Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \brief Routines for the Quickstep SCF run.
10 : !> \par History
11 : !> - Joost VandeVondele (02.2002)
12 : !> added code for: incremental (pab and gvg) update
13 : !> initialisation (init_cube, l_info)
14 : !> - Joost VandeVondele (02.2002)
15 : !> called the poisson code of the classical part
16 : !> this takes into account the spherical cutoff and allows for
17 : !> isolated systems
18 : !> - Joost VandeVondele (02.2002)
19 : !> added multiple grid feature
20 : !> changed to spherical cutoff consistently (?)
21 : !> therefore removed the gradient correct functionals
22 : !> - updated with the new QS data structures (10.04.02,MK)
23 : !> - copy_matrix replaced by transfer_matrix (11.04.02,MK)
24 : !> - nrebuild_rho and nrebuild_gvg unified (12.04.02,MK)
25 : !> - set_mo_occupation for smearing of the MO occupation numbers
26 : !> (17.04.02,MK)
27 : !> - MO level shifting added (22.04.02,MK)
28 : !> - Usage of TYPE mo_set_p_type
29 : !> - Joost VandeVondele (05.2002)
30 : !> added cholesky based diagonalisation
31 : !> - 05.2002 added pao method [fawzi]
32 : !> - parallel FFT (JGH 22.05.2002)
33 : !> - 06.2002 moved KS matrix construction to qs_build_KS_matrix.F [fawzi]
34 : !> - started to include more LSD (01.2003,Joost VandeVondele)
35 : !> - 02.2003 scf_env [fawzi]
36 : !> - got rid of nrebuild (01.2004, Joost VandeVondele)
37 : !> - 10.2004 removed pao [fawzi]
38 : !> - 03.2006 large cleaning action [Joost VandeVondele]
39 : !> - High-spin ROKS added (05.04.06,MK)
40 : !> - Mandes (10.2013)
41 : !> intermediate energy communication with external communicator added
42 : !> - kpoints (08.2014, JGH)
43 : !> - unified k-point and gamma-point code (2014.11) [Ole Schuett]
44 : !> - added extra SCF loop for CDFT constraints (12.2015) [Nico Holmberg]
45 : !> \author Matthias Krack (30.04.2001)
46 : ! **************************************************************************************************
47 : MODULE qs_scf
48 : USE atomic_kind_types, ONLY: atomic_kind_type
49 : USE cp_control_types, ONLY: dft_control_type
50 : USE cp_dbcsr_operations, ONLY: copy_dbcsr_to_fm,&
51 : dbcsr_deallocate_matrix_set
52 : USE cp_files, ONLY: close_file
53 : USE cp_fm_types, ONLY: cp_fm_create,&
54 : cp_fm_release,&
55 : cp_fm_to_fm,&
56 : cp_fm_type
57 : USE cp_log_handling, ONLY: cp_add_default_logger,&
58 : cp_get_default_logger,&
59 : cp_logger_release,&
60 : cp_logger_type,&
61 : cp_rm_default_logger,&
62 : cp_to_string
63 : USE cp_output_handling, ONLY: cp_add_iter_level,&
64 : cp_iterate,&
65 : cp_p_file,&
66 : cp_print_key_should_output,&
67 : cp_print_key_unit_nr,&
68 : cp_rm_iter_level
69 : USE cp_result_methods, ONLY: get_results,&
70 : test_for_result
71 : USE cp_result_types, ONLY: cp_result_type
72 : USE dbcsr_api, ONLY: dbcsr_copy,&
73 : dbcsr_deallocate_matrix,&
74 : dbcsr_get_info,&
75 : dbcsr_init_p,&
76 : dbcsr_p_type,&
77 : dbcsr_set,&
78 : dbcsr_type
79 : USE ec_env_types, ONLY: energy_correction_type
80 : USE input_constants, ONLY: &
81 : broyden_type_1, broyden_type_1_explicit, broyden_type_1_explicit_ls, broyden_type_1_ls, &
82 : broyden_type_2, broyden_type_2_explicit, broyden_type_2_explicit_ls, broyden_type_2_ls, &
83 : cdft2ot, history_guess, ot2cdft, ot_precond_full_all, ot_precond_full_single, &
84 : ot_precond_full_single_inverse, ot_precond_none, ot_precond_s_inverse, &
85 : outer_scf_becke_constraint, outer_scf_hirshfeld_constraint, outer_scf_optimizer_broyden, &
86 : outer_scf_optimizer_newton_ls
87 : USE input_section_types, ONLY: section_vals_get_subs_vals,&
88 : section_vals_type
89 : USE kinds, ONLY: default_path_length,&
90 : default_string_length,&
91 : dp
92 : USE kpoint_io, ONLY: write_kpoints_restart
93 : USE kpoint_types, ONLY: kpoint_type
94 : USE machine, ONLY: m_flush,&
95 : m_walltime
96 : USE mathlib, ONLY: invert_matrix
97 : USE message_passing, ONLY: mp_comm_type,&
98 : mp_para_env_type
99 : USE particle_types, ONLY: particle_type
100 : USE preconditioner, ONLY: prepare_preconditioner,&
101 : restart_preconditioner
102 : USE pw_env_types, ONLY: pw_env_get,&
103 : pw_env_type
104 : USE pw_pool_types, ONLY: pw_pool_type
105 : USE qs_block_davidson_types, ONLY: block_davidson_deallocate
106 : USE qs_cdft_scf_utils, ONLY: build_diagonal_jacobian,&
107 : create_tmp_logger,&
108 : initialize_inverse_jacobian,&
109 : prepare_jacobian_stencil,&
110 : print_inverse_jacobian,&
111 : restart_inverse_jacobian
112 : USE qs_cdft_types, ONLY: cdft_control_type
113 : USE qs_charges_types, ONLY: qs_charges_type
114 : USE qs_density_matrices, ONLY: calculate_density_matrix
115 : USE qs_density_mixing_types, ONLY: gspace_mixing_nr
116 : USE qs_diis, ONLY: qs_diis_b_clear,&
117 : qs_diis_b_clear_kp,&
118 : qs_diis_b_create,&
119 : qs_diis_b_create_kp
120 : USE qs_energy_types, ONLY: qs_energy_type
121 : USE qs_environment_types, ONLY: get_qs_env,&
122 : qs_environment_type,&
123 : set_qs_env
124 : USE qs_integrate_potential, ONLY: integrate_v_rspace
125 : USE qs_kind_types, ONLY: qs_kind_type
126 : USE qs_ks_methods, ONLY: qs_ks_update_qs_env
127 : USE qs_ks_types, ONLY: qs_ks_did_change,&
128 : qs_ks_env_type
129 : USE qs_mo_io, ONLY: write_mo_set_to_restart
130 : USE qs_mo_methods, ONLY: make_basis_simple,&
131 : make_basis_sm
132 : USE qs_mo_occupation, ONLY: set_mo_occupation
133 : USE qs_mo_types, ONLY: deallocate_mo_set,&
134 : duplicate_mo_set,&
135 : get_mo_set,&
136 : mo_set_type,&
137 : reassign_allocated_mos
138 : USE qs_ot, ONLY: qs_ot_new_preconditioner
139 : USE qs_ot_scf, ONLY: ot_scf_init,&
140 : ot_scf_read_input
141 : USE qs_outer_scf, ONLY: outer_loop_gradient,&
142 : outer_loop_optimize,&
143 : outer_loop_purge_history,&
144 : outer_loop_switch,&
145 : outer_loop_update_qs_env
146 : USE qs_rho_methods, ONLY: qs_rho_update_rho
147 : USE qs_rho_types, ONLY: qs_rho_get,&
148 : qs_rho_type
149 : USE qs_scf_initialization, ONLY: qs_scf_env_initialize
150 : USE qs_scf_loop_utils, ONLY: qs_scf_check_inner_exit,&
151 : qs_scf_check_outer_exit,&
152 : qs_scf_density_mixing,&
153 : qs_scf_inner_finalize,&
154 : qs_scf_new_mos,&
155 : qs_scf_new_mos_kp,&
156 : qs_scf_rho_update,&
157 : qs_scf_set_loop_flags
158 : USE qs_scf_output, ONLY: qs_scf_cdft_info,&
159 : qs_scf_cdft_initial_info,&
160 : qs_scf_loop_info,&
161 : qs_scf_loop_print,&
162 : qs_scf_outer_loop_info,&
163 : qs_scf_write_mos
164 : USE qs_scf_post_scf, ONLY: qs_scf_compute_properties
165 : USE qs_scf_types, ONLY: block_davidson_diag_method_nr,&
166 : block_krylov_diag_method_nr,&
167 : filter_matrix_diag_method_nr,&
168 : general_diag_method_nr,&
169 : ot_diag_method_nr,&
170 : ot_method_nr,&
171 : qs_scf_env_type,&
172 : special_diag_method_nr
173 : USE qs_wf_history_methods, ONLY: wfi_purge_history,&
174 : wfi_update
175 : USE scf_control_types, ONLY: scf_control_type
176 : #include "./base/base_uses.f90"
177 :
178 : IMPLICIT NONE
179 :
180 : PRIVATE
181 :
182 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_scf'
183 : LOGICAL, PRIVATE :: reuse_precond = .FALSE.
184 : LOGICAL, PRIVATE :: used_history = .FALSE.
185 :
186 : PUBLIC :: scf, scf_env_cleanup, scf_env_do_scf, cdft_scf
187 :
188 : CONTAINS
189 :
190 : ! **************************************************************************************************
191 : !> \brief perform an scf procedure in the given qs_env
192 : !> \param qs_env the qs_environment where to perform the scf procedure
193 : !> \param has_converged ...
194 : !> \param total_scf_steps ...
195 : !> \par History
196 : !> 02.2003 introduced scf_env, moved real work to scf_env_do_scf [fawzi]
197 : !> \author fawzi
198 : !> \note
199 : ! **************************************************************************************************
200 18135 : SUBROUTINE scf(qs_env, has_converged, total_scf_steps)
201 : TYPE(qs_environment_type), POINTER :: qs_env
202 : LOGICAL, INTENT(OUT), OPTIONAL :: has_converged
203 : INTEGER, INTENT(OUT), OPTIONAL :: total_scf_steps
204 :
205 : INTEGER :: ihistory, max_scf_tmp, tsteps
206 : LOGICAL :: converged, outer_scf_loop, should_stop
207 : LOGICAL, SAVE :: first_step_flag = .TRUE.
208 18135 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: gradient_history, variable_history
209 : TYPE(cp_logger_type), POINTER :: logger
210 : TYPE(dft_control_type), POINTER :: dft_control
211 : TYPE(qs_scf_env_type), POINTER :: scf_env
212 : TYPE(scf_control_type), POINTER :: scf_control
213 : TYPE(section_vals_type), POINTER :: dft_section, input, scf_section
214 :
215 18135 : NULLIFY (scf_env)
216 18135 : logger => cp_get_default_logger()
217 18135 : CPASSERT(ASSOCIATED(qs_env))
218 18135 : IF (PRESENT(has_converged)) THEN
219 0 : has_converged = .FALSE.
220 : END IF
221 18135 : IF (PRESENT(total_scf_steps)) THEN
222 0 : total_scf_steps = 0
223 : END IF
224 : CALL get_qs_env(qs_env, scf_env=scf_env, input=input, &
225 18135 : dft_control=dft_control, scf_control=scf_control)
226 18135 : IF (scf_control%max_scf > 0) THEN
227 :
228 17493 : dft_section => section_vals_get_subs_vals(input, "DFT")
229 17493 : scf_section => section_vals_get_subs_vals(dft_section, "SCF")
230 :
231 17493 : IF (.NOT. ASSOCIATED(scf_env)) THEN
232 5399 : CALL qs_scf_env_initialize(qs_env, scf_env)
233 : ! Moved here from qs_scf_env_initialize to be able to have more scf_env
234 5399 : CALL set_qs_env(qs_env, scf_env=scf_env)
235 : ELSE
236 12094 : CALL qs_scf_env_initialize(qs_env, scf_env)
237 : END IF
238 :
239 17493 : IF ((scf_control%density_guess .EQ. history_guess) .AND. (first_step_flag)) THEN
240 2 : max_scf_tmp = scf_control%max_scf
241 2 : scf_control%max_scf = 1
242 2 : outer_scf_loop = scf_control%outer_scf%have_scf
243 2 : scf_control%outer_scf%have_scf = .FALSE.
244 : END IF
245 :
246 17493 : IF (.NOT. dft_control%qs_control%cdft) THEN
247 : CALL scf_env_do_scf(scf_env=scf_env, scf_control=scf_control, qs_env=qs_env, &
248 17167 : converged=converged, should_stop=should_stop, total_scf_steps=tsteps)
249 : ELSE
250 : ! Third SCF loop needed for CDFT with OT to properly restart OT inner loop
251 326 : CALL cdft_scf(qs_env=qs_env, should_stop=should_stop)
252 : END IF
253 :
254 : ! If SCF has not converged, then we should not start MP2
255 17493 : IF (ASSOCIATED(qs_env%mp2_env)) qs_env%mp2_env%hf_fail = .NOT. converged
256 :
257 : ! Add the converged outer_scf SCF gradient(s)/variable(s) to history
258 17493 : IF (scf_control%outer_scf%have_scf) THEN
259 3899 : ihistory = scf_env%outer_scf%iter_count
260 : CALL get_qs_env(qs_env, gradient_history=gradient_history, &
261 3899 : variable_history=variable_history)
262 : ! We only store the latest two values
263 7828 : gradient_history(:, 1) = gradient_history(:, 2)
264 15656 : gradient_history(:, 2) = scf_env%outer_scf%gradient(:, ihistory)
265 7828 : variable_history(:, 1) = variable_history(:, 2)
266 15656 : variable_history(:, 2) = scf_env%outer_scf%variables(:, ihistory)
267 : ! Reset flag
268 3899 : IF (used_history) used_history = .FALSE.
269 : ! Update a counter and check if the Jacobian should be deallocated
270 3899 : IF (ASSOCIATED(scf_env%outer_scf%inv_jacobian)) THEN
271 64 : scf_control%outer_scf%cdft_opt_control%ijacobian(2) = scf_control%outer_scf%cdft_opt_control%ijacobian(2) + 1
272 : IF (scf_control%outer_scf%cdft_opt_control%ijacobian(2) .GE. &
273 64 : scf_control%outer_scf%cdft_opt_control%jacobian_freq(2) .AND. &
274 : scf_control%outer_scf%cdft_opt_control%jacobian_freq(2) > 0) &
275 50 : scf_env%outer_scf%deallocate_jacobian = .TRUE.
276 : END IF
277 : END IF
278 : ! *** add the converged wavefunction to the wavefunction history
279 17493 : IF ((ASSOCIATED(qs_env%wf_history)) .AND. &
280 : ((scf_control%density_guess .NE. history_guess) .OR. &
281 : (.NOT. first_step_flag))) THEN
282 17491 : IF (.NOT. dft_control%qs_control%cdft) THEN
283 17165 : CALL wfi_update(qs_env%wf_history, qs_env=qs_env, dt=1.0_dp)
284 : ELSE
285 326 : IF (dft_control%qs_control%cdft_control%should_purge) THEN
286 0 : CALL wfi_purge_history(qs_env)
287 0 : CALL outer_loop_purge_history(qs_env)
288 0 : dft_control%qs_control%cdft_control%should_purge = .FALSE.
289 : ELSE
290 326 : CALL wfi_update(qs_env%wf_history, qs_env=qs_env, dt=1.0_dp)
291 : END IF
292 : END IF
293 2 : ELSE IF ((scf_control%density_guess .EQ. history_guess) .AND. &
294 : (first_step_flag)) THEN
295 2 : scf_control%max_scf = max_scf_tmp
296 2 : scf_control%outer_scf%have_scf = outer_scf_loop
297 2 : first_step_flag = .FALSE.
298 : END IF
299 :
300 : ! *** compute properties that depend on the converged wavefunction
301 17493 : IF (.NOT. (should_stop)) CALL qs_scf_compute_properties(qs_env)
302 :
303 : ! *** cleanup
304 17493 : CALL scf_env_cleanup(scf_env)
305 17493 : IF (dft_control%qs_control%cdft) &
306 326 : CALL cdft_control_cleanup(dft_control%qs_control%cdft_control)
307 :
308 17493 : IF (PRESENT(has_converged)) THEN
309 0 : has_converged = converged
310 : END IF
311 17493 : IF (PRESENT(total_scf_steps)) THEN
312 0 : total_scf_steps = tsteps
313 : END IF
314 :
315 : END IF
316 :
317 18135 : END SUBROUTINE scf
318 :
319 : ! **************************************************************************************************
320 : !> \brief perform an scf loop
321 : !> \param scf_env the scf_env where to perform the scf procedure
322 : !> \param scf_control ...
323 : !> \param qs_env the qs_env, the scf_env lives in
324 : !> \param converged will be true / false if converged is reached
325 : !> \param should_stop ...
326 : !> \param total_scf_steps ...
327 : !> \par History
328 : !> long history, see cvs and qs_scf module history
329 : !> 02.2003 introduced scf_env [fawzi]
330 : !> 09.2005 Frozen density approximation [TdK]
331 : !> 06.2007 Check for SCF iteration count early [jgh]
332 : !> 10.2019 switch_surf_dip [SGh]
333 : !> \author Matthias Krack
334 : !> \note
335 : ! **************************************************************************************************
336 17795 : SUBROUTINE scf_env_do_scf(scf_env, scf_control, qs_env, converged, should_stop, total_scf_steps)
337 :
338 : TYPE(qs_scf_env_type), POINTER :: scf_env
339 : TYPE(scf_control_type), POINTER :: scf_control
340 : TYPE(qs_environment_type), POINTER :: qs_env
341 : LOGICAL, INTENT(OUT) :: converged, should_stop
342 : INTEGER, INTENT(OUT) :: total_scf_steps
343 :
344 : CHARACTER(LEN=*), PARAMETER :: routineN = 'scf_env_do_scf'
345 :
346 : CHARACTER(LEN=default_string_length) :: description, name
347 : INTEGER :: ext_master_id, handle, handle2, i_tmp, &
348 : ic, ispin, iter_count, output_unit, &
349 : scf_energy_message_tag, total_steps
350 : LOGICAL :: diis_step, do_kpoints, energy_only, exit_inner_loop, exit_outer_loop, &
351 : inner_loop_converged, just_energy, outer_loop_converged
352 : REAL(KIND=dp) :: t1, t2
353 : REAL(KIND=dp), DIMENSION(3) :: res_val_3
354 17795 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
355 : TYPE(cp_logger_type), POINTER :: logger
356 : TYPE(cp_result_type), POINTER :: results
357 17795 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: rho_ao_kp
358 : TYPE(dft_control_type), POINTER :: dft_control
359 : TYPE(energy_correction_type), POINTER :: ec_env
360 : TYPE(kpoint_type), POINTER :: kpoints
361 17795 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos, mos_last_converged
362 : TYPE(mp_comm_type) :: external_comm
363 : TYPE(mp_para_env_type), POINTER :: para_env
364 17795 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
365 : TYPE(pw_env_type), POINTER :: pw_env
366 : TYPE(qs_charges_type), POINTER :: qs_charges
367 : TYPE(qs_energy_type), POINTER :: energy
368 17795 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
369 : TYPE(qs_ks_env_type), POINTER :: ks_env
370 : TYPE(qs_rho_type), POINTER :: rho
371 : TYPE(section_vals_type), POINTER :: dft_section, input, scf_section
372 :
373 17795 : CALL timeset(routineN, handle)
374 :
375 17795 : NULLIFY (dft_control, rho, energy, &
376 17795 : logger, qs_charges, ks_env, mos, atomic_kind_set, qs_kind_set, &
377 17795 : particle_set, dft_section, input, &
378 17795 : scf_section, para_env, results, kpoints, pw_env, rho_ao_kp, mos_last_converged)
379 :
380 17795 : CPASSERT(ASSOCIATED(scf_env))
381 17795 : CPASSERT(ASSOCIATED(qs_env))
382 :
383 17795 : logger => cp_get_default_logger()
384 17795 : t1 = m_walltime()
385 :
386 : CALL get_qs_env(qs_env=qs_env, &
387 : energy=energy, &
388 : particle_set=particle_set, &
389 : qs_charges=qs_charges, &
390 : ks_env=ks_env, &
391 : atomic_kind_set=atomic_kind_set, &
392 : qs_kind_set=qs_kind_set, &
393 : rho=rho, &
394 : mos=mos, &
395 : input=input, &
396 : dft_control=dft_control, &
397 : do_kpoints=do_kpoints, &
398 : kpoints=kpoints, &
399 : results=results, &
400 : pw_env=pw_env, &
401 17795 : para_env=para_env)
402 :
403 17795 : CALL qs_rho_get(rho, rho_ao_kp=rho_ao_kp)
404 :
405 17795 : dft_section => section_vals_get_subs_vals(input, "DFT")
406 17795 : scf_section => section_vals_get_subs_vals(dft_section, "SCF")
407 :
408 : output_unit = cp_print_key_unit_nr(logger, scf_section, "PRINT%PROGRAM_RUN_INFO", &
409 17795 : extension=".scfLog")
410 :
411 17795 : IF (output_unit > 0) WRITE (UNIT=output_unit, FMT="(/,/,T2,A)") &
412 9080 : "SCF WAVEFUNCTION OPTIMIZATION"
413 :
414 : ! when switch_surf_dip is switched on, indicate storing mos from the last converged step
415 17795 : IF (dft_control%switch_surf_dip) THEN
416 2 : CALL get_qs_env(qs_env, mos_last_converged=mos_last_converged)
417 4 : DO ispin = 1, dft_control%nspins
418 4 : CALL reassign_allocated_mos(mos(ispin), mos_last_converged(ispin))
419 : END DO
420 2 : IF (output_unit > 0) WRITE (UNIT=output_unit, FMT="(/,/,T2,A)") &
421 1 : "COPIED mos_last_converged ---> mos"
422 : END IF
423 :
424 17795 : IF ((output_unit > 0) .AND. (.NOT. scf_control%use_ot)) THEN
425 : WRITE (UNIT=output_unit, &
426 : FMT="(/,T3,A,T12,A,T31,A,T39,A,T59,A,T75,A,/,T3,A)") &
427 6103 : "Step", "Update method", "Time", "Convergence", "Total energy", "Change", &
428 12206 : REPEAT("-", 78)
429 : END IF
430 17795 : CALL cp_add_iter_level(logger%iter_info, "QS_SCF")
431 :
432 : ! check for external communicator and if the intermediate energy should be sent
433 71180 : res_val_3(:) = -1.0_dp
434 17795 : description = "[EXT_SCF_ENER_COMM]"
435 17795 : IF (test_for_result(results, description=description)) THEN
436 : CALL get_results(results, description=description, &
437 0 : values=res_val_3, n_entries=i_tmp)
438 0 : CPASSERT(i_tmp .EQ. 3)
439 0 : IF (ALL(res_val_3(:) .LE. 0.0)) &
440 : CALL cp_abort(__LOCATION__, &
441 : " Trying to access result ("//TRIM(description)// &
442 0 : ") which is not correctly stored.")
443 0 : CALL external_comm%set_handle(NINT(res_val_3(1)))
444 : END IF
445 17795 : ext_master_id = NINT(res_val_3(2))
446 17795 : scf_energy_message_tag = NINT(res_val_3(3))
447 :
448 : ! *** outer loop of the scf, can treat other variables,
449 : ! *** such as lagrangian multipliers
450 17795 : scf_env%outer_scf%iter_count = 0
451 17795 : iter_count = 0
452 17795 : total_steps = 0
453 17795 : energy%tot_old = 0.0_dp
454 :
455 911 : scf_outer_loop: DO
456 :
457 : CALL init_scf_loop(scf_env=scf_env, qs_env=qs_env, &
458 18706 : scf_section=scf_section)
459 :
460 : CALL qs_scf_set_loop_flags(scf_env, diis_step, &
461 18706 : energy_only, just_energy, exit_inner_loop)
462 :
463 : ! decide whether to switch off dipole correction for convergence purposes
464 18706 : dft_control%surf_dip_correct_switch = dft_control%correct_surf_dip
465 18706 : IF ((dft_control%correct_surf_dip) .AND. (scf_control%outer_scf%have_scf) .AND. &
466 : (scf_env%outer_scf%iter_count > FLOOR(scf_control%outer_scf%max_scf/2.0_dp))) THEN
467 0 : IF (dft_control%switch_surf_dip) THEN
468 0 : dft_control%surf_dip_correct_switch = .FALSE.
469 0 : IF (output_unit > 0) WRITE (UNIT=output_unit, FMT="(/,/,T2,A)") &
470 0 : "SURFACE DIPOLE CORRECTION switched off"
471 : END IF
472 : END IF
473 150467 : scf_loop: DO
474 :
475 150467 : CALL timeset(routineN//"_inner_loop", handle2)
476 :
477 150467 : scf_env%iter_count = scf_env%iter_count + 1
478 150467 : iter_count = iter_count + 1
479 150467 : CALL cp_iterate(logger%iter_info, last=.FALSE., iter_nr=iter_count)
480 :
481 150467 : IF (output_unit > 0) CALL m_flush(output_unit)
482 :
483 150467 : total_steps = total_steps + 1
484 150467 : just_energy = energy_only
485 :
486 : CALL qs_ks_update_qs_env(qs_env, just_energy=just_energy, &
487 150467 : calculate_forces=.FALSE.)
488 :
489 : ! print 'heavy weight' or relatively expensive quantities
490 150467 : CALL qs_scf_loop_print(qs_env, scf_env, para_env)
491 :
492 150467 : IF (do_kpoints) THEN
493 : ! kpoints
494 5320 : CALL qs_scf_new_mos_kp(qs_env, scf_env, scf_control, diis_step)
495 : ELSE
496 : ! Gamma points only
497 145147 : CALL qs_scf_new_mos(qs_env, scf_env, scf_control, scf_section, diis_step, energy_only)
498 : END IF
499 :
500 : ! Print requested MO information (can be computationally expensive with OT)
501 150467 : CALL qs_scf_write_mos(qs_env, scf_env, final_mos=.FALSE.)
502 :
503 150467 : CALL qs_scf_density_mixing(scf_env, rho, para_env, diis_step)
504 :
505 150467 : t2 = m_walltime()
506 :
507 150467 : CALL qs_scf_loop_info(scf_env, output_unit, just_energy, t1, t2, energy)
508 :
509 150467 : IF (.NOT. just_energy) energy%tot_old = energy%total
510 :
511 : ! check for external communicator and if the intermediate energy should be sent
512 150467 : IF (scf_energy_message_tag .GT. 0) THEN
513 0 : CALL external_comm%send(energy%total, ext_master_id, scf_energy_message_tag)
514 : END IF
515 :
516 : CALL qs_scf_check_inner_exit(qs_env, scf_env, scf_control, should_stop, exit_inner_loop, &
517 150467 : inner_loop_converged, output_unit)
518 :
519 : ! In case we decide to exit we perform few more check to see if this one
520 : ! is really the last SCF step
521 150467 : IF (exit_inner_loop) THEN
522 :
523 18706 : CALL qs_scf_inner_finalize(scf_env, qs_env, diis_step, output_unit)
524 :
525 : CALL qs_scf_check_outer_exit(qs_env, scf_env, scf_control, should_stop, &
526 18706 : outer_loop_converged, exit_outer_loop)
527 :
528 : ! Let's tag the last SCF cycle so we can print informations only of the last step
529 18706 : IF (exit_outer_loop) CALL cp_iterate(logger%iter_info, last=.TRUE., iter_nr=iter_count)
530 :
531 : END IF
532 :
533 150467 : IF (do_kpoints) THEN
534 5320 : CALL write_kpoints_restart(rho_ao_kp, kpoints, scf_env, dft_section, particle_set, qs_kind_set)
535 : ELSE
536 : ! Write Wavefunction restart file
537 145147 : CALL write_mo_set_to_restart(mos, particle_set, dft_section=dft_section, qs_kind_set=qs_kind_set)
538 : END IF
539 :
540 : ! Exit if we have finished with the SCF inner loop
541 150467 : IF (exit_inner_loop) THEN
542 18706 : CALL timestop(handle2)
543 : EXIT scf_loop
544 : END IF
545 :
546 131761 : IF (.NOT. BTEST(cp_print_key_should_output(logger%iter_info, &
547 : scf_section, "PRINT%ITERATION_INFO/TIME_CUMUL"), cp_p_file)) &
548 131761 : t1 = m_walltime()
549 :
550 : ! mixing methods have the new density matrix in p_mix_new
551 131761 : IF (scf_env%mixing_method > 0) THEN
552 490656 : DO ic = 1, SIZE(rho_ao_kp, 2)
553 956283 : DO ispin = 1, dft_control%nspins
554 465627 : CALL dbcsr_get_info(rho_ao_kp(ispin, ic)%matrix, name=name) ! keep the name
555 881798 : CALL dbcsr_copy(rho_ao_kp(ispin, ic)%matrix, scf_env%p_mix_new(ispin, ic)%matrix, name=name)
556 : END DO
557 : END DO
558 : END IF
559 :
560 : CALL qs_scf_rho_update(rho, qs_env, scf_env, ks_env, &
561 131761 : mix_rho=scf_env%mixing_method >= gspace_mixing_nr)
562 :
563 131761 : CALL timestop(handle2)
564 :
565 : END DO scf_loop
566 :
567 18706 : IF (.NOT. scf_control%outer_scf%have_scf) EXIT scf_outer_loop
568 :
569 : ! In case we use the OUTER SCF loop let's print some info..
570 : CALL qs_scf_outer_loop_info(output_unit, scf_control, scf_env, &
571 5110 : energy, total_steps, should_stop, outer_loop_converged)
572 :
573 : ! save mos to converged mos if outer_loop_converged and surf_dip_correct_switch is true
574 5110 : IF (exit_outer_loop) THEN
575 4199 : IF ((dft_control%switch_surf_dip) .AND. (outer_loop_converged) .AND. &
576 : (dft_control%surf_dip_correct_switch)) THEN
577 4 : DO ispin = 1, dft_control%nspins
578 4 : CALL reassign_allocated_mos(mos_last_converged(ispin), mos(ispin))
579 : END DO
580 2 : IF (output_unit > 0) WRITE (UNIT=output_unit, FMT="(/,/,T2,A)") &
581 1 : "COPIED mos ---> mos_last_converged"
582 : END IF
583 : END IF
584 :
585 5110 : IF (exit_outer_loop) EXIT scf_outer_loop
586 :
587 : !
588 911 : CALL outer_loop_optimize(scf_env, scf_control)
589 911 : CALL outer_loop_update_qs_env(qs_env, scf_env)
590 18706 : CALL qs_ks_did_change(ks_env, potential_changed=.TRUE.)
591 :
592 : END DO scf_outer_loop
593 :
594 17795 : converged = inner_loop_converged .AND. outer_loop_converged
595 17795 : total_scf_steps = total_steps
596 :
597 17795 : IF (dft_control%qs_control%cdft) &
598 : dft_control%qs_control%cdft_control%total_steps = &
599 626 : dft_control%qs_control%cdft_control%total_steps + total_steps
600 :
601 17795 : IF (.NOT. converged) THEN
602 2098 : IF (scf_control%ignore_convergence_failure .OR. should_stop) THEN
603 2098 : CALL cp_warn(__LOCATION__, "SCF run NOT converged")
604 : ELSE
605 : CALL cp_abort(__LOCATION__, &
606 : "SCF run NOT converged. To continue the calculation "// &
607 0 : "regardless, please set the keyword IGNORE_CONVERGENCE_FAILURE.")
608 : END IF
609 : END IF
610 :
611 : ! Skip Harris functional calculation if ground-state is NOT converged
612 17795 : IF (qs_env%energy_correction) THEN
613 502 : CALL get_qs_env(qs_env, ec_env=ec_env)
614 502 : ec_env%do_skip = .FALSE.
615 502 : IF (ec_env%skip_ec .AND. .NOT. converged) ec_env%do_skip = .TRUE.
616 : END IF
617 :
618 : ! if needed copy mo_coeff dbcsr->fm for later use in post_scf!fm->dbcsr
619 38100 : DO ispin = 1, SIZE(mos) !fm -> dbcsr
620 38100 : IF (mos(ispin)%use_mo_coeff_b) THEN !fm->dbcsr
621 6919 : IF (.NOT. ASSOCIATED(mos(ispin)%mo_coeff_b)) & !fm->dbcsr
622 0 : CPABORT("mo_coeff_b is not allocated") !fm->dbcsr
623 : CALL copy_dbcsr_to_fm(mos(ispin)%mo_coeff_b, & !fm->dbcsr
624 6919 : mos(ispin)%mo_coeff) !fm -> dbcsr
625 : END IF !fm->dbcsr
626 : END DO !fm -> dbcsr
627 :
628 17795 : CALL cp_rm_iter_level(logger%iter_info, level_name="QS_SCF")
629 17795 : CALL timestop(handle)
630 :
631 17795 : END SUBROUTINE scf_env_do_scf
632 :
633 : ! **************************************************************************************************
634 : !> \brief inits those objects needed if you want to restart the scf with, say
635 : !> only a new initial guess, or different density functional or ...
636 : !> this will happen just before the scf loop starts
637 : !> \param scf_env ...
638 : !> \param qs_env ...
639 : !> \param scf_section ...
640 : !> \par History
641 : !> 03.2006 created [Joost VandeVondele]
642 : ! **************************************************************************************************
643 18706 : SUBROUTINE init_scf_loop(scf_env, qs_env, scf_section)
644 :
645 : TYPE(qs_scf_env_type), POINTER :: scf_env
646 : TYPE(qs_environment_type), POINTER :: qs_env
647 : TYPE(section_vals_type), POINTER :: scf_section
648 :
649 : CHARACTER(LEN=*), PARAMETER :: routineN = 'init_scf_loop'
650 :
651 : INTEGER :: handle, ispin, nmo, number_of_OT_envs
652 : LOGICAL :: do_kpoints, do_rotation, &
653 : has_unit_metric, is_full_all
654 : TYPE(cp_fm_type), POINTER :: mo_coeff
655 18706 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_ks, matrix_s
656 : TYPE(dbcsr_type), POINTER :: orthogonality_metric
657 : TYPE(dft_control_type), POINTER :: dft_control
658 : TYPE(kpoint_type), POINTER :: kpoints
659 18706 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
660 : TYPE(scf_control_type), POINTER :: scf_control
661 :
662 18706 : CALL timeset(routineN, handle)
663 :
664 18706 : NULLIFY (scf_control, matrix_s, matrix_ks, dft_control, mos, mo_coeff, kpoints)
665 :
666 18706 : CPASSERT(ASSOCIATED(scf_env))
667 18706 : CPASSERT(ASSOCIATED(qs_env))
668 :
669 : CALL get_qs_env(qs_env=qs_env, &
670 : scf_control=scf_control, &
671 : dft_control=dft_control, &
672 : do_kpoints=do_kpoints, &
673 : kpoints=kpoints, &
674 18706 : mos=mos)
675 :
676 : ! if using mo_coeff_b then copy to fm
677 40129 : DO ispin = 1, SIZE(mos) !fm->dbcsr
678 40129 : IF (mos(1)%use_mo_coeff_b) THEN !fm->dbcsr
679 8037 : CALL copy_dbcsr_to_fm(mos(ispin)%mo_coeff_b, mos(ispin)%mo_coeff) !fm->dbcsr
680 : END IF !fm->dbcsr
681 : END DO !fm->dbcsr
682 :
683 : ! this just guarantees that all mo_occupations match the eigenvalues, if smear
684 40129 : DO ispin = 1, dft_control%nspins
685 : ! do not reset mo_occupations if the maximum overlap method is in use
686 21423 : IF (.NOT. scf_control%diagonalization%mom) &
687 : CALL set_mo_occupation(mo_set=mos(ispin), &
688 40085 : smear=scf_control%smear)
689 : END DO
690 :
691 18706 : SELECT CASE (scf_env%method)
692 : CASE DEFAULT
693 :
694 0 : CPABORT("unknown scf method method:"//cp_to_string(scf_env%method))
695 :
696 : CASE (filter_matrix_diag_method_nr)
697 :
698 10 : IF (.NOT. scf_env%skip_diis) THEN
699 0 : IF (.NOT. ASSOCIATED(scf_env%scf_diis_buffer)) THEN
700 0 : ALLOCATE (scf_env%scf_diis_buffer)
701 0 : CALL qs_diis_b_create(scf_env%scf_diis_buffer, nbuffer=scf_control%max_diis)
702 : END IF
703 0 : CALL qs_diis_b_clear(scf_env%scf_diis_buffer)
704 : END IF
705 :
706 : CASE (general_diag_method_nr, special_diag_method_nr, block_krylov_diag_method_nr)
707 12026 : IF (.NOT. scf_env%skip_diis) THEN
708 11792 : IF (do_kpoints) THEN
709 714 : IF (.NOT. ASSOCIATED(kpoints%scf_diis_buffer)) THEN
710 120 : ALLOCATE (kpoints%scf_diis_buffer)
711 120 : CALL qs_diis_b_create_kp(kpoints%scf_diis_buffer, nbuffer=scf_control%max_diis)
712 : END IF
713 714 : CALL qs_diis_b_clear_kp(kpoints%scf_diis_buffer)
714 : ELSE
715 11078 : IF (.NOT. ASSOCIATED(scf_env%scf_diis_buffer)) THEN
716 3114 : ALLOCATE (scf_env%scf_diis_buffer)
717 3114 : CALL qs_diis_b_create(scf_env%scf_diis_buffer, nbuffer=scf_control%max_diis)
718 : END IF
719 11078 : CALL qs_diis_b_clear(scf_env%scf_diis_buffer)
720 : END IF
721 : END IF
722 :
723 : CASE (ot_diag_method_nr)
724 8 : CALL get_qs_env(qs_env, matrix_ks=matrix_ks, matrix_s=matrix_s)
725 :
726 8 : IF (.NOT. scf_env%skip_diis) THEN
727 6 : IF (.NOT. ASSOCIATED(scf_env%scf_diis_buffer)) THEN
728 6 : ALLOCATE (scf_env%scf_diis_buffer)
729 6 : CALL qs_diis_b_create(scf_env%scf_diis_buffer, nbuffer=scf_control%max_diis)
730 : END IF
731 6 : CALL qs_diis_b_clear(scf_env%scf_diis_buffer)
732 : END IF
733 :
734 : ! disable DFTB and SE for now
735 : IF (dft_control%qs_control%dftb .OR. &
736 8 : dft_control%qs_control%xtb .OR. &
737 : dft_control%qs_control%semi_empirical) THEN
738 0 : CPABORT("DFTB and SE not available with OT/DIAG")
739 : END IF
740 :
741 : ! if an old preconditioner is still around (i.e. outer SCF is active),
742 : ! remove it if this could be worthwhile
743 : CALL restart_preconditioner(qs_env, scf_env%ot_preconditioner, &
744 : scf_control%diagonalization%ot_settings%preconditioner_type, &
745 8 : dft_control%nspins)
746 :
747 : CALL prepare_preconditioner(qs_env, mos, matrix_ks, matrix_s, scf_env%ot_preconditioner, &
748 : scf_control%diagonalization%ot_settings%preconditioner_type, &
749 : scf_control%diagonalization%ot_settings%precond_solver_type, &
750 8 : scf_control%diagonalization%ot_settings%energy_gap, dft_control%nspins)
751 :
752 : CASE (block_davidson_diag_method_nr)
753 : ! Preconditioner initialized within the loop, when required
754 : CASE (ot_method_nr)
755 : CALL get_qs_env(qs_env, &
756 : has_unit_metric=has_unit_metric, &
757 : matrix_s=matrix_s, &
758 6646 : matrix_ks=matrix_ks)
759 :
760 : ! reortho the wavefunctions if we are having an outer scf and
761 : ! this is not the first iteration
762 : ! this is useful to avoid the build-up of numerical noise
763 : ! however, we can not play this trick if restricted (don't mix non-equivalent orbs)
764 6646 : IF (scf_control%do_outer_scf_reortho) THEN
765 6076 : IF (scf_control%outer_scf%have_scf .AND. .NOT. dft_control%restricted) THEN
766 4368 : IF (scf_env%outer_scf%iter_count > 0) THEN
767 1959 : DO ispin = 1, dft_control%nspins
768 1072 : CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
769 1959 : IF (has_unit_metric) THEN
770 154 : CALL make_basis_simple(mo_coeff, nmo)
771 : ELSE
772 918 : CALL make_basis_sm(mo_coeff, nmo, matrix_s(1)%matrix)
773 : END IF
774 : END DO
775 : END IF
776 : END IF
777 : ELSE
778 : ! dont need any dirty trick for the numerically stable irac algorithm.
779 : END IF
780 :
781 6646 : IF (.NOT. ASSOCIATED(scf_env%qs_ot_env)) THEN
782 :
783 : ! restricted calculations require just one set of OT orbitals
784 6646 : number_of_OT_envs = dft_control%nspins
785 6646 : IF (dft_control%restricted) number_of_OT_envs = 1
786 :
787 19938 : ALLOCATE (scf_env%qs_ot_env(number_of_OT_envs))
788 :
789 : ! XXX Joost XXX should disentangle reading input from this part
790 6646 : CALL ot_scf_read_input(scf_env%qs_ot_env, scf_section)
791 :
792 : ! keep a note that we are restricted
793 6646 : IF (dft_control%restricted) THEN
794 100 : scf_env%qs_ot_env(1)%restricted = .TRUE.
795 : ! requires rotation
796 100 : IF (.NOT. scf_env%qs_ot_env(1)%settings%do_rotation) &
797 : CALL cp_abort(__LOCATION__, &
798 : "Restricted calculation with OT requires orbital rotation. Please "// &
799 0 : "activate the OT%ROTATION keyword!")
800 : ELSE
801 14353 : scf_env%qs_ot_env(:)%restricted = .FALSE.
802 : END IF
803 :
804 : ! this will rotate the MOs to be eigen states, which is not compatible with rotation
805 : ! e.g. mo_derivs here do not yet include potentially different occupations numbers
806 6646 : do_rotation = scf_env%qs_ot_env(1)%settings%do_rotation
807 : ! only full all needs rotation
808 6646 : is_full_all = scf_env%qs_ot_env(1)%settings%preconditioner_type == ot_precond_full_all
809 6646 : IF (do_rotation .AND. is_full_all) &
810 0 : CPABORT('PRECONDITIONER FULL_ALL is not compatible with ROTATION.')
811 :
812 : ! might need the KS matrix to init properly
813 : CALL qs_ks_update_qs_env(qs_env, just_energy=.FALSE., &
814 6646 : calculate_forces=.FALSE.)
815 :
816 : ! if an old preconditioner is still around (i.e. outer SCF is active),
817 : ! remove it if this could be worthwhile
818 6646 : IF (.NOT. reuse_precond) &
819 : CALL restart_preconditioner(qs_env, scf_env%ot_preconditioner, &
820 : scf_env%qs_ot_env(1)%settings%preconditioner_type, &
821 6646 : dft_control%nspins)
822 :
823 : !
824 : ! preconditioning still needs to be done correctly with has_unit_metric
825 : ! notice that a big part of the preconditioning (S^-1) is fine anyhow
826 : !
827 6646 : IF (has_unit_metric) THEN
828 1250 : NULLIFY (orthogonality_metric)
829 : ELSE
830 5396 : orthogonality_metric => matrix_s(1)%matrix
831 : END IF
832 :
833 6646 : IF (.NOT. reuse_precond) &
834 : CALL prepare_preconditioner(qs_env, mos, matrix_ks, matrix_s, scf_env%ot_preconditioner, &
835 : scf_env%qs_ot_env(1)%settings%preconditioner_type, &
836 : scf_env%qs_ot_env(1)%settings%precond_solver_type, &
837 : scf_env%qs_ot_env(1)%settings%energy_gap, dft_control%nspins, &
838 : has_unit_metric=has_unit_metric, &
839 6646 : chol_type=scf_env%qs_ot_env(1)%settings%cholesky_type)
840 6646 : IF (reuse_precond) reuse_precond = .FALSE.
841 :
842 : CALL ot_scf_init(mo_array=mos, matrix_s=orthogonality_metric, &
843 : broyden_adaptive_sigma=qs_env%broyden_adaptive_sigma, &
844 6646 : qs_ot_env=scf_env%qs_ot_env, matrix_ks=matrix_ks(1)%matrix)
845 :
846 11289 : SELECT CASE (scf_env%qs_ot_env(1)%settings%preconditioner_type)
847 : CASE (ot_precond_none)
848 : CASE (ot_precond_full_all, ot_precond_full_single_inverse)
849 10243 : DO ispin = 1, SIZE(scf_env%qs_ot_env)
850 : CALL qs_ot_new_preconditioner(scf_env%qs_ot_env(ispin), &
851 10243 : scf_env%ot_preconditioner(ispin)%preconditioner)
852 : END DO
853 : CASE (ot_precond_s_inverse, ot_precond_full_single)
854 182 : DO ispin = 1, SIZE(scf_env%qs_ot_env)
855 : CALL qs_ot_new_preconditioner(scf_env%qs_ot_env(ispin), &
856 182 : scf_env%ot_preconditioner(1)%preconditioner)
857 : END DO
858 : CASE DEFAULT
859 8029 : DO ispin = 1, SIZE(scf_env%qs_ot_env)
860 : CALL qs_ot_new_preconditioner(scf_env%qs_ot_env(ispin), &
861 2522 : scf_env%ot_preconditioner(1)%preconditioner)
862 : END DO
863 : END SELECT
864 : END IF
865 :
866 : ! if we have non-uniform occupations we should be using rotation
867 6646 : do_rotation = scf_env%qs_ot_env(1)%settings%do_rotation
868 33359 : DO ispin = 1, SIZE(mos)
869 14653 : IF (.NOT. mos(ispin)%uniform_occupation) THEN
870 0 : CPASSERT(do_rotation)
871 : END IF
872 : END DO
873 : END SELECT
874 :
875 : ! another safety check
876 18706 : IF (dft_control%low_spin_roks) THEN
877 32 : CPASSERT(scf_env%method == ot_method_nr)
878 32 : do_rotation = scf_env%qs_ot_env(1)%settings%do_rotation
879 32 : CPASSERT(do_rotation)
880 : END IF
881 :
882 18706 : CALL timestop(handle)
883 :
884 18706 : END SUBROUTINE init_scf_loop
885 :
886 : ! **************************************************************************************************
887 : !> \brief perform cleanup operations (like releasing temporary storage)
888 : !> at the end of the scf
889 : !> \param scf_env ...
890 : !> \par History
891 : !> 02.2003 created [fawzi]
892 : !> \author fawzi
893 : ! **************************************************************************************************
894 17539 : SUBROUTINE scf_env_cleanup(scf_env)
895 : TYPE(qs_scf_env_type), INTENT(INOUT) :: scf_env
896 :
897 : CHARACTER(len=*), PARAMETER :: routineN = 'scf_env_cleanup'
898 :
899 : INTEGER :: handle
900 :
901 17539 : CALL timeset(routineN, handle)
902 :
903 : ! Release SCF work storage
904 17539 : CALL cp_fm_release(scf_env%scf_work1)
905 :
906 17539 : IF (ASSOCIATED(scf_env%scf_work2)) THEN
907 11984 : CALL cp_fm_release(scf_env%scf_work2)
908 11984 : DEALLOCATE (scf_env%scf_work2)
909 : END IF
910 17539 : IF (ASSOCIATED(scf_env%ortho)) THEN
911 9316 : CALL cp_fm_release(scf_env%ortho)
912 9316 : DEALLOCATE (scf_env%ortho)
913 : END IF
914 17539 : IF (ASSOCIATED(scf_env%ortho_m1)) THEN
915 46 : CALL cp_fm_release(scf_env%ortho_m1)
916 46 : DEALLOCATE (scf_env%ortho_m1)
917 : END IF
918 :
919 17539 : IF (ASSOCIATED(scf_env%ortho_dbcsr)) THEN
920 58 : CALL dbcsr_deallocate_matrix(scf_env%ortho_dbcsr)
921 : END IF
922 17539 : IF (ASSOCIATED(scf_env%buf1_dbcsr)) THEN
923 58 : CALL dbcsr_deallocate_matrix(scf_env%buf1_dbcsr)
924 : END IF
925 17539 : IF (ASSOCIATED(scf_env%buf2_dbcsr)) THEN
926 58 : CALL dbcsr_deallocate_matrix(scf_env%buf2_dbcsr)
927 : END IF
928 :
929 17539 : IF (ASSOCIATED(scf_env%p_mix_new)) THEN
930 11996 : CALL dbcsr_deallocate_matrix_set(scf_env%p_mix_new)
931 : END IF
932 :
933 17539 : IF (ASSOCIATED(scf_env%p_delta)) THEN
934 228 : CALL dbcsr_deallocate_matrix_set(scf_env%p_delta)
935 : END IF
936 :
937 : ! Method dependent cleanup
938 17555 : SELECT CASE (scf_env%method)
939 : CASE (ot_method_nr)
940 : !
941 : CASE (ot_diag_method_nr)
942 : !
943 : CASE (general_diag_method_nr)
944 : !
945 : CASE (special_diag_method_nr)
946 : !
947 : CASE (block_krylov_diag_method_nr)
948 : CASE (block_davidson_diag_method_nr)
949 16 : CALL block_davidson_deallocate(scf_env%block_davidson_env)
950 : CASE (filter_matrix_diag_method_nr)
951 : !
952 : CASE DEFAULT
953 17539 : CPABORT("unknown scf method method:"//cp_to_string(scf_env%method))
954 : END SELECT
955 :
956 17539 : IF (ASSOCIATED(scf_env%outer_scf%variables)) THEN
957 3903 : DEALLOCATE (scf_env%outer_scf%variables)
958 : END IF
959 17539 : IF (ASSOCIATED(scf_env%outer_scf%count)) THEN
960 3903 : DEALLOCATE (scf_env%outer_scf%count)
961 : END IF
962 17539 : IF (ASSOCIATED(scf_env%outer_scf%gradient)) THEN
963 3903 : DEALLOCATE (scf_env%outer_scf%gradient)
964 : END IF
965 17539 : IF (ASSOCIATED(scf_env%outer_scf%energy)) THEN
966 3903 : DEALLOCATE (scf_env%outer_scf%energy)
967 : END IF
968 17539 : IF (ASSOCIATED(scf_env%outer_scf%inv_jacobian) .AND. &
969 : scf_env%outer_scf%deallocate_jacobian) THEN
970 50 : DEALLOCATE (scf_env%outer_scf%inv_jacobian)
971 : END IF
972 :
973 17539 : CALL timestop(handle)
974 :
975 17539 : END SUBROUTINE scf_env_cleanup
976 :
977 : ! **************************************************************************************************
978 : !> \brief perform a CDFT scf procedure in the given qs_env
979 : !> \param qs_env the qs_environment where to perform the scf procedure
980 : !> \param should_stop flag determining if calculation should stop
981 : !> \par History
982 : !> 12.2015 Created
983 : !> \author Nico Holmberg
984 : ! **************************************************************************************************
985 652 : SUBROUTINE cdft_scf(qs_env, should_stop)
986 : TYPE(qs_environment_type), POINTER :: qs_env
987 : LOGICAL, INTENT(OUT) :: should_stop
988 :
989 : CHARACTER(len=*), PARAMETER :: routineN = 'cdft_scf'
990 :
991 : INTEGER :: handle, iatom, ispin, ivar, nmo, nvar, &
992 : output_unit, tsteps
993 : LOGICAL :: cdft_loop_converged, converged, &
994 : exit_cdft_loop, first_iteration, &
995 : my_uocc, uniform_occupation
996 326 : REAL(KIND=dp), DIMENSION(:), POINTER :: mo_occupations
997 : TYPE(cdft_control_type), POINTER :: cdft_control
998 : TYPE(cp_logger_type), POINTER :: logger
999 326 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s, rho_ao
1000 : TYPE(dft_control_type), POINTER :: dft_control
1001 326 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
1002 : TYPE(pw_env_type), POINTER :: pw_env
1003 : TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
1004 : TYPE(qs_energy_type), POINTER :: energy
1005 : TYPE(qs_ks_env_type), POINTER :: ks_env
1006 : TYPE(qs_rho_type), POINTER :: rho
1007 : TYPE(qs_scf_env_type), POINTER :: scf_env
1008 : TYPE(scf_control_type), POINTER :: scf_control
1009 : TYPE(section_vals_type), POINTER :: dft_section, input, scf_section
1010 :
1011 326 : NULLIFY (scf_env, ks_env, energy, rho, matrix_s, rho_ao, cdft_control, logger, &
1012 326 : dft_control, pw_env, auxbas_pw_pool, energy, ks_env, scf_env, dft_section, &
1013 326 : input, scf_section, scf_control, mos, mo_occupations)
1014 652 : logger => cp_get_default_logger()
1015 :
1016 326 : CPASSERT(ASSOCIATED(qs_env))
1017 : CALL get_qs_env(qs_env, scf_env=scf_env, energy=energy, &
1018 : dft_control=dft_control, scf_control=scf_control, &
1019 326 : ks_env=ks_env, input=input)
1020 :
1021 326 : CALL timeset(routineN//"_loop", handle)
1022 326 : dft_section => section_vals_get_subs_vals(input, "DFT")
1023 326 : scf_section => section_vals_get_subs_vals(dft_section, "SCF")
1024 : output_unit = cp_print_key_unit_nr(logger, scf_section, "PRINT%PROGRAM_RUN_INFO", &
1025 326 : extension=".scfLog")
1026 326 : first_iteration = .TRUE.
1027 :
1028 326 : cdft_control => dft_control%qs_control%cdft_control
1029 :
1030 326 : scf_env%outer_scf%iter_count = 0
1031 326 : cdft_control%total_steps = 0
1032 :
1033 : ! Write some info about the CDFT calculation
1034 326 : IF (output_unit > 0) THEN
1035 : WRITE (UNIT=output_unit, FMT="(/,/,T2,A)") &
1036 181 : "CDFT EXTERNAL SCF WAVEFUNCTION OPTIMIZATION"
1037 181 : CALL qs_scf_cdft_initial_info(output_unit, cdft_control)
1038 : END IF
1039 326 : IF (cdft_control%reuse_precond) THEN
1040 0 : reuse_precond = .FALSE.
1041 0 : cdft_control%nreused = 0
1042 : END IF
1043 512 : cdft_outer_loop: DO
1044 : ! Change outer_scf settings to OT settings
1045 512 : CALL outer_loop_switch(scf_env, scf_control, cdft_control, cdft2ot)
1046 : ! Solve electronic structure with fixed value of constraint
1047 : CALL scf_env_do_scf(scf_env=scf_env, scf_control=scf_control, qs_env=qs_env, &
1048 512 : converged=converged, should_stop=should_stop, total_scf_steps=tsteps)
1049 : ! Decide whether to reuse the preconditioner on the next iteration
1050 512 : IF (cdft_control%reuse_precond) THEN
1051 : ! For convergence in exactly one step, the preconditioner is always reused (assuming max_reuse > 0)
1052 : ! usually this means that the electronic structure has already converged to the correct state
1053 : ! but the constraint optimizer keeps jumping over the optimal solution
1054 : IF (scf_env%outer_scf%iter_count == 1 .AND. scf_env%iter_count == 1 &
1055 0 : .AND. cdft_control%total_steps /= 1) &
1056 0 : cdft_control%nreused = cdft_control%nreused - 1
1057 : ! SCF converged in less than precond_freq steps
1058 : IF (scf_env%outer_scf%iter_count == 1 .AND. scf_env%iter_count .LE. cdft_control%precond_freq .AND. &
1059 0 : cdft_control%total_steps /= 1 .AND. cdft_control%nreused .LT. cdft_control%max_reuse) THEN
1060 0 : reuse_precond = .TRUE.
1061 0 : cdft_control%nreused = cdft_control%nreused + 1
1062 : ELSE
1063 0 : reuse_precond = .FALSE.
1064 0 : cdft_control%nreused = 0
1065 : END IF
1066 : END IF
1067 : ! Update history purging counters
1068 512 : IF (first_iteration .AND. cdft_control%purge_history) THEN
1069 0 : cdft_control%istep = cdft_control%istep + 1
1070 0 : IF (scf_env%outer_scf%iter_count .GT. 1) THEN
1071 0 : cdft_control%nbad_conv = cdft_control%nbad_conv + 1
1072 0 : IF (cdft_control%nbad_conv .GE. cdft_control%purge_freq .AND. &
1073 : cdft_control%istep .GE. cdft_control%purge_offset) THEN
1074 0 : cdft_control%nbad_conv = 0
1075 0 : cdft_control%istep = 0
1076 0 : cdft_control%should_purge = .TRUE.
1077 : END IF
1078 : END IF
1079 : END IF
1080 512 : first_iteration = .FALSE.
1081 : ! Change outer_scf settings to CDFT settings
1082 512 : CALL outer_loop_switch(scf_env, scf_control, cdft_control, ot2cdft)
1083 : CALL qs_scf_check_outer_exit(qs_env, scf_env, scf_control, should_stop, &
1084 512 : cdft_loop_converged, exit_cdft_loop)
1085 : CALL qs_scf_cdft_info(output_unit, scf_control, scf_env, cdft_control, &
1086 : energy, cdft_control%total_steps, &
1087 512 : should_stop, cdft_loop_converged, cdft_loop=.TRUE.)
1088 512 : IF (exit_cdft_loop) EXIT cdft_outer_loop
1089 : ! Check if the inverse Jacobian needs to be calculated
1090 186 : CALL qs_calculate_inverse_jacobian(qs_env)
1091 : ! Check if a line search should be performed to find an optimal step size for the optimizer
1092 186 : CALL qs_cdft_line_search(qs_env)
1093 : ! Optimize constraint
1094 186 : CALL outer_loop_optimize(scf_env, scf_control)
1095 186 : CALL outer_loop_update_qs_env(qs_env, scf_env)
1096 512 : CALL qs_ks_did_change(ks_env, potential_changed=.TRUE.)
1097 : END DO cdft_outer_loop
1098 :
1099 326 : cdft_control%ienergy = cdft_control%ienergy + 1
1100 :
1101 : ! Store needed arrays for ET coupling calculation
1102 326 : IF (cdft_control%do_et) THEN
1103 176 : CALL get_qs_env(qs_env=qs_env, matrix_s=matrix_s, mos=mos)
1104 176 : nvar = SIZE(cdft_control%target)
1105 : ! Matrix representation of weight function
1106 708 : ALLOCATE (cdft_control%wmat(nvar))
1107 356 : DO ivar = 1, nvar
1108 180 : CALL dbcsr_init_p(cdft_control%wmat(ivar)%matrix)
1109 : CALL dbcsr_copy(cdft_control%wmat(ivar)%matrix, matrix_s(1)%matrix, &
1110 180 : name="ET_RESTRAINT_MATRIX")
1111 180 : CALL dbcsr_set(cdft_control%wmat(ivar)%matrix, 0.0_dp)
1112 : CALL integrate_v_rspace(cdft_control%group(ivar)%weight, &
1113 : hmat=cdft_control%wmat(ivar), qs_env=qs_env, &
1114 : calculate_forces=.FALSE., &
1115 356 : gapw=dft_control%qs_control%gapw)
1116 : END DO
1117 : ! Overlap matrix
1118 176 : CALL dbcsr_init_p(cdft_control%matrix_s%matrix)
1119 : CALL dbcsr_copy(cdft_control%matrix_s%matrix, matrix_s(1)%matrix, &
1120 176 : name="OVERLAP")
1121 : ! Molecular orbital coefficients
1122 176 : NULLIFY (cdft_control%mo_coeff)
1123 880 : ALLOCATE (cdft_control%mo_coeff(dft_control%nspins))
1124 528 : DO ispin = 1, dft_control%nspins
1125 : CALL cp_fm_create(matrix=cdft_control%mo_coeff(ispin), &
1126 : matrix_struct=qs_env%mos(ispin)%mo_coeff%matrix_struct, &
1127 352 : name="MO_COEFF_A"//TRIM(ADJUSTL(cp_to_string(ispin)))//"MATRIX")
1128 : CALL cp_fm_to_fm(qs_env%mos(ispin)%mo_coeff, &
1129 528 : cdft_control%mo_coeff(ispin))
1130 : END DO
1131 : ! Density matrix
1132 176 : IF (cdft_control%calculate_metric) THEN
1133 24 : CALL get_qs_env(qs_env, rho=rho)
1134 24 : CALL qs_rho_get(rho, rho_ao=rho_ao)
1135 120 : ALLOCATE (cdft_control%matrix_p(dft_control%nspins))
1136 72 : DO ispin = 1, dft_control%nspins
1137 48 : NULLIFY (cdft_control%matrix_p(ispin)%matrix)
1138 48 : CALL dbcsr_init_p(cdft_control%matrix_p(ispin)%matrix)
1139 : CALL dbcsr_copy(cdft_control%matrix_p(ispin)%matrix, rho_ao(ispin)%matrix, &
1140 72 : name="DENSITY MATRIX")
1141 : END DO
1142 : END IF
1143 : ! Copy occupation numbers if non-uniform occupation
1144 176 : uniform_occupation = .TRUE.
1145 528 : DO ispin = 1, dft_control%nspins
1146 352 : CALL get_mo_set(mo_set=mos(ispin), uniform_occupation=my_uocc)
1147 584 : uniform_occupation = uniform_occupation .AND. my_uocc
1148 : END DO
1149 176 : IF (.NOT. uniform_occupation) THEN
1150 140 : ALLOCATE (cdft_control%occupations(dft_control%nspins))
1151 84 : DO ispin = 1, dft_control%nspins
1152 : CALL get_mo_set(mo_set=mos(ispin), &
1153 : nmo=nmo, &
1154 56 : occupation_numbers=mo_occupations)
1155 168 : ALLOCATE (cdft_control%occupations(ispin)%array(nmo))
1156 588 : cdft_control%occupations(ispin)%array(1:nmo) = mo_occupations(1:nmo)
1157 : END DO
1158 : END IF
1159 : END IF
1160 :
1161 : ! Deallocate constraint storage if forces are not needed
1162 : ! In case of a simulation with multiple force_evals,
1163 : ! deallocate only if weight function should not be copied to different force_evals
1164 326 : IF (.NOT. (cdft_control%save_pot .OR. cdft_control%transfer_pot)) THEN
1165 148 : CALL get_qs_env(qs_env, pw_env=pw_env)
1166 148 : CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
1167 308 : DO iatom = 1, SIZE(cdft_control%group)
1168 160 : CALL auxbas_pw_pool%give_back_pw(cdft_control%group(iatom)%weight)
1169 308 : DEALLOCATE (cdft_control%group(iatom)%weight)
1170 : END DO
1171 148 : IF (cdft_control%atomic_charges) THEN
1172 256 : DO iatom = 1, cdft_control%natoms
1173 256 : CALL auxbas_pw_pool%give_back_pw(cdft_control%charge(iatom))
1174 : END DO
1175 84 : DEALLOCATE (cdft_control%charge)
1176 : END IF
1177 148 : IF (cdft_control%type == outer_scf_becke_constraint .AND. &
1178 : cdft_control%becke_control%cavity_confine) THEN
1179 120 : IF (.NOT. ASSOCIATED(cdft_control%becke_control%cavity_mat)) THEN
1180 110 : CALL auxbas_pw_pool%give_back_pw(cdft_control%becke_control%cavity)
1181 : ELSE
1182 10 : DEALLOCATE (cdft_control%becke_control%cavity_mat)
1183 : END IF
1184 28 : ELSE IF (cdft_control%type == outer_scf_hirshfeld_constraint) THEN
1185 20 : IF (ASSOCIATED(cdft_control%hirshfeld_control%hirshfeld_env%fnorm)) THEN
1186 0 : CALL auxbas_pw_pool%give_back_pw(cdft_control%hirshfeld_control%hirshfeld_env%fnorm)
1187 : END IF
1188 : END IF
1189 148 : IF (ASSOCIATED(cdft_control%charges_fragment)) DEALLOCATE (cdft_control%charges_fragment)
1190 148 : cdft_control%need_pot = .TRUE.
1191 148 : cdft_control%external_control = .FALSE.
1192 : END IF
1193 :
1194 326 : CALL timestop(handle)
1195 :
1196 326 : END SUBROUTINE cdft_scf
1197 :
1198 : ! **************************************************************************************************
1199 : !> \brief perform cleanup operations for cdft_control
1200 : !> \param cdft_control container for the external CDFT SCF loop variables
1201 : !> \par History
1202 : !> 12.2015 created [Nico Holmberg]
1203 : !> \author Nico Holmberg
1204 : ! **************************************************************************************************
1205 326 : SUBROUTINE cdft_control_cleanup(cdft_control)
1206 : TYPE(cdft_control_type), POINTER :: cdft_control
1207 :
1208 326 : IF (ASSOCIATED(cdft_control%constraint%variables)) &
1209 326 : DEALLOCATE (cdft_control%constraint%variables)
1210 326 : IF (ASSOCIATED(cdft_control%constraint%count)) &
1211 326 : DEALLOCATE (cdft_control%constraint%count)
1212 326 : IF (ASSOCIATED(cdft_control%constraint%gradient)) &
1213 326 : DEALLOCATE (cdft_control%constraint%gradient)
1214 326 : IF (ASSOCIATED(cdft_control%constraint%energy)) &
1215 326 : DEALLOCATE (cdft_control%constraint%energy)
1216 326 : IF (ASSOCIATED(cdft_control%constraint%inv_jacobian) .AND. &
1217 : cdft_control%constraint%deallocate_jacobian) &
1218 4 : DEALLOCATE (cdft_control%constraint%inv_jacobian)
1219 :
1220 326 : END SUBROUTINE cdft_control_cleanup
1221 :
1222 : ! **************************************************************************************************
1223 : !> \brief Calculates the finite difference inverse Jacobian
1224 : !> \param qs_env the qs_environment_type where to compute the Jacobian
1225 : !> \par History
1226 : !> 01.2017 created [Nico Holmberg]
1227 : ! **************************************************************************************************
1228 186 : SUBROUTINE qs_calculate_inverse_jacobian(qs_env)
1229 : TYPE(qs_environment_type), POINTER :: qs_env
1230 :
1231 : CHARACTER(LEN=*), PARAMETER :: routineN = 'qs_calculate_inverse_jacobian'
1232 :
1233 : CHARACTER(len=default_path_length) :: project_name
1234 : INTEGER :: counter, handle, i, ispin, iter_count, &
1235 : iwork, j, max_scf, nspins, nsteps, &
1236 : nvar, nwork, output_unit, pwork, &
1237 : tsteps, twork
1238 : LOGICAL :: converged, explicit_jacobian, &
1239 : should_build, should_stop, &
1240 : use_md_history
1241 : REAL(KIND=dp) :: inv_error, step_size
1242 186 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: coeff, dh, step_multiplier
1243 186 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: jacobian
1244 186 : REAL(KIND=dp), DIMENSION(:), POINTER :: energy
1245 186 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: gradient, inv_jacobian
1246 : TYPE(cdft_control_type), POINTER :: cdft_control
1247 : TYPE(cp_logger_type), POINTER :: logger, tmp_logger
1248 186 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: p_rmpv
1249 186 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: rho_ao_kp
1250 : TYPE(dft_control_type), POINTER :: dft_control
1251 186 : TYPE(mo_set_type), ALLOCATABLE, DIMENSION(:) :: mos_stashed
1252 186 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
1253 : TYPE(mp_para_env_type), POINTER :: para_env
1254 : TYPE(qs_energy_type), POINTER :: energy_qs
1255 : TYPE(qs_ks_env_type), POINTER :: ks_env
1256 : TYPE(qs_rho_type), POINTER :: rho
1257 : TYPE(qs_scf_env_type), POINTER :: scf_env
1258 : TYPE(scf_control_type), POINTER :: scf_control
1259 :
1260 186 : NULLIFY (energy, gradient, p_rmpv, rho_ao_kp, mos, rho, &
1261 186 : ks_env, scf_env, scf_control, dft_control, cdft_control, &
1262 186 : inv_jacobian, para_env, tmp_logger, energy_qs)
1263 372 : logger => cp_get_default_logger()
1264 :
1265 186 : CPASSERT(ASSOCIATED(qs_env))
1266 : CALL get_qs_env(qs_env, scf_env=scf_env, ks_env=ks_env, &
1267 : scf_control=scf_control, mos=mos, rho=rho, &
1268 : dft_control=dft_control, &
1269 186 : para_env=para_env, energy=energy_qs)
1270 186 : explicit_jacobian = .FALSE.
1271 186 : should_build = .FALSE.
1272 186 : use_md_history = .FALSE.
1273 186 : iter_count = scf_env%outer_scf%iter_count
1274 : ! Quick exit if optimizer does not require Jacobian
1275 186 : IF (.NOT. ASSOCIATED(scf_control%outer_scf%cdft_opt_control)) RETURN
1276 : ! Check if Jacobian should be calculated and initialize
1277 118 : CALL timeset(routineN, handle)
1278 118 : CALL initialize_inverse_jacobian(scf_control, scf_env, explicit_jacobian, should_build, used_history)
1279 118 : IF (scf_control%outer_scf%cdft_opt_control%jacobian_restart) THEN
1280 : ! Restart from previously calculated inverse Jacobian
1281 6 : should_build = .FALSE.
1282 6 : CALL restart_inverse_jacobian(qs_env)
1283 : END IF
1284 118 : IF (should_build) THEN
1285 78 : scf_env%outer_scf%deallocate_jacobian = .FALSE.
1286 : ! Actually need to (re)build the Jacobian
1287 78 : IF (explicit_jacobian) THEN
1288 : ! Build Jacobian with finite differences
1289 62 : cdft_control => dft_control%qs_control%cdft_control
1290 62 : IF (.NOT. ASSOCIATED(cdft_control)) &
1291 : CALL cp_abort(__LOCATION__, &
1292 : "Optimizers that need the explicit Jacobian can"// &
1293 0 : " only be used together with a valid CDFT constraint.")
1294 : ! Redirect output from Jacobian calculation to a new file by creating a temporary logger
1295 62 : project_name = logger%iter_info%project_name
1296 62 : CALL create_tmp_logger(para_env, project_name, "-JacobianInfo.out", output_unit, tmp_logger)
1297 : ! Save last converged state so we can roll back to it (mo_coeff and some outer_loop variables)
1298 62 : nspins = dft_control%nspins
1299 310 : ALLOCATE (mos_stashed(nspins))
1300 186 : DO ispin = 1, nspins
1301 186 : CALL duplicate_mo_set(mos_stashed(ispin), mos(ispin))
1302 : END DO
1303 62 : CALL qs_rho_get(rho, rho_ao_kp=rho_ao_kp)
1304 62 : p_rmpv => rho_ao_kp(:, 1)
1305 : ! Allocate work
1306 62 : nvar = SIZE(scf_env%outer_scf%variables, 1)
1307 62 : max_scf = scf_control%outer_scf%max_scf + 1
1308 248 : ALLOCATE (gradient(nvar, max_scf))
1309 1310 : gradient = scf_env%outer_scf%gradient
1310 186 : ALLOCATE (energy(max_scf))
1311 594 : energy = scf_env%outer_scf%energy
1312 248 : ALLOCATE (jacobian(nvar, nvar))
1313 282 : jacobian = 0.0_dp
1314 62 : nsteps = cdft_control%total_steps
1315 : ! Setup finite difference scheme
1316 62 : CALL prepare_jacobian_stencil(qs_env, output_unit, nwork, pwork, coeff, step_multiplier, dh)
1317 62 : twork = pwork - nwork
1318 148 : DO i = 1, nvar
1319 282 : jacobian(i, :) = coeff(0)*scf_env%outer_scf%gradient(i, iter_count)
1320 : END DO
1321 : ! Calculate the Jacobian by perturbing each Lagrangian and recalculating the energy self-consistently
1322 62 : CALL cp_add_default_logger(tmp_logger)
1323 148 : DO i = 1, nvar
1324 86 : IF (output_unit > 0) THEN
1325 43 : WRITE (output_unit, FMT="(A)") " "
1326 43 : WRITE (output_unit, FMT="(A)") " #####################################"
1327 : WRITE (output_unit, '(A,I3,A,I3,A)') &
1328 43 : " ### Constraint ", i, " of ", nvar, " ###"
1329 43 : WRITE (output_unit, FMT="(A)") " #####################################"
1330 : END IF
1331 86 : counter = 0
1332 332 : DO iwork = nwork, pwork
1333 184 : IF (iwork == 0) CYCLE
1334 98 : counter = counter + 1
1335 98 : IF (output_unit > 0) THEN
1336 49 : WRITE (output_unit, FMT="(A)") " #####################################"
1337 : WRITE (output_unit, '(A,I3,A,I3,A)') &
1338 49 : " ### Energy evaluation ", counter, " of ", twork, " ###"
1339 49 : WRITE (output_unit, FMT="(A)") " #####################################"
1340 : END IF
1341 98 : IF (SIZE(scf_control%outer_scf%cdft_opt_control%jacobian_step) == 1) THEN
1342 90 : step_size = scf_control%outer_scf%cdft_opt_control%jacobian_step(1)
1343 : ELSE
1344 8 : step_size = scf_control%outer_scf%cdft_opt_control%jacobian_step(i)
1345 : END IF
1346 244 : scf_env%outer_scf%variables(:, iter_count + 1) = scf_env%outer_scf%variables(:, iter_count)
1347 : scf_env%outer_scf%variables(i, iter_count + 1) = scf_env%outer_scf%variables(i, iter_count) + &
1348 98 : step_multiplier(iwork)*step_size
1349 98 : CALL outer_loop_update_qs_env(qs_env, scf_env)
1350 98 : CALL qs_ks_did_change(ks_env, potential_changed=.TRUE.)
1351 98 : CALL outer_loop_switch(scf_env, scf_control, cdft_control, cdft2ot)
1352 : CALL scf_env_do_scf(scf_env=scf_env, scf_control=scf_control, qs_env=qs_env, &
1353 98 : converged=converged, should_stop=should_stop, total_scf_steps=tsteps)
1354 98 : CALL outer_loop_switch(scf_env, scf_control, cdft_control, ot2cdft)
1355 : ! Update (iter_count + 1) element of gradient and print constraint info
1356 98 : scf_env%outer_scf%iter_count = scf_env%outer_scf%iter_count + 1
1357 98 : CALL outer_loop_gradient(qs_env, scf_env)
1358 : CALL qs_scf_cdft_info(output_unit, scf_control, scf_env, cdft_control, &
1359 : energy_qs, cdft_control%total_steps, &
1360 98 : should_stop=.FALSE., outer_loop_converged=.FALSE., cdft_loop=.FALSE.)
1361 98 : scf_env%outer_scf%iter_count = scf_env%outer_scf%iter_count - 1
1362 : ! Update Jacobian
1363 244 : DO j = 1, nvar
1364 244 : jacobian(j, i) = jacobian(j, i) + coeff(iwork)*scf_env%outer_scf%gradient(j, iter_count + 1)
1365 : END DO
1366 : ! Reset everything to last converged state
1367 244 : scf_env%outer_scf%variables(:, iter_count + 1) = 0.0_dp
1368 2026 : scf_env%outer_scf%gradient = gradient
1369 878 : scf_env%outer_scf%energy = energy
1370 98 : cdft_control%total_steps = nsteps
1371 294 : DO ispin = 1, nspins
1372 196 : CALL deallocate_mo_set(mos(ispin))
1373 196 : CALL duplicate_mo_set(mos(ispin), mos_stashed(ispin))
1374 : CALL calculate_density_matrix(mos(ispin), &
1375 294 : p_rmpv(ispin)%matrix)
1376 : END DO
1377 98 : CALL qs_rho_update_rho(rho, qs_env=qs_env)
1378 368 : CALL qs_ks_did_change(qs_env%ks_env, rho_changed=.TRUE.)
1379 : END DO
1380 : END DO
1381 62 : CALL cp_rm_default_logger()
1382 62 : CALL cp_logger_release(tmp_logger)
1383 : ! Finalize and invert Jacobian
1384 148 : DO j = 1, nvar
1385 282 : DO i = 1, nvar
1386 220 : jacobian(i, j) = jacobian(i, j)/dh(j)
1387 : END DO
1388 : END DO
1389 62 : IF (.NOT. ASSOCIATED(scf_env%outer_scf%inv_jacobian)) &
1390 136 : ALLOCATE (scf_env%outer_scf%inv_jacobian(nvar, nvar))
1391 62 : inv_jacobian => scf_env%outer_scf%inv_jacobian
1392 62 : CALL invert_matrix(jacobian, inv_jacobian, inv_error)
1393 62 : scf_control%outer_scf%cdft_opt_control%broyden_update = .FALSE.
1394 : ! Release temporary storage
1395 186 : DO ispin = 1, nspins
1396 186 : CALL deallocate_mo_set(mos_stashed(ispin))
1397 : END DO
1398 62 : DEALLOCATE (mos_stashed, jacobian, gradient, energy, coeff, step_multiplier, dh)
1399 186 : IF (output_unit > 0) THEN
1400 : WRITE (output_unit, FMT="(/,A)") &
1401 31 : " ================================== JACOBIAN CALCULATED =================================="
1402 31 : CALL close_file(unit_number=output_unit)
1403 : END IF
1404 : ELSE
1405 : ! Build a strictly diagonal Jacobian from history and invert it
1406 16 : CALL build_diagonal_jacobian(qs_env, used_history)
1407 : END IF
1408 : END IF
1409 118 : IF (ASSOCIATED(scf_env%outer_scf%inv_jacobian) .AND. para_env%is_source()) THEN
1410 : ! Write restart file for inverse Jacobian
1411 55 : CALL print_inverse_jacobian(logger, scf_env%outer_scf%inv_jacobian, iter_count)
1412 : END IF
1413 : ! Update counter
1414 118 : scf_control%outer_scf%cdft_opt_control%ijacobian(1) = scf_control%outer_scf%cdft_opt_control%ijacobian(1) + 1
1415 118 : CALL timestop(handle)
1416 :
1417 372 : END SUBROUTINE qs_calculate_inverse_jacobian
1418 :
1419 : ! **************************************************************************************************
1420 : !> \brief Perform backtracking line search to find the optimal step size for the CDFT constraint
1421 : !> optimizer. Assumes that the CDFT gradient function is a smooth function of the constraint
1422 : !> variables.
1423 : !> \param qs_env the qs_environment_type where to perform the line search
1424 : !> \par History
1425 : !> 02.2017 created [Nico Holmberg]
1426 : ! **************************************************************************************************
1427 186 : SUBROUTINE qs_cdft_line_search(qs_env)
1428 : TYPE(qs_environment_type), POINTER :: qs_env
1429 :
1430 : CHARACTER(LEN=*), PARAMETER :: routineN = 'qs_cdft_line_search'
1431 :
1432 : CHARACTER(len=default_path_length) :: project_name
1433 : INTEGER :: handle, i, ispin, iter_count, &
1434 : max_linesearch, max_scf, nspins, &
1435 : nsteps, nvar, output_unit, tsteps
1436 : LOGICAL :: continue_ls, continue_ls_exit, converged, do_linesearch, found_solution, &
1437 : reached_maxls, should_exit, should_stop, sign_changed
1438 186 : LOGICAL, ALLOCATABLE, DIMENSION(:) :: positive_sign
1439 : REAL(KIND=dp) :: alpha, alpha_ls, factor, norm_ls
1440 186 : REAL(KIND=dp), DIMENSION(:), POINTER :: energy
1441 186 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: gradient, inv_jacobian
1442 : REAL(KIND=dp), EXTERNAL :: dnrm2
1443 : TYPE(cdft_control_type), POINTER :: cdft_control
1444 : TYPE(cp_logger_type), POINTER :: logger, tmp_logger
1445 186 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: p_rmpv
1446 186 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: rho_ao_kp
1447 : TYPE(dft_control_type), POINTER :: dft_control
1448 186 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
1449 : TYPE(mp_para_env_type), POINTER :: para_env
1450 : TYPE(qs_energy_type), POINTER :: energy_qs
1451 : TYPE(qs_ks_env_type), POINTER :: ks_env
1452 : TYPE(qs_rho_type), POINTER :: rho
1453 : TYPE(qs_scf_env_type), POINTER :: scf_env
1454 : TYPE(scf_control_type), POINTER :: scf_control
1455 :
1456 186 : CALL timeset(routineN, handle)
1457 :
1458 186 : NULLIFY (energy, gradient, p_rmpv, rho_ao_kp, mos, rho, &
1459 186 : ks_env, scf_env, scf_control, dft_control, &
1460 186 : cdft_control, inv_jacobian, para_env, &
1461 186 : tmp_logger, energy_qs)
1462 186 : logger => cp_get_default_logger()
1463 :
1464 186 : CPASSERT(ASSOCIATED(qs_env))
1465 : CALL get_qs_env(qs_env, scf_env=scf_env, ks_env=ks_env, &
1466 : scf_control=scf_control, mos=mos, rho=rho, &
1467 : dft_control=dft_control, &
1468 186 : para_env=para_env, energy=energy_qs)
1469 186 : do_linesearch = .FALSE.
1470 186 : SELECT CASE (scf_control%outer_scf%optimizer)
1471 : CASE DEFAULT
1472 : do_linesearch = .FALSE.
1473 : CASE (outer_scf_optimizer_newton_ls)
1474 24 : do_linesearch = .TRUE.
1475 : CASE (outer_scf_optimizer_broyden)
1476 186 : SELECT CASE (scf_control%outer_scf%cdft_opt_control%broyden_type)
1477 : CASE (broyden_type_1, broyden_type_2, broyden_type_1_explicit, broyden_type_2_explicit)
1478 0 : do_linesearch = .FALSE.
1479 : CASE (broyden_type_1_ls, broyden_type_1_explicit_ls, broyden_type_2_ls, broyden_type_2_explicit_ls)
1480 0 : cdft_control => dft_control%qs_control%cdft_control
1481 0 : IF (.NOT. ASSOCIATED(cdft_control)) &
1482 : CALL cp_abort(__LOCATION__, &
1483 : "Optimizers that perform a line search can"// &
1484 0 : " only be used together with a valid CDFT constraint")
1485 0 : IF (ASSOCIATED(scf_env%outer_scf%inv_jacobian)) &
1486 : do_linesearch = .TRUE.
1487 : END SELECT
1488 : END SELECT
1489 : IF (do_linesearch) THEN
1490 8 : BLOCK
1491 8 : TYPE(mo_set_type), DIMENSION(:), ALLOCATABLE :: mos_ls, mos_stashed
1492 8 : cdft_control => dft_control%qs_control%cdft_control
1493 8 : IF (.NOT. ASSOCIATED(cdft_control)) &
1494 : CALL cp_abort(__LOCATION__, &
1495 : "Optimizers that perform a line search can"// &
1496 0 : " only be used together with a valid CDFT constraint")
1497 8 : CPASSERT(ASSOCIATED(scf_env%outer_scf%inv_jacobian))
1498 8 : CPASSERT(ASSOCIATED(scf_control%outer_scf%cdft_opt_control))
1499 8 : alpha = scf_control%outer_scf%cdft_opt_control%newton_step_save
1500 8 : iter_count = scf_env%outer_scf%iter_count
1501 : ! Redirect output from line search procedure to a new file by creating a temporary logger
1502 8 : project_name = logger%iter_info%project_name
1503 8 : CALL create_tmp_logger(para_env, project_name, "-LineSearch.out", output_unit, tmp_logger)
1504 : ! Save last converged state so we can roll back to it (mo_coeff and some outer_loop variables)
1505 8 : nspins = dft_control%nspins
1506 40 : ALLOCATE (mos_stashed(nspins))
1507 24 : DO ispin = 1, nspins
1508 24 : CALL duplicate_mo_set(mos_stashed(ispin), mos(ispin))
1509 : END DO
1510 8 : CALL qs_rho_get(rho, rho_ao_kp=rho_ao_kp)
1511 8 : p_rmpv => rho_ao_kp(:, 1)
1512 8 : nsteps = cdft_control%total_steps
1513 : ! Allocate work
1514 8 : nvar = SIZE(scf_env%outer_scf%variables, 1)
1515 8 : max_scf = scf_control%outer_scf%max_scf + 1
1516 8 : max_linesearch = scf_control%outer_scf%cdft_opt_control%max_ls
1517 8 : continue_ls = scf_control%outer_scf%cdft_opt_control%continue_ls
1518 8 : factor = scf_control%outer_scf%cdft_opt_control%factor_ls
1519 8 : continue_ls_exit = .FALSE.
1520 8 : found_solution = .FALSE.
1521 32 : ALLOCATE (gradient(nvar, max_scf))
1522 104 : gradient = scf_env%outer_scf%gradient
1523 24 : ALLOCATE (energy(max_scf))
1524 56 : energy = scf_env%outer_scf%energy
1525 8 : reached_maxls = .FALSE.
1526 : ! Broyden optimizers: perform update of inv_jacobian if necessary
1527 8 : IF (scf_control%outer_scf%cdft_opt_control%broyden_update) THEN
1528 0 : CALL outer_loop_optimize(scf_env, scf_control)
1529 : ! Reset the variables and prevent a reupdate of inv_jacobian
1530 0 : scf_env%outer_scf%variables(:, iter_count + 1) = 0
1531 0 : scf_control%outer_scf%cdft_opt_control%broyden_update = .FALSE.
1532 : END IF
1533 : ! Print some info
1534 8 : IF (output_unit > 0) THEN
1535 : WRITE (output_unit, FMT="(/,A)") &
1536 4 : " ================================== LINE SEARCH STARTED =================================="
1537 : WRITE (output_unit, FMT="(A,I5,A)") &
1538 4 : " Evaluating optimal step size for optimizer using a maximum of", max_linesearch, " steps"
1539 4 : IF (continue_ls) THEN
1540 : WRITE (output_unit, FMT="(A)") &
1541 2 : " Line search continues until best step size is found or max steps are reached"
1542 : END IF
1543 : WRITE (output_unit, '(/,A,F5.3)') &
1544 4 : " Initial step size: ", alpha
1545 : WRITE (output_unit, '(/,A,F5.3)') &
1546 4 : " Step size update factor: ", factor
1547 : WRITE (output_unit, '(/,A,I10,A,I10)') &
1548 4 : " Energy evaluation: ", cdft_control%ienergy, ", CDFT SCF iteration: ", iter_count
1549 : END IF
1550 : ! Perform backtracking line search
1551 8 : CALL cp_add_default_logger(tmp_logger)
1552 16 : DO i = 1, max_linesearch
1553 16 : IF (output_unit > 0) THEN
1554 8 : WRITE (output_unit, FMT="(A)") " "
1555 8 : WRITE (output_unit, FMT="(A)") " #####################################"
1556 : WRITE (output_unit, '(A,I10,A)') &
1557 8 : " ### Line search step: ", i, " ###"
1558 8 : WRITE (output_unit, FMT="(A)") " #####################################"
1559 : END IF
1560 16 : inv_jacobian => scf_env%outer_scf%inv_jacobian
1561 : ! Newton update of CDFT variables with a step size of alpha
1562 : scf_env%outer_scf%variables(:, iter_count + 1) = scf_env%outer_scf%variables(:, iter_count) - alpha* &
1563 144 : MATMUL(inv_jacobian, scf_env%outer_scf%gradient(:, iter_count))
1564 : ! With updated CDFT variables, perform SCF
1565 16 : CALL outer_loop_update_qs_env(qs_env, scf_env)
1566 16 : CALL qs_ks_did_change(ks_env, potential_changed=.TRUE.)
1567 16 : CALL outer_loop_switch(scf_env, scf_control, cdft_control, cdft2ot)
1568 : CALL scf_env_do_scf(scf_env=scf_env, scf_control=scf_control, qs_env=qs_env, &
1569 16 : converged=converged, should_stop=should_stop, total_scf_steps=tsteps)
1570 16 : CALL outer_loop_switch(scf_env, scf_control, cdft_control, ot2cdft)
1571 : ! Update (iter_count + 1) element of gradient and print constraint info
1572 16 : scf_env%outer_scf%iter_count = scf_env%outer_scf%iter_count + 1
1573 16 : CALL outer_loop_gradient(qs_env, scf_env)
1574 : CALL qs_scf_cdft_info(output_unit, scf_control, scf_env, cdft_control, &
1575 : energy_qs, cdft_control%total_steps, &
1576 16 : should_stop=.FALSE., outer_loop_converged=.FALSE., cdft_loop=.FALSE.)
1577 16 : scf_env%outer_scf%iter_count = scf_env%outer_scf%iter_count - 1
1578 : ! Store sign of initial gradient for each variable for continue_ls
1579 16 : IF (continue_ls .AND. .NOT. ALLOCATED(positive_sign)) THEN
1580 12 : ALLOCATE (positive_sign(nvar))
1581 8 : DO ispin = 1, nvar
1582 8 : positive_sign(ispin) = scf_env%outer_scf%gradient(ispin, iter_count + 1) .GE. 0.0_dp
1583 : END DO
1584 : END IF
1585 : ! Check if the L2 norm of the gradient decreased
1586 16 : inv_jacobian => scf_env%outer_scf%inv_jacobian
1587 16 : IF (dnrm2(nvar, scf_env%outer_scf%gradient(:, iter_count + 1), 1) .LT. &
1588 : dnrm2(nvar, scf_env%outer_scf%gradient(:, iter_count), 1)) THEN
1589 : ! Optimal step size found
1590 14 : IF (.NOT. continue_ls) THEN
1591 : should_exit = .TRUE.
1592 : ELSE
1593 : ! But line search continues for at least one more iteration in an attempt to find a better solution
1594 : ! if max number of steps is not exceeded
1595 10 : IF (found_solution) THEN
1596 : ! Check if the norm also decreased w.r.t. to previously found solution
1597 6 : IF (dnrm2(nvar, scf_env%outer_scf%gradient(:, iter_count + 1), 1) .GT. norm_ls) THEN
1598 : ! Norm increased => accept previous solution and exit
1599 : continue_ls_exit = .TRUE.
1600 : END IF
1601 : END IF
1602 : ! Store current state and the value of alpha
1603 10 : IF (.NOT. continue_ls_exit) THEN
1604 10 : should_exit = .FALSE.
1605 10 : alpha_ls = alpha
1606 10 : found_solution = .TRUE.
1607 10 : norm_ls = dnrm2(nvar, scf_env%outer_scf%gradient(:, iter_count + 1), 1)
1608 : ! Check if the sign of the gradient has changed for all variables (w.r.t initial gradient)
1609 : ! In this case we should exit because further line search steps will just increase the norm
1610 10 : sign_changed = .TRUE.
1611 20 : DO ispin = 1, nvar
1612 : sign_changed = sign_changed .AND. (positive_sign(ispin) .NEQV. &
1613 28 : scf_env%outer_scf%gradient(ispin, iter_count + 1) .GE. 0.0_dp)
1614 : END DO
1615 10 : IF (.NOT. ALLOCATED(mos_ls)) THEN
1616 20 : ALLOCATE (mos_ls(nspins))
1617 : ELSE
1618 18 : DO ispin = 1, nspins
1619 18 : CALL deallocate_mo_set(mos_ls(ispin))
1620 : END DO
1621 : END IF
1622 30 : DO ispin = 1, nspins
1623 30 : CALL duplicate_mo_set(mos_ls(ispin), mos(ispin))
1624 : END DO
1625 10 : alpha = alpha*factor
1626 : ! Exit on last iteration
1627 10 : IF (i == max_linesearch) continue_ls_exit = .TRUE.
1628 : ! Exit if constraint target is satisfied to requested tolerance
1629 30 : IF (SQRT(MAXVAL(scf_env%outer_scf%gradient(:, scf_env%outer_scf%iter_count + 1)**2)) .LT. &
1630 : scf_control%outer_scf%eps_scf) &
1631 2 : continue_ls_exit = .TRUE.
1632 : ! Exit if line search jumped over the optimal step length
1633 10 : IF (sign_changed) continue_ls_exit = .TRUE.
1634 : END IF
1635 : END IF
1636 : ELSE
1637 : ! Gradient increased => alpha is too large (if the gradient function is smooth)
1638 2 : should_exit = .FALSE.
1639 : ! Update alpha using Armijo's scheme
1640 2 : alpha = alpha*factor
1641 : END IF
1642 16 : IF (continue_ls_exit) THEN
1643 : ! Continuation of line search did not yield a better alpha, use previously located solution and exit
1644 4 : alpha = alpha_ls
1645 12 : DO ispin = 1, nspins
1646 8 : CALL deallocate_mo_set(mos(ispin))
1647 8 : CALL duplicate_mo_set(mos(ispin), mos_ls(ispin))
1648 : CALL calculate_density_matrix(mos(ispin), &
1649 8 : p_rmpv(ispin)%matrix)
1650 12 : CALL deallocate_mo_set(mos_ls(ispin))
1651 : END DO
1652 4 : CALL qs_rho_update_rho(rho, qs_env=qs_env)
1653 4 : CALL qs_ks_did_change(qs_env%ks_env, rho_changed=.TRUE.)
1654 4 : DEALLOCATE (mos_ls)
1655 : should_exit = .TRUE.
1656 : END IF
1657 : ! Reached max steps and SCF converged: continue with last iterated step size
1658 12 : IF (.NOT. should_exit .AND. &
1659 : (i == max_linesearch .AND. converged .AND. .NOT. found_solution)) THEN
1660 0 : should_exit = .TRUE.
1661 0 : reached_maxls = .TRUE.
1662 0 : alpha = alpha*(1.0_dp/factor)
1663 : END IF
1664 : ! Reset outer SCF environment to last converged state
1665 32 : scf_env%outer_scf%variables(:, iter_count + 1) = 0.0_dp
1666 208 : scf_env%outer_scf%gradient = gradient
1667 112 : scf_env%outer_scf%energy = energy
1668 : ! Exit line search if a suitable step size was found
1669 16 : IF (should_exit) EXIT
1670 : ! Reset the electronic structure
1671 8 : cdft_control%total_steps = nsteps
1672 24 : DO ispin = 1, nspins
1673 16 : CALL deallocate_mo_set(mos(ispin))
1674 16 : CALL duplicate_mo_set(mos(ispin), mos_stashed(ispin))
1675 : CALL calculate_density_matrix(mos(ispin), &
1676 24 : p_rmpv(ispin)%matrix)
1677 : END DO
1678 8 : CALL qs_rho_update_rho(rho, qs_env=qs_env)
1679 24 : CALL qs_ks_did_change(qs_env%ks_env, rho_changed=.TRUE.)
1680 : END DO
1681 8 : scf_control%outer_scf%cdft_opt_control%newton_step = alpha
1682 8 : IF (.NOT. should_exit) THEN
1683 : CALL cp_warn(__LOCATION__, &
1684 0 : "Line search did not converge. CDFT SCF proceeds with fixed step size.")
1685 0 : scf_control%outer_scf%cdft_opt_control%newton_step = scf_control%outer_scf%cdft_opt_control%newton_step_save
1686 : END IF
1687 8 : IF (reached_maxls) &
1688 : CALL cp_warn(__LOCATION__, &
1689 0 : "Line search did not converge. CDFT SCF proceeds with lasted iterated step size.")
1690 8 : CALL cp_rm_default_logger()
1691 8 : CALL cp_logger_release(tmp_logger)
1692 : ! Release temporary storage
1693 24 : DO ispin = 1, nspins
1694 24 : CALL deallocate_mo_set(mos_stashed(ispin))
1695 : END DO
1696 8 : DEALLOCATE (mos_stashed, gradient, energy)
1697 8 : IF (ALLOCATED(positive_sign)) DEALLOCATE (positive_sign)
1698 20 : IF (output_unit > 0) THEN
1699 : WRITE (output_unit, FMT="(/,A)") &
1700 4 : " ================================== LINE SEARCH COMPLETE =================================="
1701 4 : CALL close_file(unit_number=output_unit)
1702 : END IF
1703 : END BLOCK
1704 : END IF
1705 :
1706 186 : CALL timestop(handle)
1707 :
1708 186 : END SUBROUTINE qs_cdft_line_search
1709 :
1710 16 : END MODULE qs_scf
|