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 : MODULE qs_scf_output
9 : USE admm_types, ONLY: admm_type
10 : USE admm_utils, ONLY: admm_correct_for_eigenvalues,&
11 : admm_uncorrect_for_eigenvalues
12 : USE cp_blacs_env, ONLY: cp_blacs_env_type
13 : USE cp_control_types, ONLY: dft_control_type
14 : USE cp_dbcsr_api, ONLY: dbcsr_p_type,&
15 : dbcsr_type
16 : USE cp_dbcsr_output, ONLY: cp_dbcsr_write_sparse_matrix
17 : USE cp_fm_struct, ONLY: cp_fm_struct_create,&
18 : cp_fm_struct_release,&
19 : cp_fm_struct_type
20 : USE cp_fm_types, ONLY: cp_fm_init_random,&
21 : cp_fm_type
22 : USE cp_log_handling, ONLY: cp_get_default_logger,&
23 : cp_logger_type
24 : USE cp_output_handling, ONLY: cp_p_file,&
25 : cp_print_key_finished_output,&
26 : cp_print_key_should_output,&
27 : cp_print_key_unit_nr
28 : USE cp_units, ONLY: cp_unit_from_cp2k
29 : USE input_constants, ONLY: &
30 : becke_cutoff_element, becke_cutoff_global, cdft_alpha_constraint, cdft_beta_constraint, &
31 : cdft_charge_constraint, cdft_magnetization_constraint, ot_precond_full_all, &
32 : outer_scf_becke_constraint, outer_scf_hirshfeld_constraint, outer_scf_optimizer_bisect, &
33 : outer_scf_optimizer_broyden, outer_scf_optimizer_diis, outer_scf_optimizer_newton, &
34 : outer_scf_optimizer_newton_ls, outer_scf_optimizer_sd, outer_scf_optimizer_secant, &
35 : radius_covalent, radius_default, radius_single, radius_user, radius_vdw, &
36 : shape_function_density, shape_function_gaussian, smear_fermi_dirac, smear_gaussian, &
37 : smear_mp, smear_mv
38 : USE input_section_types, ONLY: section_get_ivals,&
39 : section_vals_get_subs_vals,&
40 : section_vals_type,&
41 : section_vals_val_get
42 : USE kahan_sum, ONLY: accurate_sum
43 : USE kinds, ONLY: default_string_length,&
44 : dp
45 : USE kpoint_types, ONLY: kpoint_type
46 : USE machine, ONLY: m_flush
47 : USE message_passing, ONLY: mp_para_env_type
48 : USE particle_types, ONLY: particle_type
49 : USE physcon, ONLY: evolt,&
50 : kcalmol
51 : USE preconditioner_types, ONLY: preconditioner_type
52 : USE ps_implicit_types, ONLY: MIXED_BC,&
53 : MIXED_PERIODIC_BC,&
54 : NEUMANN_BC,&
55 : PERIODIC_BC
56 : USE pw_env_types, ONLY: pw_env_type
57 : USE pw_poisson_types, ONLY: pw_poisson_implicit
58 : USE qmmm_image_charge, ONLY: print_image_coefficients
59 : USE qs_cdft_opt_types, ONLY: cdft_opt_type_write
60 : USE qs_cdft_types, ONLY: cdft_control_type
61 : USE qs_charges_types, ONLY: qs_charges_type
62 : USE qs_energy_types, ONLY: qs_energy_type
63 : USE qs_environment_types, ONLY: get_qs_env,&
64 : qs_environment_type
65 : USE qs_kind_types, ONLY: qs_kind_type
66 : USE qs_mo_io, ONLY: write_mo_set_to_output_unit
67 : USE qs_mo_methods, ONLY: calculate_magnitude,&
68 : calculate_orthonormality,&
69 : calculate_subspace_eigenvalues
70 : USE qs_mo_occupation, ONLY: set_mo_occupation
71 : USE qs_mo_types, ONLY: allocate_mo_set,&
72 : deallocate_mo_set,&
73 : get_mo_set,&
74 : init_mo_set,&
75 : mo_set_type
76 : USE qs_ot_eigensolver, ONLY: ot_eigensolver
77 : USE qs_rho_types, ONLY: qs_rho_get,&
78 : qs_rho_type
79 : USE qs_sccs, ONLY: print_sccs_results
80 : USE qs_scf_types, ONLY: ot_method_nr,&
81 : qs_scf_env_type,&
82 : special_diag_method_nr
83 : USE scf_control_types, ONLY: scf_control_type
84 : #include "./base/base_uses.f90"
85 :
86 : IMPLICIT NONE
87 :
88 : PRIVATE
89 :
90 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_scf_output'
91 :
92 : PUBLIC :: qs_scf_loop_info, &
93 : qs_scf_print_summary, &
94 : qs_scf_loop_print, &
95 : qs_scf_outer_loop_info, &
96 : qs_scf_initial_info, &
97 : qs_scf_write_mos, &
98 : qs_scf_cdft_info, &
99 : qs_scf_cdft_initial_info, &
100 : qs_scf_cdft_constraint_info
101 :
102 : CONTAINS
103 :
104 : ! **************************************************************************************************
105 : !> \brief writes a summary of information after scf
106 : !> \param output_unit ...
107 : !> \param qs_env ...
108 : ! **************************************************************************************************
109 23699 : SUBROUTINE qs_scf_print_summary(output_unit, qs_env)
110 : INTEGER, INTENT(IN) :: output_unit
111 : TYPE(qs_environment_type), POINTER :: qs_env
112 :
113 : INTEGER :: nelectron_total
114 : LOGICAL :: gapw, gapw_xc, qmmm
115 : TYPE(dft_control_type), POINTER :: dft_control
116 : TYPE(qs_charges_type), POINTER :: qs_charges
117 : TYPE(qs_energy_type), POINTER :: energy
118 : TYPE(qs_rho_type), POINTER :: rho
119 : TYPE(qs_scf_env_type), POINTER :: scf_env
120 :
121 23699 : NULLIFY (rho, energy, dft_control, scf_env, qs_charges)
122 : CALL get_qs_env(qs_env=qs_env, rho=rho, energy=energy, dft_control=dft_control, &
123 23699 : scf_env=scf_env, qs_charges=qs_charges)
124 :
125 23699 : gapw = dft_control%qs_control%gapw
126 23699 : gapw_xc = dft_control%qs_control%gapw_xc
127 23699 : qmmm = qs_env%qmmm
128 23699 : nelectron_total = scf_env%nelectron
129 :
130 : CALL qs_scf_print_scf_summary(output_unit, rho, qs_charges, energy, nelectron_total, &
131 23699 : dft_control, qmmm, qs_env, gapw, gapw_xc)
132 :
133 23699 : END SUBROUTINE qs_scf_print_summary
134 :
135 : ! **************************************************************************************************
136 : !> \brief writes basic information at the beginning of an scf run
137 : !> \param output_unit ...
138 : !> \param mos ...
139 : !> \param dft_control ...
140 : !> \param ndep ...
141 : ! **************************************************************************************************
142 25345 : SUBROUTINE qs_scf_initial_info(output_unit, mos, dft_control, ndep)
143 : INTEGER :: output_unit
144 : TYPE(mo_set_type), DIMENSION(:), INTENT(IN) :: mos
145 : TYPE(dft_control_type), POINTER :: dft_control
146 : INTEGER, INTENT(IN) :: ndep
147 :
148 : CHARACTER(LEN=*), PARAMETER :: routineN = 'qs_scf_initial_info'
149 :
150 : INTEGER :: handle, homo, ispin, nao, &
151 : nelectron_spin, nmo
152 :
153 25345 : CALL timeset(routineN, handle)
154 :
155 25345 : IF (output_unit > 0) THEN
156 27356 : DO ispin = 1, dft_control%nspins
157 : CALL get_mo_set(mo_set=mos(ispin), &
158 : homo=homo, &
159 : nelectron=nelectron_spin, &
160 : nao=nao, &
161 14506 : nmo=nmo)
162 14506 : IF (dft_control%nspins > 1) THEN
163 3312 : WRITE (UNIT=output_unit, FMT="(/,T2,A,I2)") "Spin", ispin
164 : END IF
165 : WRITE (UNIT=output_unit, FMT="(/,(T2,A,T71,I10))") &
166 14506 : "Number of electrons:", nelectron_spin, &
167 14506 : "Number of occupied orbitals:", homo, &
168 56368 : "Number of molecular orbitals:", nmo
169 : END DO
170 : WRITE (UNIT=output_unit, FMT="(/,(T2,A,T71,I10))") &
171 12850 : "Number of orbital functions:", nao, &
172 25700 : "Number of independent orbital functions:", nao - ndep
173 : END IF
174 :
175 25345 : CALL timestop(handle)
176 :
177 25345 : END SUBROUTINE qs_scf_initial_info
178 :
179 : ! **************************************************************************************************
180 : !> \brief Write the MO eigenvector, eigenvalues, and occupation numbers to the output unit
181 : !> \param qs_env ...
182 : !> \param scf_env ...
183 : !> \param final_mos ...
184 : !> \par History
185 : !> - Revise MO printout to enable eigenvalues with OT (05.05.2021, MK)
186 : ! **************************************************************************************************
187 935568 : SUBROUTINE qs_scf_write_mos(qs_env, scf_env, final_mos)
188 : TYPE(qs_environment_type), POINTER :: qs_env
189 : TYPE(qs_scf_env_type), POINTER :: scf_env
190 : LOGICAL, INTENT(IN) :: final_mos
191 :
192 : CHARACTER(LEN=*), PARAMETER :: routineN = 'qs_scf_write_mos'
193 :
194 : CHARACTER(LEN=2) :: solver_method
195 : CHARACTER(LEN=3*default_string_length) :: message
196 : CHARACTER(LEN=5) :: spin
197 : CHARACTER(LEN=default_string_length), &
198 233892 : DIMENSION(:), POINTER :: tmpstringlist
199 : INTEGER :: handle, homo, ikp, ispin, iw, kpoint, &
200 : nao, nelectron, nkp, nmo, nspin, numo
201 : INTEGER, DIMENSION(2) :: nmos_occ
202 233892 : INTEGER, DIMENSION(:), POINTER :: mo_index_range
203 : LOGICAL :: do_kpoints, do_printout, print_eigvals, &
204 : print_eigvecs, print_mo_info, &
205 : print_occup, print_occup_stats
206 : REAL(KIND=dp) :: flexible_electron_count, maxocc, n_el_f, &
207 : occup_stats_occ_threshold
208 233892 : REAL(KIND=dp), DIMENSION(:), POINTER :: mo_eigenvalues, umo_eigenvalues
209 : TYPE(admm_type), POINTER :: admm_env
210 : TYPE(cp_blacs_env_type), POINTER :: blacs_env
211 : TYPE(cp_fm_struct_type), POINTER :: fm_struct_tmp
212 : TYPE(cp_fm_type), POINTER :: mo_coeff, umo_coeff
213 : TYPE(cp_logger_type), POINTER :: logger
214 233892 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: ks, s
215 : TYPE(dbcsr_type), POINTER :: matrix_ks, matrix_s, mo_coeff_deriv
216 : TYPE(dft_control_type), POINTER :: dft_control
217 : TYPE(kpoint_type), POINTER :: kpoints
218 233892 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
219 : TYPE(mo_set_type), POINTER :: mo_set, umo_set
220 : TYPE(mp_para_env_type), POINTER :: para_env
221 233892 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
222 : TYPE(preconditioner_type), POINTER :: local_preconditioner
223 : TYPE(qs_environment_type), POINTER :: cart_overlap_qs_env
224 233892 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
225 : TYPE(scf_control_type), POINTER :: scf_control
226 : TYPE(section_vals_type), POINTER :: dft_section, input
227 :
228 233892 : CALL timeset(routineN, handle)
229 :
230 233892 : CPASSERT(ASSOCIATED(qs_env))
231 :
232 : ! Retrieve the required information for the requested print output
233 : CALL get_qs_env(qs_env, &
234 : blacs_env=blacs_env, &
235 : dft_control=dft_control, &
236 : do_kpoints=do_kpoints, &
237 : input=input, &
238 : qs_kind_set=qs_kind_set, &
239 : para_env=para_env, &
240 : particle_set=particle_set, &
241 233892 : scf_control=scf_control)
242 :
243 : ! Quick return, if no printout of MO information is requested
244 233892 : dft_section => section_vals_get_subs_vals(input, "DFT")
245 233892 : CALL section_vals_val_get(dft_section, "PRINT%MO%EIGENVALUES", l_val=print_eigvals)
246 233892 : CALL section_vals_val_get(dft_section, "PRINT%MO%EIGENVECTORS", l_val=print_eigvecs)
247 233892 : CALL section_vals_val_get(dft_section, "PRINT%MO%OCCUPATION_NUMBERS", l_val=print_occup)
248 233892 : CALL section_vals_val_get(dft_section, "PRINT%MO%OCCUPATION_NUMBERS_STATS", c_vals=tmpstringlist)
249 :
250 233892 : print_occup_stats = .FALSE.
251 233892 : occup_stats_occ_threshold = 1e-6_dp
252 233892 : IF (SIZE(tmpstringlist) > 0) THEN ! the lone_keyword_c_vals doesn't work as advertised, handle it manually
253 233892 : print_occup_stats = .TRUE.
254 233892 : IF (LEN_TRIM(tmpstringlist(1)) > 0) &
255 233892 : READ (tmpstringlist(1), *) print_occup_stats
256 : END IF
257 233892 : IF (SIZE(tmpstringlist) > 1) &
258 233892 : READ (tmpstringlist(2), *) occup_stats_occ_threshold
259 :
260 233892 : logger => cp_get_default_logger()
261 233892 : print_mo_info = (cp_print_key_should_output(logger%iter_info, dft_section, "PRINT%MO") /= 0)
262 :
263 233892 : IF ((.NOT. print_mo_info) .OR. (.NOT. (print_eigvals .OR. print_eigvecs .OR. print_occup .OR. print_occup_stats))) THEN
264 226448 : CALL timestop(handle)
265 226448 : RETURN
266 : END IF
267 :
268 7444 : NULLIFY (fm_struct_tmp)
269 7444 : NULLIFY (mo_coeff)
270 7444 : NULLIFY (mo_coeff_deriv)
271 7444 : NULLIFY (mo_eigenvalues)
272 7444 : NULLIFY (mo_set)
273 7444 : NULLIFY (umo_coeff)
274 7444 : NULLIFY (umo_eigenvalues)
275 7444 : NULLIFY (umo_set)
276 :
277 7444 : do_printout = .TRUE.
278 7444 : nspin = dft_control%nspins
279 7444 : nmos_occ = 0
280 :
281 : ! Check, if we have k points
282 7444 : IF (do_kpoints) THEN
283 22 : CALL get_qs_env(qs_env, kpoints=kpoints)
284 22 : nkp = SIZE(kpoints%kp_env)
285 : ELSE
286 7422 : CALL get_qs_env(qs_env, matrix_ks=ks, matrix_s=s)
287 7422 : CPASSERT(ASSOCIATED(ks))
288 7422 : CPASSERT(ASSOCIATED(s))
289 : nkp = 1
290 : END IF
291 :
292 12484 : kp_loop: DO ikp = 1, nkp
293 :
294 8008 : IF (do_kpoints) THEN
295 586 : mos => kpoints%kp_env(ikp)%kpoint_env%mos(1, :)
296 586 : kpoint = ikp
297 : ELSE
298 7422 : CALL get_qs_env(qs_env, matrix_ks=ks, mos=mos)
299 7422 : kpoint = 0 ! Gamma point only
300 : END IF
301 8008 : CPASSERT(ASSOCIATED(mos))
302 :
303 : ! Prepare MO information for printout
304 17818 : DO ispin = 1, nspin
305 :
306 : ! Calculate MO eigenvalues and eigenvector when OT is used
307 8302 : IF (scf_env%method == ot_method_nr) THEN
308 :
309 3190 : solver_method = "OT"
310 :
311 3190 : IF (do_kpoints) THEN
312 0 : CPABORT("The OT method is not implemented for k points")
313 : END IF
314 :
315 3190 : IF (final_mos) THEN
316 :
317 222 : matrix_ks => ks(ispin)%matrix
318 222 : matrix_s => s(1)%matrix
319 :
320 : ! With ADMM, we have to modify the Kohn-Sham matrix
321 222 : IF (dft_control%do_admm) THEN
322 0 : CALL get_qs_env(qs_env, admm_env=admm_env)
323 0 : CALL admm_correct_for_eigenvalues(ispin, admm_env, matrix_ks)
324 : END IF
325 :
326 222 : mo_set => mos(ispin)
327 : CALL get_mo_set(mo_set=mo_set, &
328 : mo_coeff=mo_coeff, &
329 : eigenvalues=mo_eigenvalues, &
330 : homo=homo, &
331 : maxocc=maxocc, &
332 : nelectron=nelectron, &
333 : n_el_f=n_el_f, &
334 : nao=nao, &
335 : nmo=nmo, &
336 222 : flexible_electron_count=flexible_electron_count)
337 :
338 222 : IF (ASSOCIATED(qs_env%mo_derivs)) THEN
339 222 : mo_coeff_deriv => qs_env%mo_derivs(ispin)%matrix
340 : ELSE
341 0 : mo_coeff_deriv => NULL()
342 : END IF
343 :
344 : ! Update the eigenvalues of the occupied orbitals
345 : CALL calculate_subspace_eigenvalues(orbitals=mo_coeff, &
346 : ks_matrix=matrix_ks, &
347 : evals_arg=mo_eigenvalues, &
348 222 : co_rotate_dbcsr=mo_coeff_deriv)
349 222 : CALL set_mo_occupation(mo_set=mo_set)
350 :
351 : ! Retrieve the index of the last MO for which a printout is requested
352 222 : mo_index_range => section_get_ivals(dft_section, "PRINT%MO%MO_INDEX_RANGE")
353 222 : CPASSERT(ASSOCIATED(mo_index_range))
354 222 : IF (mo_index_range(2) < 0) THEN
355 0 : numo = nao - homo
356 : ELSE
357 222 : numo = MIN(mo_index_range(2) - homo, nao - homo)
358 : END IF
359 :
360 : ! Calculate the unoccupied MO set (umo_set) with OT if needed
361 222 : IF (numo > 0) THEN
362 :
363 : ! Create temporary virtual MO set for printout
364 : CALL cp_fm_struct_create(fm_struct_tmp, &
365 : context=blacs_env, &
366 : para_env=para_env, &
367 : nrow_global=nao, &
368 20 : ncol_global=numo)
369 20 : ALLOCATE (umo_set)
370 : CALL allocate_mo_set(mo_set=umo_set, &
371 : nao=nao, &
372 : nmo=numo, &
373 : nelectron=0, &
374 : n_el_f=n_el_f, &
375 : maxocc=maxocc, &
376 20 : flexible_electron_count=flexible_electron_count)
377 : CALL init_mo_set(mo_set=umo_set, &
378 : fm_struct=fm_struct_tmp, &
379 20 : name="Temporary MO set (unoccupied MOs only) for printout")
380 20 : CALL cp_fm_struct_release(fm_struct_tmp)
381 : CALL get_mo_set(mo_set=umo_set, &
382 : mo_coeff=umo_coeff, &
383 20 : eigenvalues=umo_eigenvalues)
384 :
385 : ! Prepare printout of the additional unoccupied MOs when OT is being employed
386 20 : CALL cp_fm_init_random(umo_coeff)
387 :
388 : ! The FULL_ALL preconditioner makes not much sense for the unoccupied orbitals
389 20 : NULLIFY (local_preconditioner)
390 20 : IF (ASSOCIATED(scf_env%ot_preconditioner)) THEN
391 20 : local_preconditioner => scf_env%ot_preconditioner(1)%preconditioner
392 20 : IF (local_preconditioner%in_use == ot_precond_full_all) THEN
393 0 : NULLIFY (local_preconditioner)
394 : END IF
395 : END IF
396 :
397 : ! Calculate the MO information for the request MO index range
398 : CALL ot_eigensolver(matrix_h=matrix_ks, &
399 : matrix_s=matrix_s, &
400 : matrix_c_fm=umo_coeff, &
401 : matrix_orthogonal_space_fm=mo_coeff, &
402 : eps_gradient=scf_control%eps_lumos, &
403 : preconditioner=local_preconditioner, &
404 : iter_max=scf_control%max_iter_lumos, &
405 20 : size_ortho_space=nmo)
406 :
407 : CALL calculate_subspace_eigenvalues(orbitals=umo_coeff, &
408 : ks_matrix=matrix_ks, &
409 20 : evals_arg=umo_eigenvalues)
410 20 : CALL set_mo_occupation(mo_set=umo_set)
411 :
412 : END IF ! numo > 0
413 :
414 : ! With ADMM, we have to undo the modification of the Kohn-Sham matrix
415 222 : IF (dft_control%do_admm) THEN
416 0 : CALL admm_uncorrect_for_eigenvalues(ispin, admm_env, matrix_ks)
417 : END IF
418 :
419 : ELSE
420 :
421 : message = "The MO information is only calculated after SCF convergence "// &
422 2968 : "is achieved when the orbital transformation (OT) method is used"
423 2968 : CPWARN(TRIM(message))
424 2968 : do_printout = .FALSE.
425 2968 : EXIT kp_loop
426 :
427 : END IF ! final MOs
428 :
429 : ELSE
430 :
431 5112 : solver_method = "TD"
432 5112 : mo_set => mos(ispin)
433 5112 : NULLIFY (umo_set)
434 :
435 : END IF ! OT is used
436 :
437 : ! Print MO information
438 5334 : NULLIFY (cart_overlap_qs_env)
439 5334 : IF ((ikp == 1) .AND. (ispin == 1)) cart_overlap_qs_env => qs_env
440 5334 : IF (nspin > 1) THEN
441 294 : SELECT CASE (ispin)
442 : CASE (1)
443 294 : spin = "ALPHA"
444 : CASE (2)
445 294 : spin = "BETA"
446 : CASE DEFAULT
447 588 : CPABORT("Invalid spin")
448 : END SELECT
449 588 : IF (ASSOCIATED(umo_set)) THEN
450 : CALL write_mo_set_to_output_unit(mo_set, qs_kind_set, particle_set, dft_section, 4, kpoint, &
451 : final_mos=final_mos, spin=TRIM(spin), solver_method=solver_method, &
452 12 : umo_set=umo_set, qs_env=cart_overlap_qs_env)
453 : ELSE
454 : CALL write_mo_set_to_output_unit(mo_set, qs_kind_set, particle_set, dft_section, 4, kpoint, &
455 : final_mos=final_mos, spin=TRIM(spin), solver_method=solver_method, &
456 576 : qs_env=cart_overlap_qs_env)
457 : END IF
458 : ELSE
459 4746 : IF (ASSOCIATED(umo_set)) THEN
460 : CALL write_mo_set_to_output_unit(mo_set, qs_kind_set, particle_set, dft_section, 4, kpoint, &
461 : final_mos=final_mos, solver_method=solver_method, &
462 8 : umo_set=umo_set, qs_env=cart_overlap_qs_env)
463 : ELSE
464 : CALL write_mo_set_to_output_unit(mo_set, qs_kind_set, particle_set, dft_section, 4, kpoint, &
465 : final_mos=final_mos, solver_method=solver_method, &
466 4738 : qs_env=cart_overlap_qs_env)
467 : END IF
468 : END IF
469 :
470 46544 : nmos_occ(ispin) = MAX(nmos_occ(ispin), COUNT(mo_set%occupation_numbers > occup_stats_occ_threshold))
471 :
472 : ! Deallocate temporary objects needed for OT
473 5334 : IF (scf_env%method == ot_method_nr) THEN
474 222 : IF (ASSOCIATED(umo_set)) THEN
475 20 : CALL deallocate_mo_set(umo_set)
476 20 : DEALLOCATE (umo_set)
477 : END IF
478 222 : NULLIFY (matrix_ks)
479 222 : NULLIFY (matrix_s)
480 : END IF
481 10374 : NULLIFY (mo_set)
482 :
483 : END DO ! ispin
484 :
485 : END DO kp_loop
486 :
487 7444 : IF (do_printout .AND. print_mo_info .AND. print_occup_stats) THEN
488 : iw = cp_print_key_unit_nr(logger, dft_section, "PRINT%MO", &
489 : ignore_should_output=print_mo_info, &
490 0 : extension=".MOLog")
491 0 : IF (iw > 0) THEN
492 0 : IF (SIZE(mos) > 1) THEN
493 0 : WRITE (UNIT=iw, FMT="(A,I4)") " MO| Total occupied (ALPHA):", nmos_occ(1)
494 0 : WRITE (UNIT=iw, FMT="(A,I4)") " MO| Total occupied (BETA): ", nmos_occ(2)
495 : ELSE
496 0 : WRITE (UNIT=iw, FMT="(A,I4)") " MO| Total occupied: ", nmos_occ(1)
497 : END IF
498 0 : WRITE (UNIT=iw, FMT="(A)") ""
499 : END IF
500 : CALL cp_print_key_finished_output(iw, logger, dft_section, "PRINT%MO", &
501 0 : ignore_should_output=print_mo_info)
502 : END IF
503 :
504 7444 : CALL timestop(handle)
505 :
506 233892 : END SUBROUTINE qs_scf_write_mos
507 :
508 : ! **************************************************************************************************
509 : !> \brief writes basic information obtained in a scf outer loop step
510 : !> \param output_unit ...
511 : !> \param scf_control ...
512 : !> \param scf_env ...
513 : !> \param energy ...
514 : !> \param total_steps ...
515 : !> \param should_stop ...
516 : !> \param outer_loop_converged ...
517 : ! **************************************************************************************************
518 5399 : SUBROUTINE qs_scf_outer_loop_info(output_unit, scf_control, scf_env, &
519 : energy, total_steps, should_stop, outer_loop_converged)
520 : INTEGER :: output_unit
521 : TYPE(scf_control_type), POINTER :: scf_control
522 : TYPE(qs_scf_env_type), POINTER :: scf_env
523 : TYPE(qs_energy_type), POINTER :: energy
524 : INTEGER :: total_steps
525 : LOGICAL, INTENT(IN) :: should_stop, outer_loop_converged
526 :
527 : REAL(KIND=dp) :: outer_loop_eps
528 :
529 10798 : outer_loop_eps = SQRT(MAXVAL(scf_env%outer_scf%gradient(:, scf_env%outer_scf%iter_count)**2))
530 5399 : IF (output_unit > 0) WRITE (output_unit, '(/,T3,A,I4,A,E10.2,A,F22.10)') &
531 2820 : "outer SCF iter = ", scf_env%outer_scf%iter_count, &
532 5640 : " RMS gradient = ", outer_loop_eps, " energy =", energy%total
533 :
534 5399 : IF (outer_loop_converged) THEN
535 4443 : IF (output_unit > 0) WRITE (output_unit, '(T3,A,I4,A,I4,A,/)') &
536 2335 : "outer SCF loop converged in", scf_env%outer_scf%iter_count, &
537 4670 : " iterations or ", total_steps, " steps"
538 : ELSE IF (scf_env%outer_scf%iter_count > scf_control%outer_scf%max_scf &
539 956 : .OR. should_stop) THEN
540 108 : IF (output_unit > 0) WRITE (output_unit, '(T3,A,I4,A,I4,A,/)') &
541 54 : "outer SCF loop FAILED to converge after ", &
542 108 : scf_env%outer_scf%iter_count, " iterations or ", total_steps, " steps"
543 : END IF
544 :
545 5399 : END SUBROUTINE qs_scf_outer_loop_info
546 :
547 : ! **************************************************************************************************
548 : !> \brief writes basic information obtained in a scf step
549 : !> \param scf_env ...
550 : !> \param output_unit ...
551 : !> \param just_energy ...
552 : !> \param t1 ...
553 : !> \param t2 ...
554 : !> \param energy ...
555 : ! **************************************************************************************************
556 215881 : SUBROUTINE qs_scf_loop_info(scf_env, output_unit, just_energy, t1, t2, energy)
557 :
558 : TYPE(qs_scf_env_type), POINTER :: scf_env
559 : INTEGER :: output_unit
560 : LOGICAL :: just_energy
561 : REAL(KIND=dp) :: t1, t2
562 : TYPE(qs_energy_type), POINTER :: energy
563 :
564 215881 : IF ((output_unit > 0) .AND. scf_env%print_iter_line) THEN
565 109112 : IF (just_energy) THEN
566 : WRITE (UNIT=output_unit, &
567 : FMT="(T2,A,1X,A,T20,E8.2,1X,F6.1,16X,F20.10)") &
568 7217 : " -", TRIM(scf_env%iter_method), scf_env%iter_param, t2 - t1, energy%total
569 : ELSE
570 92363 : IF ((ABS(scf_env%iter_delta) < 1.0E-8_dp) .OR. &
571 101895 : (ABS(scf_env%iter_delta) >= 1.0E5_dp)) THEN
572 : WRITE (UNIT=output_unit, &
573 : FMT="(T2,I5,1X,A,T20,E8.2,1X,F6.1,1X,ES14.4,1X,F20.10,1X,ES9.2)") &
574 9532 : scf_env%iter_count, TRIM(scf_env%iter_method), scf_env%iter_param, &
575 19064 : t2 - t1, scf_env%iter_delta, energy%total, energy%total - energy%tot_old
576 : ELSE
577 : WRITE (UNIT=output_unit, &
578 : FMT="(T2,I5,1X,A,T20,E8.2,1X,F6.1,1X,F14.8,1X,F20.10,1X,ES9.2)") &
579 92363 : scf_env%iter_count, TRIM(scf_env%iter_method), scf_env%iter_param, &
580 184726 : t2 - t1, scf_env%iter_delta, energy%total, energy%total - energy%tot_old
581 : END IF
582 : END IF
583 : END IF
584 :
585 215881 : END SUBROUTINE qs_scf_loop_info
586 :
587 : ! **************************************************************************************************
588 : !> \brief writes rather detailed summary of densities and energies
589 : !> after the SCF
590 : !> \param output_unit ...
591 : !> \param rho ...
592 : !> \param qs_charges ...
593 : !> \param energy ...
594 : !> \param nelectron_total ...
595 : !> \param dft_control ...
596 : !> \param qmmm ...
597 : !> \param qs_env ...
598 : !> \param gapw ...
599 : !> \param gapw_xc ...
600 : !> \par History
601 : !> 03.2006 created [Joost VandeVondele]
602 : !> 10.2019 print dipole moment [SGh]
603 : !> 11.2022 print SCCS results [MK]
604 : ! **************************************************************************************************
605 23699 : SUBROUTINE qs_scf_print_scf_summary(output_unit, rho, qs_charges, energy, nelectron_total, &
606 : dft_control, qmmm, qs_env, gapw, gapw_xc)
607 : INTEGER, INTENT(IN) :: output_unit
608 : TYPE(qs_rho_type), POINTER :: rho
609 : TYPE(qs_charges_type), POINTER :: qs_charges
610 : TYPE(qs_energy_type), POINTER :: energy
611 : INTEGER, INTENT(IN) :: nelectron_total
612 : TYPE(dft_control_type), POINTER :: dft_control
613 : LOGICAL, INTENT(IN) :: qmmm
614 : TYPE(qs_environment_type), POINTER :: qs_env
615 : LOGICAL, INTENT(IN) :: gapw, gapw_xc
616 :
617 : CHARACTER(LEN=*), PARAMETER :: routineN = 'qs_scf_print_scf_summary'
618 :
619 : INTEGER :: bc, handle, ispin, psolver
620 : REAL(kind=dp) :: e_extrapolated, exc1_energy, exc_energy, &
621 : implicit_ps_ehartree, tot1_h, tot1_s
622 23699 : REAL(KIND=dp), DIMENSION(:), POINTER :: tot_rho_r
623 : TYPE(pw_env_type), POINTER :: pw_env
624 : TYPE(scf_control_type), POINTER :: scf_control
625 :
626 23699 : NULLIFY (tot_rho_r, pw_env)
627 23699 : CALL timeset(routineN, handle)
628 :
629 23699 : CALL get_qs_env(qs_env=qs_env, pw_env=pw_env, scf_control=scf_control)
630 23699 : psolver = pw_env%poisson_env%parameters%solver
631 :
632 23699 : IF (output_unit > 0) THEN
633 12037 : CALL qs_rho_get(rho, tot_rho_r=tot_rho_r)
634 12037 : IF (.NOT. (dft_control%qs_control%semi_empirical .OR. &
635 : dft_control%qs_control%xtb .OR. &
636 : dft_control%qs_control%dftb)) THEN
637 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T41,2F20.10))") &
638 6111 : "Electronic density on regular grids: ", &
639 6111 : accurate_sum(tot_rho_r), &
640 6111 : accurate_sum(tot_rho_r) + nelectron_total, &
641 6111 : "Core density on regular grids:", &
642 6111 : qs_charges%total_rho_core_rspace, &
643 : qs_charges%total_rho_core_rspace + &
644 : qs_charges%total_rho1_hard_nuc - &
645 12222 : REAL(nelectron_total + dft_control%charge, dp)
646 :
647 6111 : IF (dft_control%correct_surf_dip) THEN
648 : WRITE (UNIT=output_unit, FMT="((T3,A,/,T3,A,T41,F20.10))") &
649 5 : "Total dipole moment perpendicular to ", &
650 5 : "the slab [electrons-Angstroem]: ", &
651 10 : qs_env%surface_dipole_moment
652 : END IF
653 :
654 6111 : IF (gapw) THEN
655 1095 : tot1_h = qs_charges%total_rho1_hard(1)
656 1095 : tot1_s = qs_charges%total_rho1_soft(1)
657 1307 : DO ispin = 2, dft_control%nspins
658 212 : tot1_h = tot1_h + qs_charges%total_rho1_hard(ispin)
659 1307 : tot1_s = tot1_s + qs_charges%total_rho1_soft(ispin)
660 : END DO
661 : WRITE (UNIT=output_unit, FMT="((T3,A,T41,2F20.10))") &
662 1095 : "Hard and soft densities (Lebedev):", &
663 2190 : tot1_h, tot1_s
664 : WRITE (UNIT=output_unit, FMT="(T3,A,T41,F20.10)") &
665 1095 : "Total Rho_soft + Rho1_hard - Rho1_soft (r-space): ", &
666 1095 : accurate_sum(tot_rho_r) + tot1_h - tot1_s, &
667 1095 : "Total charge density (r-space): ", &
668 : accurate_sum(tot_rho_r) + tot1_h - tot1_s &
669 : + qs_charges%total_rho_core_rspace &
670 2190 : + qs_charges%total_rho1_hard_nuc
671 1095 : IF (qs_charges%total_rho1_hard_nuc /= 0.0_dp) THEN
672 : WRITE (UNIT=output_unit, FMT="(T3,A,T41,F20.10)") &
673 4 : "Total CNEO nuc. char. den. (Lebedev): ", &
674 4 : qs_charges%total_rho1_hard_nuc, &
675 4 : "Total CNEO soft char. den. (Lebedev): ", &
676 4 : qs_charges%total_rho1_soft_nuc_lebedev, &
677 4 : "Total CNEO soft char. den. (r-space): ", &
678 4 : qs_charges%total_rho1_soft_nuc_rspace, &
679 4 : "Total soft Rho_e+n+0 (g-space):", &
680 8 : qs_charges%total_rho_gspace
681 : ELSE
682 : WRITE (UNIT=output_unit, FMT="(T3,A,T41,F20.10)") &
683 1091 : "Total Rho_soft + Rho0_soft (g-space):", &
684 2182 : qs_charges%total_rho_gspace
685 : END IF
686 : ! only add total_rho1_hard_nuc for gapw as cneo requires gapw
687 : ELSE
688 : WRITE (UNIT=output_unit, FMT="(T3,A,T41,F20.10)") &
689 5016 : "Total charge density on r-space grids: ", &
690 : accurate_sum(tot_rho_r) + &
691 5016 : qs_charges%total_rho_core_rspace, &
692 5016 : "Total charge density g-space grids: ", &
693 10032 : qs_charges%total_rho_gspace
694 : END IF
695 : END IF
696 12037 : IF (dft_control%qs_control%semi_empirical) THEN
697 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
698 1917 : "Core-core repulsion energy [eV]: ", energy%core_overlap*evolt, &
699 1917 : "Core Hamiltonian energy [eV]: ", energy%core*evolt, &
700 1917 : "Two-electron integral energy [eV]: ", energy%hartree*evolt, &
701 1917 : "Electronic energy [eV]: ", &
702 3834 : (energy%core + 0.5_dp*energy%hartree)*evolt
703 1917 : IF (energy%dispersion /= 0.0_dp) &
704 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
705 8 : "Dispersion energy [eV]: ", energy%dispersion*evolt
706 10120 : ELSEIF (dft_control%qs_control%dftb) THEN
707 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
708 1129 : "Core Hamiltonian energy: ", energy%core, &
709 1129 : "Repulsive potential energy: ", energy%repulsive, &
710 1129 : "Electronic energy: ", energy%hartree, &
711 2258 : "Dispersion energy: ", energy%dispersion
712 1129 : IF (energy%dftb3 /= 0.0_dp) &
713 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
714 376 : "DFTB3 3rd order energy: ", energy%dftb3
715 1129 : IF (energy%efield /= 0.0_dp) &
716 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
717 16 : "Electric field interaction energy: ", energy%efield
718 8991 : ELSEIF (dft_control%qs_control%xtb) THEN
719 2880 : IF (dft_control%qs_control%xtb_control%do_tblite) THEN
720 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
721 1137 : "Core Hamiltonian energy: ", energy%core, &
722 1137 : "Repulsive potential energy: ", energy%repulsive, &
723 1137 : "Electrostatic energy: ", energy%el_stat, &
724 1137 : "Self-consistent dispersion energy: ", energy%dispersion_sc, &
725 1137 : "Non-self consistent dispersion energy: ", energy%dispersion, &
726 2274 : "Correction for halogen bonding: ", energy%xtb_xb_inter
727 : ELSE
728 1743 : IF (dft_control%qs_control%xtb_control%gfn_type == 0) THEN
729 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
730 0 : "Core Hamiltonian energy: ", energy%core, &
731 0 : "Repulsive potential energy: ", energy%repulsive, &
732 0 : "SRB Correction energy: ", energy%srb, &
733 0 : "Charge equilibration energy: ", energy%eeq, &
734 0 : "Dispersion energy: ", energy%dispersion
735 1743 : ELSEIF (dft_control%qs_control%xtb_control%gfn_type == 1) THEN
736 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
737 1743 : "Core Hamiltonian energy: ", energy%core, &
738 1743 : "Repulsive potential energy: ", energy%repulsive, &
739 1743 : "Electronic energy: ", energy%hartree, &
740 1743 : "DFTB3 3rd order energy: ", energy%dftb3, &
741 3486 : "Dispersion energy: ", energy%dispersion
742 1743 : IF (dft_control%qs_control%xtb_control%xb_interaction) &
743 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
744 1723 : "Correction for halogen bonding: ", energy%xtb_xb_inter
745 0 : ELSEIF (dft_control%qs_control%xtb_control%gfn_type == 2) THEN
746 0 : CPABORT("gfn_typ 2 NYA")
747 : ELSE
748 0 : CPABORT("invalid gfn_typ")
749 : END IF
750 : END IF
751 2880 : IF (dft_control%qs_control%xtb_control%do_nonbonded) &
752 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
753 12 : "Correction for nonbonded interactions: ", energy%xtb_nonbonded
754 2880 : IF (energy%efield /= 0.0_dp) &
755 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
756 394 : "Electric field interaction energy: ", energy%efield
757 : ELSE
758 6111 : IF (dft_control%do_admm) THEN
759 530 : exc_energy = energy%exc + energy%exc_aux_fit
760 530 : IF (gapw .OR. gapw_xc) exc1_energy = energy%exc1 + energy%exc1_aux_fit
761 : ELSE
762 5581 : exc_energy = energy%exc
763 5581 : IF (gapw .OR. gapw_xc) exc1_energy = energy%exc1
764 : END IF
765 :
766 6111 : IF (psolver == pw_poisson_implicit) THEN
767 60 : implicit_ps_ehartree = pw_env%poisson_env%implicit_env%ehartree
768 60 : bc = pw_env%poisson_env%parameters%ps_implicit_params%boundary_condition
769 41 : SELECT CASE (bc)
770 : CASE (MIXED_PERIODIC_BC, MIXED_BC)
771 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
772 41 : "Overlap energy of the core charge distribution:", energy%core_overlap, &
773 41 : "Self energy of the core charge distribution: ", energy%core_self, &
774 41 : "Core Hamiltonian energy: ", energy%core, &
775 41 : "Hartree energy: ", implicit_ps_ehartree, &
776 41 : "Electric enthalpy: ", energy%hartree, &
777 82 : "Exchange-correlation energy: ", exc_energy
778 : CASE (PERIODIC_BC, NEUMANN_BC)
779 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
780 19 : "Overlap energy of the core charge distribution:", energy%core_overlap, &
781 19 : "Self energy of the core charge distribution: ", energy%core_self, &
782 19 : "Core Hamiltonian energy: ", energy%core, &
783 19 : "Hartree energy: ", energy%hartree, &
784 79 : "Exchange-correlation energy: ", exc_energy
785 : END SELECT
786 : ELSE
787 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
788 6051 : "Overlap energy of the core charge distribution:", energy%core_overlap, &
789 6051 : "Self energy of the core charge distribution: ", energy%core_self, &
790 6051 : "Core Hamiltonian energy: ", energy%core, &
791 6051 : "Hartree energy: ", energy%hartree, &
792 12102 : "Exchange-correlation energy: ", exc_energy
793 : END IF
794 6111 : IF (energy%e_hartree /= 0.0_dp) &
795 : WRITE (UNIT=output_unit, FMT="(T3,A,/,T3,A,T56,F25.14)") &
796 44 : "Coulomb Electron-Electron Interaction Energy ", &
797 88 : "- Already included in the total Hartree term ", energy%e_hartree
798 6111 : IF (energy%ex /= 0.0_dp) &
799 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
800 1183 : "Hartree-Fock Exchange energy: ", energy%ex
801 6111 : IF (energy%dispersion /= 0.0_dp) &
802 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
803 178 : "Dispersion energy: ", energy%dispersion
804 6111 : IF (energy%gcp /= 0.0_dp) &
805 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
806 2 : "gCP energy: ", energy%gcp
807 6111 : IF (energy%efield /= 0.0_dp) &
808 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
809 447 : "Electric field interaction energy: ", energy%efield
810 6111 : IF (gapw) THEN
811 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
812 1095 : "GAPW| Exc from hard and soft atomic rho1: ", exc1_energy, &
813 2190 : "GAPW| local Eh = 1 center integrals: ", energy%hartree_1c
814 : END IF
815 6111 : IF (gapw_xc) THEN
816 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
817 203 : "GAPW_XC| Exc from hard and soft atomic rho1: ", exc1_energy
818 : END IF
819 6111 : IF (energy%core_cneo /= 0.0_dp) THEN
820 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
821 4 : "CNEO| quantum nuclear core energy: ", energy%core_cneo
822 : END IF
823 : END IF
824 12037 : IF (dft_control%hairy_probes .EQV. .TRUE.) THEN
825 : WRITE (UNIT=output_unit, FMT="((T3,A,T56,F25.14))") &
826 2 : "Electronic entropic energy:", energy%kTS
827 : WRITE (UNIT=output_unit, FMT="((T3,A,T56,F25.14))") &
828 2 : "Fermi energy:", energy%efermi
829 : END IF
830 12037 : IF (dft_control%smear) THEN
831 1467 : SELECT CASE (scf_control%smear%method)
832 : CASE (smear_gaussian, smear_mp, smear_mv)
833 : ! kTS does not have physical meaning in these smearing methods
834 : WRITE (UNIT=output_unit, FMT="((T3,A,T56,F25.14))") &
835 63 : "Smearing free energy correction:", energy%kTS
836 : CASE DEFAULT
837 : WRITE (UNIT=output_unit, FMT="((T3,A,T56,F25.14))") &
838 1404 : "Electronic entropic energy:", energy%kTS
839 : END SELECT
840 : WRITE (UNIT=output_unit, FMT="((T3,A,T56,F25.14))") &
841 1404 : "Fermi energy:", energy%efermi
842 : END IF
843 12037 : IF (dft_control%dft_plus_u) THEN
844 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
845 50 : "DFT+U energy:", energy%dft_plus_u
846 : END IF
847 12037 : IF (dft_control%do_sccs) THEN
848 6 : WRITE (UNIT=output_unit, FMT="(A)") ""
849 6 : CALL print_sccs_results(energy, dft_control%sccs_control, output_unit)
850 : END IF
851 12037 : IF (qmmm) THEN
852 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
853 1746 : "QM/MM Electrostatic energy: ", energy%qmmm_el
854 1746 : IF (qs_env%qmmm_env_qm%image_charge) THEN
855 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
856 10 : "QM/MM image charge energy: ", energy%image_charge
857 : END IF
858 : END IF
859 12037 : IF (dft_control%qs_control%mulliken_restraint) THEN
860 : WRITE (UNIT=output_unit, FMT="(T3,A,T56,F25.14)") &
861 3 : "Mulliken restraint energy: ", energy%mulliken
862 : END IF
863 12037 : IF (dft_control%qs_control%semi_empirical) THEN
864 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
865 1917 : "Total energy [eV]: ", energy%total*evolt
866 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
867 1917 : "Atomic reference energy [eV]: ", energy%core_self*evolt, &
868 1917 : "Heat of formation [kcal/mol]: ", &
869 3834 : (energy%total + energy%core_self)*kcalmol
870 : ELSE
871 : WRITE (UNIT=output_unit, FMT="(/,(T3,A,T56,F25.14))") &
872 10120 : "Total energy: ", energy%total
873 10120 : IF (dft_control%smear) THEN
874 2738 : SELECT CASE (scf_control%smear%method)
875 : CASE (smear_fermi_dirac)
876 1334 : e_extrapolated = energy%total - 0.5_dp*energy%kTS
877 : WRITE (UNIT=output_unit, FMT="((T3,A,T56,F25.14))") &
878 1334 : "Total energy (extrapolated to T->0): ", e_extrapolated
879 : CASE (smear_gaussian)
880 61 : e_extrapolated = energy%total - 0.5_dp*energy%kTS
881 : WRITE (UNIT=output_unit, FMT="((T3,A,T56,F25.14))") &
882 1404 : "Total energy (extrapolated to sigma->0): ", e_extrapolated
883 : CASE (smear_mp, smear_mv)
884 : ! Sigma->0 extrapolation does not apply to MP or MV method.
885 : END SELECT
886 : END IF
887 : END IF
888 12037 : IF (qmmm) THEN
889 1746 : IF (qs_env%qmmm_env_qm%image_charge) THEN
890 10 : CALL print_image_coefficients(qs_env%image_coeff, qs_env)
891 : END IF
892 : END IF
893 12037 : CALL m_flush(output_unit)
894 : END IF
895 :
896 23699 : CALL timestop(handle)
897 :
898 23699 : END SUBROUTINE qs_scf_print_scf_summary
899 :
900 : ! **************************************************************************************************
901 : !> \brief collects the 'heavy duty' printing tasks out of the SCF loop
902 : !> \param qs_env ...
903 : !> \param scf_env ...
904 : !> \param para_env ...
905 : !> \par History
906 : !> 03.2006 created [Joost VandeVondele]
907 : ! **************************************************************************************************
908 656133 : SUBROUTINE qs_scf_loop_print(qs_env, scf_env, para_env)
909 : TYPE(qs_environment_type), POINTER :: qs_env
910 : TYPE(qs_scf_env_type), POINTER :: scf_env
911 : TYPE(mp_para_env_type), POINTER :: para_env
912 :
913 : CHARACTER(LEN=*), PARAMETER :: routineN = 'qs_scf_loop_print'
914 :
915 : INTEGER :: after, handle, ic, ispin, iw
916 : LOGICAL :: do_kpoints, omit_headers
917 : REAL(KIND=dp) :: mo_mag_max, mo_mag_min, orthonormality
918 : TYPE(cp_logger_type), POINTER :: logger
919 218711 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_ks, matrix_p, matrix_s
920 : TYPE(dft_control_type), POINTER :: dft_control
921 218711 : TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
922 : TYPE(qs_rho_type), POINTER :: rho
923 : TYPE(section_vals_type), POINTER :: dft_section, input, scf_section
924 :
925 437422 : logger => cp_get_default_logger()
926 218711 : CALL timeset(routineN, handle)
927 :
928 : CALL get_qs_env(qs_env=qs_env, input=input, dft_control=dft_control, &
929 218711 : do_kpoints=do_kpoints)
930 :
931 218711 : dft_section => section_vals_get_subs_vals(input, "DFT")
932 218711 : scf_section => section_vals_get_subs_vals(dft_section, "SCF")
933 :
934 218711 : CALL section_vals_val_get(input, "DFT%PRINT%AO_MATRICES%OMIT_HEADERS", l_val=omit_headers)
935 467933 : DO ispin = 1, dft_control%nspins
936 :
937 249222 : IF (BTEST(cp_print_key_should_output(logger%iter_info, &
938 : dft_section, "PRINT%AO_MATRICES/DENSITY"), cp_p_file)) THEN
939 6888 : CALL get_qs_env(qs_env, rho=rho)
940 6888 : CALL qs_rho_get(rho, rho_ao_kp=matrix_p)
941 : iw = cp_print_key_unit_nr(logger, dft_section, "PRINT%AO_MATRICES/DENSITY", &
942 6888 : extension=".Log")
943 6888 : CALL section_vals_val_get(dft_section, "PRINT%AO_MATRICES%NDIGITS", i_val=after)
944 6888 : after = MIN(MAX(after, 1), 16)
945 13776 : DO ic = 1, SIZE(matrix_p, 2)
946 : CALL cp_dbcsr_write_sparse_matrix(matrix_p(ispin, ic)%matrix, 4, after, qs_env, para_env, &
947 13776 : output_unit=iw, omit_headers=omit_headers)
948 : END DO
949 : CALL cp_print_key_finished_output(iw, logger, dft_section, &
950 6888 : "PRINT%AO_MATRICES/DENSITY")
951 : END IF
952 :
953 249222 : IF (BTEST(cp_print_key_should_output(logger%iter_info, &
954 218711 : dft_section, "PRINT%AO_MATRICES/KOHN_SHAM_MATRIX"), cp_p_file)) THEN
955 : iw = cp_print_key_unit_nr(logger, dft_section, "PRINT%AO_MATRICES/KOHN_SHAM_MATRIX", &
956 5762 : extension=".Log")
957 5762 : CALL section_vals_val_get(dft_section, "PRINT%AO_MATRICES%NDIGITS", i_val=after)
958 5762 : after = MIN(MAX(after, 1), 16)
959 5762 : CALL get_qs_env(qs_env=qs_env, matrix_ks_kp=matrix_ks)
960 11524 : DO ic = 1, SIZE(matrix_ks, 2)
961 11524 : IF (dft_control%qs_control%semi_empirical) THEN
962 : CALL cp_dbcsr_write_sparse_matrix(matrix_ks(ispin, ic)%matrix, 4, after, qs_env, para_env, &
963 5758 : scale=evolt, output_unit=iw, omit_headers=omit_headers)
964 : ELSE
965 : CALL cp_dbcsr_write_sparse_matrix(matrix_ks(ispin, ic)%matrix, 4, after, qs_env, para_env, &
966 4 : output_unit=iw, omit_headers=omit_headers)
967 : END IF
968 : END DO
969 : CALL cp_print_key_finished_output(iw, logger, dft_section, &
970 5762 : "PRINT%AO_MATRICES/KOHN_SHAM_MATRIX")
971 : END IF
972 :
973 : END DO
974 :
975 218711 : IF (BTEST(cp_print_key_should_output(logger%iter_info, &
976 : scf_section, "PRINT%MO_ORTHONORMALITY"), cp_p_file)) THEN
977 1172 : IF (do_kpoints) THEN
978 : iw = cp_print_key_unit_nr(logger, scf_section, "PRINT%MO_ORTHONORMALITY", &
979 16 : extension=".scfLog")
980 16 : IF (iw > 0) THEN
981 : WRITE (iw, '(T8,A)') &
982 8 : " K-points: Maximum deviation from MO S-orthonormality not determined"
983 : END IF
984 : CALL cp_print_key_finished_output(iw, logger, scf_section, &
985 16 : "PRINT%MO_ORTHONORMALITY")
986 : ELSE
987 1156 : CALL get_qs_env(qs_env, mos=mos)
988 1156 : IF (scf_env%method == special_diag_method_nr) THEN
989 58 : CALL calculate_orthonormality(orthonormality, mos)
990 : ELSE
991 1098 : CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrix_s)
992 1098 : CALL calculate_orthonormality(orthonormality, mos, matrix_s(1, 1)%matrix)
993 : END IF
994 : iw = cp_print_key_unit_nr(logger, scf_section, "PRINT%MO_ORTHONORMALITY", &
995 1156 : extension=".scfLog")
996 1156 : IF (iw > 0) THEN
997 : WRITE (iw, '(T8,A,T61,E20.4)') &
998 578 : " Maximum deviation from MO S-orthonormality", orthonormality
999 : END IF
1000 : CALL cp_print_key_finished_output(iw, logger, scf_section, &
1001 1156 : "PRINT%MO_ORTHONORMALITY")
1002 : END IF
1003 : END IF
1004 218711 : IF (BTEST(cp_print_key_should_output(logger%iter_info, &
1005 : scf_section, "PRINT%MO_MAGNITUDE"), cp_p_file)) THEN
1006 1172 : IF (do_kpoints) THEN
1007 : iw = cp_print_key_unit_nr(logger, scf_section, "PRINT%MO_MAGNITUDE", &
1008 16 : extension=".scfLog")
1009 16 : IF (iw > 0) THEN
1010 : WRITE (iw, '(T8,A)') &
1011 8 : " K-points: Minimum/Maximum MO magnitude not determined"
1012 : END IF
1013 : CALL cp_print_key_finished_output(iw, logger, scf_section, &
1014 16 : "PRINT%MO_MAGNITUDE")
1015 : ELSE
1016 1156 : CALL get_qs_env(qs_env, mos=mos)
1017 1156 : CALL calculate_magnitude(mos, mo_mag_min, mo_mag_max)
1018 : iw = cp_print_key_unit_nr(logger, scf_section, "PRINT%MO_MAGNITUDE", &
1019 1156 : extension=".scfLog")
1020 1156 : IF (iw > 0) THEN
1021 : WRITE (iw, '(T8,A,T41,2E20.4)') &
1022 578 : " Minimum/Maximum MO magnitude ", mo_mag_min, mo_mag_max
1023 : END IF
1024 : CALL cp_print_key_finished_output(iw, logger, scf_section, &
1025 1156 : "PRINT%MO_MAGNITUDE")
1026 : END IF
1027 : END IF
1028 :
1029 218711 : CALL timestop(handle)
1030 :
1031 218711 : END SUBROUTINE qs_scf_loop_print
1032 :
1033 : ! **************************************************************************************************
1034 : !> \brief writes CDFT constraint information and optionally CDFT scf loop info
1035 : !> \param output_unit where to write the information
1036 : !> \param scf_control settings of the SCF loop
1037 : !> \param scf_env the env which holds convergence data
1038 : !> \param cdft_control the env which holds information about the constraint
1039 : !> \param energy the total energy
1040 : !> \param total_steps the total number of performed SCF iterations
1041 : !> \param should_stop if the calculation should stop
1042 : !> \param outer_loop_converged logical which determines if the CDFT SCF loop converged
1043 : !> \param cdft_loop logical which determines a CDFT SCF loop is active
1044 : !> \par History
1045 : !> 12.2015 created [Nico Holmberg]
1046 : ! **************************************************************************************************
1047 682 : SUBROUTINE qs_scf_cdft_info(output_unit, scf_control, scf_env, cdft_control, &
1048 : energy, total_steps, should_stop, outer_loop_converged, &
1049 : cdft_loop)
1050 : INTEGER :: output_unit
1051 : TYPE(scf_control_type), POINTER :: scf_control
1052 : TYPE(qs_scf_env_type), POINTER :: scf_env
1053 : TYPE(cdft_control_type), POINTER :: cdft_control
1054 : TYPE(qs_energy_type), POINTER :: energy
1055 : INTEGER :: total_steps
1056 : LOGICAL, INTENT(IN) :: should_stop, outer_loop_converged, &
1057 : cdft_loop
1058 :
1059 : REAL(KIND=dp) :: outer_loop_eps
1060 :
1061 682 : IF (cdft_loop) THEN
1062 1222 : outer_loop_eps = SQRT(MAXVAL(scf_env%outer_scf%gradient(:, scf_env%outer_scf%iter_count)**2))
1063 568 : IF (output_unit > 0) WRITE (output_unit, '(/,T3,A,I4,A,E10.2,A,F22.10)') &
1064 306 : "CDFT SCF iter = ", scf_env%outer_scf%iter_count, &
1065 612 : " RMS gradient = ", outer_loop_eps, " energy =", energy%total
1066 568 : IF (outer_loop_converged) THEN
1067 280 : IF (output_unit > 0) WRITE (output_unit, '(T3,A,I4,A,I4,A,/)') &
1068 159 : "CDFT SCF loop converged in", scf_env%outer_scf%iter_count, &
1069 318 : " iterations or ", total_steps, " steps"
1070 : END IF
1071 : IF ((scf_env%outer_scf%iter_count > scf_control%outer_scf%max_scf .OR. should_stop) &
1072 568 : .AND. .NOT. outer_loop_converged) THEN
1073 64 : IF (output_unit > 0) WRITE (output_unit, '(T3,A,I4,A,I4,A,/)') &
1074 32 : "CDFT SCF loop FAILED to converge after ", &
1075 64 : scf_env%outer_scf%iter_count, " iterations or ", total_steps, " steps"
1076 : END IF
1077 : END IF
1078 682 : CALL qs_scf_cdft_constraint_info(output_unit, cdft_control)
1079 :
1080 682 : END SUBROUTINE qs_scf_cdft_info
1081 :
1082 : ! **************************************************************************************************
1083 : !> \brief writes information about the CDFT env
1084 : !> \param output_unit where to write the information
1085 : !> \param cdft_control the CDFT env that stores information about the constraint calculation
1086 : !> \par History
1087 : !> 12.2015 created [Nico Holmberg]
1088 : ! **************************************************************************************************
1089 191 : SUBROUTINE qs_scf_cdft_initial_info(output_unit, cdft_control)
1090 : INTEGER :: output_unit
1091 : TYPE(cdft_control_type), POINTER :: cdft_control
1092 :
1093 191 : IF (output_unit > 0) THEN
1094 : WRITE (output_unit, '(/,A)') &
1095 191 : " ---------------------------------- CDFT --------------------------------------"
1096 : WRITE (output_unit, '(A)') &
1097 191 : " Optimizing a density constraint in an external SCF loop "
1098 191 : WRITE (output_unit, '(A)') " "
1099 206 : SELECT CASE (cdft_control%type)
1100 : CASE (outer_scf_hirshfeld_constraint)
1101 15 : WRITE (output_unit, '(A)') " Type of constraint: Hirshfeld"
1102 : CASE (outer_scf_becke_constraint)
1103 191 : WRITE (output_unit, '(A)') " Type of constraint: Becke"
1104 : END SELECT
1105 191 : WRITE (output_unit, '(A,I8)') " Number of constraints: ", SIZE(cdft_control%group)
1106 191 : WRITE (output_unit, '(A,L8)') " Using fragment densities:", cdft_control%fragment_density
1107 191 : WRITE (output_unit, '(A)') " "
1108 191 : IF (cdft_control%atomic_charges) WRITE (output_unit, '(A,/)') " Calculating atomic CDFT charges"
1109 191 : SELECT CASE (cdft_control%constraint_control%optimizer)
1110 : CASE (outer_scf_optimizer_sd)
1111 : WRITE (output_unit, '(A)') &
1112 0 : " Minimizer : SD : steepest descent"
1113 : CASE (outer_scf_optimizer_diis)
1114 : WRITE (output_unit, '(A)') &
1115 15 : " Minimizer : DIIS : direct inversion"
1116 : WRITE (output_unit, '(A)') &
1117 15 : " in the iterative subspace"
1118 : WRITE (output_unit, '(A,I3,A)') &
1119 15 : " using ", &
1120 30 : cdft_control%constraint_control%diis_buffer_length, " DIIS vectors"
1121 : CASE (outer_scf_optimizer_bisect)
1122 : WRITE (output_unit, '(A)') &
1123 115 : " Minimizer : BISECT : gradient bisection"
1124 : WRITE (output_unit, '(A,I3)') &
1125 115 : " using a trust count of", &
1126 230 : cdft_control%constraint_control%bisect_trust_count
1127 : CASE (outer_scf_optimizer_broyden, outer_scf_optimizer_newton, &
1128 : outer_scf_optimizer_newton_ls)
1129 : CALL cdft_opt_type_write(cdft_control%constraint_control%cdft_opt_control, &
1130 60 : cdft_control%constraint_control%optimizer, output_unit)
1131 : CASE (outer_scf_optimizer_secant)
1132 1 : WRITE (output_unit, '(A)') " Minimizer : Secant"
1133 : CASE DEFAULT
1134 191 : CPABORT("")
1135 : END SELECT
1136 : WRITE (output_unit, '(/,A,L7)') &
1137 191 : " Reusing OT preconditioner: ", cdft_control%reuse_precond
1138 191 : IF (cdft_control%reuse_precond) THEN
1139 : WRITE (output_unit, '(A,I3,A,I3,A)') &
1140 0 : " using old preconditioner for up to ", &
1141 0 : cdft_control%max_reuse, " subsequent CDFT SCF"
1142 : WRITE (output_unit, '(A,I3,A,I3,A)') &
1143 0 : " iterations if the relevant loop converged in less than ", &
1144 0 : cdft_control%precond_freq, " steps"
1145 : END IF
1146 206 : SELECT CASE (cdft_control%type)
1147 : CASE (outer_scf_hirshfeld_constraint)
1148 15 : WRITE (output_unit, '(/,A)') " Hirshfeld constraint settings"
1149 15 : WRITE (output_unit, '(A)') " "
1150 204 : SELECT CASE (cdft_control%hirshfeld_control%shape_function)
1151 : CASE (shape_function_gaussian)
1152 : WRITE (output_unit, '(A, A8)') &
1153 13 : " Shape function type: ", "Gaussian"
1154 : WRITE (output_unit, '(A)', ADVANCE='NO') &
1155 13 : " Type of Gaussian: "
1156 17 : SELECT CASE (cdft_control%hirshfeld_control%gaussian_shape)
1157 : CASE (radius_default)
1158 2 : WRITE (output_unit, '(A13)') "Default"
1159 : CASE (radius_covalent)
1160 11 : WRITE (output_unit, '(A13)') "Covalent"
1161 : CASE (radius_single)
1162 0 : WRITE (output_unit, '(A13)') "Fixed radius"
1163 : CASE (radius_vdw)
1164 0 : WRITE (output_unit, '(A13)') "Van der Waals"
1165 : CASE (radius_user)
1166 13 : WRITE (output_unit, '(A13)') "User-defined"
1167 :
1168 : END SELECT
1169 : CASE (shape_function_density)
1170 : WRITE (output_unit, '(A, A8)') &
1171 15 : " Shape function type: ", "Density"
1172 : END SELECT
1173 : CASE (outer_scf_becke_constraint)
1174 176 : WRITE (output_unit, '(/, A)') " Becke constraint settings"
1175 176 : WRITE (output_unit, '(A)') " "
1176 283 : SELECT CASE (cdft_control%becke_control%cutoff_type)
1177 : CASE (becke_cutoff_global)
1178 : WRITE (output_unit, '(A,F8.3,A)') &
1179 107 : " Cutoff for partitioning :", cp_unit_from_cp2k(cdft_control%becke_control%rglobal, &
1180 214 : "angstrom"), " angstrom"
1181 : CASE (becke_cutoff_element)
1182 : WRITE (output_unit, '(A)') &
1183 176 : " Using element specific cutoffs for partitioning"
1184 : END SELECT
1185 : WRITE (output_unit, '(A,L7)') &
1186 176 : " Skipping distant gpoints: ", cdft_control%becke_control%should_skip
1187 : WRITE (output_unit, '(A,L7)') &
1188 176 : " Precompute gradients : ", cdft_control%becke_control%in_memory
1189 176 : WRITE (output_unit, '(A)') " "
1190 176 : IF (cdft_control%becke_control%adjust) &
1191 : WRITE (output_unit, '(A)') &
1192 110 : " Using atomic radii to generate a heteronuclear charge partitioning"
1193 176 : WRITE (output_unit, '(A)') " "
1194 367 : IF (.NOT. cdft_control%becke_control%cavity_confine) THEN
1195 : WRITE (output_unit, '(A)') &
1196 19 : " No confinement is active"
1197 : ELSE
1198 157 : WRITE (output_unit, '(A)') " Confinement using a Gaussian shaped cavity is active"
1199 158 : SELECT CASE (cdft_control%becke_control%cavity_shape)
1200 : CASE (radius_single)
1201 : WRITE (output_unit, '(A,F8.4, A)') &
1202 1 : " Type of Gaussian : Fixed radius: ", &
1203 2 : cp_unit_from_cp2k(cdft_control%becke_control%rcavity, "angstrom"), " angstrom"
1204 : CASE (radius_covalent)
1205 : WRITE (output_unit, '(A)') &
1206 1 : " Type of Gaussian : Covalent radius "
1207 : CASE (radius_vdw)
1208 : WRITE (output_unit, '(A)') &
1209 154 : " Type of Gaussian : vdW radius "
1210 : CASE (radius_user)
1211 : WRITE (output_unit, '(A)') &
1212 157 : " Type of Gaussian : User radius "
1213 : END SELECT
1214 : WRITE (output_unit, '(A,ES12.4)') &
1215 157 : " Cavity threshold : ", cdft_control%becke_control%eps_cavity
1216 : END IF
1217 : END SELECT
1218 : WRITE (output_unit, '(/,A)') &
1219 191 : " ---------------------------------- CDFT --------------------------------------"
1220 : END IF
1221 :
1222 191 : END SUBROUTINE qs_scf_cdft_initial_info
1223 :
1224 : ! **************************************************************************************************
1225 : !> \brief writes CDFT constraint information
1226 : !> \param output_unit where to write the information
1227 : !> \param cdft_control the env which holds information about the constraint
1228 : !> \par History
1229 : !> 08.2018 separated from qs_scf_cdft_info to make code callable elsewhere [Nico Holmberg]
1230 : ! **************************************************************************************************
1231 4080 : SUBROUTINE qs_scf_cdft_constraint_info(output_unit, cdft_control)
1232 : INTEGER :: output_unit
1233 : TYPE(cdft_control_type), POINTER :: cdft_control
1234 :
1235 : INTEGER :: igroup
1236 :
1237 4080 : IF (output_unit > 0) THEN
1238 2203 : SELECT CASE (cdft_control%type)
1239 : CASE (outer_scf_hirshfeld_constraint)
1240 : WRITE (output_unit, '(/,T3,A,T60)') &
1241 61 : '------------------- Hirshfeld constraint information -------------------'
1242 : CASE (outer_scf_becke_constraint)
1243 : WRITE (output_unit, '(/,T3,A,T60)') &
1244 2081 : '--------------------- Becke constraint information ---------------------'
1245 : CASE DEFAULT
1246 2142 : CPABORT("Unknown CDFT constraint.")
1247 : END SELECT
1248 4839 : DO igroup = 1, SIZE(cdft_control%target)
1249 2697 : IF (igroup > 1) WRITE (output_unit, '(T3,A)') ' '
1250 : WRITE (output_unit, '(T3,A,T54,(3X,I18))') &
1251 2697 : 'Atomic group :', igroup
1252 4284 : SELECT CASE (cdft_control%group(igroup)%constraint_type)
1253 : CASE (cdft_charge_constraint)
1254 1587 : IF (cdft_control%group(igroup)%is_fragment_constraint) THEN
1255 : WRITE (output_unit, '(T3,A,T42,A)') &
1256 6 : 'Type of constraint :', ADJUSTR('Charge density constraint (frag.)')
1257 : ELSE
1258 : WRITE (output_unit, '(T3,A,T50,A)') &
1259 1581 : 'Type of constraint :', ADJUSTR('Charge density constraint')
1260 : END IF
1261 : CASE (cdft_magnetization_constraint)
1262 8 : IF (cdft_control%group(igroup)%is_fragment_constraint) THEN
1263 : WRITE (output_unit, '(T3,A,T35,A)') &
1264 6 : 'Type of constraint :', ADJUSTR('Magnetization density constraint (frag.)')
1265 : ELSE
1266 : WRITE (output_unit, '(T3,A,T43,A)') &
1267 2 : 'Type of constraint :', ADJUSTR('Magnetization density constraint')
1268 : END IF
1269 : CASE (cdft_alpha_constraint)
1270 551 : IF (cdft_control%group(igroup)%is_fragment_constraint) THEN
1271 : WRITE (output_unit, '(T3,A,T38,A)') &
1272 0 : 'Type of constraint :', ADJUSTR('Alpha spin density constraint (frag.)')
1273 : ELSE
1274 : WRITE (output_unit, '(T3,A,T46,A)') &
1275 551 : 'Type of constraint :', ADJUSTR('Alpha spin density constraint')
1276 : END IF
1277 : CASE (cdft_beta_constraint)
1278 551 : IF (cdft_control%group(igroup)%is_fragment_constraint) THEN
1279 : WRITE (output_unit, '(T3,A,T39,A)') &
1280 0 : 'Type of constraint :', ADJUSTR('Beta spin density constraint (frag.)')
1281 : ELSE
1282 : WRITE (output_unit, '(T3,A,T47,A)') &
1283 551 : 'Type of constraint :', ADJUSTR('Beta spin density constraint')
1284 : END IF
1285 : CASE DEFAULT
1286 2697 : CPABORT("Unknown constraint type.")
1287 : END SELECT
1288 : WRITE (output_unit, '(T3,A,T54,(3X,F18.12))') &
1289 2697 : 'Target value of constraint :', cdft_control%target(igroup)
1290 : WRITE (output_unit, '(T3,A,T54,(3X,F18.12))') &
1291 2697 : 'Current value of constraint :', cdft_control%value(igroup)
1292 : WRITE (output_unit, '(T3,A,T59,(3X,ES13.3))') &
1293 2697 : 'Deviation from target :', cdft_control%value(igroup) - cdft_control%target(igroup)
1294 : WRITE (output_unit, '(T3,A,T54,(3X,F18.12))') &
1295 4839 : 'Strength of constraint :', cdft_control%strength(igroup)
1296 : END DO
1297 : WRITE (output_unit, '(T3,A)') &
1298 2142 : '------------------------------------------------------------------------'
1299 : END IF
1300 :
1301 4080 : END SUBROUTINE qs_scf_cdft_constraint_info
1302 :
1303 : END MODULE qs_scf_output
|