Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \brief Utility routines for qs_scf
10 : ! **************************************************************************************************
11 : MODULE qs_scf_initialization
12 : USE atomic_kind_types, ONLY: atomic_kind_type
13 : USE cp_control_types, ONLY: dft_control_type
14 : USE cp_dbcsr_api, ONLY: dbcsr_create,&
15 : dbcsr_init_p,&
16 : dbcsr_p_type,&
17 : dbcsr_type,&
18 : dbcsr_type_no_symmetry
19 : USE cp_dbcsr_operations, ONLY: copy_dbcsr_to_fm,&
20 : copy_fm_to_dbcsr,&
21 : cp_dbcsr_m_by_n_from_row_template,&
22 : cp_dbcsr_sm_fm_multiply
23 : USE cp_dbcsr_output, ONLY: write_fm_with_basis_info
24 : USE cp_fm_basic_linalg, ONLY: cp_fm_column_scale,&
25 : cp_fm_row_scale,&
26 : cp_fm_transpose,&
27 : cp_fm_triangular_invert
28 : USE cp_fm_cholesky, ONLY: cp_fm_cholesky_decompose
29 : USE cp_fm_diag, ONLY: FM_DIAG_TYPE_CUSOLVER,&
30 : choose_eigv_solver,&
31 : cp_fm_power,&
32 : cusolver_generalized,&
33 : diag_type
34 : USE cp_fm_pool_types, ONLY: cp_fm_pool_p_type,&
35 : fm_pool_get_el_struct
36 : USE cp_fm_struct, ONLY: cp_fm_struct_create,&
37 : cp_fm_struct_get,&
38 : cp_fm_struct_release,&
39 : cp_fm_struct_type
40 : USE cp_fm_types, ONLY: cp_fm_create,&
41 : cp_fm_get_info,&
42 : cp_fm_release,&
43 : cp_fm_set_all,&
44 : cp_fm_to_fm,&
45 : cp_fm_to_fm_triangular,&
46 : cp_fm_type
47 : USE cp_log_handling, ONLY: cp_get_default_logger,&
48 : cp_logger_type,&
49 : cp_to_string
50 : USE cp_output_handling, ONLY: cp_p_file,&
51 : cp_print_key_finished_output,&
52 : cp_print_key_should_output,&
53 : cp_print_key_unit_nr
54 : USE hairy_probes, ONLY: AO_boundaries
55 : USE input_constants, ONLY: &
56 : broy_mix, cholesky_dbcsr, cholesky_inverse, cholesky_off, diag_block_davidson, &
57 : diag_block_krylov, diag_filter_matrix, diag_ot, diag_standard, direct_p_mix, kerker_mix, &
58 : multisec_mix, no_mix, ot2cdft, outer_scf_none, plus_u_lowdin, pulay_mix, &
59 : smeagol_runtype_emtransport, wfi_frozen_method_nr, wfi_ps_method_nr, &
60 : wfi_use_guess_method_nr
61 : USE input_section_types, ONLY: section_vals_get_subs_vals,&
62 : section_vals_type,&
63 : section_vals_val_get
64 : USE kinds, ONLY: dp
65 : USE kpoint_types, ONLY: kpoint_type
66 : USE message_passing, ONLY: mp_para_env_type
67 : USE parallel_gemm_api, ONLY: parallel_gemm
68 : USE particle_types, ONLY: particle_type
69 : USE pw_types, ONLY: pw_c1d_gs_type
70 : USE qmmm_image_charge, ONLY: conditional_calc_image_matrix
71 : USE qs_block_davidson_types, ONLY: block_davidson_allocate,&
72 : block_davidson_env_create
73 : USE qs_cdft_opt_types, ONLY: cdft_opt_type_copy
74 : USE qs_density_mixing_types, ONLY: direct_mixing_nr,&
75 : mixing_storage_create,&
76 : mixing_storage_release,&
77 : no_mixing_nr
78 : USE qs_environment_types, ONLY: get_qs_env,&
79 : qs_environment_type,&
80 : set_qs_env
81 : USE qs_fb_distribution_methods, ONLY: fb_distribution_build
82 : USE qs_fb_env_methods, ONLY: fb_env_build_atomic_halos,&
83 : fb_env_build_rcut_auto,&
84 : fb_env_read_input,&
85 : fb_env_write_info
86 : USE qs_fb_env_types, ONLY: fb_env_create,&
87 : fb_env_has_data
88 : USE qs_harris_types, ONLY: harris_type
89 : USE qs_harris_utils, ONLY: harris_density_update
90 : USE qs_initial_guess, ONLY: calculate_first_density_matrix
91 : USE qs_kind_types, ONLY: get_qs_kind,&
92 : qs_kind_type,&
93 : set_qs_kind
94 : USE qs_ks_types, ONLY: qs_ks_did_change
95 : USE qs_matrix_pools, ONLY: mpools_get
96 : USE qs_mixing_utils, ONLY: charge_mixing_init,&
97 : mixing_allocate,&
98 : mixing_init
99 : USE qs_mo_occupation, ONLY: set_mo_occupation
100 : USE qs_mo_types, ONLY: get_mo_set,&
101 : init_mo_set,&
102 : mo_set_type,&
103 : set_mo_set
104 : USE qs_outer_scf, ONLY: outer_loop_extrapolate,&
105 : outer_loop_switch,&
106 : outer_loop_variables_count
107 : USE qs_rho_atom_types, ONLY: rho_atom_type
108 : USE qs_rho_methods, ONLY: duplicate_rho_type,&
109 : qs_rho_update_rho
110 : USE qs_rho_types, ONLY: qs_rho_create,&
111 : qs_rho_get,&
112 : qs_rho_type
113 : USE qs_scf_diagonalization, ONLY: diag_kp_smat,&
114 : diag_subspace_allocate
115 : USE qs_scf_lanczos, ONLY: krylov_space_allocate
116 : USE qs_scf_output, ONLY: qs_scf_initial_info
117 : USE qs_scf_types, ONLY: &
118 : block_davidson_diag_method_nr, block_krylov_diag_method_nr, diag_subspace_env_create, &
119 : filter_matrix_diag_method_nr, general_diag_method_nr, krylov_space_create, &
120 : ot_diag_method_nr, ot_method_nr, qs_scf_env_type, scf_env_create, smeagol_method_nr, &
121 : special_diag_method_nr
122 : USE qs_wf_history_methods, ONLY: reorthogonalize_vectors,&
123 : wfi_extrapolate,&
124 : wfi_get_method_label,&
125 : wfi_update
126 : USE scf_control_types, ONLY: scf_control_type
127 : USE xas_env_types, ONLY: xas_environment_type
128 : USE xas_restart, ONLY: xas_initialize_rho
129 : #include "./base/base_uses.f90"
130 :
131 : IMPLICIT NONE
132 :
133 : PRIVATE
134 :
135 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_scf_initialization'
136 :
137 : PUBLIC:: qs_scf_env_initialize, qs_scf_env_init_basic
138 :
139 : CONTAINS
140 :
141 : ! **************************************************************************************************
142 : !> \brief initializes input parameters if needed or restores values from
143 : !> previous runs to fill scf_env with the values required for scf
144 : !> \param qs_env the qs_environment where to perform the scf procedure
145 : !> \param scf_env ...
146 : !> \param scf_control ...
147 : !> \param scf_section ...
148 : ! **************************************************************************************************
149 21205 : SUBROUTINE qs_scf_env_initialize(qs_env, scf_env, scf_control, scf_section)
150 : TYPE(qs_environment_type), POINTER :: qs_env
151 : TYPE(qs_scf_env_type), POINTER :: scf_env
152 : TYPE(scf_control_type), OPTIONAL, POINTER :: scf_control
153 : TYPE(section_vals_type), OPTIONAL, POINTER :: scf_section
154 :
155 : INTEGER :: ip, np
156 21205 : TYPE(atomic_kind_type), POINTER :: atomic_kind_set(:)
157 : TYPE(dft_control_type), POINTER :: dft_control
158 21205 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
159 21205 : TYPE(particle_type), POINTER :: particle_set(:)
160 21205 : TYPE(qs_kind_type), POINTER :: qs_kind_set(:)
161 : TYPE(scf_control_type), POINTER :: my_scf_control
162 : TYPE(section_vals_type), POINTER :: dft_section, input, my_scf_section
163 :
164 21205 : CALL get_qs_env(qs_env, input=input, dft_control=dft_control)
165 :
166 : !Initialize Hairy Probe calculation
167 21205 : IF (dft_control%hairy_probes .EQV. .TRUE.) THEN
168 : CALL get_qs_env(qs_env, mos=mos, &
169 4 : atomic_kind_set=atomic_kind_set, qs_kind_set=qs_kind_set, particle_set=particle_set)
170 4 : np = SIZE(dft_control%probe)
171 12 : DO ip = 1, np
172 : CALL AO_boundaries(probe=dft_control%probe(ip), atomic_kind_set=atomic_kind_set, qs_kind_set=qs_kind_set, &
173 12 : particle_set=particle_set, nAO=mos(1)%nao) !FIX THIS!
174 : END DO
175 : END IF
176 :
177 21205 : IF (PRESENT(scf_control)) THEN
178 82 : my_scf_control => scf_control
179 : ELSE
180 21123 : CALL get_qs_env(qs_env, scf_control=my_scf_control)
181 : END IF
182 :
183 21205 : dft_section => section_vals_get_subs_vals(input, "DFT")
184 21205 : IF (PRESENT(scf_section)) THEN
185 82 : my_scf_section => scf_section
186 : ELSE
187 21123 : my_scf_section => section_vals_get_subs_vals(dft_section, "SCF")
188 : END IF
189 :
190 21205 : CALL qs_scf_ensure_scf_env(qs_env, scf_env)
191 :
192 21205 : CALL section_vals_val_get(my_scf_section, "CHOLESKY", i_val=scf_env%cholesky_method)
193 :
194 21205 : CALL qs_scf_ensure_mos(qs_env)
195 :
196 : ! set flags for diagonalization
197 : CALL qs_scf_ensure_diagonalization(scf_env, my_scf_section, qs_env, &
198 21205 : my_scf_control, qs_env%has_unit_metric)
199 : ! set parameters for mixing/DIIS during scf
200 21205 : CALL qs_scf_ensure_mixing(my_scf_control, my_scf_section, scf_env, dft_control)
201 :
202 21205 : CALL qs_scf_ensure_work_matrices(qs_env, scf_env)
203 :
204 21205 : CALL qs_scf_ensure_mixing_store(qs_env, scf_env)
205 :
206 : ! Initialize outer loop variables: handle CDFT and regular outer loop separately
207 21205 : IF (dft_control%qs_control%cdft) THEN
208 : CALL qs_scf_ensure_cdft_loop_vars(qs_env, scf_env, dft_control, &
209 326 : scf_control=my_scf_control)
210 : ELSE
211 20879 : CALL qs_scf_ensure_outer_loop_vars(scf_env, my_scf_control)
212 : END IF
213 :
214 21205 : CALL init_scf_run(scf_env, qs_env, my_scf_section, my_scf_control)
215 :
216 21205 : END SUBROUTINE qs_scf_env_initialize
217 :
218 : ! **************************************************************************************************
219 : !> \brief initializes input parameters if needed for non-scf calclulations using diagonalization
220 : !> \param qs_env the qs_environment where to perform the scf procedure
221 : !> \param scf_env ...
222 : ! **************************************************************************************************
223 2 : SUBROUTINE qs_scf_env_init_basic(qs_env, scf_env)
224 : TYPE(qs_environment_type), POINTER :: qs_env
225 : TYPE(qs_scf_env_type), POINTER :: scf_env
226 :
227 : TYPE(dft_control_type), POINTER :: dft_control
228 : TYPE(scf_control_type), POINTER :: scf_control
229 : TYPE(section_vals_type), POINTER :: dft_section, input, scf_section
230 :
231 2 : CALL get_qs_env(qs_env, input=input, dft_control=dft_control)
232 :
233 2 : CALL get_qs_env(qs_env, scf_control=scf_control)
234 2 : dft_section => section_vals_get_subs_vals(input, "DFT")
235 2 : scf_section => section_vals_get_subs_vals(dft_section, "SCF")
236 :
237 2 : CALL qs_scf_ensure_scf_env(qs_env, scf_env)
238 :
239 2 : CALL section_vals_val_get(scf_section, "CHOLESKY", i_val=scf_env%cholesky_method)
240 2 : scf_control%use_diag = .TRUE.
241 2 : scf_control%diagonalization%method = diag_standard
242 :
243 2 : CALL qs_scf_ensure_mos(qs_env)
244 :
245 : ! set flags for diagonalization
246 : CALL qs_scf_ensure_diagonalization(scf_env, scf_section, qs_env, &
247 2 : scf_control, qs_env%has_unit_metric)
248 2 : CALL qs_scf_ensure_work_matrices(qs_env, scf_env)
249 :
250 2 : CALL init_scf_run(scf_env, qs_env, scf_section, scf_control)
251 :
252 2 : END SUBROUTINE qs_scf_env_init_basic
253 :
254 : ! **************************************************************************************************
255 : !> \brief makes sure scf_env is allocated (might already be from before)
256 : !> in case it is present the g-space mixing storage is reset
257 : !> \param qs_env ...
258 : !> \param scf_env ...
259 : ! **************************************************************************************************
260 21207 : SUBROUTINE qs_scf_ensure_scf_env(qs_env, scf_env)
261 : TYPE(qs_environment_type), POINTER :: qs_env
262 : TYPE(qs_scf_env_type), POINTER :: scf_env
263 :
264 21207 : TYPE(pw_c1d_gs_type), DIMENSION(:), POINTER :: rho_g
265 : TYPE(qs_rho_type), POINTER :: rho
266 :
267 21207 : NULLIFY (rho_g)
268 :
269 27776 : IF (.NOT. ASSOCIATED(scf_env)) THEN ! i.e. for MD this is associated on the second step (it so seems)
270 6569 : ALLOCATE (scf_env)
271 6569 : CALL scf_env_create(scf_env)
272 : ELSE
273 : ! Reallocate mixing store, if the g space grid (cell) has changed
274 14706 : SELECT CASE (scf_env%mixing_method)
275 : CASE (kerker_mix, pulay_mix, broy_mix, multisec_mix)
276 14638 : IF (ASSOCIATED(scf_env%mixing_store)) THEN
277 : ! The current mixing_store data structure does not allow for an unique
278 : ! grid comparison, but the probability that the 1d lengths of the old and
279 : ! the new grid are accidentily equal is rather low
280 68 : CALL get_qs_env(qs_env, rho=rho)
281 68 : CALL qs_rho_get(rho, rho_g=rho_g)
282 68 : IF (ASSOCIATED(scf_env%mixing_store%rhoin)) THEN
283 32 : IF (SIZE(rho_g(1)%pw_grid%gsq) /= SIZE(scf_env%mixing_store%rhoin(1)%cc)) THEN
284 0 : CALL mixing_storage_release(scf_env%mixing_store)
285 0 : DEALLOCATE (scf_env%mixing_store)
286 : END IF
287 : END IF
288 : END IF
289 : END SELECT
290 : END IF
291 :
292 21207 : END SUBROUTINE qs_scf_ensure_scf_env
293 :
294 : ! **************************************************************************************************
295 : !> \brief performs allocation of outer SCF variables
296 : !> \param scf_env the SCF environment which contains the outer SCF variables
297 : !> \param scf_control control settings for the outer SCF loop
298 : !> \param nvar (optional) set number of outer SCF variables externally if CDFT SCF is active
299 : ! **************************************************************************************************
300 21205 : SUBROUTINE qs_scf_ensure_outer_loop_vars(scf_env, scf_control, nvar)
301 : TYPE(qs_scf_env_type), POINTER :: scf_env
302 : TYPE(scf_control_type), POINTER :: scf_control
303 : INTEGER, OPTIONAL :: nvar
304 :
305 : INTEGER :: nhistory, nvariables
306 :
307 21205 : IF (scf_control%outer_scf%have_scf) THEN
308 4139 : nhistory = scf_control%outer_scf%max_scf + 1
309 4139 : IF (PRESENT(nvar)) THEN
310 326 : IF (nvar > 0) THEN
311 : nvariables = nvar
312 : ELSE
313 0 : nvariables = outer_loop_variables_count(scf_control)
314 : END IF
315 : ELSE
316 3813 : nvariables = outer_loop_variables_count(scf_control)
317 : END IF
318 16556 : ALLOCATE (scf_env%outer_scf%variables(nvariables, nhistory))
319 12417 : ALLOCATE (scf_env%outer_scf%count(nhistory))
320 77593 : scf_env%outer_scf%count = 0
321 12417 : ALLOCATE (scf_env%outer_scf%gradient(nvariables, nhistory))
322 12417 : ALLOCATE (scf_env%outer_scf%energy(nhistory))
323 : END IF
324 :
325 21205 : END SUBROUTINE qs_scf_ensure_outer_loop_vars
326 :
327 : ! **************************************************************************************************
328 : !> \brief performs allocation of CDFT SCF variables
329 : !> \param qs_env the qs_env where to perform the allocation
330 : !> \param scf_env the currently active scf_env
331 : !> \param dft_control the dft_control that holds the cdft_control type
332 : !> \param scf_control the currently active scf_control
333 : ! **************************************************************************************************
334 326 : SUBROUTINE qs_scf_ensure_cdft_loop_vars(qs_env, scf_env, dft_control, scf_control)
335 : TYPE(qs_environment_type), POINTER :: qs_env
336 : TYPE(qs_scf_env_type), POINTER :: scf_env
337 : TYPE(dft_control_type), POINTER :: dft_control
338 : TYPE(scf_control_type), POINTER :: scf_control
339 :
340 : INTEGER :: nhistory, nvariables
341 : LOGICAL :: do_kpoints
342 326 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: gradient_history, outer_scf_history, &
343 326 : variable_history
344 :
345 326 : NULLIFY (outer_scf_history, gradient_history, variable_history)
346 326 : CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints)
347 : ! Test kpoints
348 326 : IF (do_kpoints) &
349 0 : CPABORT("CDFT calculation not possible with kpoints")
350 : ! Check that OUTER_SCF section in DFT&SCF is active
351 : ! This section must always be active to facilitate
352 : ! switching of the CDFT and SCF control parameters in outer_loop_switch
353 326 : IF (.NOT. scf_control%outer_scf%have_scf) &
354 0 : CPABORT("Section SCF&OUTER_SCF must be active for CDFT calculations.")
355 : ! Initialize CDFT and outer_loop variables (constraint settings active in scf_control)
356 326 : IF (dft_control%qs_control%cdft_control%constraint_control%have_scf) THEN
357 326 : nhistory = dft_control%qs_control%cdft_control%constraint_control%max_scf + 1
358 326 : IF (scf_control%outer_scf%type /= outer_scf_none) THEN
359 : nvariables = outer_loop_variables_count(scf_control, &
360 62 : dft_control%qs_control%cdft_control)
361 : ELSE
362 : ! First iteration: scf_control has not yet been updated
363 264 : nvariables = SIZE(dft_control%qs_control%cdft_control%target)
364 : END IF
365 1304 : ALLOCATE (dft_control%qs_control%cdft_control%constraint%variables(nvariables, nhistory))
366 978 : ALLOCATE (dft_control%qs_control%cdft_control%constraint%count(nhistory))
367 2246 : dft_control%qs_control%cdft_control%constraint%count = 0
368 978 : ALLOCATE (dft_control%qs_control%cdft_control%constraint%gradient(nvariables, nhistory))
369 978 : ALLOCATE (dft_control%qs_control%cdft_control%constraint%energy(nhistory))
370 326 : CALL qs_scf_ensure_outer_loop_vars(scf_env, scf_control, nvariables)
371 : END IF
372 : ! Executed only on first call (OT settings active in scf_control)
373 : ! Save OT settings and constraint initial values in CDFT control
374 : ! Then switch to constraint outer_scf settings for proper initialization of history
375 326 : IF (scf_control%outer_scf%have_scf) THEN
376 326 : IF (scf_control%outer_scf%type == outer_scf_none) THEN
377 264 : dft_control%qs_control%cdft_control%ot_control%have_scf = .TRUE.
378 264 : dft_control%qs_control%cdft_control%ot_control%max_scf = scf_control%outer_scf%max_scf
379 264 : dft_control%qs_control%cdft_control%ot_control%eps_scf = scf_control%outer_scf%eps_scf
380 264 : dft_control%qs_control%cdft_control%ot_control%step_size = scf_control%outer_scf%step_size
381 264 : dft_control%qs_control%cdft_control%ot_control%type = scf_control%outer_scf%type
382 264 : dft_control%qs_control%cdft_control%ot_control%optimizer = scf_control%outer_scf%optimizer
383 264 : dft_control%qs_control%cdft_control%ot_control%diis_buffer_length = scf_control%outer_scf%diis_buffer_length
384 264 : dft_control%qs_control%cdft_control%ot_control%bisect_trust_count = scf_control%outer_scf%bisect_trust_count
385 : CALL cdft_opt_type_copy(dft_control%qs_control%cdft_control%ot_control%cdft_opt_control, &
386 264 : scf_control%outer_scf%cdft_opt_control)
387 : ! In case constraint and OT extrapolation orders are different, make sure to use former
388 264 : nvariables = SIZE(dft_control%qs_control%cdft_control%target)
389 : IF (scf_control%outer_scf%extrapolation_order /= &
390 : dft_control%qs_control%cdft_control%constraint_control%extrapolation_order &
391 264 : .OR. nvariables /= 1) THEN
392 256 : DEALLOCATE (qs_env%outer_scf_history)
393 256 : DEALLOCATE (qs_env%gradient_history)
394 256 : DEALLOCATE (qs_env%variable_history)
395 256 : nhistory = dft_control%qs_control%cdft_control%constraint_control%extrapolation_order
396 1024 : ALLOCATE (outer_scf_history(nvariables, nhistory))
397 768 : ALLOCATE (gradient_history(nvariables, 2))
398 1324 : gradient_history = 0.0_dp
399 512 : ALLOCATE (variable_history(nvariables, 2))
400 1324 : variable_history = 0.0_dp
401 : CALL set_qs_env(qs_env, outer_scf_history=outer_scf_history, &
402 256 : gradient_history=gradient_history, variable_history=variable_history)
403 : END IF
404 264 : CALL outer_loop_switch(scf_env, scf_control, dft_control%qs_control%cdft_control, ot2cdft)
405 : END IF
406 : END IF
407 :
408 326 : END SUBROUTINE qs_scf_ensure_cdft_loop_vars
409 :
410 : ! **************************************************************************************************
411 : !> \brief performs allocation of the mixing storage
412 : !> \param qs_env ...
413 : !> \param scf_env ...
414 : ! **************************************************************************************************
415 21205 : SUBROUTINE qs_scf_ensure_mixing_store(qs_env, scf_env)
416 : TYPE(qs_environment_type), POINTER :: qs_env
417 : TYPE(qs_scf_env_type), POINTER :: scf_env
418 :
419 : TYPE(dft_control_type), POINTER :: dft_control
420 :
421 21205 : NULLIFY (dft_control)
422 21205 : CALL get_qs_env(qs_env=qs_env, dft_control=dft_control)
423 :
424 21205 : IF (scf_env%mixing_method > 0) THEN
425 : CALL mixing_allocate(qs_env, scf_env%mixing_method, scf_env%p_mix_new, &
426 : scf_env%p_delta, dft_control%nspins, &
427 15088 : scf_env%mixing_store)
428 : ELSE
429 6117 : NULLIFY (scf_env%p_mix_new)
430 : END IF
431 :
432 21205 : END SUBROUTINE qs_scf_ensure_mixing_store
433 :
434 : ! **************************************************************************************************
435 : !> \brief Performs allocation of the SCF work matrices
436 : !> In case of kpoints we probably don't need most of these matrices,
437 : !> maybe we have to initialize some matrices in the fm_pool in kpoints
438 : !> \param qs_env ...
439 : !> \param scf_env ...
440 : ! **************************************************************************************************
441 63621 : SUBROUTINE qs_scf_ensure_work_matrices(qs_env, scf_env)
442 :
443 : TYPE(qs_environment_type), POINTER :: qs_env
444 : TYPE(qs_scf_env_type), POINTER :: scf_env
445 :
446 : CHARACTER(LEN=*), PARAMETER :: routineN = 'qs_scf_ensure_work_matrices'
447 :
448 : INTEGER :: handle, is, nao, nrow_block, nw
449 : LOGICAL :: do_kpoints
450 21207 : TYPE(cp_fm_pool_p_type), DIMENSION(:), POINTER :: ao_mo_fm_pools
451 : TYPE(cp_fm_struct_type), POINTER :: ao_ao_fmstruct, ao_mo_fmstruct
452 21207 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s
453 : TYPE(dbcsr_type), POINTER :: ref_matrix
454 : TYPE(dft_control_type), POINTER :: dft_control
455 21207 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
456 : TYPE(scf_control_type), POINTER :: scf_control
457 :
458 21207 : CALL timeset(routineN, handle)
459 :
460 21207 : NULLIFY (ao_mo_fm_pools, ao_mo_fmstruct, ao_ao_fmstruct, dft_control, matrix_s, mos)
461 :
462 : CALL get_qs_env(qs_env=qs_env, &
463 : dft_control=dft_control, &
464 : matrix_s_kp=matrix_s, &
465 : mos=mos, &
466 : scf_control=scf_control, &
467 21207 : do_kpoints=do_kpoints)
468 21207 : CALL mpools_get(qs_env%mpools, ao_mo_fm_pools=ao_mo_fm_pools)
469 :
470 : ! create an ao_ao parallel matrix structure
471 21207 : ao_mo_fmstruct => fm_pool_get_el_struct(ao_mo_fm_pools(1)%pool)
472 21207 : CALL cp_fm_struct_get(ao_mo_fmstruct, nrow_block=nrow_block)
473 21207 : CALL get_mo_set(mos(1), nao=nao)
474 : CALL cp_fm_struct_create(fmstruct=ao_ao_fmstruct, &
475 : nrow_block=nrow_block, &
476 : ncol_block=nrow_block, &
477 : nrow_global=nao, &
478 : ncol_global=nao, &
479 21207 : template_fmstruct=ao_mo_fmstruct)
480 :
481 21207 : IF ((scf_env%method /= ot_method_nr) .AND. &
482 : (scf_env%method /= block_davidson_diag_method_nr)) THEN
483 15074 : IF (.NOT. ASSOCIATED(scf_env%scf_work1)) THEN
484 13626 : nw = dft_control%nspins
485 13626 : IF (do_kpoints) nw = 4
486 57966 : ALLOCATE (scf_env%scf_work1(nw))
487 30714 : DO is = 1, SIZE(scf_env%scf_work1)
488 : CALL cp_fm_create(scf_env%scf_work1(is), &
489 : matrix_struct=ao_ao_fmstruct, &
490 30714 : name="SCF-WORK_MATRIX-1-"//TRIM(ADJUSTL(cp_to_string(is))))
491 : END DO
492 : END IF
493 : IF ((.NOT. ASSOCIATED(scf_env%ortho)) .AND. &
494 15074 : (scf_env%method /= ot_diag_method_nr) .AND. &
495 : (scf_env%method /= special_diag_method_nr)) THEN
496 : ! Initialize fm matrix to store the Cholesky decomposition
497 10962 : ALLOCATE (scf_env%ortho)
498 : CALL cp_fm_create(scf_env%ortho, &
499 : matrix_struct=ao_ao_fmstruct, &
500 10962 : name="SCF-ORTHO_MATRIX")
501 : ! Initialize dbcsr matrix to store the Cholesky decomposition
502 10962 : IF (scf_env%cholesky_method == cholesky_dbcsr) THEN
503 58 : ref_matrix => matrix_s(1, 1)%matrix
504 58 : CALL dbcsr_init_p(scf_env%ortho_dbcsr)
505 : CALL dbcsr_create(scf_env%ortho_dbcsr, template=ref_matrix, &
506 58 : matrix_type=dbcsr_type_no_symmetry)
507 58 : CALL dbcsr_init_p(scf_env%buf1_dbcsr)
508 : CALL dbcsr_create(scf_env%buf1_dbcsr, template=ref_matrix, &
509 58 : matrix_type=dbcsr_type_no_symmetry)
510 58 : CALL dbcsr_init_p(scf_env%buf2_dbcsr)
511 : CALL dbcsr_create(scf_env%buf2_dbcsr, template=ref_matrix, &
512 58 : matrix_type=dbcsr_type_no_symmetry)
513 10904 : ELSE IF (scf_env%cholesky_method == cholesky_inverse .OR. &
514 : (scf_control%level_shift /= 0.0_dp .AND. &
515 : scf_env%cholesky_method == cholesky_off)) THEN
516 48 : ALLOCATE (scf_env%ortho_m1)
517 : CALL cp_fm_create(scf_env%ortho_m1, &
518 : matrix_struct=ao_ao_fmstruct, &
519 48 : name="SCF-ORTHO_MATRIX-1")
520 : END IF
521 : END IF
522 15074 : IF (.NOT. ASSOCIATED(scf_env%scf_work2)) THEN
523 13626 : ALLOCATE (scf_env%scf_work2)
524 : CALL cp_fm_create(scf_env%scf_work2, &
525 : matrix_struct=ao_ao_fmstruct, &
526 13626 : name="SCF-WORK_MATRIX-2")
527 : END IF
528 : END IF
529 :
530 21207 : IF (dft_control%dft_plus_u) THEN
531 92 : IF (dft_control%plus_u_method_id == plus_u_lowdin) THEN
532 14 : IF (.NOT. ASSOCIATED(scf_env%s_half)) THEN
533 10 : ALLOCATE (scf_env%s_half)
534 : CALL cp_fm_create(scf_env%s_half, &
535 : matrix_struct=ao_ao_fmstruct, &
536 10 : name="S**(1/2) MATRIX")
537 : END IF
538 : END IF
539 : END IF
540 :
541 21207 : IF (do_kpoints) THEN
542 940 : IF (.NOT. ASSOCIATED(scf_env%scf_work1)) THEN
543 0 : nw = 4
544 0 : ALLOCATE (scf_env%scf_work1(nw))
545 0 : DO is = 1, SIZE(scf_env%scf_work1)
546 : CALL cp_fm_create(scf_env%scf_work1(is), &
547 : matrix_struct=ao_ao_fmstruct, &
548 0 : name="SCF-WORK_MATRIX-1-"//TRIM(ADJUSTL(cp_to_string(is))))
549 : END DO
550 : END IF
551 : END IF
552 :
553 21207 : CALL cp_fm_struct_release(ao_ao_fmstruct)
554 :
555 21207 : CALL timestop(handle)
556 :
557 21207 : END SUBROUTINE qs_scf_ensure_work_matrices
558 :
559 : ! **************************************************************************************************
560 : !> \brief performs allocation of the MO matrices
561 : !> \param qs_env ...
562 : ! **************************************************************************************************
563 21207 : SUBROUTINE qs_scf_ensure_mos(qs_env)
564 : TYPE(qs_environment_type), POINTER :: qs_env
565 :
566 : CHARACTER(len=*), PARAMETER :: routineN = 'qs_scf_ensure_mos'
567 :
568 : INTEGER :: handle, ic, ik, ikk, ispin, nmo, nmo_mat
569 21207 : TYPE(cp_fm_pool_p_type), DIMENSION(:), POINTER :: ao_mo_fm_pools
570 : TYPE(cp_fm_type), POINTER :: mo_coeff, mo_coeff_last
571 21207 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: mo_derivs
572 21207 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s
573 : TYPE(dbcsr_type), POINTER :: mo_coeff_b
574 : TYPE(dft_control_type), POINTER :: dft_control
575 : TYPE(kpoint_type), POINTER :: kpoints
576 21207 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos, mos_last_converged
577 21207 : TYPE(mo_set_type), DIMENSION(:, :), POINTER :: mos_k
578 : TYPE(xas_environment_type), POINTER :: xas_env
579 :
580 21207 : CALL timeset(routineN, handle)
581 :
582 21207 : NULLIFY (ao_mo_fm_pools, dft_control, mos, xas_env, matrix_s, mos_last_converged, mo_coeff_last)
583 :
584 : CALL get_qs_env(qs_env=qs_env, &
585 : dft_control=dft_control, &
586 : mos=mos, &
587 : matrix_s_kp=matrix_s, &
588 21207 : xas_env=xas_env)
589 21207 : CALL mpools_get(qs_env%mpools, ao_mo_fm_pools=ao_mo_fm_pools)
590 21207 : IF (dft_control%switch_surf_dip) THEN
591 2 : CALL get_qs_env(qs_env, mos_last_converged=mos_last_converged)
592 : END IF
593 :
594 21207 : nmo_mat = dft_control%nspins
595 21207 : IF (dft_control%restricted) nmo_mat = 1 ! right now, there might be more mos than needed derivs
596 :
597 : ! Finish initialization of the MOs
598 21207 : CPASSERT(ASSOCIATED(mos))
599 44982 : DO ispin = 1, SIZE(mos)
600 23775 : CALL get_mo_set(mos(ispin), mo_coeff=mo_coeff, mo_coeff_b=mo_coeff_b)
601 23775 : IF (.NOT. ASSOCIATED(mo_coeff)) THEN
602 : CALL init_mo_set(mos(ispin), &
603 : fm_pool=ao_mo_fm_pools(ispin)%pool, &
604 8060 : name="qs_env%mo"//TRIM(ADJUSTL(cp_to_string(ispin))))
605 : END IF
606 44982 : IF (.NOT. ASSOCIATED(mo_coeff_b)) THEN
607 8060 : CALL cp_fm_get_info(mos(ispin)%mo_coeff, ncol_global=nmo)
608 8060 : CALL dbcsr_init_p(mos(ispin)%mo_coeff_b)
609 : CALL cp_dbcsr_m_by_n_from_row_template(mos(ispin)%mo_coeff_b, template=matrix_s(1, 1)%matrix, n=nmo, &
610 8060 : sym=dbcsr_type_no_symmetry)
611 : END IF
612 : END DO
613 : ! Get the mo_derivs OK if needed
614 21207 : IF (qs_env%requires_mo_derivs) THEN
615 6123 : CALL get_qs_env(qs_env, mo_derivs=mo_derivs)
616 6123 : IF (.NOT. ASSOCIATED(mo_derivs)) THEN
617 9067 : ALLOCATE (mo_derivs(nmo_mat))
618 4833 : DO ispin = 1, nmo_mat
619 2716 : CALL get_mo_set(mos(ispin), mo_coeff_b=mo_coeff_b)
620 2716 : NULLIFY (mo_derivs(ispin)%matrix)
621 2716 : CALL dbcsr_init_p(mo_derivs(ispin)%matrix)
622 : CALL dbcsr_create(mo_derivs(ispin)%matrix, template=mo_coeff_b, &
623 4833 : name="mo_derivs", matrix_type=dbcsr_type_no_symmetry)
624 : END DO
625 2117 : CALL set_qs_env(qs_env, mo_derivs=mo_derivs)
626 : END IF
627 :
628 : ELSE
629 : ! nothing should be done
630 : END IF
631 :
632 : ! Finish initialization of the MOs for ADMM and derivs if needed ***
633 21207 : IF (dft_control%do_admm) THEN
634 910 : IF (dft_control%restricted) CPABORT("ROKS with ADMM is not implemented")
635 : END IF
636 :
637 : ! Finish initialization of mos_last_converged [SGh]
638 21207 : IF (dft_control%switch_surf_dip) THEN
639 2 : CPASSERT(ASSOCIATED(mos_last_converged))
640 4 : DO ispin = 1, SIZE(mos_last_converged)
641 2 : CALL get_mo_set(mos_last_converged(ispin), mo_coeff=mo_coeff_last)
642 4 : IF (.NOT. ASSOCIATED(mo_coeff_last)) THEN
643 : CALL init_mo_set(mos_last_converged(ispin), &
644 : fm_ref=mos(ispin)%mo_coeff, &
645 2 : name="qs_env%mos_last_converged"//TRIM(ADJUSTL(cp_to_string(ispin))))
646 : END IF
647 : END DO
648 : END IF
649 : ! kpoints: we have to initialize all the k-point MOs
650 21207 : CALL get_qs_env(qs_env=qs_env, kpoints=kpoints)
651 21207 : IF (kpoints%nkp /= 0) THEN
652 : ! check for some incompatible options
653 940 : IF (qs_env%requires_mo_derivs) THEN
654 2 : CPWARN("MO derivative methods flag has been switched off for kpoint calculation")
655 : ! we switch it off to make band structure calculations
656 : ! possible for OT gamma point calculations
657 2 : qs_env%requires_mo_derivs = .FALSE.
658 : END IF
659 940 : IF (dft_control%do_xas_calculation) &
660 0 : CPABORT("No XAS implemented with kpoints")
661 940 : IF (qs_env%do_rixs) &
662 0 : CPABORT("RIXS not implemented with kpoints")
663 3852 : DO ik = 1, SIZE(kpoints%kp_env)
664 2912 : CALL mpools_get(kpoints%mpools, ao_mo_fm_pools=ao_mo_fm_pools)
665 2912 : mos_k => kpoints%kp_env(ik)%kpoint_env%mos
666 2912 : ikk = kpoints%kp_range(1) + ik - 1
667 2912 : CPASSERT(ASSOCIATED(mos_k))
668 7246 : DO ispin = 1, SIZE(mos_k, 2)
669 13080 : DO ic = 1, SIZE(mos_k, 1)
670 6774 : CALL get_mo_set(mos_k(ic, ispin), mo_coeff=mo_coeff, mo_coeff_b=mo_coeff_b)
671 6774 : IF (.NOT. ASSOCIATED(mo_coeff)) THEN
672 : CALL init_mo_set(mos_k(ic, ispin), &
673 : fm_pool=ao_mo_fm_pools(ispin)%pool, &
674 : name="kpoints_"//TRIM(ADJUSTL(cp_to_string(ikk)))// &
675 2754 : "%mo"//TRIM(ADJUSTL(cp_to_string(ispin))))
676 : END IF
677 : ! no sparse matrix representation of kpoint MO vectors
678 10168 : CPASSERT(.NOT. ASSOCIATED(mo_coeff_b))
679 : END DO
680 : END DO
681 : END DO
682 : END IF
683 :
684 21207 : CALL timestop(handle)
685 :
686 21207 : END SUBROUTINE qs_scf_ensure_mos
687 :
688 : ! **************************************************************************************************
689 : !> \brief sets flag for mixing/DIIS during scf
690 : !> \param scf_control ...
691 : !> \param scf_section ...
692 : !> \param scf_env ...
693 : !> \param dft_control ...
694 : ! **************************************************************************************************
695 21205 : SUBROUTINE qs_scf_ensure_mixing(scf_control, scf_section, scf_env, dft_control)
696 : TYPE(scf_control_type), POINTER :: scf_control
697 : TYPE(section_vals_type), POINTER :: scf_section
698 : TYPE(qs_scf_env_type), POINTER :: scf_env
699 : TYPE(dft_control_type), POINTER :: dft_control
700 :
701 : TYPE(section_vals_type), POINTER :: mixing_section
702 :
703 21205 : SELECT CASE (scf_control%mixing_method)
704 : CASE (no_mix)
705 0 : scf_env%mixing_method = no_mixing_nr
706 0 : scf_env%p_mix_alpha = 1.0_dp
707 : CASE (direct_p_mix, kerker_mix, pulay_mix, broy_mix, multisec_mix)
708 21205 : scf_env%mixing_method = scf_control%mixing_method
709 21205 : mixing_section => section_vals_get_subs_vals(scf_section, "MIXING")
710 21205 : IF (.NOT. ASSOCIATED(scf_env%mixing_store)) THEN
711 26268 : ALLOCATE (scf_env%mixing_store)
712 : CALL mixing_storage_create(scf_env%mixing_store, mixing_section, scf_env%mixing_method, &
713 6567 : dft_control%qs_control%cutoff)
714 : END IF
715 : CASE DEFAULT
716 21205 : CPABORT("Unknown mixing method")
717 : END SELECT
718 :
719 : ! Disable DIIS for OT and g-space density mixing methods
720 21205 : IF (scf_env%method == ot_method_nr) THEN
721 : ! No mixing is used with OT
722 6117 : scf_env%mixing_method = no_mixing_nr
723 6117 : scf_env%p_mix_alpha = 1.0_dp
724 6117 : scf_env%skip_diis = .TRUE.
725 : END IF
726 :
727 21205 : IF (scf_control%use_diag .AND. scf_env%mixing_method == no_mixing_nr) THEN
728 0 : CPABORT("Diagonalization procedures without mixing are not recommendable")
729 : END IF
730 :
731 21205 : IF (scf_env%mixing_method > direct_mixing_nr) THEN
732 280 : scf_env%skip_diis = .TRUE.
733 280 : scf_env%p_mix_alpha = scf_env%mixing_store%alpha
734 280 : IF (scf_env%mixing_store%beta == 0.0_dp) THEN
735 0 : CPABORT("Mixing employing the Kerker damping factor needs BETA /= 0.0")
736 : END IF
737 : END IF
738 :
739 21205 : IF (scf_env%mixing_method == direct_mixing_nr) THEN
740 14808 : scf_env%p_mix_alpha = scf_env%mixing_store%alpha
741 14808 : IF (scf_control%eps_diis < scf_control%eps_scf) THEN
742 42 : scf_env%skip_diis = .TRUE.
743 42 : CPWARN("the DIIS scheme is disabled, since EPS_DIIS < EPS_SCF")
744 : END IF
745 : END IF
746 :
747 21205 : END SUBROUTINE qs_scf_ensure_mixing
748 :
749 : ! **************************************************************************************************
750 : !> \brief sets flags for diagonalization and ensure that everything is
751 : !> allocated
752 : !> \param scf_env ...
753 : !> \param scf_section ...
754 : !> \param qs_env ...
755 : !> \param scf_control ...
756 : !> \param has_unit_metric ...
757 : ! **************************************************************************************************
758 21207 : SUBROUTINE qs_scf_ensure_diagonalization(scf_env, scf_section, qs_env, &
759 : scf_control, has_unit_metric)
760 : TYPE(qs_scf_env_type), POINTER :: scf_env
761 : TYPE(section_vals_type), POINTER :: scf_section
762 : TYPE(qs_environment_type), POINTER :: qs_env
763 : TYPE(scf_control_type), POINTER :: scf_control
764 : LOGICAL :: has_unit_metric
765 :
766 : INTEGER :: ispin, nao, nmo
767 : LOGICAL :: do_kpoints, need_coeff_b, not_se_or_tb
768 : TYPE(cp_fm_type), POINTER :: mo_coeff
769 : TYPE(dft_control_type), POINTER :: dft_control
770 21207 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
771 :
772 21207 : CALL get_qs_env(qs_env=qs_env, do_kpoints=do_kpoints, dft_control=dft_control, mos=mos)
773 : not_se_or_tb = .NOT. (dft_control%qs_control%dftb .OR. dft_control%qs_control%xtb .OR. &
774 21207 : dft_control%qs_control%semi_empirical)
775 21207 : need_coeff_b = .FALSE.
776 21207 : scf_env%needs_ortho = .FALSE.
777 :
778 21207 : IF (dft_control%smeagol_control%smeagol_enabled .AND. &
779 : dft_control%smeagol_control%run_type == smeagol_runtype_emtransport) THEN
780 0 : scf_env%method = smeagol_method_nr
781 0 : scf_env%skip_diis = .TRUE.
782 0 : scf_control%use_diag = .FALSE.
783 :
784 0 : IF (.NOT. do_kpoints) THEN
785 0 : CPABORT("SMEAGOL requires kpoint calculations")
786 : END IF
787 0 : CPWARN_IF(scf_control%use_ot, "OT is irrelevant to NEGF method")
788 : END IF
789 :
790 21207 : IF (scf_control%use_diag) THEN
791 : ! sanity check whether combinations are allowed
792 15090 : IF (dft_control%restricted) &
793 0 : CPABORT("OT only for restricted (ROKS)")
794 15122 : SELECT CASE (scf_control%diagonalization%method)
795 : CASE (diag_ot, diag_block_krylov, diag_block_davidson)
796 32 : IF (.NOT. not_se_or_tb) &
797 15090 : CPABORT("TB and SE not possible with OT diagonalization")
798 : END SELECT
799 30138 : SELECT CASE (scf_control%diagonalization%method)
800 : ! Diagonalization: additional check whether we are in an orthonormal basis
801 : CASE (diag_standard)
802 15048 : scf_env%method = general_diag_method_nr
803 15048 : scf_env%needs_ortho = (.NOT. has_unit_metric) .AND. (.NOT. do_kpoints)
804 : IF (diag_type == FM_DIAG_TYPE_CUSOLVER .AND. &
805 : cusolver_generalized .AND. &
806 15048 : scf_control%level_shift == 0.0_dp .AND. &
807 : scf_env%cholesky_method /= cholesky_off) THEN
808 0 : scf_env%needs_ortho = .FALSE.
809 : END IF
810 15048 : IF (has_unit_metric) THEN
811 2656 : scf_env%method = special_diag_method_nr
812 : END IF
813 : ! OT Diagonalization: not possible with ROKS
814 : CASE (diag_ot)
815 8 : IF (dft_control%roks) &
816 0 : CPABORT("ROKS with OT diagonalization not possible")
817 8 : IF (do_kpoints) &
818 0 : CPABORT("OT diagonalization not possible with kpoint calculations")
819 8 : scf_env%method = ot_diag_method_nr
820 8 : need_coeff_b = .TRUE.
821 : ! Block Krylov diagonlization: not possible with ROKS,
822 : ! allocation of additional matrices is needed
823 : CASE (diag_block_krylov)
824 8 : IF (dft_control%roks) &
825 0 : CPABORT("ROKS with block PF diagonalization not possible")
826 8 : IF (do_kpoints) &
827 0 : CPABORT("Block Krylov diagonalization not possible with kpoint calculations")
828 8 : scf_env%method = block_krylov_diag_method_nr
829 8 : scf_env%needs_ortho = .TRUE.
830 8 : IF (.NOT. ASSOCIATED(scf_env%krylov_space)) &
831 4 : CALL krylov_space_create(scf_env%krylov_space, scf_section)
832 8 : CALL krylov_space_allocate(scf_env%krylov_space, scf_control, mos)
833 : ! Block davidson diagonlization: allocation of additional matrices is needed
834 : CASE (diag_block_davidson)
835 16 : IF (do_kpoints) &
836 0 : CPABORT("Block Davidson diagonalization not possible with kpoint calculations")
837 16 : scf_env%method = block_davidson_diag_method_nr
838 16 : IF (.NOT. ASSOCIATED(scf_env%block_davidson_env)) &
839 : CALL block_davidson_env_create(scf_env%block_davidson_env, dft_control%nspins, &
840 12 : scf_section)
841 34 : DO ispin = 1, dft_control%nspins
842 18 : CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff, nao=nao, nmo=nmo)
843 34 : CALL block_davidson_allocate(scf_env%block_davidson_env(ispin), mo_coeff, nao, nmo)
844 : END DO
845 10 : need_coeff_b = .TRUE.
846 : ! Filter matrix diagonalisation method
847 : CASE (diag_filter_matrix)
848 10 : scf_env%method = filter_matrix_diag_method_nr
849 10 : IF (.NOT. fb_env_has_data(scf_env%filter_matrix_env)) THEN
850 10 : CALL fb_env_create(scf_env%filter_matrix_env)
851 : END IF
852 10 : CALL fb_env_read_input(scf_env%filter_matrix_env, scf_section)
853 10 : CALL fb_env_build_rcut_auto(scf_env%filter_matrix_env, qs_env)
854 10 : CALL fb_env_write_info(scf_env%filter_matrix_env, qs_env, scf_section)
855 10 : CALL fb_distribution_build(scf_env%filter_matrix_env, qs_env, scf_section)
856 10 : CALL fb_env_build_atomic_halos(scf_env%filter_matrix_env, qs_env, scf_section)
857 : CASE DEFAULT
858 15090 : CPABORT("Unknown diagonalization method")
859 : END SELECT
860 : ! Check if subspace diagonlization is requested: allocation of additional matrices is needed
861 15090 : IF (scf_control%do_diag_sub) THEN
862 2 : scf_env%needs_ortho = .TRUE.
863 2 : IF (.NOT. ASSOCIATED(scf_env%subspace_env)) &
864 : CALL diag_subspace_env_create(scf_env%subspace_env, scf_section, &
865 2 : dft_control%qs_control%cutoff)
866 2 : CALL diag_subspace_allocate(scf_env%subspace_env, qs_env, mos)
867 2 : IF (do_kpoints) &
868 0 : CPABORT("No subspace diagonlization with kpoint calculation")
869 : END IF
870 : ! OT: check if OT is used instead of diagonalization. Not possible with added MOS at the moment
871 6117 : ELSEIF (scf_control%use_ot) THEN
872 6117 : scf_env%method = ot_method_nr
873 6117 : need_coeff_b = .TRUE.
874 18351 : IF (SUM(ABS(scf_control%added_mos)) > 0) &
875 0 : CPABORT("OT with ADDED_MOS/=0 not implemented")
876 6117 : IF (dft_control%restricted .AND. dft_control%nspins /= 2) &
877 0 : CPABORT("nspin must be 2 for restricted (ROKS)")
878 6117 : IF (do_kpoints) &
879 0 : CPABORT("OT not possible with kpoint calculations")
880 0 : ELSEIF (scf_env%method /= smeagol_method_nr) THEN
881 0 : CPABORT("OT or DIAGONALIZATION have to be set")
882 : END IF
883 44982 : DO ispin = 1, dft_control%nspins
884 44982 : mos(ispin)%use_mo_coeff_b = need_coeff_b
885 : END DO
886 :
887 21207 : END SUBROUTINE qs_scf_ensure_diagonalization
888 :
889 : ! **************************************************************************************************
890 : !> \brief performs those initialisations that need to be done only once
891 : !> (e.g. that only depend on the atomic positions)
892 : !> this will be called in scf
893 : !> \param scf_env ...
894 : !> \param qs_env ...
895 : !> \param scf_section ...
896 : !> \param scf_control ...
897 : !> \par History
898 : !> 03.2006 created [Joost VandeVondele]
899 : ! **************************************************************************************************
900 21207 : SUBROUTINE init_scf_run(scf_env, qs_env, scf_section, scf_control)
901 :
902 : TYPE(qs_scf_env_type), POINTER :: scf_env
903 : TYPE(qs_environment_type), POINTER :: qs_env
904 : TYPE(section_vals_type), POINTER :: scf_section
905 : TYPE(scf_control_type), POINTER :: scf_control
906 :
907 : CHARACTER(LEN=*), PARAMETER :: routineN = 'init_scf_run'
908 :
909 : INTEGER :: after, handle, homo, ii, ikind, ispin, &
910 : iw, nao, ndep, needed_evals, nmo, &
911 : output_unit
912 : LOGICAL :: dft_plus_u_atom, do_kpoints, &
913 : init_u_ramping_each_scf, omit_headers, &
914 : s_minus_half_available
915 : REAL(KIND=dp) :: u_ramping
916 21207 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: evals
917 21207 : REAL(KIND=dp), DIMENSION(:), POINTER :: eigenvalues
918 : TYPE(cp_fm_struct_type), POINTER :: fm_struct
919 : TYPE(cp_fm_type) :: evecs, fm_w
920 : TYPE(cp_fm_type), POINTER :: mo_coeff
921 : TYPE(cp_logger_type), POINTER :: logger
922 21207 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s
923 21207 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s_kp
924 : TYPE(dft_control_type), POINTER :: dft_control
925 : TYPE(kpoint_type), POINTER :: kpoints
926 21207 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
927 : TYPE(mp_para_env_type), POINTER :: para_env
928 21207 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
929 : TYPE(qs_kind_type), POINTER :: qs_kind
930 : TYPE(qs_rho_type), POINTER :: rho
931 : TYPE(xas_environment_type), POINTER :: xas_env
932 :
933 21207 : CALL timeset(routineN, handle)
934 :
935 21207 : NULLIFY (qs_kind_set, matrix_s, dft_control, mos, qs_kind, rho, xas_env, mo_coeff)
936 :
937 21207 : logger => cp_get_default_logger()
938 :
939 21207 : CPASSERT(ASSOCIATED(scf_env))
940 21207 : CPASSERT(ASSOCIATED(qs_env))
941 21207 : NULLIFY (para_env)
942 :
943 21207 : s_minus_half_available = .FALSE.
944 : CALL get_qs_env(qs_env, &
945 : dft_control=dft_control, &
946 : qs_kind_set=qs_kind_set, &
947 : mos=mos, &
948 : rho=rho, &
949 : nelectron_total=scf_env%nelectron, &
950 : do_kpoints=do_kpoints, &
951 : para_env=para_env, &
952 21207 : xas_env=xas_env)
953 :
954 : !Check restricted optimizers available for tblite library
955 21207 : IF (dft_control%qs_control%xtb_control%do_tblite) THEN
956 894 : IF (dft_control%lsd) THEN
957 0 : CPABORT("LSD option not compatible with tblite library.")
958 : END IF
959 894 : IF (scf_env%method == ot_method_nr) THEN
960 0 : CPABORT("OT SCF option not compatible with tblite library.")
961 : END IF
962 : END IF
963 :
964 : ! Calculate ortho matrix
965 21207 : ndep = 0
966 21207 : IF (scf_env%needs_ortho) THEN
967 11460 : CALL get_qs_env(qs_env, matrix_s=matrix_s)
968 11460 : CALL copy_dbcsr_to_fm(matrix_s(1)%matrix, scf_env%ortho)
969 11460 : IF (scf_env%cholesky_method > cholesky_off) THEN
970 11412 : CALL cp_fm_cholesky_decompose(scf_env%ortho)
971 11412 : IF (scf_env%cholesky_method == cholesky_dbcsr) THEN
972 58 : CALL cp_fm_triangular_invert(scf_env%ortho)
973 58 : CALL cp_fm_set_all(scf_env%scf_work2, 0.0_dp)
974 58 : CALL cp_fm_to_fm_triangular(scf_env%ortho, scf_env%scf_work2, "U")
975 58 : CALL copy_fm_to_dbcsr(scf_env%scf_work2, scf_env%ortho_dbcsr)
976 11354 : ELSE IF (scf_env%cholesky_method == cholesky_inverse) THEN
977 34 : CALL cp_fm_to_fm(scf_env%ortho, scf_env%ortho_m1)
978 34 : CALL cp_fm_triangular_invert(scf_env%ortho_m1)
979 : END IF
980 : ELSE
981 48 : CALL cp_fm_get_info(scf_env%ortho, ncol_global=nao)
982 144 : ALLOCATE (evals(nao))
983 1908 : evals = 0
984 :
985 48 : CALL cp_fm_create(evecs, scf_env%ortho%matrix_struct)
986 :
987 : ! Perform an EVD
988 48 : CALL choose_eigv_solver(scf_env%ortho, evecs, evals)
989 :
990 : ! Determine the number of neglectable eigenvalues assuming that the eigenvalues are in ascending order
991 : ! (Required by Lapack)
992 : ndep = 0
993 112 : DO ii = 1, nao
994 112 : IF (evals(ii) > scf_control%eps_eigval) THEN
995 48 : ndep = ii - 1
996 48 : EXIT
997 : END IF
998 : END DO
999 48 : needed_evals = nao - ndep
1000 :
1001 : ! Set the eigenvalue of the eigenvectors belonging to the linear subspace to zero
1002 112 : evals(1:ndep) = 0.0_dp
1003 : ! Determine the eigenvalues of the inverse square root
1004 1844 : evals(ndep + 1:nao) = 1.0_dp/SQRT(evals(ndep + 1:nao))
1005 :
1006 : ! Create reduced matrices
1007 48 : NULLIFY (fm_struct)
1008 : CALL cp_fm_struct_create(fm_struct, template_fmstruct=scf_env%ortho%matrix_struct, &
1009 48 : nrow_global=nao, ncol_global=needed_evals)
1010 :
1011 48 : ALLOCATE (scf_env%ortho_red, scf_env%scf_work2_red)
1012 48 : CALL cp_fm_create(scf_env%ortho_red, fm_struct)
1013 48 : CALL cp_fm_create(scf_env%scf_work2_red, fm_struct)
1014 48 : CALL cp_fm_struct_release(fm_struct)
1015 :
1016 48 : IF (scf_control%level_shift /= 0.0_dp) THEN
1017 : CALL cp_fm_struct_create(fm_struct, template_fmstruct=scf_env%ortho%matrix_struct, &
1018 6 : nrow_global=needed_evals, ncol_global=nao)
1019 :
1020 6 : ALLOCATE (scf_env%ortho_m1_red)
1021 6 : CALL cp_fm_create(scf_env%ortho_m1_red, fm_struct)
1022 6 : CALL cp_fm_struct_release(fm_struct)
1023 : END IF
1024 :
1025 206 : ALLOCATE (scf_env%scf_work1_red(SIZE(scf_env%scf_work1)))
1026 110 : DO ispin = 1, SIZE(scf_env%scf_work1)
1027 : CALL cp_fm_struct_create(fm_struct, template_fmstruct=scf_env%ortho%matrix_struct, &
1028 62 : nrow_global=needed_evals, ncol_global=needed_evals)
1029 62 : CALL cp_fm_create(scf_env%scf_work1_red(ispin), fm_struct)
1030 110 : CALL cp_fm_struct_release(fm_struct)
1031 : END DO
1032 :
1033 : ! Scale the eigenvalues and copy them to
1034 48 : CALL cp_fm_to_fm(evecs, scf_env%ortho_red, needed_evals, ndep + 1, 1)
1035 :
1036 48 : IF (scf_control%level_shift /= 0.0_dp) THEN
1037 6 : CALL cp_fm_transpose(scf_env%ortho_red, scf_env%ortho_m1_red)
1038 : END IF
1039 :
1040 48 : CALL cp_fm_column_scale(scf_env%ortho_red, evals(ndep + 1:))
1041 :
1042 : ! Copy the linear dependent columns to the MO sets and set their orbital energies
1043 : ! to a very large value to reduce the probability of occupying them
1044 110 : DO ispin = 1, SIZE(mos)
1045 62 : CALL get_mo_set(mos(ispin), nmo=nmo, mo_coeff=mo_coeff, homo=homo, eigenvalues=eigenvalues)
1046 62 : IF (needed_evals < nmo) THEN
1047 2 : IF (needed_evals < homo) THEN
1048 : CALL cp_abort(__LOCATION__, &
1049 : "The numerical rank of the overlap matrix is lower than the "// &
1050 : "number of orbitals to be occupied! Check the geometry or increase "// &
1051 0 : "EPS_DEFAULT or EPS_PGF_ORB!")
1052 : END IF
1053 : CALL cp_warn(__LOCATION__, &
1054 : "The numerical rank of the overlap matrix is lower than the number of requested MOs! "// &
1055 : "Reduce the number of MOs to the number of available MOs. If necessary, "// &
1056 2 : "request a lower number of MOs or increase EPS_DEFAULT or EPS_PGF_ORB.")
1057 2 : CALL set_mo_set(mos(ispin), nmo=needed_evals)
1058 : END IF
1059 : ! Copy the last columns to mo_coeff if the container is large enough
1060 62 : CALL cp_fm_to_fm(evecs, mo_coeff, MIN(ndep, MAX(0, nmo - needed_evals)), 1, needed_evals + 1)
1061 : ! Set the corresponding eigenvalues to a large value
1062 : ! This prevents their occupation but still keeps the information on them
1063 182 : eigenvalues(needed_evals + 1:MIN(nao, nmo)) = 1.0_dp/scf_control%eps_eigval
1064 : END DO
1065 :
1066 : ! Obtain ortho from (P)DGEMM, skip the linear dependent columns
1067 : CALL parallel_gemm("N", "T", nao, nao, needed_evals, 1.0_dp, scf_env%ortho_red, evecs, &
1068 48 : 0.0_dp, scf_env%ortho, b_first_col=ndep + 1)
1069 :
1070 48 : IF (scf_control%level_shift /= 0.0_dp) THEN
1071 : ! We need SQRT(evals) of the eigenvalues of H, so 1/SQRT(evals) of ortho_red
1072 168 : evals(ndep + 1:nao) = 1.0_dp/evals(ndep + 1:nao)
1073 6 : CALL cp_fm_row_scale(scf_env%ortho_m1_red, evals(ndep + 1:))
1074 :
1075 : CALL parallel_gemm("T", "T", nao, nao, needed_evals, 1.0_dp, scf_env%ortho_m1_red, evecs, &
1076 6 : 0.0_dp, scf_env%ortho_m1, b_first_col=ndep + 1)
1077 : END IF
1078 :
1079 48 : CALL cp_fm_release(evecs)
1080 :
1081 144 : s_minus_half_available = .TRUE.
1082 : END IF
1083 :
1084 11460 : IF (BTEST(cp_print_key_should_output(logger%iter_info, &
1085 : qs_env%input, "DFT%PRINT%AO_MATRICES/ORTHO"), cp_p_file)) THEN
1086 : iw = cp_print_key_unit_nr(logger, qs_env%input, "DFT%PRINT%AO_MATRICES/ORTHO", &
1087 4 : extension=".Log")
1088 4 : CALL section_vals_val_get(qs_env%input, "DFT%PRINT%AO_MATRICES%NDIGITS", i_val=after)
1089 4 : CALL section_vals_val_get(qs_env%input, "DFT%PRINT%AO_MATRICES%OMIT_HEADERS", l_val=omit_headers)
1090 4 : after = MIN(MAX(after, 1), 16)
1091 : CALL write_fm_with_basis_info(scf_env%ortho, 4, after, qs_env, &
1092 4 : para_env, output_unit=iw, omit_headers=omit_headers)
1093 : CALL cp_print_key_finished_output(iw, logger, qs_env%input, &
1094 4 : "DFT%PRINT%AO_MATRICES/ORTHO")
1095 : END IF
1096 : END IF
1097 :
1098 21207 : CALL get_mo_set(mo_set=mos(1), nao=nao)
1099 :
1100 : ! DFT+U methods based on Lowdin charges need S^(1/2)
1101 21207 : IF (dft_control%dft_plus_u) THEN
1102 92 : IF (dft_control%plus_u_method_id == plus_u_lowdin) THEN
1103 14 : IF (do_kpoints) THEN
1104 0 : CALL get_qs_env(qs_env, kpoints=kpoints, matrix_s_kp=matrix_s_kp)
1105 0 : CALL diag_kp_smat(matrix_s_kp, kpoints, scf_env%scf_work1)
1106 : ELSE
1107 14 : CALL get_qs_env(qs_env, matrix_s=matrix_s)
1108 14 : IF (s_minus_half_available) THEN
1109 : CALL cp_dbcsr_sm_fm_multiply(matrix_s(1)%matrix, scf_env%ortho, &
1110 0 : scf_env%s_half, nao)
1111 : ELSE
1112 14 : CALL copy_dbcsr_to_fm(matrix_s(1)%matrix, scf_env%s_half)
1113 14 : CALL cp_fm_create(fm_w, scf_env%s_half%matrix_struct)
1114 14 : CALL cp_fm_power(scf_env%s_half, fm_w, 0.5_dp, scf_control%eps_eigval, ndep)
1115 14 : CALL cp_fm_release(fm_w)
1116 : END IF
1117 : END IF
1118 : END IF
1119 276 : DO ikind = 1, SIZE(qs_kind_set)
1120 184 : qs_kind => qs_kind_set(ikind)
1121 : CALL get_qs_kind(qs_kind=qs_kind, &
1122 : dft_plus_u_atom=dft_plus_u_atom, &
1123 : u_ramping=u_ramping, &
1124 184 : init_u_ramping_each_scf=init_u_ramping_each_scf)
1125 276 : IF (dft_plus_u_atom .AND. (u_ramping /= 0.0_dp)) THEN
1126 24 : IF (init_u_ramping_each_scf) THEN
1127 12 : CALL set_qs_kind(qs_kind=qs_kind, u_minus_j=0.0_dp)
1128 : END IF
1129 : END IF
1130 : END DO
1131 : END IF
1132 :
1133 : ! extrapolate outer loop variables
1134 21207 : IF (scf_control%outer_scf%have_scf) THEN
1135 4141 : CALL outer_loop_extrapolate(qs_env)
1136 : END IF
1137 :
1138 : ! initializes rho and the mos
1139 21207 : IF (ASSOCIATED(qs_env%xas_env)) THEN
1140 : ! if just optimized wfn, e.g. ground state
1141 : ! changes come from a perturbation, e.g., the occupation numbers
1142 : ! it could be generalized for other cases, at the moment used only for core level spectroscopy
1143 : ! initialize the density with the localized mos
1144 82 : CALL xas_initialize_rho(qs_env, scf_env, scf_control)
1145 : ELSE
1146 : CALL scf_env_initial_rho_setup(scf_env, qs_env=qs_env, &
1147 21125 : scf_section=scf_section, scf_control=scf_control)
1148 : END IF
1149 :
1150 : ! Frozen density approximation
1151 21207 : IF (ASSOCIATED(qs_env%wf_history)) THEN
1152 21207 : IF (qs_env%wf_history%interpolation_method_nr == wfi_frozen_method_nr) THEN
1153 12 : IF (.NOT. ASSOCIATED(qs_env%wf_history%past_states(1)%snapshot)) THEN
1154 4 : CALL wfi_update(qs_env%wf_history, qs_env=qs_env, dt=1.0_dp)
1155 4 : ALLOCATE (qs_env%wf_history%past_states(1)%snapshot%rho_frozen)
1156 4 : CALL qs_rho_create(qs_env%wf_history%past_states(1)%snapshot%rho_frozen)
1157 : CALL duplicate_rho_type(rho_input=rho, &
1158 : rho_output=qs_env%wf_history%past_states(1)%snapshot%rho_frozen, &
1159 4 : qs_env=qs_env)
1160 : END IF
1161 : END IF
1162 : END IF
1163 :
1164 : !image charge method, calculate image_matrix if required
1165 21207 : IF (qs_env%qmmm) THEN
1166 3802 : IF (qs_env%qmmm .AND. qs_env%qmmm_env_qm%image_charge) THEN
1167 : CALL conditional_calc_image_matrix(qs_env=qs_env, &
1168 20 : qmmm_env=qs_env%qmmm_env_qm)
1169 : END IF
1170 : END IF
1171 :
1172 : output_unit = cp_print_key_unit_nr(logger, scf_section, "PRINT%PROGRAM_RUN_INFO", &
1173 21207 : extension=".scfLog")
1174 21207 : CALL qs_scf_initial_info(output_unit, mos, dft_control, ndep)
1175 : CALL cp_print_key_finished_output(output_unit, logger, scf_section, &
1176 21207 : "PRINT%PROGRAM_RUN_INFO")
1177 :
1178 21207 : CALL timestop(handle)
1179 :
1180 42414 : END SUBROUTINE init_scf_run
1181 :
1182 : ! **************************************************************************************************
1183 : !> \brief Initializes rho and the mos, so that an scf cycle can start
1184 : !> \param scf_env the scf env in which to do the scf
1185 : !> \param qs_env the qs env the scf_env lives in
1186 : !> \param scf_section ...
1187 : !> \param scf_control ...
1188 : !> \par History
1189 : !> 02.2003 created [fawzi]
1190 : !> \author fawzi
1191 : ! **************************************************************************************************
1192 21125 : SUBROUTINE scf_env_initial_rho_setup(scf_env, qs_env, scf_section, scf_control)
1193 : TYPE(qs_scf_env_type), POINTER :: scf_env
1194 : TYPE(qs_environment_type), POINTER :: qs_env
1195 : TYPE(section_vals_type), POINTER :: scf_section
1196 : TYPE(scf_control_type), POINTER :: scf_control
1197 :
1198 : CHARACTER(len=*), PARAMETER :: routineN = 'scf_env_initial_rho_setup'
1199 :
1200 : INTEGER :: extrapolation_method_nr, handle, ispin, &
1201 : nmo, output_unit
1202 : LOGICAL :: do_harris, orthogonal_wf
1203 : TYPE(cp_fm_type), POINTER :: mo_coeff
1204 : TYPE(cp_logger_type), POINTER :: logger
1205 : TYPE(dft_control_type), POINTER :: dft_control
1206 : TYPE(harris_type), POINTER :: harris_env
1207 21125 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
1208 : TYPE(mp_para_env_type), POINTER :: para_env
1209 : TYPE(qs_rho_type), POINTER :: rho
1210 21125 : TYPE(rho_atom_type), DIMENSION(:), POINTER :: rho_atom
1211 :
1212 21125 : CALL timeset(routineN, handle)
1213 21125 : NULLIFY (mo_coeff, rho, dft_control, para_env, mos)
1214 21125 : logger => cp_get_default_logger()
1215 21125 : CPASSERT(ASSOCIATED(scf_env))
1216 21125 : CPASSERT(ASSOCIATED(qs_env))
1217 :
1218 : CALL get_qs_env(qs_env, &
1219 : rho=rho, &
1220 : mos=mos, &
1221 : dft_control=dft_control, &
1222 21125 : para_env=para_env)
1223 :
1224 21125 : do_harris = qs_env%harris_method
1225 :
1226 21125 : extrapolation_method_nr = wfi_use_guess_method_nr
1227 21125 : IF (ASSOCIATED(qs_env%wf_history)) THEN
1228 : CALL wfi_extrapolate(qs_env%wf_history, &
1229 : qs_env=qs_env, dt=1.0_dp, &
1230 : extrapolation_method_nr=extrapolation_method_nr, &
1231 21125 : orthogonal_wf=orthogonal_wf)
1232 : ! wfi_use_guess_method_nr the wavefunctions are not yet initialized
1233 : IF ((.NOT. orthogonal_wf) .AND. &
1234 21125 : (scf_env%method == ot_method_nr) .AND. &
1235 : (.NOT. (extrapolation_method_nr == wfi_use_guess_method_nr))) THEN
1236 0 : DO ispin = 1, SIZE(mos)
1237 0 : CALL get_mo_set(mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
1238 0 : CALL reorthogonalize_vectors(qs_env, v_matrix=mo_coeff, n_col=nmo)
1239 0 : IF (dft_control%hairy_probes .EQV. .TRUE.) THEN
1240 0 : scf_control%smear%do_smear = .FALSE.
1241 : CALL set_mo_occupation(mo_set=mos(ispin), &
1242 0 : smear=scf_control%smear, probe=dft_control%probe)
1243 : ELSE
1244 : CALL set_mo_occupation(mo_set=mos(ispin), &
1245 0 : smear=scf_control%smear)
1246 : END IF
1247 : END DO
1248 : END IF
1249 : END IF
1250 :
1251 21125 : IF (.NOT. do_harris) THEN
1252 : output_unit = cp_print_key_unit_nr(logger, scf_section, "PRINT%PROGRAM_RUN_INFO", &
1253 21109 : extension=".scfLog")
1254 21109 : IF (output_unit > 0) THEN
1255 : WRITE (UNIT=output_unit, FMT="(/,T2,A,I0)") &
1256 : "Extrapolation method: "// &
1257 10737 : TRIM(wfi_get_method_label(extrapolation_method_nr))
1258 10737 : IF (extrapolation_method_nr == wfi_ps_method_nr) THEN
1259 : WRITE (UNIT=output_unit, FMT="(T2,A,I0,A)") &
1260 188 : "Extrapolation order: ", &
1261 376 : MAX((MIN(qs_env%wf_history%memory_depth, qs_env%wf_history%snapshot_count) - 1), 0)
1262 : END IF
1263 : END IF
1264 : CALL cp_print_key_finished_output(output_unit, logger, scf_section, &
1265 21109 : "PRINT%PROGRAM_RUN_INFO")
1266 : END IF
1267 :
1268 : IF (do_harris) THEN
1269 16 : CALL get_qs_env(qs_env, harris_env=harris_env)
1270 16 : CALL harris_density_update(qs_env, harris_env)
1271 16 : CALL qs_rho_update_rho(rho, qs_env=qs_env)
1272 16 : CALL qs_ks_did_change(qs_env%ks_env, rho_changed=.TRUE.)
1273 21109 : ELSE IF (extrapolation_method_nr == wfi_use_guess_method_nr) THEN
1274 7121 : CALL calculate_first_density_matrix(scf_env=scf_env, qs_env=qs_env)
1275 7121 : CALL qs_rho_update_rho(rho, qs_env=qs_env)
1276 7121 : CALL qs_ks_did_change(qs_env%ks_env, rho_changed=.TRUE.)
1277 : END IF
1278 :
1279 : ! Some preparation for the mixing
1280 21125 : IF (scf_env%mixing_method > 1) THEN
1281 274 : IF (dft_control%qs_control%gapw) THEN
1282 40 : CALL get_qs_env(qs_env=qs_env, rho_atom_set=rho_atom)
1283 : CALL mixing_init(scf_env%mixing_method, rho, scf_env%mixing_store, &
1284 40 : para_env, rho_atom=rho_atom)
1285 234 : ELSEIF (dft_control%qs_control%dftb .OR. dft_control%qs_control%xtb) THEN
1286 54 : CALL charge_mixing_init(scf_env%mixing_store)
1287 180 : ELSEIF (dft_control%qs_control%semi_empirical) THEN
1288 0 : CPABORT('SE Code not possible')
1289 : ELSE
1290 : CALL mixing_init(scf_env%mixing_method, rho, scf_env%mixing_store, &
1291 180 : para_env)
1292 : END IF
1293 : END IF
1294 :
1295 44736 : DO ispin = 1, SIZE(mos) !fm->dbcsr
1296 44736 : IF (mos(ispin)%use_mo_coeff_b) THEN
1297 : CALL copy_fm_to_dbcsr(mos(ispin)%mo_coeff, &
1298 7127 : mos(ispin)%mo_coeff_b) !fm->dbcsr
1299 : END IF
1300 : END DO !fm->dbcsr
1301 :
1302 21125 : CALL timestop(handle)
1303 :
1304 21125 : END SUBROUTINE scf_env_initial_rho_setup
1305 :
1306 : END MODULE qs_scf_initialization
|