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