178#include "./base/base_uses.f90"
184 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_scf'
185 LOGICAL,
PRIVATE :: reuse_precond = .false.
186 LOGICAL,
PRIVATE :: used_history = .false.
202 SUBROUTINE scf(qs_env, has_converged, total_scf_steps)
204 LOGICAL,
INTENT(OUT),
OPTIONAL :: has_converged
205 INTEGER,
INTENT(OUT),
OPTIONAL :: total_scf_steps
207 INTEGER :: ihistory, max_scf_tmp, tsteps
208 LOGICAL :: converged, outer_scf_loop, should_stop
209 LOGICAL,
SAVE :: first_step_flag = .true.
210 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: gradient_history, variable_history
219 cpassert(
ASSOCIATED(qs_env))
220 IF (
PRESENT(has_converged))
THEN
221 has_converged = .false.
223 IF (
PRESENT(total_scf_steps))
THEN
226 CALL get_qs_env(qs_env, scf_env=scf_env, input=input, &
227 dft_control=dft_control, scf_control=scf_control)
228 IF (scf_control%max_scf > 0)
THEN
233 IF (.NOT.
ASSOCIATED(scf_env))
THEN
241 IF ((scf_control%density_guess ==
history_guess) .AND. (first_step_flag))
THEN
242 max_scf_tmp = scf_control%max_scf
243 scf_control%max_scf = 1
244 outer_scf_loop = scf_control%outer_scf%have_scf
245 scf_control%outer_scf%have_scf = .false.
248 IF (.NOT. dft_control%qs_control%cdft)
THEN
249 CALL scf_env_do_scf(scf_env=scf_env, scf_control=scf_control, qs_env=qs_env, &
250 converged=converged, should_stop=should_stop, total_scf_steps=tsteps)
253 CALL cdft_scf(qs_env=qs_env, should_stop=should_stop)
257 IF (
ASSOCIATED(qs_env%mp2_env)) qs_env%mp2_env%hf_fail = .NOT. converged
260 IF (scf_control%outer_scf%have_scf)
THEN
261 ihistory = scf_env%outer_scf%iter_count
262 CALL get_qs_env(qs_env, gradient_history=gradient_history, &
263 variable_history=variable_history)
265 gradient_history(:, 1) = gradient_history(:, 2)
266 gradient_history(:, 2) = scf_env%outer_scf%gradient(:, ihistory)
267 variable_history(:, 1) = variable_history(:, 2)
268 variable_history(:, 2) = scf_env%outer_scf%variables(:, ihistory)
270 IF (used_history) used_history = .false.
272 IF (
ASSOCIATED(scf_env%outer_scf%inv_jacobian))
THEN
273 scf_control%outer_scf%cdft_opt_control%ijacobian(2) = scf_control%outer_scf%cdft_opt_control%ijacobian(2) + 1
274 IF (scf_control%outer_scf%cdft_opt_control%ijacobian(2) >= &
275 scf_control%outer_scf%cdft_opt_control%jacobian_freq(2) .AND. &
276 scf_control%outer_scf%cdft_opt_control%jacobian_freq(2) > 0) &
277 scf_env%outer_scf%deallocate_jacobian = .true.
281 IF ((
ASSOCIATED(qs_env%wf_history)) .AND. &
283 (.NOT. first_step_flag)))
THEN
284 IF (.NOT. dft_control%qs_control%cdft)
THEN
285 CALL wfi_update(qs_env%wf_history, qs_env=qs_env, dt=1.0_dp)
287 IF (dft_control%qs_control%cdft_control%should_purge)
THEN
290 dft_control%qs_control%cdft_control%should_purge = .false.
292 CALL wfi_update(qs_env%wf_history, qs_env=qs_env, dt=1.0_dp)
295 ELSE IF ((scf_control%density_guess ==
history_guess) .AND. &
296 (first_step_flag))
THEN
297 scf_control%max_scf = max_scf_tmp
298 scf_control%outer_scf%have_scf = outer_scf_loop
299 first_step_flag = .false.
306 IF (.NOT. (should_stop))
THEN
315 IF (dft_control%qs_control%cdft) &
316 CALL cdft_control_cleanup(dft_control%qs_control%cdft_control)
318 IF (
PRESENT(has_converged))
THEN
319 has_converged = converged
321 IF (
PRESENT(total_scf_steps))
THEN
322 total_scf_steps = tsteps
346 SUBROUTINE scf_env_do_scf(scf_env, scf_control, qs_env, converged, should_stop, total_scf_steps)
351 LOGICAL,
INTENT(OUT) :: converged, should_stop
352 INTEGER,
INTENT(OUT) :: total_scf_steps
354 CHARACTER(LEN=*),
PARAMETER :: routinen =
'scf_env_do_scf'
356 CHARACTER(LEN=default_string_length) :: description, name
357 INTEGER :: ext_master_id, handle, handle2, i_tmp, &
358 ic, ispin, iter_count, output_unit, &
359 scf_energy_message_tag, total_steps
360 LOGICAL :: diis_step, do_kpoints, energy_only, exit_inner_loop, exit_outer_loop, &
361 inner_loop_converged, just_energy, outer_loop_converged
362 REAL(kind=
dp) :: t1, t2
363 REAL(kind=
dp),
DIMENSION(3) :: res_val_3
368 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: rho_ao_kp
372 TYPE(
mo_set_type),
DIMENSION(:),
POINTER :: mos, mos_last_converged
379 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
384 CALL timeset(routinen, handle)
386 NULLIFY (dft_control, rho, energy, &
387 logger, qs_charges, ks_env, mos, atomic_kind_set, qs_kind_set, &
388 particle_set, dft_section, input, &
389 scf_section, para_env, results, kpoints, pw_env, rho_ao_kp, mos_last_converged)
391 cpassert(
ASSOCIATED(scf_env))
392 cpassert(
ASSOCIATED(qs_env))
399 particle_set=particle_set, &
400 qs_charges=qs_charges, &
402 atomic_kind_set=atomic_kind_set, &
403 qs_kind_set=qs_kind_set, &
407 dft_control=dft_control, &
408 do_kpoints=do_kpoints, &
422 IF (output_unit > 0)
WRITE (unit=output_unit, fmt=
"(/,/,T2,A)") &
423 "SCF WAVEFUNCTION OPTIMIZATION"
426 IF (dft_control%switch_surf_dip)
THEN
427 CALL get_qs_env(qs_env, mos_last_converged=mos_last_converged)
428 DO ispin = 1, dft_control%nspins
431 IF (output_unit > 0)
WRITE (unit=output_unit, fmt=
"(/,/,T2,A)") &
432 "COPIED mos_last_converged ---> mos"
435 IF ((output_unit > 0) .AND. (.NOT. scf_control%use_ot))
THEN
436 WRITE (unit=output_unit, &
437 fmt=
"(/,T3,A,T12,A,T31,A,T39,A,T59,A,T75,A,/,T3,A)") &
438 "Step",
"Update method",
"Time",
"Convergence",
"Total energy",
"Change", &
444 res_val_3(:) = -1.0_dp
445 description =
"[EXT_SCF_ENER_COMM]"
447 CALL get_results(results, description=description, &
448 values=res_val_3, n_entries=i_tmp)
450 IF (all(res_val_3(:) <= 0.0)) &
451 CALL cp_abort(__location__, &
452 " Trying to access result ("//trim(description)// &
453 ") which is not correctly stored.")
454 CALL external_comm%set_handle(nint(res_val_3(1)))
456 ext_master_id = nint(res_val_3(2))
457 scf_energy_message_tag = nint(res_val_3(3))
461 scf_env%outer_scf%iter_count = 0
464 energy%tot_old = 0.0_dp
469 scf_section=scf_section)
472 energy_only, just_energy, exit_inner_loop)
475 dft_control%surf_dip_correct_switch = dft_control%correct_surf_dip
476 IF ((dft_control%correct_surf_dip) .AND. (scf_control%outer_scf%have_scf) .AND. &
477 (scf_env%outer_scf%iter_count > floor(scf_control%outer_scf%max_scf/2.0_dp)))
THEN
478 IF (dft_control%switch_surf_dip)
THEN
479 dft_control%surf_dip_correct_switch = .false.
480 IF (output_unit > 0)
WRITE (unit=output_unit, fmt=
"(/,/,T2,A)") &
481 "SURFACE DIPOLE CORRECTION switched off"
487 CALL timeset(routinen//
"_inner_loop", handle2)
489 IF (.NOT. just_energy) scf_env%iter_count = scf_env%iter_count + 1
490 iter_count = iter_count + 1
491 CALL cp_iterate(logger%iter_info, last=.false., iter_nr=iter_count)
493 IF (output_unit > 0)
CALL m_flush(output_unit)
495 total_steps = total_steps + 1
496 just_energy = energy_only
499 calculate_forces=.false.)
506 IF (dft_control%hairy_probes .EQV. .true.)
THEN
507 scf_control%smear%do_smear = .false.
508 CALL qs_scf_new_mos_kp(qs_env, scf_env, scf_control, diis_step, dft_control%probe)
514 IF (dft_control%hairy_probes .EQV. .true.)
THEN
515 scf_control%smear%do_smear = .false.
516 CALL qs_scf_new_mos(qs_env, scf_env, scf_control, scf_section, diis_step, energy_only, &
519 CALL qs_scf_new_mos(qs_env, scf_env, scf_control, scf_section, diis_step, energy_only)
526 IF (dft_control%qs_control%xtb_control%do_tblite)
THEN
527 cpassert(scf_env%mixing_method > 0)
539 IF (.NOT. just_energy) energy%tot_old = energy%total
542 IF (scf_energy_message_tag > 0)
THEN
543 CALL external_comm%send(energy%total, ext_master_id, scf_energy_message_tag)
547 exit_inner_loop, inner_loop_converged, output_unit)
551 IF (exit_inner_loop)
THEN
556 outer_loop_converged, exit_outer_loop)
559 IF (exit_outer_loop)
CALL cp_iterate(logger%iter_info, last=.true., iter_nr=iter_count)
569 CALL get_ks_env(ks_env=ks_env, matrix_ks=matrix_ks)
579 IF (exit_inner_loop)
THEN
580 CALL timestop(handle2)
585 scf_section,
"PRINT%ITERATION_INFO/TIME_CUMUL"),
cp_p_file)) &
589 IF (scf_env%mixing_method > 0)
THEN
590 DO ic = 1,
SIZE(rho_ao_kp, 2)
591 DO ispin = 1, dft_control%nspins
593 CALL dbcsr_copy(rho_ao_kp(ispin, ic)%matrix, scf_env%p_mix_new(ispin, ic)%matrix, name=name)
601 CALL timestop(handle2)
605 IF (.NOT. scf_control%outer_scf%have_scf)
EXIT scf_outer_loop
609 energy, total_steps, should_stop, outer_loop_converged)
612 IF (exit_outer_loop)
THEN
613 IF ((dft_control%switch_surf_dip) .AND. (outer_loop_converged) .AND. &
614 (dft_control%surf_dip_correct_switch))
THEN
615 DO ispin = 1, dft_control%nspins
618 IF (output_unit > 0)
WRITE (unit=output_unit, fmt=
"(/,/,T2,A)") &
619 "COPIED mos ---> mos_last_converged"
623 IF (exit_outer_loop)
EXIT scf_outer_loop
630 END DO scf_outer_loop
632 converged = inner_loop_converged .AND. outer_loop_converged
633 total_scf_steps = total_steps
635 IF (dft_control%qs_control%cdft) &
636 dft_control%qs_control%cdft_control%total_steps = &
637 dft_control%qs_control%cdft_control%total_steps + total_steps
639 IF (.NOT. converged)
THEN
640 IF (scf_control%ignore_convergence_failure .OR. should_stop)
THEN
641 CALL cp_warn(__location__,
"SCF run NOT converged")
643 CALL cp_abort(__location__, &
644 "SCF run NOT converged. To continue the calculation "// &
645 "regardless, please set the keyword IGNORE_CONVERGENCE_FAILURE.")
650 IF (qs_env%energy_correction)
THEN
652 ec_env%do_skip = .false.
653 IF (ec_env%skip_ec .AND. .NOT. converged) ec_env%do_skip = .true.
657 DO ispin = 1,
SIZE(mos)
658 IF (mos(ispin)%use_mo_coeff_b)
THEN
659 IF (.NOT.
ASSOCIATED(mos(ispin)%mo_coeff_b)) &
660 cpabort(
"mo_coeff_b is not allocated")
667 CALL timestop(handle)
687 CHARACTER(LEN=*),
PARAMETER :: routinen =
'init_scf_loop'
689 INTEGER :: handle, ispin, nmo, number_of_ot_envs
690 LOGICAL :: do_kpoints, do_rotation, &
691 has_unit_metric, is_full_all
693 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_ks, matrix_s
694 TYPE(
dbcsr_type),
POINTER :: orthogonality_metric
700 CALL timeset(routinen, handle)
702 NULLIFY (scf_control, matrix_s, matrix_ks, dft_control, mos, mo_coeff, kpoints)
704 cpassert(
ASSOCIATED(scf_env))
705 cpassert(
ASSOCIATED(qs_env))
708 scf_control=scf_control, &
709 dft_control=dft_control, &
710 do_kpoints=do_kpoints, &
715 DO ispin = 1,
SIZE(mos)
716 IF (mos(1)%use_mo_coeff_b)
THEN
722 DO ispin = 1, dft_control%nspins
724 IF (.NOT. scf_control%diagonalization%mom)
THEN
727 IF (dft_control%hairy_probes .EQV. .true.)
THEN
728 IF (scf_env%outer_scf%iter_count > 0)
THEN
729 scf_control%smear%do_smear = .false.
731 smear=scf_control%smear, &
732 probe=dft_control%probe)
736 smear=scf_control%smear)
741 SELECT CASE (scf_env%method)
744 cpabort(
"unknown scf method method:"//
cp_to_string(scf_env%method))
748 IF (.NOT. scf_env%skip_diis)
THEN
749 IF (.NOT.
ASSOCIATED(scf_env%scf_diis_buffer))
THEN
750 ALLOCATE (scf_env%scf_diis_buffer)
751 CALL qs_diis_b_create(scf_env%scf_diis_buffer, nbuffer=scf_control%max_diis)
757 IF (.NOT. scf_env%skip_diis)
THEN
759 IF (.NOT.
ASSOCIATED(kpoints%scf_diis_buffer))
THEN
760 ALLOCATE (kpoints%scf_diis_buffer)
765 IF (.NOT.
ASSOCIATED(scf_env%scf_diis_buffer))
THEN
766 ALLOCATE (scf_env%scf_diis_buffer)
767 CALL qs_diis_b_create(scf_env%scf_diis_buffer, nbuffer=scf_control%max_diis)
774 CALL get_qs_env(qs_env, matrix_ks=matrix_ks, matrix_s=matrix_s)
776 IF (.NOT. scf_env%skip_diis)
THEN
777 IF (.NOT.
ASSOCIATED(scf_env%scf_diis_buffer))
THEN
778 ALLOCATE (scf_env%scf_diis_buffer)
779 CALL qs_diis_b_create(scf_env%scf_diis_buffer, nbuffer=scf_control%max_diis)
785 IF (dft_control%qs_control%dftb .OR. &
786 dft_control%qs_control%xtb .OR. &
787 dft_control%qs_control%semi_empirical)
THEN
788 cpabort(
"DFTB and SE not available with OT/DIAG")
794 scf_control%diagonalization%ot_settings%preconditioner_type, &
798 scf_control%diagonalization%ot_settings%preconditioner_type, &
799 scf_control%diagonalization%ot_settings%precond_solver_type, &
800 scf_control%diagonalization%ot_settings%energy_gap, dft_control%nspins)
806 has_unit_metric=has_unit_metric, &
814 IF (scf_control%do_outer_scf_reortho)
THEN
815 IF (scf_control%outer_scf%have_scf .AND. .NOT. dft_control%restricted)
THEN
816 IF (scf_env%outer_scf%iter_count > 0)
THEN
817 DO ispin = 1, dft_control%nspins
818 CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
819 IF (has_unit_metric)
THEN
831 IF (.NOT.
ASSOCIATED(scf_env%qs_ot_env))
THEN
834 number_of_ot_envs = dft_control%nspins
835 IF (dft_control%restricted) number_of_ot_envs = 1
837 ALLOCATE (scf_env%qs_ot_env(number_of_ot_envs))
840 IF (scf_env%outer_scf%iter_count > 0)
THEN
841 IF (scf_env%iter_delta < scf_control%eps_diis)
THEN
842 scf_env%qs_ot_env(1)%settings%ot_state = 1
848 IF (scf_env%outer_scf%iter_count > 0)
THEN
849 IF (scf_env%qs_ot_env(1)%settings%ot_state == 1)
THEN
850 scf_control%max_scf = max(scf_env%qs_ot_env(1)%settings%max_scf_diis, &
856 IF (dft_control%restricted)
THEN
857 scf_env%qs_ot_env(1)%restricted = .true.
859 IF (.NOT. scf_env%qs_ot_env(1)%settings%do_rotation) &
860 CALL cp_abort(__location__, &
861 "Restricted calculation with OT requires orbital rotation. Please "// &
862 "activate the OT%ROTATION keyword!")
864 scf_env%qs_ot_env(:)%restricted = .false.
869 do_rotation = scf_env%qs_ot_env(1)%settings%do_rotation
872 IF (do_rotation .AND. is_full_all) &
873 cpabort(
'PRECONDITIONER FULL_ALL is not compatible with ROTATION.')
877 calculate_forces=.false.)
881 IF (.NOT. reuse_precond) &
883 scf_env%qs_ot_env(1)%settings%preconditioner_type, &
890 IF (has_unit_metric)
THEN
891 NULLIFY (orthogonality_metric)
893 orthogonality_metric => matrix_s(1)%matrix
896 IF (.NOT. reuse_precond) &
898 scf_env%qs_ot_env(1)%settings%preconditioner_type, &
899 scf_env%qs_ot_env(1)%settings%precond_solver_type, &
900 scf_env%qs_ot_env(1)%settings%energy_gap, dft_control%nspins, &
901 has_unit_metric=has_unit_metric, &
902 chol_type=scf_env%qs_ot_env(1)%settings%cholesky_type)
903 IF (reuse_precond) reuse_precond = .false.
905 CALL ot_scf_init(mo_array=mos, matrix_s=orthogonality_metric, &
906 broyden_adaptive_sigma=qs_env%broyden_adaptive_sigma, &
907 qs_ot_env=scf_env%qs_ot_env, matrix_ks=matrix_ks(1)%matrix)
909 SELECT CASE (scf_env%qs_ot_env(1)%settings%preconditioner_type)
912 DO ispin = 1,
SIZE(scf_env%qs_ot_env)
914 scf_env%ot_preconditioner(ispin)%preconditioner)
917 DO ispin = 1,
SIZE(scf_env%qs_ot_env)
919 scf_env%ot_preconditioner(1)%preconditioner)
922 DO ispin = 1,
SIZE(scf_env%qs_ot_env)
924 scf_env%ot_preconditioner(1)%preconditioner)
930 do_rotation = scf_env%qs_ot_env(1)%settings%do_rotation
931 DO ispin = 1,
SIZE(mos)
932 IF (.NOT. mos(ispin)%uniform_occupation)
THEN
933 cpassert(do_rotation)
939 IF (dft_control%low_spin_roks)
THEN
941 do_rotation = scf_env%qs_ot_env(1)%settings%do_rotation
942 cpassert(do_rotation)
945 CALL timestop(handle)
960 CHARACTER(len=*),
PARAMETER :: routinen =
'scf_env_cleanup'
964 CALL timeset(routinen, handle)
969 IF (
ASSOCIATED(scf_env%scf_work1_red))
THEN
972 IF (
ASSOCIATED(scf_env%scf_work2))
THEN
974 DEALLOCATE (scf_env%scf_work2)
975 NULLIFY (scf_env%scf_work2)
977 IF (
ASSOCIATED(scf_env%scf_work2_red))
THEN
979 DEALLOCATE (scf_env%scf_work2_red)
980 NULLIFY (scf_env%scf_work2_red)
982 IF (
ASSOCIATED(scf_env%ortho))
THEN
984 DEALLOCATE (scf_env%ortho)
985 NULLIFY (scf_env%ortho)
987 IF (
ASSOCIATED(scf_env%ortho_red))
THEN
989 DEALLOCATE (scf_env%ortho_red)
990 NULLIFY (scf_env%ortho_red)
992 IF (
ASSOCIATED(scf_env%ortho_m1))
THEN
994 DEALLOCATE (scf_env%ortho_m1)
995 NULLIFY (scf_env%ortho_m1)
997 IF (
ASSOCIATED(scf_env%ortho_m1_red))
THEN
999 DEALLOCATE (scf_env%ortho_m1_red)
1000 NULLIFY (scf_env%ortho_m1_red)
1003 IF (
ASSOCIATED(scf_env%ortho_dbcsr))
THEN
1006 IF (
ASSOCIATED(scf_env%buf1_dbcsr))
THEN
1009 IF (
ASSOCIATED(scf_env%buf2_dbcsr))
THEN
1013 IF (
ASSOCIATED(scf_env%p_mix_new))
THEN
1017 IF (
ASSOCIATED(scf_env%p_delta))
THEN
1022 SELECT CASE (scf_env%method)
1039 cpabort(
"unknown scf method method:"//
cp_to_string(scf_env%method))
1042 IF (
ASSOCIATED(scf_env%outer_scf%variables))
THEN
1043 DEALLOCATE (scf_env%outer_scf%variables)
1045 IF (
ASSOCIATED(scf_env%outer_scf%count))
THEN
1046 DEALLOCATE (scf_env%outer_scf%count)
1048 IF (
ASSOCIATED(scf_env%outer_scf%gradient))
THEN
1049 DEALLOCATE (scf_env%outer_scf%gradient)
1051 IF (
ASSOCIATED(scf_env%outer_scf%energy))
THEN
1052 DEALLOCATE (scf_env%outer_scf%energy)
1054 IF (
ASSOCIATED(scf_env%outer_scf%inv_jacobian) .AND. &
1055 scf_env%outer_scf%deallocate_jacobian)
THEN
1056 DEALLOCATE (scf_env%outer_scf%inv_jacobian)
1059 CALL timestop(handle)
1073 LOGICAL,
INTENT(OUT) :: should_stop
1075 CHARACTER(len=*),
PARAMETER :: routinen =
'cdft_scf'
1077 INTEGER :: handle, iatom, ispin, ivar, nmo, nvar, &
1079 LOGICAL :: cdft_loop_converged, converged, &
1080 exit_cdft_loop, first_iteration, &
1081 my_uocc, uniform_occupation
1082 REAL(kind=
dp),
DIMENSION(:),
POINTER :: mo_occupations
1085 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_s, rho_ao
1097 NULLIFY (scf_env, ks_env, energy, rho, matrix_s, rho_ao, cdft_control, logger, &
1098 dft_control, pw_env, auxbas_pw_pool, energy, ks_env, scf_env, dft_section, &
1099 input, scf_section, scf_control, mos, mo_occupations)
1102 cpassert(
ASSOCIATED(qs_env))
1103 CALL get_qs_env(qs_env, scf_env=scf_env, energy=energy, &
1104 dft_control=dft_control, scf_control=scf_control, &
1105 ks_env=ks_env, input=input)
1107 CALL timeset(routinen//
"_loop", handle)
1111 extension=
".scfLog")
1112 first_iteration = .true.
1114 cdft_control => dft_control%qs_control%cdft_control
1116 scf_env%outer_scf%iter_count = 0
1117 cdft_control%total_steps = 0
1120 IF (output_unit > 0)
THEN
1121 WRITE (unit=output_unit, fmt=
"(/,/,T2,A)") &
1122 "CDFT EXTERNAL SCF WAVEFUNCTION OPTIMIZATION"
1125 IF (cdft_control%reuse_precond)
THEN
1126 reuse_precond = .false.
1127 cdft_control%nreused = 0
1133 CALL scf_env_do_scf(scf_env=scf_env, scf_control=scf_control, qs_env=qs_env, &
1134 converged=converged, should_stop=should_stop, total_scf_steps=tsteps)
1136 IF (cdft_control%reuse_precond)
THEN
1140 IF (scf_env%outer_scf%iter_count == 1 .AND. scf_env%iter_count == 1 &
1141 .AND. cdft_control%total_steps /= 1) &
1142 cdft_control%nreused = cdft_control%nreused - 1
1144 IF (scf_env%outer_scf%iter_count == 1 .AND. scf_env%iter_count <= cdft_control%precond_freq .AND. &
1145 cdft_control%total_steps /= 1 .AND. cdft_control%nreused < cdft_control%max_reuse)
THEN
1146 reuse_precond = .true.
1147 cdft_control%nreused = cdft_control%nreused + 1
1149 reuse_precond = .false.
1150 cdft_control%nreused = 0
1154 IF (first_iteration .AND. cdft_control%purge_history)
THEN
1155 cdft_control%istep = cdft_control%istep + 1
1156 IF (scf_env%outer_scf%iter_count > 1)
THEN
1157 cdft_control%nbad_conv = cdft_control%nbad_conv + 1
1158 IF (cdft_control%nbad_conv >= cdft_control%purge_freq .AND. &
1159 cdft_control%istep >= cdft_control%purge_offset)
THEN
1160 cdft_control%nbad_conv = 0
1161 cdft_control%istep = 0
1162 cdft_control%should_purge = .true.
1166 first_iteration = .false.
1170 cdft_loop_converged, exit_cdft_loop)
1172 energy, cdft_control%total_steps, &
1173 should_stop, cdft_loop_converged, cdft_loop=.true.)
1174 IF (exit_cdft_loop)
EXIT cdft_outer_loop
1176 CALL qs_calculate_inverse_jacobian(qs_env)
1178 CALL qs_cdft_line_search(qs_env)
1183 END DO cdft_outer_loop
1185 cdft_control%ienergy = cdft_control%ienergy + 1
1188 IF (cdft_control%do_et)
THEN
1189 CALL get_qs_env(qs_env=qs_env, matrix_s=matrix_s, mos=mos)
1190 nvar =
SIZE(cdft_control%target)
1192 ALLOCATE (cdft_control%wmat(nvar))
1195 CALL dbcsr_copy(cdft_control%wmat(ivar)%matrix, matrix_s(1)%matrix, &
1196 name=
"ET_RESTRAINT_MATRIX")
1197 CALL dbcsr_set(cdft_control%wmat(ivar)%matrix, 0.0_dp)
1198 CALL integrate_v_rspace(cdft_control%group(ivar)%weight, &
1199 hmat=cdft_control%wmat(ivar), qs_env=qs_env, &
1200 calculate_forces=.false., &
1201 gapw=dft_control%qs_control%gapw)
1205 CALL dbcsr_copy(cdft_control%matrix_s%matrix, matrix_s(1)%matrix, &
1208 NULLIFY (cdft_control%mo_coeff)
1209 ALLOCATE (cdft_control%mo_coeff(dft_control%nspins))
1210 DO ispin = 1, dft_control%nspins
1211 CALL cp_fm_create(matrix=cdft_control%mo_coeff(ispin), &
1212 matrix_struct=qs_env%mos(ispin)%mo_coeff%matrix_struct, &
1213 name=
"MO_COEFF_A"//trim(adjustl(
cp_to_string(ispin)))//
"MATRIX")
1215 cdft_control%mo_coeff(ispin))
1218 IF (cdft_control%calculate_metric)
THEN
1221 ALLOCATE (cdft_control%matrix_p(dft_control%nspins))
1222 DO ispin = 1, dft_control%nspins
1223 NULLIFY (cdft_control%matrix_p(ispin)%matrix)
1225 CALL dbcsr_copy(cdft_control%matrix_p(ispin)%matrix, rho_ao(ispin)%matrix, &
1226 name=
"DENSITY MATRIX")
1230 uniform_occupation = .true.
1231 DO ispin = 1, dft_control%nspins
1232 CALL get_mo_set(mo_set=mos(ispin), uniform_occupation=my_uocc)
1233 uniform_occupation = uniform_occupation .AND. my_uocc
1235 IF (.NOT. uniform_occupation)
THEN
1236 ALLOCATE (cdft_control%occupations(dft_control%nspins))
1237 DO ispin = 1, dft_control%nspins
1240 occupation_numbers=mo_occupations)
1241 ALLOCATE (cdft_control%occupations(ispin)%array(nmo))
1242 cdft_control%occupations(ispin)%array(1:nmo) = mo_occupations(1:nmo)
1250 IF (.NOT. (cdft_control%save_pot .OR. cdft_control%transfer_pot))
THEN
1252 CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
1253 DO iatom = 1,
SIZE(cdft_control%group)
1254 CALL auxbas_pw_pool%give_back_pw(cdft_control%group(iatom)%weight)
1255 DEALLOCATE (cdft_control%group(iatom)%weight)
1257 IF (cdft_control%atomic_charges)
THEN
1258 DO iatom = 1, cdft_control%natoms
1259 CALL auxbas_pw_pool%give_back_pw(cdft_control%charge(iatom))
1261 DEALLOCATE (cdft_control%charge)
1264 cdft_control%becke_control%cavity_confine)
THEN
1265 IF (.NOT.
ASSOCIATED(cdft_control%becke_control%cavity_mat))
THEN
1266 CALL auxbas_pw_pool%give_back_pw(cdft_control%becke_control%cavity)
1268 DEALLOCATE (cdft_control%becke_control%cavity_mat)
1271 IF (
ASSOCIATED(cdft_control%hirshfeld_control%hirshfeld_env%fnorm))
THEN
1272 CALL auxbas_pw_pool%give_back_pw(cdft_control%hirshfeld_control%hirshfeld_env%fnorm)
1275 IF (
ASSOCIATED(cdft_control%charges_fragment))
DEALLOCATE (cdft_control%charges_fragment)
1276 cdft_control%need_pot = .true.
1277 cdft_control%external_control = .false.
1280 CALL timestop(handle)
1291 SUBROUTINE cdft_control_cleanup(cdft_control)
1294 IF (
ASSOCIATED(cdft_control%constraint%variables)) &
1295 DEALLOCATE (cdft_control%constraint%variables)
1296 IF (
ASSOCIATED(cdft_control%constraint%count)) &
1297 DEALLOCATE (cdft_control%constraint%count)
1298 IF (
ASSOCIATED(cdft_control%constraint%gradient)) &
1299 DEALLOCATE (cdft_control%constraint%gradient)
1300 IF (
ASSOCIATED(cdft_control%constraint%energy)) &
1301 DEALLOCATE (cdft_control%constraint%energy)
1302 IF (
ASSOCIATED(cdft_control%constraint%inv_jacobian) .AND. &
1303 cdft_control%constraint%deallocate_jacobian) &
1304 DEALLOCATE (cdft_control%constraint%inv_jacobian)
1306 END SUBROUTINE cdft_control_cleanup
1314 SUBROUTINE qs_calculate_inverse_jacobian(qs_env)
1317 CHARACTER(LEN=*),
PARAMETER :: routinen =
'qs_calculate_inverse_jacobian'
1319 CHARACTER(len=default_path_length) :: project_name
1320 INTEGER :: counter, handle, i, ispin, iter_count, &
1321 iwork, j, max_scf, nspins, nsteps, &
1322 nvar, nwork, output_unit, pwork, &
1324 LOGICAL :: converged, explicit_jacobian, &
1325 should_build, should_stop, &
1327 REAL(kind=
dp) :: inv_error, step_size
1328 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: coeff, dh, step_multiplier
1329 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: jacobian
1330 REAL(kind=
dp),
DIMENSION(:),
POINTER :: energy
1331 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: gradient, inv_jacobian
1335 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: rho_ao_kp
1337 TYPE(
mo_set_type),
ALLOCATABLE,
DIMENSION(:) :: mos_stashed
1346 NULLIFY (energy, gradient, p_rmpv, rho_ao_kp, mos, rho, &
1347 ks_env, scf_env, scf_control, dft_control, cdft_control, &
1348 inv_jacobian, para_env, tmp_logger, energy_qs)
1351 cpassert(
ASSOCIATED(qs_env))
1352 CALL get_qs_env(qs_env, scf_env=scf_env, ks_env=ks_env, &
1353 scf_control=scf_control, mos=mos, rho=rho, &
1354 dft_control=dft_control, &
1355 para_env=para_env, energy=energy_qs)
1356 explicit_jacobian = .false.
1357 should_build = .false.
1358 use_md_history = .false.
1359 iter_count = scf_env%outer_scf%iter_count
1361 IF (.NOT.
ASSOCIATED(scf_control%outer_scf%cdft_opt_control))
RETURN
1363 CALL timeset(routinen, handle)
1365 IF (scf_control%outer_scf%cdft_opt_control%jacobian_restart)
THEN
1367 should_build = .false.
1370 IF (should_build)
THEN
1371 scf_env%outer_scf%deallocate_jacobian = .false.
1373 IF (explicit_jacobian)
THEN
1375 cdft_control => dft_control%qs_control%cdft_control
1376 IF (.NOT.
ASSOCIATED(cdft_control)) &
1377 CALL cp_abort(__location__, &
1378 "Optimizers that need the explicit Jacobian can"// &
1379 " only be used together with a valid CDFT constraint.")
1381 project_name = logger%iter_info%project_name
1382 CALL create_tmp_logger(para_env, project_name,
"-JacobianInfo.out", output_unit, tmp_logger)
1384 nspins = dft_control%nspins
1385 ALLOCATE (mos_stashed(nspins))
1386 DO ispin = 1, nspins
1390 p_rmpv => rho_ao_kp(:, 1)
1392 nvar =
SIZE(scf_env%outer_scf%variables, 1)
1393 max_scf = scf_control%outer_scf%max_scf + 1
1394 ALLOCATE (gradient(nvar, max_scf))
1395 gradient = scf_env%outer_scf%gradient
1396 ALLOCATE (energy(max_scf))
1397 energy = scf_env%outer_scf%energy
1398 ALLOCATE (jacobian(nvar, nvar))
1400 nsteps = cdft_control%total_steps
1403 twork = pwork - nwork
1405 jacobian(i, :) = coeff(0)*scf_env%outer_scf%gradient(i, iter_count)
1410 IF (output_unit > 0)
THEN
1411 WRITE (output_unit, fmt=
"(A)")
" "
1412 WRITE (output_unit, fmt=
"(A)")
" #####################################"
1413 WRITE (output_unit,
'(A,I3,A,I3,A)') &
1414 " ### Constraint ", i,
" of ", nvar,
" ###"
1415 WRITE (output_unit, fmt=
"(A)")
" #####################################"
1418 DO iwork = nwork, pwork
1419 IF (iwork == 0) cycle
1420 counter = counter + 1
1421 IF (output_unit > 0)
THEN
1422 WRITE (output_unit, fmt=
"(A)")
" #####################################"
1423 WRITE (output_unit,
'(A,I3,A,I3,A)') &
1424 " ### Energy evaluation ", counter,
" of ", twork,
" ###"
1425 WRITE (output_unit, fmt=
"(A)")
" #####################################"
1427 IF (
SIZE(scf_control%outer_scf%cdft_opt_control%jacobian_step) == 1)
THEN
1428 step_size = scf_control%outer_scf%cdft_opt_control%jacobian_step(1)
1430 step_size = scf_control%outer_scf%cdft_opt_control%jacobian_step(i)
1432 scf_env%outer_scf%variables(:, iter_count + 1) = scf_env%outer_scf%variables(:, iter_count)
1433 scf_env%outer_scf%variables(i, iter_count + 1) = scf_env%outer_scf%variables(i, iter_count) + &
1434 step_multiplier(iwork)*step_size
1438 CALL scf_env_do_scf(scf_env=scf_env, scf_control=scf_control, qs_env=qs_env, &
1439 converged=converged, should_stop=should_stop, total_scf_steps=tsteps)
1442 scf_env%outer_scf%iter_count = scf_env%outer_scf%iter_count + 1
1445 energy_qs, cdft_control%total_steps, &
1446 should_stop=.false., outer_loop_converged=.false., cdft_loop=.false.)
1447 scf_env%outer_scf%iter_count = scf_env%outer_scf%iter_count - 1
1450 jacobian(j, i) = jacobian(j, i) + coeff(iwork)*scf_env%outer_scf%gradient(j, iter_count + 1)
1453 scf_env%outer_scf%variables(:, iter_count + 1) = 0.0_dp
1454 scf_env%outer_scf%gradient = gradient
1455 scf_env%outer_scf%energy = energy
1456 cdft_control%total_steps = nsteps
1457 DO ispin = 1, nspins
1461 p_rmpv(ispin)%matrix)
1472 jacobian(i, j) = jacobian(i, j)/dh(j)
1475 IF (.NOT.
ASSOCIATED(scf_env%outer_scf%inv_jacobian)) &
1476 ALLOCATE (scf_env%outer_scf%inv_jacobian(nvar, nvar))
1477 inv_jacobian => scf_env%outer_scf%inv_jacobian
1479 scf_control%outer_scf%cdft_opt_control%broyden_update = .false.
1481 DO ispin = 1, nspins
1484 DEALLOCATE (mos_stashed, jacobian, gradient, energy, coeff, step_multiplier, dh)
1485 IF (output_unit > 0)
THEN
1486 WRITE (output_unit, fmt=
"(/,A)") &
1487 " ================================== JACOBIAN CALCULATED =================================="
1495 IF (
ASSOCIATED(scf_env%outer_scf%inv_jacobian) .AND. para_env%is_source())
THEN
1500 scf_control%outer_scf%cdft_opt_control%ijacobian(1) = scf_control%outer_scf%cdft_opt_control%ijacobian(1) + 1
1501 CALL timestop(handle)
1503 END SUBROUTINE qs_calculate_inverse_jacobian
1513 SUBROUTINE qs_cdft_line_search(qs_env)
1516 CHARACTER(LEN=*),
PARAMETER :: routinen =
'qs_cdft_line_search'
1518 CHARACTER(len=default_path_length) :: project_name
1519 INTEGER :: handle, i, ispin, iter_count, &
1520 max_linesearch, max_scf, nspins, &
1521 nsteps, nvar, output_unit, tsteps
1522 LOGICAL :: continue_ls, continue_ls_exit, converged, do_linesearch, found_solution, &
1523 reached_maxls, should_exit, should_stop, sign_changed
1524 LOGICAL,
ALLOCATABLE,
DIMENSION(:) :: positive_sign
1525 REAL(kind=
dp) :: alpha, alpha_ls, factor, norm_ls
1526 REAL(kind=
dp),
DIMENSION(:),
POINTER :: energy
1527 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: gradient, inv_jacobian
1528 REAL(kind=
dp),
EXTERNAL :: dnrm2
1532 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: rho_ao_kp
1542 CALL timeset(routinen, handle)
1544 NULLIFY (energy, gradient, p_rmpv, rho_ao_kp, mos, rho, &
1545 ks_env, scf_env, scf_control, dft_control, &
1546 cdft_control, inv_jacobian, para_env, &
1547 tmp_logger, energy_qs)
1550 cpassert(
ASSOCIATED(qs_env))
1551 CALL get_qs_env(qs_env, scf_env=scf_env, ks_env=ks_env, &
1552 scf_control=scf_control, mos=mos, rho=rho, &
1553 dft_control=dft_control, &
1554 para_env=para_env, energy=energy_qs)
1555 do_linesearch = .false.
1556 SELECT CASE (scf_control%outer_scf%optimizer)
1558 do_linesearch = .false.
1560 do_linesearch = .true.
1562 SELECT CASE (scf_control%outer_scf%cdft_opt_control%broyden_type)
1564 do_linesearch = .false.
1566 cdft_control => dft_control%qs_control%cdft_control
1567 IF (.NOT.
ASSOCIATED(cdft_control)) &
1568 CALL cp_abort(__location__, &
1569 "Optimizers that perform a line search can"// &
1570 " only be used together with a valid CDFT constraint")
1571 IF (
ASSOCIATED(scf_env%outer_scf%inv_jacobian)) &
1572 do_linesearch = .true.
1575 IF (do_linesearch)
THEN
1577 TYPE(
mo_set_type),
DIMENSION(:),
ALLOCATABLE :: mos_ls, mos_stashed
1578 cdft_control => dft_control%qs_control%cdft_control
1579 IF (.NOT.
ASSOCIATED(cdft_control)) &
1580 CALL cp_abort(__location__, &
1581 "Optimizers that perform a line search can"// &
1582 " only be used together with a valid CDFT constraint")
1583 cpassert(
ASSOCIATED(scf_env%outer_scf%inv_jacobian))
1584 cpassert(
ASSOCIATED(scf_control%outer_scf%cdft_opt_control))
1585 alpha = scf_control%outer_scf%cdft_opt_control%newton_step_save
1586 iter_count = scf_env%outer_scf%iter_count
1588 project_name = logger%iter_info%project_name
1589 CALL create_tmp_logger(para_env, project_name,
"-LineSearch.out", output_unit, tmp_logger)
1591 nspins = dft_control%nspins
1592 ALLOCATE (mos_stashed(nspins))
1593 DO ispin = 1, nspins
1597 p_rmpv => rho_ao_kp(:, 1)
1598 nsteps = cdft_control%total_steps
1600 nvar =
SIZE(scf_env%outer_scf%variables, 1)
1601 max_scf = scf_control%outer_scf%max_scf + 1
1602 max_linesearch = scf_control%outer_scf%cdft_opt_control%max_ls
1603 continue_ls = scf_control%outer_scf%cdft_opt_control%continue_ls
1604 factor = scf_control%outer_scf%cdft_opt_control%factor_ls
1605 continue_ls_exit = .false.
1606 found_solution = .false.
1607 ALLOCATE (gradient(nvar, max_scf))
1608 gradient = scf_env%outer_scf%gradient
1609 ALLOCATE (energy(max_scf))
1610 energy = scf_env%outer_scf%energy
1611 reached_maxls = .false.
1613 IF (scf_control%outer_scf%cdft_opt_control%broyden_update)
THEN
1616 scf_env%outer_scf%variables(:, iter_count + 1) = 0
1617 scf_control%outer_scf%cdft_opt_control%broyden_update = .false.
1620 IF (output_unit > 0)
THEN
1621 WRITE (output_unit, fmt=
"(/,A)") &
1622 " ================================== LINE SEARCH STARTED =================================="
1623 WRITE (output_unit, fmt=
"(A,I5,A)") &
1624 " Evaluating optimal step size for optimizer using a maximum of", max_linesearch,
" steps"
1625 IF (continue_ls)
THEN
1626 WRITE (output_unit, fmt=
"(A)") &
1627 " Line search continues until best step size is found or max steps are reached"
1629 WRITE (output_unit,
'(/,A,F5.3)') &
1630 " Initial step size: ", alpha
1631 WRITE (output_unit,
'(/,A,F5.3)') &
1632 " Step size update factor: ", factor
1633 WRITE (output_unit,
'(/,A,I10,A,I10)') &
1634 " Energy evaluation: ", cdft_control%ienergy,
", CDFT SCF iteration: ", iter_count
1638 DO i = 1, max_linesearch
1639 IF (output_unit > 0)
THEN
1640 WRITE (output_unit, fmt=
"(A)")
" "
1641 WRITE (output_unit, fmt=
"(A)")
" #####################################"
1642 WRITE (output_unit,
'(A,I10,A)') &
1643 " ### Line search step: ", i,
" ###"
1644 WRITE (output_unit, fmt=
"(A)")
" #####################################"
1646 inv_jacobian => scf_env%outer_scf%inv_jacobian
1648 scf_env%outer_scf%variables(:, iter_count + 1) = scf_env%outer_scf%variables(:, iter_count) - alpha* &
1649 matmul(inv_jacobian, scf_env%outer_scf%gradient(:, iter_count))
1654 CALL scf_env_do_scf(scf_env=scf_env, scf_control=scf_control, qs_env=qs_env, &
1655 converged=converged, should_stop=should_stop, total_scf_steps=tsteps)
1658 scf_env%outer_scf%iter_count = scf_env%outer_scf%iter_count + 1
1661 energy_qs, cdft_control%total_steps, &
1662 should_stop=.false., outer_loop_converged=.false., cdft_loop=.false.)
1663 scf_env%outer_scf%iter_count = scf_env%outer_scf%iter_count - 1
1665 IF (continue_ls .AND. .NOT.
ALLOCATED(positive_sign))
THEN
1666 ALLOCATE (positive_sign(nvar))
1668 positive_sign(ispin) = scf_env%outer_scf%gradient(ispin, iter_count + 1) >= 0.0_dp
1672 inv_jacobian => scf_env%outer_scf%inv_jacobian
1673 IF (dnrm2(nvar, scf_env%outer_scf%gradient(:, iter_count + 1), 1) < &
1674 dnrm2(nvar, scf_env%outer_scf%gradient(:, iter_count), 1))
THEN
1676 IF (.NOT. continue_ls)
THEN
1677 should_exit = .true.
1681 IF (found_solution)
THEN
1683 IF (dnrm2(nvar, scf_env%outer_scf%gradient(:, iter_count + 1), 1) > norm_ls)
THEN
1685 continue_ls_exit = .true.
1689 IF (.NOT. continue_ls_exit)
THEN
1690 should_exit = .false.
1692 found_solution = .true.
1693 norm_ls = dnrm2(nvar, scf_env%outer_scf%gradient(:, iter_count + 1), 1)
1696 sign_changed = .true.
1698 sign_changed = sign_changed .AND. (positive_sign(ispin) .NEQV. &
1699 scf_env%outer_scf%gradient(ispin, iter_count + 1) >= 0.0_dp)
1701 IF (.NOT.
ALLOCATED(mos_ls))
THEN
1702 ALLOCATE (mos_ls(nspins))
1704 DO ispin = 1, nspins
1708 DO ispin = 1, nspins
1711 alpha = alpha*factor
1713 IF (i == max_linesearch) continue_ls_exit = .true.
1715 IF (sqrt(maxval(scf_env%outer_scf%gradient(:, scf_env%outer_scf%iter_count + 1)**2)) < &
1716 scf_control%outer_scf%eps_scf) &
1717 continue_ls_exit = .true.
1719 IF (sign_changed) continue_ls_exit = .true.
1724 should_exit = .false.
1726 alpha = alpha*factor
1728 IF (continue_ls_exit)
THEN
1731 DO ispin = 1, nspins
1735 p_rmpv(ispin)%matrix)
1741 should_exit = .true.
1744 IF (.NOT. should_exit .AND. &
1745 (i == max_linesearch .AND. converged .AND. .NOT. found_solution))
THEN
1746 should_exit = .true.
1747 reached_maxls = .true.
1748 alpha = alpha*(1.0_dp/factor)
1751 scf_env%outer_scf%variables(:, iter_count + 1) = 0.0_dp
1752 scf_env%outer_scf%gradient = gradient
1753 scf_env%outer_scf%energy = energy
1755 IF (should_exit)
EXIT
1757 cdft_control%total_steps = nsteps
1758 DO ispin = 1, nspins
1762 p_rmpv(ispin)%matrix)
1767 scf_control%outer_scf%cdft_opt_control%newton_step = alpha
1768 IF (.NOT. should_exit)
THEN
1769 CALL cp_warn(__location__, &
1770 "Line search did not converge. CDFT SCF proceeds with fixed step size.")
1771 scf_control%outer_scf%cdft_opt_control%newton_step = scf_control%outer_scf%cdft_opt_control%newton_step_save
1773 IF (reached_maxls) &
1774 CALL cp_warn(__location__, &
1775 "Line search did not converge. CDFT SCF proceeds with lasted iterated step size.")
1779 DO ispin = 1, nspins
1782 DEALLOCATE (mos_stashed, gradient, energy)
1783 IF (
ALLOCATED(positive_sign))
DEALLOCATE (positive_sign)
1784 IF (output_unit > 0)
THEN
1785 WRITE (output_unit, fmt=
"(/,A)") &
1786 " ================================== LINE SEARCH COMPLETE =================================="
1792 CALL timestop(handle)
1794 END SUBROUTINE qs_cdft_line_search
Define the atomic kind types and their sub types.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_deallocate_matrix(matrix)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
subroutine, public dbcsr_init_p(matrix)
...
subroutine, public dbcsr_set(matrix, alpha)
...
DBCSR operations in CP2K.
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
Utility routines to open and close files. Tracking of preconnections.
subroutine, public close_file(unit_number, file_status, keep_preconnection)
Close an open file given by its logical unit number. Optionally, keep the file and unit preconnected.
represent a full matrix distributed on many processors
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
various routines to log and control the output. The idea is that decisions about where to log should ...
subroutine, public cp_rm_default_logger()
the cousin of cp_add_default_logger, decrements the stack, so that the default logger is what it has ...
subroutine, public cp_logger_release(logger)
releases this logger
subroutine, public cp_add_default_logger(logger)
adds a default logger. MUST be called before logging occours
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)
...
integer, parameter, public cp_p_file
subroutine, public cp_iterate(iteration_info, last, iter_nr, increment, iter_nr_out)
adds one to the actual iteration
subroutine, public cp_rm_iter_level(iteration_info, level_name, n_rlevel_att)
Removes an iteration level.
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...
subroutine, public cp_add_iter_level(iteration_info, level_name, n_rlevel_new)
Adds an iteration level.
set of type/routines to handle the storage of results in force_envs
logical function, public test_for_result(results, description)
test for a certain result in the result_list
set of type/routines to handle the storage of results in force_envs
Types needed for a for a Energy Correction.
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
integer, parameter, public default_path_length
Restart file for k point calculations.
subroutine, public write_kpoints_restart(denmat, kpoints, scf_env, dft_section, particle_set, qs_kind_set)
...
Types and basic routines needed for a kpoint calculation.
Machine interface based on Fortran 2003 and POSIX.
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
Collection of simple mathematical functions and subroutines.
Interface to the message passing library MPI.
Define the data structure for the particle information.
computes preconditioners, and implements methods to apply them currently used in qs_ot
subroutine, public restart_preconditioner(qs_env, preconditioner, prec_type, nspins)
Allows for a restart of the preconditioner depending on the method it purges all arrays or keeps them...
subroutine, public prepare_preconditioner(qs_env, mos, matrix_ks, matrix_s, ot_preconditioner, prec_type, solver_type, energy_gap, nspins, has_unit_metric, convert_to_dbcsr, chol_type, full_mo_set)
...
container for various plainwaves related things
subroutine, public pw_env_get(pw_env, pw_pools, cube_info, gridlevel_info, auxbas_pw_pool, auxbas_grid, auxbas_rs_desc, auxbas_rs_grid, rs_descs, rs_grids, xc_pw_pool, vdw_pw_pool, poisson_env, interp_section)
returns the various attributes of the pw env
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
module that contains the algorithms to perform an iterative diagonalization by the block-Davidson app...
subroutine, public block_davidson_deallocate(bdav_env)
...
Auxiliary routines for performing a constrained DFT SCF run with Quickstep.
subroutine, public prepare_jacobian_stencil(qs_env, output_unit, nwork, pwork, coeff, step_multiplier, dh)
Prepares the finite difference stencil for computing the Jacobian. The constraints are re-evaluated b...
subroutine, public build_diagonal_jacobian(qs_env, used_history)
Builds a strictly diagonal inverse Jacobian from MD/SCF history.
subroutine, public print_inverse_jacobian(logger, inv_jacobian, iter_count)
Prints the finite difference inverse Jacobian to file.
subroutine, public create_tmp_logger(para_env, project_name, suffix, output_unit, tmp_logger)
Creates a temporary logger for redirecting output to a new file.
subroutine, public restart_inverse_jacobian(qs_env)
Restarts the finite difference inverse Jacobian.
subroutine, public initialize_inverse_jacobian(scf_control, scf_env, explicit_jacobian, should_build, used_history)
Checks if the inverse Jacobian should be calculated and initializes the calculation.
Defines CDFT control structures.
container for information about total charges on the grids
collects routines that calculate density matrices
module that contains the definitions of the scf types
integer, parameter, public gspace_mixing_nr
Apply the direct inversion in the iterative subspace (DIIS) of Pulay in the framework of an SCF itera...
pure subroutine, public qs_diis_b_clear(diis_buffer)
clears the buffer
subroutine, public qs_diis_b_create(diis_buffer, nbuffer)
Allocates an SCF DIIS buffer.
subroutine, public qs_diis_b_create_kp(diis_buffer, nbuffer)
Allocates an SCF DIIS buffer for k-points.
pure subroutine, public qs_diis_b_clear_kp(diis_buffer)
clears the buffer
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.
subroutine, public set_qs_env(qs_env, super_cell, mos, qmmm, qmmm_periodic, mimic, ewald_env, ewald_pw, mpools, rho_external, external_vxc, mask, scf_control, rel_control, qs_charges, ks_env, ks_qmmm_env, wf_history, scf_env, active_space, input, oce, rho_atom_set, rho0_atom_set, rho0_mpole, run_rtp, rtp, rhoz_set, rhoz_tot, ecoul_1c, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, efield, rhoz_cneo_set, linres_control, xas_env, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, ls_scf_env, do_transport, transport_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, harris_env, gcp_env, mp2_env, bs_env, kg_env, force, kpoints, wanniercentres, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Set the QUICKSTEP environment.
Integrate single or product functions over a potential on a RS grid.
Define the quickstep kind type and their sub types.
routines that build the Kohn-Sham matrix (i.e calculate the coulomb and xc parts
subroutine, public evaluate_core_matrix_traces(qs_env, rho_ao_ext)
Calculates the traces of the core matrices and the density matrix.
subroutine, public qs_ks_update_qs_env(qs_env, calculate_forces, just_energy, print_active)
updates the Kohn Sham matrix of the given qs_env (facility method)
subroutine, public qs_ks_did_change(ks_env, s_mstruct_changed, rho_changed, potential_changed, full_reset)
tells that some of the things relevant to the ks calculation did change. has to be called when change...
subroutine, public get_ks_env(ks_env, v_hartree_rspace, s_mstruct_changed, rho_changed, exc_accint, potential_changed, forces_up_to_date, complex_ks, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, kinetic, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_ks_im_kp, rho, rho_xc, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, vee, neighbor_list_id, sab_orb, sab_all, sac_ae, sac_ppl, sac_lri, sap_ppnl, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_vdw, sab_scp, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, task_list, task_list_soft, kpoints, do_kpoints, atomic_kind_set, qs_kind_set, cell, cell_ref, use_ref_cell, particle_set, energy, force, local_particles, local_molecules, molecule_kind_set, molecule_set, subsys, cp_subsys, virial, results, atprop, nkind, natom, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env, nelectron_total, nelectron_spin)
...
Definition and initialisation of the mo data type.
subroutine, public write_mo_set_to_restart(mo_array, particle_set, dft_section, qs_kind_set, matrix_ks)
...
collects routines that perform operations directly related to MOs
subroutine, public make_basis_simple(vmatrix, ncol)
given a set of vectors, return an orthogonal (C^T C == 1) set spanning the same space (notice,...
subroutine, public make_basis_sm(vmatrix, ncol, matrix_s)
returns an S-orthonormal basis v (v^T S v ==1)
Set occupation of molecular orbitals.
Definition and initialisation of the mo data type.
subroutine, public duplicate_mo_set(mo_set_new, mo_set_old)
allocate a new mo_set, and copy the old data
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 reassign_allocated_mos(mo_set_new, mo_set_old)
reassign an already allocated mo_set
basic functionality for using ot in the scf routines.
subroutine, public ot_scf_read_input(qs_ot_env, scf_section)
...
subroutine, public ot_scf_init(mo_array, matrix_s, qs_ot_env, matrix_ks, broyden_adaptive_sigma)
...
subroutine, public qs_ot_new_preconditioner(qs_ot_env, preconditioner)
gets ready to use the preconditioner/ or renew the preconditioner only keeps a pointer to the precond...
Routines for performing an outer scf loop.
subroutine, public outer_loop_purge_history(qs_env)
purges outer_scf_history zeroing everything except the latest value of the outer_scf variable stored ...
subroutine, public outer_loop_switch(scf_env, scf_control, cdft_control, dir)
switch between two outer_scf envs stored in cdft_control
subroutine, public outer_loop_optimize(scf_env, scf_control)
optimizes the parameters of the outer_scf
subroutine, public outer_loop_update_qs_env(qs_env, scf_env)
propagates the updated variables to wherever they need to be set in qs_env
subroutine, public outer_loop_gradient(qs_env, scf_env)
computes the gradient wrt to the outer loop variables
methods of the rho structure (defined in qs_rho_types)
subroutine, public qs_rho_update_rho(rho_struct, qs_env, rho_xc_external, local_rho_set, task_list_external, task_list_external_soft, pw_env_external, para_env_external)
updates rho_r and rho_g to the rhorho_ao. if use_kinetic_energy_density also computes tau_r and tau_g...
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...
Utility routines for qs_scf.
subroutine, public qs_scf_env_initialize(qs_env, scf_env, scf_control, scf_section)
initializes input parameters if needed or restores values from previous runs to fill scf_env with the...
Utility routines for qs_scf.
subroutine, public qs_scf_check_inner_exit(qs_env, scf_env, scf_control, should_stop, just_energy, exit_inner_loop, inner_loop_converged, output_unit)
checks whether exit conditions for inner loop are satisfied
subroutine, public qs_scf_inner_finalize(scf_env, qs_env, diis_step, output_unit)
Performs the necessary steps before leaving innner scf loop.
subroutine, public qs_scf_set_loop_flags(scf_env, diis_step, energy_only, just_energy, exit_inner_loop)
computes properties for a given hamiltonian using the current wfn
subroutine, public qs_scf_rho_update(rho, qs_env, scf_env, ks_env, mix_rho)
Performs the updates rho (takes care of mixing as well)
subroutine, public qs_scf_new_mos_kp(qs_env, scf_env, scf_control, diis_step, probe)
Updates MOs and density matrix using diagonalization Kpoint code.
subroutine, public qs_scf_check_outer_exit(qs_env, scf_env, scf_control, should_stop, outer_loop_converged, exit_outer_loop)
checks whether exit conditions for outer loop are satisfied
subroutine, public qs_scf_density_mixing(scf_env, rho, para_env, diis_step)
Performs the requested density mixing if any needed.
subroutine, public qs_scf_new_mos(qs_env, scf_env, scf_control, scf_section, diis_step, energy_only, probe)
takes known energy and derivatives and produces new wfns and or density matrix
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_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
Utility routines for qs_scf.
subroutine, public qs_scf_compute_properties(qs_env, wf_type, do_mp2)
computes properties for a given hamilonian using the current wfn
module that contains the definitions of the scf types
integer, parameter, public ot_diag_method_nr
integer, parameter, public filter_matrix_diag_method_nr
integer, parameter, public block_davidson_diag_method_nr
integer, parameter, public smeagol_method_nr
integer, parameter, public ot_method_nr
integer, parameter, public special_diag_method_nr
integer, parameter, public block_krylov_diag_method_nr
integer, parameter, public general_diag_method_nr
Routines for the Quickstep SCF run.
subroutine, public init_scf_loop(scf_env, qs_env, scf_section)
inits those objects needed if you want to restart the scf with, say only a new initial guess,...
subroutine, public scf(qs_env, has_converged, total_scf_steps)
perform an scf procedure in the given qs_env
subroutine, public scf_env_cleanup(scf_env)
perform cleanup operations (like releasing temporary storage) at the end of the scf
subroutine, public scf_env_do_scf(scf_env, scf_control, qs_env, converged, should_stop, total_scf_steps)
perform an scf loop
subroutine, public cdft_scf(qs_env, should_stop)
perform a CDFT scf procedure in the given qs_env
Storage of past states of the qs_env. Methods to interpolate (or actually normally extrapolate) the n...
subroutine, public wfi_purge_history(qs_env)
purges wf_history retaining only the latest snapshot
subroutine, public wfi_update(wf_history, qs_env, dt)
updates the snapshot buffer, taking a new snapshot
parameters that control an scf iteration
subroutine, public run_smeagol_emtrans(qs_env, last, iter, rho_ao_kp)
Run NEGF/SMEAGOL transport calculation.
subroutine, public run_smeagol_bulktrans(qs_env)
Save overlap, Kohn-Sham, electron density, and energy-density matrices of semi-infinite electrodes in...
subroutine, public tb_update_charges(qs_env, dft_control, tb, calculate_forces, use_rho)
...
subroutine, public tb_get_energy(qs_env, tb, energy)
...
Provides all information about an atomic kind.
type of a logger, at the moment it contains just a print level starting at which level it should be l...
contains arbitrary information which need to be stored
Contains information on the energy correction functional for KG.
Contains information about kpoints.
stores all the informations relevant to an mpi environment
contained for different pw related things
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Container for information about total charges on the grids.
Provides all information about a quickstep kind.
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...
keeps the density in various representations, keeping track of which ones are valid.