(git:1155b05)
Loading...
Searching...
No Matches
qs_scf_output.F
Go to the documentation of this file.
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
9 USE admm_types, ONLY: admm_type
14 USE cp_dbcsr_api, ONLY: dbcsr_p_type,&
24 USE cp_output_handling, ONLY: cp_p_file,&
29 USE input_constants, ONLY: &
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
49 USE physcon, ONLY: evolt,&
52 USE ps_implicit_types, ONLY: mixed_bc,&
56 USE pw_env_types, ONLY: pw_env_type
71 USE qs_mo_types, ONLY: allocate_mo_set,&
77 USE qs_rho_types, ONLY: qs_rho_get,&
80 USE qs_scf_types, ONLY: ot_method_nr,&
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, &
101
102CONTAINS
103
104! **************************************************************************************************
105!> \brief writes a summary of information after scf
106!> \param output_unit ...
107!> \param qs_env ...
108! **************************************************************************************************
109 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 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 scf_env=scf_env, qs_charges=qs_charges)
124
125 gapw = dft_control%qs_control%gapw
126 gapw_xc = dft_control%qs_control%gapw_xc
127 qmmm = qs_env%qmmm
128 nelectron_total = scf_env%nelectron
129
130 CALL qs_scf_print_scf_summary(output_unit, rho, qs_charges, energy, nelectron_total, &
131 dft_control, qmmm, qs_env, gapw, gapw_xc)
132
133 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 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 CALL timeset(routinen, handle)
154
155 IF (output_unit > 0) THEN
156 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 nmo=nmo)
162 IF (dft_control%nspins > 1) THEN
163 WRITE (unit=output_unit, fmt="(/,T2,A,I2)") "Spin", ispin
164 END IF
165 WRITE (unit=output_unit, fmt="(/,(T2,A,T71,I10))") &
166 "Number of electrons:", nelectron_spin, &
167 "Number of occupied orbitals:", homo, &
168 "Number of molecular orbitals:", nmo
169 END DO
170 WRITE (unit=output_unit, fmt="(/,(T2,A,T71,I10))") &
171 "Number of orbital functions:", nao, &
172 "Number of independent orbital functions:", nao - ndep
173 END IF
174
175 CALL timestop(handle)
176
177 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 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 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 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 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 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 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 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
222 TYPE(preconditioner_type), POINTER :: local_preconditioner
223 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
224 TYPE(scf_control_type), POINTER :: scf_control
225 TYPE(section_vals_type), POINTER :: dft_section, input
226
227 CALL timeset(routinen, handle)
228
229 cpassert(ASSOCIATED(qs_env))
230
231 ! Retrieve the required information for the requested print output
232 CALL get_qs_env(qs_env, &
233 blacs_env=blacs_env, &
234 dft_control=dft_control, &
235 do_kpoints=do_kpoints, &
236 input=input, &
237 qs_kind_set=qs_kind_set, &
238 para_env=para_env, &
239 particle_set=particle_set, &
240 scf_control=scf_control)
241
242 ! Quick return, if no printout of MO information is requested
243 dft_section => section_vals_get_subs_vals(input, "DFT")
244 CALL section_vals_val_get(dft_section, "PRINT%MO%EIGENVALUES", l_val=print_eigvals)
245 CALL section_vals_val_get(dft_section, "PRINT%MO%EIGENVECTORS", l_val=print_eigvecs)
246 CALL section_vals_val_get(dft_section, "PRINT%MO%OCCUPATION_NUMBERS", l_val=print_occup)
247 CALL section_vals_val_get(dft_section, "PRINT%MO%OCCUPATION_NUMBERS_STATS", c_vals=tmpstringlist)
248
249 print_occup_stats = .false.
250 occup_stats_occ_threshold = 1e-6_dp
251 IF (SIZE(tmpstringlist) > 0) THEN ! the lone_keyword_c_vals doesn't work as advertised, handle it manually
252 print_occup_stats = .true.
253 IF (len_trim(tmpstringlist(1)) > 0) &
254 READ (tmpstringlist(1), *) print_occup_stats
255 END IF
256 IF (SIZE(tmpstringlist) > 1) &
257 READ (tmpstringlist(2), *) occup_stats_occ_threshold
258
259 logger => cp_get_default_logger()
260 print_mo_info = (cp_print_key_should_output(logger%iter_info, dft_section, "PRINT%MO") /= 0)
261
262 IF ((.NOT. print_mo_info) .OR. (.NOT. (print_eigvals .OR. print_eigvecs .OR. print_occup .OR. print_occup_stats))) THEN
263 CALL timestop(handle)
264 RETURN
265 END IF
266
267 NULLIFY (fm_struct_tmp)
268 NULLIFY (mo_coeff)
269 NULLIFY (mo_coeff_deriv)
270 NULLIFY (mo_eigenvalues)
271 NULLIFY (mo_set)
272 NULLIFY (umo_coeff)
273 NULLIFY (umo_eigenvalues)
274 NULLIFY (umo_set)
275
276 do_printout = .true.
277 nspin = dft_control%nspins
278 nmos_occ = 0
279
280 ! Check, if we have k points
281 IF (do_kpoints) THEN
282 CALL get_qs_env(qs_env, kpoints=kpoints)
283 nkp = SIZE(kpoints%kp_env)
284 ELSE
285 CALL get_qs_env(qs_env, matrix_ks=ks, matrix_s=s)
286 cpassert(ASSOCIATED(ks))
287 cpassert(ASSOCIATED(s))
288 nkp = 1
289 END IF
290
291 kp_loop: DO ikp = 1, nkp
292
293 IF (do_kpoints) THEN
294 mos => kpoints%kp_env(ikp)%kpoint_env%mos(1, :)
295 kpoint = ikp
296 ELSE
297 CALL get_qs_env(qs_env, matrix_ks=ks, mos=mos)
298 kpoint = 0 ! Gamma point only
299 END IF
300 cpassert(ASSOCIATED(mos))
301
302 ! Prepare MO information for printout
303 DO ispin = 1, nspin
304
305 ! Calculate MO eigenvalues and eigenvector when OT is used
306 IF (scf_env%method == ot_method_nr) THEN
307
308 solver_method = "OT"
309
310 IF (do_kpoints) THEN
311 cpabort("The OT method is not implemented for k points")
312 END IF
313
314 IF (final_mos) THEN
315
316 matrix_ks => ks(ispin)%matrix
317 matrix_s => s(1)%matrix
318
319 ! With ADMM, we have to modify the Kohn-Sham matrix
320 IF (dft_control%do_admm) THEN
321 CALL get_qs_env(qs_env, admm_env=admm_env)
322 CALL admm_correct_for_eigenvalues(ispin, admm_env, matrix_ks)
323 END IF
324
325 mo_set => mos(ispin)
326 CALL get_mo_set(mo_set=mo_set, &
327 mo_coeff=mo_coeff, &
328 eigenvalues=mo_eigenvalues, &
329 homo=homo, &
330 maxocc=maxocc, &
331 nelectron=nelectron, &
332 n_el_f=n_el_f, &
333 nao=nao, &
334 nmo=nmo, &
335 flexible_electron_count=flexible_electron_count)
336
337 IF (ASSOCIATED(qs_env%mo_derivs)) THEN
338 mo_coeff_deriv => qs_env%mo_derivs(ispin)%matrix
339 ELSE
340 mo_coeff_deriv => null()
341 END IF
342
343 ! Update the eigenvalues of the occupied orbitals
344 CALL calculate_subspace_eigenvalues(orbitals=mo_coeff, &
345 ks_matrix=matrix_ks, &
346 evals_arg=mo_eigenvalues, &
347 co_rotate_dbcsr=mo_coeff_deriv)
348 CALL set_mo_occupation(mo_set=mo_set)
349
350 ! Retrieve the index of the last MO for which a printout is requested
351 mo_index_range => section_get_ivals(dft_section, "PRINT%MO%MO_INDEX_RANGE")
352 cpassert(ASSOCIATED(mo_index_range))
353 IF (mo_index_range(2) < 0) THEN
354 numo = nao - homo
355 ELSE
356 numo = min(mo_index_range(2) - homo, nao - homo)
357 END IF
358
359 ! Calculate the unoccupied MO set (umo_set) with OT if needed
360 IF (numo > 0) THEN
361
362 ! Create temporary virtual MO set for printout
363 CALL cp_fm_struct_create(fm_struct_tmp, &
364 context=blacs_env, &
365 para_env=para_env, &
366 nrow_global=nao, &
367 ncol_global=numo)
368 ALLOCATE (umo_set)
369 CALL allocate_mo_set(mo_set=umo_set, &
370 nao=nao, &
371 nmo=numo, &
372 nelectron=0, &
373 n_el_f=n_el_f, &
374 maxocc=maxocc, &
375 flexible_electron_count=flexible_electron_count)
376 CALL init_mo_set(mo_set=umo_set, &
377 fm_struct=fm_struct_tmp, &
378 name="Temporary MO set (unoccupied MOs only) for printout")
379 CALL cp_fm_struct_release(fm_struct_tmp)
380 CALL get_mo_set(mo_set=umo_set, &
381 mo_coeff=umo_coeff, &
382 eigenvalues=umo_eigenvalues)
383
384 ! Prepare printout of the additional unoccupied MOs when OT is being employed
385 CALL cp_fm_init_random(umo_coeff)
386
387 ! The FULL_ALL preconditioner makes not much sense for the unoccupied orbitals
388 NULLIFY (local_preconditioner)
389 IF (ASSOCIATED(scf_env%ot_preconditioner)) THEN
390 local_preconditioner => scf_env%ot_preconditioner(1)%preconditioner
391 IF (local_preconditioner%in_use == ot_precond_full_all) THEN
392 NULLIFY (local_preconditioner)
393 END IF
394 END IF
395
396 ! Calculate the MO information for the request MO index range
397 CALL ot_eigensolver(matrix_h=matrix_ks, &
398 matrix_s=matrix_s, &
399 matrix_c_fm=umo_coeff, &
400 matrix_orthogonal_space_fm=mo_coeff, &
401 eps_gradient=scf_control%eps_lumos, &
402 preconditioner=local_preconditioner, &
403 iter_max=scf_control%max_iter_lumos, &
404 size_ortho_space=nmo)
405
406 CALL calculate_subspace_eigenvalues(orbitals=umo_coeff, &
407 ks_matrix=matrix_ks, &
408 evals_arg=umo_eigenvalues)
409 CALL set_mo_occupation(mo_set=umo_set)
410
411 END IF ! numo > 0
412
413 ! With ADMM, we have to undo the modification of the Kohn-Sham matrix
414 IF (dft_control%do_admm) THEN
415 CALL admm_uncorrect_for_eigenvalues(ispin, admm_env, matrix_ks)
416 END IF
417
418 ELSE
419
420 message = "The MO information is only calculated after SCF convergence "// &
421 "is achieved when the orbital transformation (OT) method is used"
422 cpwarn(trim(message))
423 do_printout = .false.
424 EXIT kp_loop
425
426 END IF ! final MOs
427
428 ELSE
429
430 solver_method = "TD"
431 mo_set => mos(ispin)
432 NULLIFY (umo_set)
433
434 END IF ! OT is used
435
436 ! Print MO information
437 IF (nspin > 1) THEN
438 SELECT CASE (ispin)
439 CASE (1)
440 spin = "ALPHA"
441 CASE (2)
442 spin = "BETA"
443 CASE DEFAULT
444 cpabort("Invalid spin")
445 END SELECT
446 IF (ASSOCIATED(umo_set)) THEN
447 CALL write_mo_set_to_output_unit(mo_set, qs_kind_set, particle_set, dft_section, 4, kpoint, &
448 final_mos=final_mos, spin=trim(spin), solver_method=solver_method, &
449 umo_set=umo_set)
450 ELSE
451 CALL write_mo_set_to_output_unit(mo_set, qs_kind_set, particle_set, dft_section, 4, kpoint, &
452 final_mos=final_mos, spin=trim(spin), solver_method=solver_method)
453 END IF
454 ELSE
455 IF (ASSOCIATED(umo_set)) THEN
456 CALL write_mo_set_to_output_unit(mo_set, qs_kind_set, particle_set, dft_section, 4, kpoint, &
457 final_mos=final_mos, solver_method=solver_method, &
458 umo_set=umo_set)
459 ELSE
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 END IF
463 END IF
464
465 nmos_occ(ispin) = max(nmos_occ(ispin), count(mo_set%occupation_numbers > occup_stats_occ_threshold))
466
467 ! Deallocate temporary objects needed for OT
468 IF (scf_env%method == ot_method_nr) THEN
469 IF (ASSOCIATED(umo_set)) THEN
470 CALL deallocate_mo_set(umo_set)
471 DEALLOCATE (umo_set)
472 END IF
473 NULLIFY (matrix_ks)
474 NULLIFY (matrix_s)
475 END IF
476 NULLIFY (mo_set)
477
478 END DO ! ispin
479
480 END DO kp_loop
481
482 IF (do_printout .AND. print_mo_info .AND. print_occup_stats) THEN
483 iw = cp_print_key_unit_nr(logger, dft_section, "PRINT%MO", &
484 ignore_should_output=print_mo_info, &
485 extension=".MOLog")
486 IF (iw > 0) THEN
487 IF (SIZE(mos) > 1) THEN
488 WRITE (unit=iw, fmt="(A,I4)") " MO| Total occupied (ALPHA):", nmos_occ(1)
489 WRITE (unit=iw, fmt="(A,I4)") " MO| Total occupied (BETA): ", nmos_occ(2)
490 ELSE
491 WRITE (unit=iw, fmt="(A,I4)") " MO| Total occupied: ", nmos_occ(1)
492 END IF
493 WRITE (unit=iw, fmt="(A)") ""
494 END IF
495 CALL cp_print_key_finished_output(iw, logger, dft_section, "PRINT%MO", &
496 ignore_should_output=print_mo_info)
497 END IF
498
499 CALL timestop(handle)
500
501 END SUBROUTINE qs_scf_write_mos
502
503! **************************************************************************************************
504!> \brief writes basic information obtained in a scf outer loop step
505!> \param output_unit ...
506!> \param scf_control ...
507!> \param scf_env ...
508!> \param energy ...
509!> \param total_steps ...
510!> \param should_stop ...
511!> \param outer_loop_converged ...
512! **************************************************************************************************
513 SUBROUTINE qs_scf_outer_loop_info(output_unit, scf_control, scf_env, &
514 energy, total_steps, should_stop, outer_loop_converged)
515 INTEGER :: output_unit
516 TYPE(scf_control_type), POINTER :: scf_control
517 TYPE(qs_scf_env_type), POINTER :: scf_env
518 TYPE(qs_energy_type), POINTER :: energy
519 INTEGER :: total_steps
520 LOGICAL, INTENT(IN) :: should_stop, outer_loop_converged
521
522 REAL(kind=dp) :: outer_loop_eps
523
524 outer_loop_eps = sqrt(maxval(scf_env%outer_scf%gradient(:, scf_env%outer_scf%iter_count)**2))
525 IF (output_unit > 0) WRITE (output_unit, '(/,T3,A,I4,A,E10.2,A,F22.10)') &
526 "outer SCF iter = ", scf_env%outer_scf%iter_count, &
527 " RMS gradient = ", outer_loop_eps, " energy =", energy%total
528
529 IF (outer_loop_converged) THEN
530 IF (output_unit > 0) WRITE (output_unit, '(T3,A,I4,A,I4,A,/)') &
531 "outer SCF loop converged in", scf_env%outer_scf%iter_count, &
532 " iterations or ", total_steps, " steps"
533 ELSE IF (scf_env%outer_scf%iter_count > scf_control%outer_scf%max_scf &
534 .OR. should_stop) THEN
535 IF (output_unit > 0) WRITE (output_unit, '(T3,A,I4,A,I4,A,/)') &
536 "outer SCF loop FAILED to converge after ", &
537 scf_env%outer_scf%iter_count, " iterations or ", total_steps, " steps"
538 END IF
539
540 END SUBROUTINE qs_scf_outer_loop_info
541
542! **************************************************************************************************
543!> \brief writes basic information obtained in a scf step
544!> \param scf_env ...
545!> \param output_unit ...
546!> \param just_energy ...
547!> \param t1 ...
548!> \param t2 ...
549!> \param energy ...
550! **************************************************************************************************
551 SUBROUTINE qs_scf_loop_info(scf_env, output_unit, just_energy, t1, t2, energy)
552
553 TYPE(qs_scf_env_type), POINTER :: scf_env
554 INTEGER :: output_unit
555 LOGICAL :: just_energy
556 REAL(kind=dp) :: t1, t2
557 TYPE(qs_energy_type), POINTER :: energy
558
559 IF ((output_unit > 0) .AND. scf_env%print_iter_line) THEN
560 IF (just_energy) THEN
561 WRITE (unit=output_unit, &
562 fmt="(T2,A,1X,A,T20,E8.2,1X,F6.1,16X,F20.10)") &
563 " -", trim(scf_env%iter_method), scf_env%iter_param, t2 - t1, energy%total
564 ELSE
565 IF ((abs(scf_env%iter_delta) < 1.0e-8_dp) .OR. &
566 (abs(scf_env%iter_delta) >= 1.0e5_dp)) THEN
567 WRITE (unit=output_unit, &
568 fmt="(T2,I5,1X,A,T20,E8.2,1X,F6.1,1X,ES14.4,1X,F20.10,1X,ES9.2)") &
569 scf_env%iter_count, trim(scf_env%iter_method), scf_env%iter_param, &
570 t2 - t1, scf_env%iter_delta, energy%total, energy%total - energy%tot_old
571 ELSE
572 WRITE (unit=output_unit, &
573 fmt="(T2,I5,1X,A,T20,E8.2,1X,F6.1,1X,F14.8,1X,F20.10,1X,ES9.2)") &
574 scf_env%iter_count, trim(scf_env%iter_method), scf_env%iter_param, &
575 t2 - t1, scf_env%iter_delta, energy%total, energy%total - energy%tot_old
576 END IF
577 END IF
578 END IF
579
580 END SUBROUTINE qs_scf_loop_info
581
582! **************************************************************************************************
583!> \brief writes rather detailed summary of densities and energies
584!> after the SCF
585!> \param output_unit ...
586!> \param rho ...
587!> \param qs_charges ...
588!> \param energy ...
589!> \param nelectron_total ...
590!> \param dft_control ...
591!> \param qmmm ...
592!> \param qs_env ...
593!> \param gapw ...
594!> \param gapw_xc ...
595!> \par History
596!> 03.2006 created [Joost VandeVondele]
597!> 10.2019 print dipole moment [SGh]
598!> 11.2022 print SCCS results [MK]
599! **************************************************************************************************
600 SUBROUTINE qs_scf_print_scf_summary(output_unit, rho, qs_charges, energy, nelectron_total, &
601 dft_control, qmmm, qs_env, gapw, gapw_xc)
602 INTEGER, INTENT(IN) :: output_unit
603 TYPE(qs_rho_type), POINTER :: rho
604 TYPE(qs_charges_type), POINTER :: qs_charges
605 TYPE(qs_energy_type), POINTER :: energy
606 INTEGER, INTENT(IN) :: nelectron_total
607 TYPE(dft_control_type), POINTER :: dft_control
608 LOGICAL, INTENT(IN) :: qmmm
609 TYPE(qs_environment_type), POINTER :: qs_env
610 LOGICAL, INTENT(IN) :: gapw, gapw_xc
611
612 CHARACTER(LEN=*), PARAMETER :: routinen = 'qs_scf_print_scf_summary'
613
614 INTEGER :: bc, handle, ispin, psolver
615 REAL(kind=dp) :: e_extrapolated, exc1_energy, exc_energy, &
616 implicit_ps_ehartree, tot1_h, tot1_s
617 REAL(kind=dp), DIMENSION(:), POINTER :: tot_rho_r
618 TYPE(pw_env_type), POINTER :: pw_env
619 TYPE(scf_control_type), POINTER :: scf_control
620
621 NULLIFY (tot_rho_r, pw_env)
622 CALL timeset(routinen, handle)
623
624 CALL get_qs_env(qs_env=qs_env, pw_env=pw_env, scf_control=scf_control)
625 psolver = pw_env%poisson_env%parameters%solver
626
627 IF (output_unit > 0) THEN
628 CALL qs_rho_get(rho, tot_rho_r=tot_rho_r)
629 IF (.NOT. (dft_control%qs_control%semi_empirical .OR. &
630 dft_control%qs_control%xtb .OR. &
631 dft_control%qs_control%dftb)) THEN
632 WRITE (unit=output_unit, fmt="(/,(T3,A,T41,2F20.10))") &
633 "Electronic density on regular grids: ", &
634 accurate_sum(tot_rho_r), &
635 accurate_sum(tot_rho_r) + nelectron_total, &
636 "Core density on regular grids:", &
637 qs_charges%total_rho_core_rspace, &
638 qs_charges%total_rho_core_rspace + &
639 qs_charges%total_rho1_hard_nuc - &
640 REAL(nelectron_total + dft_control%charge, dp)
641
642 IF (dft_control%correct_surf_dip) THEN
643 WRITE (unit=output_unit, fmt="((T3,A,/,T3,A,T41,F20.10))") &
644 "Total dipole moment perpendicular to ", &
645 "the slab [electrons-Angstroem]: ", &
646 qs_env%surface_dipole_moment
647 END IF
648
649 IF (gapw) THEN
650 tot1_h = qs_charges%total_rho1_hard(1)
651 tot1_s = qs_charges%total_rho1_soft(1)
652 DO ispin = 2, dft_control%nspins
653 tot1_h = tot1_h + qs_charges%total_rho1_hard(ispin)
654 tot1_s = tot1_s + qs_charges%total_rho1_soft(ispin)
655 END DO
656 WRITE (unit=output_unit, fmt="((T3,A,T41,2F20.10))") &
657 "Hard and soft densities (Lebedev):", &
658 tot1_h, tot1_s
659 WRITE (unit=output_unit, fmt="(T3,A,T41,F20.10)") &
660 "Total Rho_soft + Rho1_hard - Rho1_soft (r-space): ", &
661 accurate_sum(tot_rho_r) + tot1_h - tot1_s, &
662 "Total charge density (r-space): ", &
663 accurate_sum(tot_rho_r) + tot1_h - tot1_s &
664 + qs_charges%total_rho_core_rspace &
665 + qs_charges%total_rho1_hard_nuc
666 IF (qs_charges%total_rho1_hard_nuc /= 0.0_dp) THEN
667 WRITE (unit=output_unit, fmt="(T3,A,T41,F20.10)") &
668 "Total CNEO nuc. char. den. (Lebedev): ", &
669 qs_charges%total_rho1_hard_nuc, &
670 "Total CNEO soft char. den. (Lebedev): ", &
671 qs_charges%total_rho1_soft_nuc_lebedev, &
672 "Total CNEO soft char. den. (r-space): ", &
673 qs_charges%total_rho1_soft_nuc_rspace, &
674 "Total soft Rho_e+n+0 (g-space):", &
675 qs_charges%total_rho_gspace
676 ELSE
677 WRITE (unit=output_unit, fmt="(T3,A,T41,F20.10)") &
678 "Total Rho_soft + Rho0_soft (g-space):", &
679 qs_charges%total_rho_gspace
680 END IF
681 ! only add total_rho1_hard_nuc for gapw as cneo requires gapw
682 ELSE
683 WRITE (unit=output_unit, fmt="(T3,A,T41,F20.10)") &
684 "Total charge density on r-space grids: ", &
685 accurate_sum(tot_rho_r) + &
686 qs_charges%total_rho_core_rspace, &
687 "Total charge density g-space grids: ", &
688 qs_charges%total_rho_gspace
689 END IF
690 END IF
691 IF (dft_control%qs_control%semi_empirical) THEN
692 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
693 "Core-core repulsion energy [eV]: ", energy%core_overlap*evolt, &
694 "Core Hamiltonian energy [eV]: ", energy%core*evolt, &
695 "Two-electron integral energy [eV]: ", energy%hartree*evolt, &
696 "Electronic energy [eV]: ", &
697 (energy%core + 0.5_dp*energy%hartree)*evolt
698 IF (energy%dispersion /= 0.0_dp) &
699 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
700 "Dispersion energy [eV]: ", energy%dispersion*evolt
701 ELSEIF (dft_control%qs_control%dftb) THEN
702 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
703 "Core Hamiltonian energy: ", energy%core, &
704 "Repulsive potential energy: ", energy%repulsive, &
705 "Electronic energy: ", energy%hartree, &
706 "Dispersion energy: ", energy%dispersion
707 IF (energy%dftb3 /= 0.0_dp) &
708 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
709 "DFTB3 3rd order energy: ", energy%dftb3
710 IF (energy%efield /= 0.0_dp) &
711 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
712 "Electric field interaction energy: ", energy%efield
713 ELSEIF (dft_control%qs_control%xtb) THEN
714 IF (dft_control%qs_control%xtb_control%do_tblite) THEN
715 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
716 "Core Hamiltonian energy: ", energy%core, &
717 "Repulsive potential energy: ", energy%repulsive, &
718 "Electrostatic energy: ", energy%el_stat, &
719 "Self-consistent dispersion energy: ", energy%dispersion_sc, &
720 "Non-self consistent dispersion energy: ", energy%dispersion, &
721 "Correction for halogen bonding: ", energy%xtb_xb_inter
722 ELSE
723 IF (dft_control%qs_control%xtb_control%gfn_type == 0) THEN
724 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
725 "Core Hamiltonian energy: ", energy%core, &
726 "Repulsive potential energy: ", energy%repulsive, &
727 "SRB Correction energy: ", energy%srb, &
728 "Charge equilibration energy: ", energy%eeq, &
729 "Dispersion energy: ", energy%dispersion
730 ELSEIF (dft_control%qs_control%xtb_control%gfn_type == 1) THEN
731 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
732 "Core Hamiltonian energy: ", energy%core, &
733 "Repulsive potential energy: ", energy%repulsive, &
734 "Electronic energy: ", energy%hartree, &
735 "DFTB3 3rd order energy: ", energy%dftb3, &
736 "Dispersion energy: ", energy%dispersion
737 IF (dft_control%qs_control%xtb_control%xb_interaction) &
738 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
739 "Correction for halogen bonding: ", energy%xtb_xb_inter
740 ELSEIF (dft_control%qs_control%xtb_control%gfn_type == 2) THEN
741 cpabort("gfn_typ 2 NYA")
742 ELSE
743 cpabort("invalid gfn_typ")
744 END IF
745 END IF
746 IF (dft_control%qs_control%xtb_control%do_nonbonded) &
747 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
748 "Correction for nonbonded interactions: ", energy%xtb_nonbonded
749 IF (energy%efield /= 0.0_dp) &
750 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
751 "Electric field interaction energy: ", energy%efield
752 ELSE
753 IF (dft_control%do_admm) THEN
754 exc_energy = energy%exc + energy%exc_aux_fit
755 IF (gapw .OR. gapw_xc) exc1_energy = energy%exc1 + energy%exc1_aux_fit
756 ELSE
757 exc_energy = energy%exc
758 IF (gapw .OR. gapw_xc) exc1_energy = energy%exc1
759 END IF
760
761 IF (psolver == pw_poisson_implicit) THEN
762 implicit_ps_ehartree = pw_env%poisson_env%implicit_env%ehartree
763 bc = pw_env%poisson_env%parameters%ps_implicit_params%boundary_condition
764 SELECT CASE (bc)
766 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
767 "Overlap energy of the core charge distribution:", energy%core_overlap, &
768 "Self energy of the core charge distribution: ", energy%core_self, &
769 "Core Hamiltonian energy: ", energy%core, &
770 "Hartree energy: ", implicit_ps_ehartree, &
771 "Electric enthalpy: ", energy%hartree, &
772 "Exchange-correlation energy: ", exc_energy
773 CASE (periodic_bc, neumann_bc)
774 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
775 "Overlap energy of the core charge distribution:", energy%core_overlap, &
776 "Self energy of the core charge distribution: ", energy%core_self, &
777 "Core Hamiltonian energy: ", energy%core, &
778 "Hartree energy: ", energy%hartree, &
779 "Exchange-correlation energy: ", exc_energy
780 END SELECT
781 ELSE
782 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
783 "Overlap energy of the core charge distribution:", energy%core_overlap, &
784 "Self energy of the core charge distribution: ", energy%core_self, &
785 "Core Hamiltonian energy: ", energy%core, &
786 "Hartree energy: ", energy%hartree, &
787 "Exchange-correlation energy: ", exc_energy
788 END IF
789 IF (energy%e_hartree /= 0.0_dp) &
790 WRITE (unit=output_unit, fmt="(T3,A,/,T3,A,T56,F25.14)") &
791 "Coulomb Electron-Electron Interaction Energy ", &
792 "- Already included in the total Hartree term ", energy%e_hartree
793 IF (energy%ex /= 0.0_dp) &
794 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
795 "Hartree-Fock Exchange energy: ", energy%ex
796 IF (energy%dispersion /= 0.0_dp) &
797 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
798 "Dispersion energy: ", energy%dispersion
799 IF (energy%gcp /= 0.0_dp) &
800 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
801 "gCP energy: ", energy%gcp
802 IF (energy%efield /= 0.0_dp) &
803 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
804 "Electric field interaction energy: ", energy%efield
805 IF (gapw) THEN
806 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
807 "GAPW| Exc from hard and soft atomic rho1: ", exc1_energy, &
808 "GAPW| local Eh = 1 center integrals: ", energy%hartree_1c
809 END IF
810 IF (gapw_xc) THEN
811 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
812 "GAPW_XC| Exc from hard and soft atomic rho1: ", exc1_energy
813 END IF
814 IF (energy%core_cneo /= 0.0_dp) THEN
815 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
816 "CNEO| quantum nuclear core energy: ", energy%core_cneo
817 END IF
818 END IF
819 IF (dft_control%hairy_probes .EQV. .true.) THEN
820 WRITE (unit=output_unit, fmt="((T3,A,T56,F25.14))") &
821 "Electronic entropic energy:", energy%kTS
822 WRITE (unit=output_unit, fmt="((T3,A,T56,F25.14))") &
823 "Fermi energy:", energy%efermi
824 END IF
825 IF (dft_control%smear) THEN
826 SELECT CASE (scf_control%smear%method)
828 ! kTS does not have physical meaning in these smearing methods
829 WRITE (unit=output_unit, fmt="((T3,A,T56,F25.14))") &
830 "Smearing free energy correction:", energy%kTS
831 CASE DEFAULT
832 WRITE (unit=output_unit, fmt="((T3,A,T56,F25.14))") &
833 "Electronic entropic energy:", energy%kTS
834 END SELECT
835 WRITE (unit=output_unit, fmt="((T3,A,T56,F25.14))") &
836 "Fermi energy:", energy%efermi
837 END IF
838 IF (dft_control%dft_plus_u) THEN
839 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
840 "DFT+U energy:", energy%dft_plus_u
841 END IF
842 IF (dft_control%do_sccs) THEN
843 WRITE (unit=output_unit, fmt="(A)") ""
844 CALL print_sccs_results(energy, dft_control%sccs_control, output_unit)
845 END IF
846 IF (qmmm) THEN
847 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
848 "QM/MM Electrostatic energy: ", energy%qmmm_el
849 IF (qs_env%qmmm_env_qm%image_charge) THEN
850 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
851 "QM/MM image charge energy: ", energy%image_charge
852 END IF
853 END IF
854 IF (dft_control%qs_control%mulliken_restraint) THEN
855 WRITE (unit=output_unit, fmt="(T3,A,T56,F25.14)") &
856 "Mulliken restraint energy: ", energy%mulliken
857 END IF
858 IF (dft_control%qs_control%semi_empirical) THEN
859 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
860 "Total energy [eV]: ", energy%total*evolt
861 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
862 "Atomic reference energy [eV]: ", energy%core_self*evolt, &
863 "Heat of formation [kcal/mol]: ", &
864 (energy%total + energy%core_self)*kcalmol
865 ELSE
866 WRITE (unit=output_unit, fmt="(/,(T3,A,T56,F25.14))") &
867 "Total energy: ", energy%total
868 IF (dft_control%smear) THEN
869 SELECT CASE (scf_control%smear%method)
870 CASE (smear_fermi_dirac)
871 e_extrapolated = energy%total - 0.5_dp*energy%kTS
872 WRITE (unit=output_unit, fmt="((T3,A,T56,F25.14))") &
873 "Total energy (extrapolated to T->0): ", e_extrapolated
874 CASE (smear_gaussian)
875 e_extrapolated = energy%total - 0.5_dp*energy%kTS
876 WRITE (unit=output_unit, fmt="((T3,A,T56,F25.14))") &
877 "Total energy (extrapolated to sigma->0): ", e_extrapolated
878 CASE (smear_mp, smear_mv)
879 ! Sigma->0 extrapolation does not apply to MP or MV method.
880 END SELECT
881 END IF
882 END IF
883 IF (qmmm) THEN
884 IF (qs_env%qmmm_env_qm%image_charge) THEN
885 CALL print_image_coefficients(qs_env%image_coeff, qs_env)
886 END IF
887 END IF
888 CALL m_flush(output_unit)
889 END IF
890
891 CALL timestop(handle)
892
893 END SUBROUTINE qs_scf_print_scf_summary
894
895! **************************************************************************************************
896!> \brief collects the 'heavy duty' printing tasks out of the SCF loop
897!> \param qs_env ...
898!> \param scf_env ...
899!> \param para_env ...
900!> \par History
901!> 03.2006 created [Joost VandeVondele]
902! **************************************************************************************************
903 SUBROUTINE qs_scf_loop_print(qs_env, scf_env, para_env)
904 TYPE(qs_environment_type), POINTER :: qs_env
905 TYPE(qs_scf_env_type), POINTER :: scf_env
906 TYPE(mp_para_env_type), POINTER :: para_env
907
908 CHARACTER(LEN=*), PARAMETER :: routinen = 'qs_scf_loop_print'
909
910 INTEGER :: after, handle, ic, ispin, iw
911 LOGICAL :: do_kpoints, omit_headers
912 REAL(kind=dp) :: mo_mag_max, mo_mag_min, orthonormality
913 TYPE(cp_logger_type), POINTER :: logger
914 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_ks, matrix_p, matrix_s
915 TYPE(dft_control_type), POINTER :: dft_control
916 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
917 TYPE(qs_rho_type), POINTER :: rho
918 TYPE(section_vals_type), POINTER :: dft_section, input, scf_section
919
920 logger => cp_get_default_logger()
921 CALL timeset(routinen, handle)
922
923 CALL get_qs_env(qs_env=qs_env, input=input, dft_control=dft_control, &
924 do_kpoints=do_kpoints)
925
926 dft_section => section_vals_get_subs_vals(input, "DFT")
927 scf_section => section_vals_get_subs_vals(dft_section, "SCF")
928
929 CALL section_vals_val_get(input, "DFT%PRINT%AO_MATRICES%OMIT_HEADERS", l_val=omit_headers)
930 DO ispin = 1, dft_control%nspins
931
932 IF (btest(cp_print_key_should_output(logger%iter_info, &
933 dft_section, "PRINT%AO_MATRICES/DENSITY"), cp_p_file)) THEN
934 CALL get_qs_env(qs_env, rho=rho)
935 CALL qs_rho_get(rho, rho_ao_kp=matrix_p)
936 iw = cp_print_key_unit_nr(logger, dft_section, "PRINT%AO_MATRICES/DENSITY", &
937 extension=".Log")
938 CALL section_vals_val_get(dft_section, "PRINT%AO_MATRICES%NDIGITS", i_val=after)
939 after = min(max(after, 1), 16)
940 DO ic = 1, SIZE(matrix_p, 2)
941 CALL cp_dbcsr_write_sparse_matrix(matrix_p(ispin, ic)%matrix, 4, after, qs_env, para_env, &
942 output_unit=iw, omit_headers=omit_headers)
943 END DO
944 CALL cp_print_key_finished_output(iw, logger, dft_section, &
945 "PRINT%AO_MATRICES/DENSITY")
946 END IF
947
948 IF (btest(cp_print_key_should_output(logger%iter_info, &
949 dft_section, "PRINT%AO_MATRICES/KOHN_SHAM_MATRIX"), cp_p_file)) THEN
950 iw = cp_print_key_unit_nr(logger, dft_section, "PRINT%AO_MATRICES/KOHN_SHAM_MATRIX", &
951 extension=".Log")
952 CALL section_vals_val_get(dft_section, "PRINT%AO_MATRICES%NDIGITS", i_val=after)
953 after = min(max(after, 1), 16)
954 CALL get_qs_env(qs_env=qs_env, matrix_ks_kp=matrix_ks)
955 DO ic = 1, SIZE(matrix_ks, 2)
956 IF (dft_control%qs_control%semi_empirical) THEN
957 CALL cp_dbcsr_write_sparse_matrix(matrix_ks(ispin, ic)%matrix, 4, after, qs_env, para_env, &
958 scale=evolt, output_unit=iw, omit_headers=omit_headers)
959 ELSE
960 CALL cp_dbcsr_write_sparse_matrix(matrix_ks(ispin, ic)%matrix, 4, after, qs_env, para_env, &
961 output_unit=iw, omit_headers=omit_headers)
962 END IF
963 END DO
964 CALL cp_print_key_finished_output(iw, logger, dft_section, &
965 "PRINT%AO_MATRICES/KOHN_SHAM_MATRIX")
966 END IF
967
968 END DO
969
970 IF (btest(cp_print_key_should_output(logger%iter_info, &
971 scf_section, "PRINT%MO_ORTHONORMALITY"), cp_p_file)) THEN
972 IF (do_kpoints) THEN
973 iw = cp_print_key_unit_nr(logger, scf_section, "PRINT%MO_ORTHONORMALITY", &
974 extension=".scfLog")
975 IF (iw > 0) THEN
976 WRITE (iw, '(T8,A)') &
977 " K-points: Maximum deviation from MO S-orthonormality not determined"
978 END IF
979 CALL cp_print_key_finished_output(iw, logger, scf_section, &
980 "PRINT%MO_ORTHONORMALITY")
981 ELSE
982 CALL get_qs_env(qs_env, mos=mos)
983 IF (scf_env%method == special_diag_method_nr) THEN
984 CALL calculate_orthonormality(orthonormality, mos)
985 ELSE
986 CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrix_s)
987 CALL calculate_orthonormality(orthonormality, mos, matrix_s(1, 1)%matrix)
988 END IF
989 iw = cp_print_key_unit_nr(logger, scf_section, "PRINT%MO_ORTHONORMALITY", &
990 extension=".scfLog")
991 IF (iw > 0) THEN
992 WRITE (iw, '(T8,A,T61,E20.4)') &
993 " Maximum deviation from MO S-orthonormality", orthonormality
994 END IF
995 CALL cp_print_key_finished_output(iw, logger, scf_section, &
996 "PRINT%MO_ORTHONORMALITY")
997 END IF
998 END IF
999 IF (btest(cp_print_key_should_output(logger%iter_info, &
1000 scf_section, "PRINT%MO_MAGNITUDE"), cp_p_file)) THEN
1001 IF (do_kpoints) THEN
1002 iw = cp_print_key_unit_nr(logger, scf_section, "PRINT%MO_MAGNITUDE", &
1003 extension=".scfLog")
1004 IF (iw > 0) THEN
1005 WRITE (iw, '(T8,A)') &
1006 " K-points: Minimum/Maximum MO magnitude not determined"
1007 END IF
1008 CALL cp_print_key_finished_output(iw, logger, scf_section, &
1009 "PRINT%MO_MAGNITUDE")
1010 ELSE
1011 CALL get_qs_env(qs_env, mos=mos)
1012 CALL calculate_magnitude(mos, mo_mag_min, mo_mag_max)
1013 iw = cp_print_key_unit_nr(logger, scf_section, "PRINT%MO_MAGNITUDE", &
1014 extension=".scfLog")
1015 IF (iw > 0) THEN
1016 WRITE (iw, '(T8,A,T41,2E20.4)') &
1017 " Minimum/Maximum MO magnitude ", mo_mag_min, mo_mag_max
1018 END IF
1019 CALL cp_print_key_finished_output(iw, logger, scf_section, &
1020 "PRINT%MO_MAGNITUDE")
1021 END IF
1022 END IF
1023
1024 CALL timestop(handle)
1025
1026 END SUBROUTINE qs_scf_loop_print
1027
1028! **************************************************************************************************
1029!> \brief writes CDFT constraint information and optionally CDFT scf loop info
1030!> \param output_unit where to write the information
1031!> \param scf_control settings of the SCF loop
1032!> \param scf_env the env which holds convergence data
1033!> \param cdft_control the env which holds information about the constraint
1034!> \param energy the total energy
1035!> \param total_steps the total number of performed SCF iterations
1036!> \param should_stop if the calculation should stop
1037!> \param outer_loop_converged logical which determines if the CDFT SCF loop converged
1038!> \param cdft_loop logical which determines a CDFT SCF loop is active
1039!> \par History
1040!> 12.2015 created [Nico Holmberg]
1041! **************************************************************************************************
1042 SUBROUTINE qs_scf_cdft_info(output_unit, scf_control, scf_env, cdft_control, &
1043 energy, total_steps, should_stop, outer_loop_converged, &
1044 cdft_loop)
1045 INTEGER :: output_unit
1046 TYPE(scf_control_type), POINTER :: scf_control
1047 TYPE(qs_scf_env_type), POINTER :: scf_env
1048 TYPE(cdft_control_type), POINTER :: cdft_control
1049 TYPE(qs_energy_type), POINTER :: energy
1050 INTEGER :: total_steps
1051 LOGICAL, INTENT(IN) :: should_stop, outer_loop_converged, &
1052 cdft_loop
1053
1054 REAL(kind=dp) :: outer_loop_eps
1055
1056 IF (cdft_loop) THEN
1057 outer_loop_eps = sqrt(maxval(scf_env%outer_scf%gradient(:, scf_env%outer_scf%iter_count)**2))
1058 IF (output_unit > 0) WRITE (output_unit, '(/,T3,A,I4,A,E10.2,A,F22.10)') &
1059 "CDFT SCF iter = ", scf_env%outer_scf%iter_count, &
1060 " RMS gradient = ", outer_loop_eps, " energy =", energy%total
1061 IF (outer_loop_converged) THEN
1062 IF (output_unit > 0) WRITE (output_unit, '(T3,A,I4,A,I4,A,/)') &
1063 "CDFT SCF loop converged in", scf_env%outer_scf%iter_count, &
1064 " iterations or ", total_steps, " steps"
1065 END IF
1066 IF ((scf_env%outer_scf%iter_count > scf_control%outer_scf%max_scf .OR. should_stop) &
1067 .AND. .NOT. outer_loop_converged) THEN
1068 IF (output_unit > 0) WRITE (output_unit, '(T3,A,I4,A,I4,A,/)') &
1069 "CDFT SCF loop FAILED to converge after ", &
1070 scf_env%outer_scf%iter_count, " iterations or ", total_steps, " steps"
1071 END IF
1072 END IF
1073 CALL qs_scf_cdft_constraint_info(output_unit, cdft_control)
1074
1075 END SUBROUTINE qs_scf_cdft_info
1076
1077! **************************************************************************************************
1078!> \brief writes information about the CDFT env
1079!> \param output_unit where to write the information
1080!> \param cdft_control the CDFT env that stores information about the constraint calculation
1081!> \par History
1082!> 12.2015 created [Nico Holmberg]
1083! **************************************************************************************************
1084 SUBROUTINE qs_scf_cdft_initial_info(output_unit, cdft_control)
1085 INTEGER :: output_unit
1086 TYPE(cdft_control_type), POINTER :: cdft_control
1087
1088 IF (output_unit > 0) THEN
1089 WRITE (output_unit, '(/,A)') &
1090 " ---------------------------------- CDFT --------------------------------------"
1091 WRITE (output_unit, '(A)') &
1092 " Optimizing a density constraint in an external SCF loop "
1093 WRITE (output_unit, '(A)') " "
1094 SELECT CASE (cdft_control%type)
1096 WRITE (output_unit, '(A)') " Type of constraint: Hirshfeld"
1098 WRITE (output_unit, '(A)') " Type of constraint: Becke"
1099 END SELECT
1100 WRITE (output_unit, '(A,I8)') " Number of constraints: ", SIZE(cdft_control%group)
1101 WRITE (output_unit, '(A,L8)') " Using fragment densities:", cdft_control%fragment_density
1102 WRITE (output_unit, '(A)') " "
1103 IF (cdft_control%atomic_charges) WRITE (output_unit, '(A,/)') " Calculating atomic CDFT charges"
1104 SELECT CASE (cdft_control%constraint_control%optimizer)
1106 WRITE (output_unit, '(A)') &
1107 " Minimizer : SD : steepest descent"
1109 WRITE (output_unit, '(A)') &
1110 " Minimizer : DIIS : direct inversion"
1111 WRITE (output_unit, '(A)') &
1112 " in the iterative subspace"
1113 WRITE (output_unit, '(A,I3,A)') &
1114 " using ", &
1115 cdft_control%constraint_control%diis_buffer_length, " DIIS vectors"
1117 WRITE (output_unit, '(A)') &
1118 " Minimizer : BISECT : gradient bisection"
1119 WRITE (output_unit, '(A,I3)') &
1120 " using a trust count of", &
1121 cdft_control%constraint_control%bisect_trust_count
1124 CALL cdft_opt_type_write(cdft_control%constraint_control%cdft_opt_control, &
1125 cdft_control%constraint_control%optimizer, output_unit)
1127 WRITE (output_unit, '(A)') " Minimizer : Secant"
1128 CASE DEFAULT
1129 cpabort("")
1130 END SELECT
1131 WRITE (output_unit, '(/,A,L7)') &
1132 " Reusing OT preconditioner: ", cdft_control%reuse_precond
1133 IF (cdft_control%reuse_precond) THEN
1134 WRITE (output_unit, '(A,I3,A,I3,A)') &
1135 " using old preconditioner for up to ", &
1136 cdft_control%max_reuse, " subsequent CDFT SCF"
1137 WRITE (output_unit, '(A,I3,A,I3,A)') &
1138 " iterations if the relevant loop converged in less than ", &
1139 cdft_control%precond_freq, " steps"
1140 END IF
1141 SELECT CASE (cdft_control%type)
1143 WRITE (output_unit, '(/,A)') " Hirshfeld constraint settings"
1144 WRITE (output_unit, '(A)') " "
1145 SELECT CASE (cdft_control%hirshfeld_control%shape_function)
1147 WRITE (output_unit, '(A, A8)') &
1148 " Shape function type: ", "Gaussian"
1149 WRITE (output_unit, '(A)', advance='NO') &
1150 " Type of Gaussian: "
1151 SELECT CASE (cdft_control%hirshfeld_control%gaussian_shape)
1152 CASE (radius_default)
1153 WRITE (output_unit, '(A13)') "Default"
1154 CASE (radius_covalent)
1155 WRITE (output_unit, '(A13)') "Covalent"
1156 CASE (radius_single)
1157 WRITE (output_unit, '(A13)') "Fixed radius"
1158 CASE (radius_vdw)
1159 WRITE (output_unit, '(A13)') "Van der Waals"
1160 CASE (radius_user)
1161 WRITE (output_unit, '(A13)') "User-defined"
1162
1163 END SELECT
1165 WRITE (output_unit, '(A, A8)') &
1166 " Shape function type: ", "Density"
1167 END SELECT
1169 WRITE (output_unit, '(/, A)') " Becke constraint settings"
1170 WRITE (output_unit, '(A)') " "
1171 SELECT CASE (cdft_control%becke_control%cutoff_type)
1172 CASE (becke_cutoff_global)
1173 WRITE (output_unit, '(A,F8.3,A)') &
1174 " Cutoff for partitioning :", cp_unit_from_cp2k(cdft_control%becke_control%rglobal, &
1175 "angstrom"), " angstrom"
1177 WRITE (output_unit, '(A)') &
1178 " Using element specific cutoffs for partitioning"
1179 END SELECT
1180 WRITE (output_unit, '(A,L7)') &
1181 " Skipping distant gpoints: ", cdft_control%becke_control%should_skip
1182 WRITE (output_unit, '(A,L7)') &
1183 " Precompute gradients : ", cdft_control%becke_control%in_memory
1184 WRITE (output_unit, '(A)') " "
1185 IF (cdft_control%becke_control%adjust) &
1186 WRITE (output_unit, '(A)') &
1187 " Using atomic radii to generate a heteronuclear charge partitioning"
1188 WRITE (output_unit, '(A)') " "
1189 IF (.NOT. cdft_control%becke_control%cavity_confine) THEN
1190 WRITE (output_unit, '(A)') &
1191 " No confinement is active"
1192 ELSE
1193 WRITE (output_unit, '(A)') " Confinement using a Gaussian shaped cavity is active"
1194 SELECT CASE (cdft_control%becke_control%cavity_shape)
1195 CASE (radius_single)
1196 WRITE (output_unit, '(A,F8.4, A)') &
1197 " Type of Gaussian : Fixed radius: ", &
1198 cp_unit_from_cp2k(cdft_control%becke_control%rcavity, "angstrom"), " angstrom"
1199 CASE (radius_covalent)
1200 WRITE (output_unit, '(A)') &
1201 " Type of Gaussian : Covalent radius "
1202 CASE (radius_vdw)
1203 WRITE (output_unit, '(A)') &
1204 " Type of Gaussian : vdW radius "
1205 CASE (radius_user)
1206 WRITE (output_unit, '(A)') &
1207 " Type of Gaussian : User radius "
1208 END SELECT
1209 WRITE (output_unit, '(A,ES12.4)') &
1210 " Cavity threshold : ", cdft_control%becke_control%eps_cavity
1211 END IF
1212 END SELECT
1213 WRITE (output_unit, '(/,A)') &
1214 " ---------------------------------- CDFT --------------------------------------"
1215 END IF
1216
1217 END SUBROUTINE qs_scf_cdft_initial_info
1218
1219! **************************************************************************************************
1220!> \brief writes CDFT constraint information
1221!> \param output_unit where to write the information
1222!> \param cdft_control the env which holds information about the constraint
1223!> \par History
1224!> 08.2018 separated from qs_scf_cdft_info to make code callable elsewhere [Nico Holmberg]
1225! **************************************************************************************************
1226 SUBROUTINE qs_scf_cdft_constraint_info(output_unit, cdft_control)
1227 INTEGER :: output_unit
1228 TYPE(cdft_control_type), POINTER :: cdft_control
1229
1230 INTEGER :: igroup
1231
1232 IF (output_unit > 0) THEN
1233 SELECT CASE (cdft_control%type)
1235 WRITE (output_unit, '(/,T3,A,T60)') &
1236 '------------------- Hirshfeld constraint information -------------------'
1238 WRITE (output_unit, '(/,T3,A,T60)') &
1239 '--------------------- Becke constraint information ---------------------'
1240 CASE DEFAULT
1241 cpabort("Unknown CDFT constraint.")
1242 END SELECT
1243 DO igroup = 1, SIZE(cdft_control%target)
1244 IF (igroup > 1) WRITE (output_unit, '(T3,A)') ' '
1245 WRITE (output_unit, '(T3,A,T54,(3X,I18))') &
1246 'Atomic group :', igroup
1247 SELECT CASE (cdft_control%group(igroup)%constraint_type)
1249 IF (cdft_control%group(igroup)%is_fragment_constraint) THEN
1250 WRITE (output_unit, '(T3,A,T42,A)') &
1251 'Type of constraint :', adjustr('Charge density constraint (frag.)')
1252 ELSE
1253 WRITE (output_unit, '(T3,A,T50,A)') &
1254 'Type of constraint :', adjustr('Charge density constraint')
1255 END IF
1257 IF (cdft_control%group(igroup)%is_fragment_constraint) THEN
1258 WRITE (output_unit, '(T3,A,T35,A)') &
1259 'Type of constraint :', adjustr('Magnetization density constraint (frag.)')
1260 ELSE
1261 WRITE (output_unit, '(T3,A,T43,A)') &
1262 'Type of constraint :', adjustr('Magnetization density constraint')
1263 END IF
1265 IF (cdft_control%group(igroup)%is_fragment_constraint) THEN
1266 WRITE (output_unit, '(T3,A,T38,A)') &
1267 'Type of constraint :', adjustr('Alpha spin density constraint (frag.)')
1268 ELSE
1269 WRITE (output_unit, '(T3,A,T46,A)') &
1270 'Type of constraint :', adjustr('Alpha spin density constraint')
1271 END IF
1273 IF (cdft_control%group(igroup)%is_fragment_constraint) THEN
1274 WRITE (output_unit, '(T3,A,T39,A)') &
1275 'Type of constraint :', adjustr('Beta spin density constraint (frag.)')
1276 ELSE
1277 WRITE (output_unit, '(T3,A,T47,A)') &
1278 'Type of constraint :', adjustr('Beta spin density constraint')
1279 END IF
1280 CASE DEFAULT
1281 cpabort("Unknown constraint type.")
1282 END SELECT
1283 WRITE (output_unit, '(T3,A,T54,(3X,F18.12))') &
1284 'Target value of constraint :', cdft_control%target(igroup)
1285 WRITE (output_unit, '(T3,A,T54,(3X,F18.12))') &
1286 'Current value of constraint :', cdft_control%value(igroup)
1287 WRITE (output_unit, '(T3,A,T59,(3X,ES13.3))') &
1288 'Deviation from target :', cdft_control%value(igroup) - cdft_control%target(igroup)
1289 WRITE (output_unit, '(T3,A,T54,(3X,F18.12))') &
1290 'Strength of constraint :', cdft_control%strength(igroup)
1291 END DO
1292 WRITE (output_unit, '(T3,A)') &
1293 '------------------------------------------------------------------------'
1294 END IF
1295
1296 END SUBROUTINE qs_scf_cdft_constraint_info
1297
1298END MODULE qs_scf_output
Types and set/get functions for auxiliary density matrix methods.
Definition admm_types.F:15
Contains methods used in the context of density fitting.
Definition admm_utils.F:15
subroutine, public admm_uncorrect_for_eigenvalues(ispin, admm_env, ks_matrix)
...
Definition admm_utils.F:127
subroutine, public admm_correct_for_eigenvalues(ispin, admm_env, ks_matrix)
...
Definition admm_utils.F:53
methods related to the blacs parallel environment
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
DBCSR output in CP2K.
subroutine, public cp_dbcsr_write_sparse_matrix(sparse_matrix, before, after, qs_env, para_env, first_row, last_row, first_col, last_col, scale, output_unit, omit_headers)
...
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_init_random(matrix, ncol, start_col)
fills a matrix with random numbers
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
integer, parameter, public cp_p_file
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
unit conversion facility
Definition cp_units.F:30
real(kind=dp) function, public cp_unit_from_cp2k(value, unit_str, defaults, power)
converts from the internal cp2k units to the given unit
Definition cp_units.F:1178
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public smear_fermi_dirac
integer, parameter, public radius_vdw
integer, parameter, public outer_scf_optimizer_sd
integer, parameter, public cdft_beta_constraint
integer, parameter, public cdft_magnetization_constraint
integer, parameter, public outer_scf_optimizer_bisect
integer, parameter, public outer_scf_optimizer_secant
integer, parameter, public becke_cutoff_element
integer, parameter, public radius_default
integer, parameter, public outer_scf_optimizer_broyden
integer, parameter, public cdft_charge_constraint
integer, parameter, public outer_scf_becke_constraint
integer, parameter, public radius_user
integer, parameter, public smear_gaussian
integer, parameter, public shape_function_density
integer, parameter, public outer_scf_hirshfeld_constraint
integer, parameter, public radius_covalent
integer, parameter, public smear_mv
integer, parameter, public shape_function_gaussian
integer, parameter, public radius_single
integer, parameter, public outer_scf_optimizer_newton_ls
integer, parameter, public outer_scf_optimizer_newton
integer, parameter, public smear_mp
integer, parameter, public becke_cutoff_global
integer, parameter, public cdft_alpha_constraint
integer, parameter, public outer_scf_optimizer_diis
integer, parameter, public ot_precond_full_all
objects that represent the structure of input sections and the data contained in an input section
integer function, dimension(:), pointer, public section_get_ivals(section_vals, keyword_name)
...
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
sums arrays of real/complex numbers with much reduced round-off as compared to a naive implementation...
Definition kahan_sum.F:29
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
Types and basic routines needed for a kpoint calculation.
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:136
Interface to the message passing library MPI.
Define the data structure for the particle information.
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public kcalmol
Definition physcon.F:171
real(kind=dp), parameter, public evolt
Definition physcon.F:183
types of preconditioners
computes preconditioners, and implements methods to apply them currently used in qs_ot
Types containing essential information for running implicit (iterative) Poisson solver.
integer, parameter, public neumann_bc
integer, parameter, public mixed_bc
integer, parameter, public mixed_periodic_bc
integer, parameter, public periodic_bc
container for various plainwaves related things
functions related to the poisson solver on regular grids
integer, parameter, public pw_poisson_implicit
Routines for image charge calculation within QM/MM.
subroutine, public print_image_coefficients(image_coeff, qs_env)
Print image coefficients.
Control parameters for optimizers that work with CDFT constraints.
subroutine, public cdft_opt_type_write(cdft_opt_control, optimizer, output_unit)
writes information about the CDFT optimizer object
Defines CDFT control structures.
container for information about total charges on the grids
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, mimic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
Definition and initialisation of the mo data type.
Definition qs_mo_io.F:21
subroutine, public write_mo_set_to_output_unit(mo_set, qs_kind_set, particle_set, dft_section, before, kpoint, final_mos, spin, solver_method, rtp, cpart, sim_step, umo_set)
Write MO information to output file (eigenvalues, occupation numbers, coefficients)
Definition qs_mo_io.F:1011
collects routines that perform operations directly related to MOs
subroutine, public calculate_magnitude(mo_array, mo_mag_min, mo_mag_max)
...
subroutine, public calculate_orthonormality(orthonormality, mo_array, matrix_s)
...
Set occupation of molecular orbitals.
Definition and initialisation of the mo data type.
Definition qs_mo_types.F:22
subroutine, public allocate_mo_set(mo_set, nao, nmo, nelectron, n_el_f, maxocc, flexible_electron_count)
Allocates a mo set and partially initializes it (nao,nmo,nelectron, and flexible_electron_count are v...
subroutine, public deallocate_mo_set(mo_set)
Deallocate a wavefunction data structure.
subroutine, public get_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, mo_coeff, mo_coeff_b, uniform_occupation, kts, mu, flexible_electron_count)
Get the components of a MO set data structure.
subroutine, public init_mo_set(mo_set, fm_pool, fm_ref, fm_struct, name)
initializes an allocated mo_set. eigenvalues, mo_coeff, occupation_numbers are valid only after this ...
an eigen-space solver for the generalised symmetric eigenvalue problem for sparse matrices,...
subroutine, public ot_eigensolver(matrix_h, matrix_s, matrix_orthogonal_space_fm, matrix_c_fm, preconditioner, eps_gradient, iter_max, size_ortho_space, silent, ot_settings)
...
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public qs_rho_get(rho_struct, rho_ao, rho_ao_im, rho_ao_kp, rho_ao_im_kp, rho_r, drho_r, rho_g, drho_g, tau_r, tau_g, rho_r_valid, drho_r_valid, rho_g_valid, drho_g_valid, tau_r_valid, tau_g_valid, tot_rho_r, tot_rho_g, rho_r_sccs, soft_valid, complex_rho_ao)
returns info about the density described by this object. If some representation is not available an e...
Self-consistent continuum solvation (SCCS) model implementation.
Definition qs_sccs.F:29
subroutine, public print_sccs_results(energy, sccs_control, output_unit)
Print SCCS results.
Definition qs_sccs.F:1211
subroutine, public qs_scf_print_summary(output_unit, qs_env)
writes a summary of information after scf
subroutine, public qs_scf_cdft_constraint_info(output_unit, cdft_control)
writes CDFT constraint information
subroutine, public qs_scf_loop_print(qs_env, scf_env, para_env)
collects the 'heavy duty' printing tasks out of the SCF loop
subroutine, public qs_scf_outer_loop_info(output_unit, scf_control, scf_env, energy, total_steps, should_stop, outer_loop_converged)
writes basic information obtained in a scf outer loop step
subroutine, public qs_scf_initial_info(output_unit, mos, dft_control, ndep)
writes basic information at the beginning of an scf run
subroutine, public qs_scf_write_mos(qs_env, scf_env, final_mos)
Write the MO eigenvector, eigenvalues, and occupation numbers to the output unit.
subroutine, public qs_scf_loop_info(scf_env, output_unit, just_energy, t1, t2, energy)
writes basic information obtained in a scf step
subroutine, public qs_scf_cdft_info(output_unit, scf_control, scf_env, cdft_control, energy, total_steps, should_stop, outer_loop_converged, cdft_loop)
writes CDFT constraint information and optionally CDFT scf loop info
subroutine, public qs_scf_cdft_initial_info(output_unit, cdft_control)
writes information about the CDFT env
module that contains the definitions of the scf types
integer, parameter, public ot_method_nr
integer, parameter, public special_diag_method_nr
parameters that control an scf iteration
stores some data used in wavefunction fitting
Definition admm_types.F:120
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
keeps the information about the structure of a full matrix
represent a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
Contains information about kpoints.
stores all the informations relevant to an mpi environment
contained for different pw related things
Container for information about total charges on the grids.
Provides all information about a quickstep kind.
keeps the density in various representations, keeping track of which ones are valid.