119#include "./base/base_uses.f90"
124 LOGICAL,
PRIVATE,
PARAMETER :: debug_this_module = .true.
125 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_wf_history_methods'
141 SUBROUTINE wfs_create(snapshot)
142 TYPE(qs_wf_snapshot_type),
INTENT(OUT) :: snapshot
144 NULLIFY (snapshot%wf, snapshot%rho_r, &
145 snapshot%rho_g, snapshot%rho_ao, snapshot%rho_ao_kp, &
146 snapshot%overlap, snapshot%wf_kp, snapshot%overlap_cfm_kp, &
149 END SUBROUTINE wfs_create
162 SUBROUTINE wfs_update(snapshot, wf_history, qs_env, dt)
163 TYPE(qs_wf_snapshot_type),
POINTER :: snapshot
164 TYPE(qs_wf_history_type),
POINTER :: wf_history
165 TYPE(qs_environment_type),
POINTER :: qs_env
166 REAL(KIND=
dp),
INTENT(in),
OPTIONAL :: dt
168 CHARACTER(len=*),
PARAMETER :: routineN =
'wfs_update'
170 INTEGER :: handle, ic, igroup, ik, ikp, img, &
171 indx_ft, ispin, kplocal, nc, nimg, &
172 nkp_all, nkp_grps, nspin_kp, nspins
173 INTEGER,
DIMENSION(2) :: kp_range
174 INTEGER,
DIMENSION(:, :),
POINTER :: kp_dist
175 INTEGER,
DIMENSION(:, :, :),
POINTER :: cell_to_index
177 REAL(KIND=
dp),
DIMENSION(:, :),
POINTER :: xkp
178 TYPE(copy_info_type),
ALLOCATABLE,
DIMENSION(:, :) :: info_ft
179 TYPE(cp_fm_pool_p_type),
DIMENSION(:),
POINTER :: ao_ao_fm_pools, ao_mo_pools
180 TYPE(cp_fm_struct_type),
POINTER :: ao_ao_struct_ft
181 TYPE(cp_fm_type) :: fmdummy_ft, fmlocal_ft
182 TYPE(cp_fm_type),
POINTER :: mo_coeff
183 TYPE(dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_s, rho_ao
184 TYPE(dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_s_kp, rho_ao_kp
185 TYPE(dbcsr_type),
POINTER :: cmat_ft, rmat_ft, tmpmat_ft
186 TYPE(dft_control_type),
POINTER :: dft_control
187 TYPE(kpoint_env_type),
POINTER :: kp
188 TYPE(kpoint_type),
POINTER :: kpoints
189 TYPE(mo_set_type),
DIMENSION(:),
POINTER :: mos
190 TYPE(mp_para_env_type),
POINTER :: para_env_ft
191 TYPE(neighbor_list_set_p_type),
DIMENSION(:), &
193 TYPE(pw_c1d_gs_type),
DIMENSION(:),
POINTER :: rho_g
194 TYPE(pw_env_type),
POINTER :: pw_env
195 TYPE(pw_pool_type),
POINTER :: auxbas_pw_pool
196 TYPE(pw_r3d_rs_type),
DIMENSION(:),
POINTER :: rho_r
197 TYPE(qs_matrix_pools_type),
POINTER :: mpools_kp
198 TYPE(qs_rho_type),
POINTER :: rho
199 TYPE(qs_scf_env_type),
POINTER :: scf_env
201 CALL timeset(routinen, handle)
203 NULLIFY (pw_env, auxbas_pw_pool, ao_mo_pools, ao_ao_fm_pools, dft_control, mos, mo_coeff, &
204 rho, rho_r, rho_g, rho_ao, matrix_s, matrix_s_kp, kpoints, kp, &
205 kp_dist, cell_to_index, xkp, sab_nl, scf_env, mpools_kp, para_env_ft, &
206 rmat_ft, cmat_ft, tmpmat_ft, ao_ao_struct_ft)
208 dft_control=dft_control, rho=rho)
209 CALL mpools_get(qs_env%mpools, ao_mo_fm_pools=ao_mo_pools)
210 CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
212 cpassert(
ASSOCIATED(wf_history))
213 cpassert(
ASSOCIATED(dft_control))
214 IF (.NOT.
ASSOCIATED(snapshot))
THEN
216 CALL wfs_create(snapshot)
218 cpassert(wf_history%ref_count > 0)
220 nspins = dft_control%nspins
222 IF (
PRESENT(dt)) snapshot%dt = dt
223 IF (wf_history%store_wf)
THEN
225 IF (.NOT.
ASSOCIATED(snapshot%wf))
THEN
228 cpassert(nspins ==
SIZE(snapshot%wf))
231 CALL get_mo_set(mos(ispin), mo_coeff=mo_coeff)
238 IF (wf_history%store_rho_r)
THEN
240 cpassert(
ASSOCIATED(rho_r))
241 IF (.NOT.
ASSOCIATED(snapshot%rho_r))
THEN
242 ALLOCATE (snapshot%rho_r(nspins))
244 CALL auxbas_pw_pool%create_pw(snapshot%rho_r(ispin))
248 CALL pw_copy(rho_r(ispin), snapshot%rho_r(ispin))
250 ELSE IF (
ASSOCIATED(snapshot%rho_r))
THEN
251 DO ispin = 1,
SIZE(snapshot%rho_r)
252 CALL auxbas_pw_pool%give_back_pw(snapshot%rho_r(ispin))
254 DEALLOCATE (snapshot%rho_r)
257 IF (wf_history%store_rho_g)
THEN
259 cpassert(
ASSOCIATED(rho_g))
260 IF (.NOT.
ASSOCIATED(snapshot%rho_g))
THEN
261 ALLOCATE (snapshot%rho_g(nspins))
263 CALL auxbas_pw_pool%create_pw(snapshot%rho_g(ispin))
267 CALL pw_copy(rho_g(ispin), snapshot%rho_g(ispin))
269 ELSE IF (
ASSOCIATED(snapshot%rho_g))
THEN
270 DO ispin = 1,
SIZE(snapshot%rho_g)
271 CALL auxbas_pw_pool%give_back_pw(snapshot%rho_g(ispin))
273 DEALLOCATE (snapshot%rho_g)
276 IF (
ASSOCIATED(snapshot%rho_ao))
THEN
280 IF (wf_history%store_rho_ao)
THEN
282 cpassert(
ASSOCIATED(rho_ao))
286 ALLOCATE (snapshot%rho_ao(ispin)%matrix)
287 CALL dbcsr_copy(snapshot%rho_ao(ispin)%matrix, rho_ao(ispin)%matrix)
291 IF (
ASSOCIATED(snapshot%rho_ao_kp))
THEN
295 IF (wf_history%store_rho_ao_kp)
THEN
297 cpassert(
ASSOCIATED(rho_ao_kp))
299 nimg = dft_control%nimages
303 ALLOCATE (snapshot%rho_ao_kp(ispin, img)%matrix)
304 CALL dbcsr_copy(snapshot%rho_ao_kp(ispin, img)%matrix, &
305 rho_ao_kp(ispin, img)%matrix)
310 IF (
ASSOCIATED(snapshot%overlap))
THEN
314 IF (wf_history%store_overlap)
THEN
316 cpassert(
ASSOCIATED(matrix_s))
317 cpassert(
ASSOCIATED(matrix_s(1)%matrix))
318 ALLOCATE (snapshot%overlap)
319 CALL dbcsr_copy(snapshot%overlap, matrix_s(1)%matrix)
323 IF (
ASSOCIATED(kpoints))
THEN
324 IF (
ASSOCIATED(kpoints%kp_env))
THEN
326 IF (wf_history%store_wf_kp)
THEN
328 kplocal = kp_range(2) - kp_range(1) + 1
329 nspin_kp =
SIZE(kpoints%kp_env(1)%kpoint_env%mos, 2)
330 nc =
SIZE(kpoints%kp_env(1)%kpoint_env%mos, 1)
332 IF (
ASSOCIATED(snapshot%wf_kp))
THEN
333 DO ikp = 1,
SIZE(snapshot%wf_kp, 1)
334 DO ic = 1,
SIZE(snapshot%wf_kp, 2)
335 DO ispin = 1,
SIZE(snapshot%wf_kp, 3)
340 DEALLOCATE (snapshot%wf_kp)
343 ALLOCATE (snapshot%wf_kp(kplocal, nc, nspin_kp))
345 kp => kpoints%kp_env(ikp)%kpoint_env
346 DO ispin = 1, nspin_kp
348 CALL get_mo_set(kp%mos(ic, ispin), mo_coeff=mo_coeff)
350 mo_coeff%matrix_struct, &
352 CALL cp_fm_to_fm(mo_coeff, snapshot%wf_kp(ikp, ic, ispin))
362 IF (wf_history%store_overlap_kp)
THEN
363 CALL get_qs_env(qs_env, matrix_s_kp=matrix_s_kp, scf_env=scf_env)
364 CALL get_kpoint_info(kpoints, nkp=nkp_all, xkp=xkp, kp_range=kp_range, &
365 nkp_groups=nkp_grps, kp_dist=kp_dist, &
366 sab_nl=sab_nl, cell_to_index=cell_to_index)
367 kplocal = kp_range(2) - kp_range(1) + 1
368 para_env_ft => kpoints%blacs_env_all%para_env
371 ALLOCATE (rmat_ft, cmat_ft, tmpmat_ft)
372 CALL dbcsr_create(rmat_ft, template=matrix_s_kp(1, 1)%matrix, &
373 matrix_type=dbcsr_type_symmetric)
374 CALL dbcsr_create(cmat_ft, template=matrix_s_kp(1, 1)%matrix, &
375 matrix_type=dbcsr_type_antisymmetric)
376 CALL dbcsr_create(tmpmat_ft, template=matrix_s_kp(1, 1)%matrix, &
377 matrix_type=dbcsr_type_no_symmetry)
383 CALL mpools_get(mpools_kp, ao_ao_fm_pools=ao_ao_fm_pools)
387 IF (
ASSOCIATED(snapshot%overlap_cfm_kp))
THEN
388 DO ikp = 1,
SIZE(snapshot%overlap_cfm_kp)
391 DEALLOCATE (snapshot%overlap_cfm_kp)
393 ALLOCATE (snapshot%overlap_cfm_kp(kplocal))
398 ALLOCATE (info_ft(kplocal*nkp_grps, 2))
403 DO igroup = 1, nkp_grps
404 ik = kp_dist(1, igroup) + ikp - 1
405 my_kpgrp = (ik >= kp_range(1) .AND. ik <= kp_range(2))
406 indx_ft = indx_ft + 1
410 CALL rskp_transform(rmatrix=rmat_ft, cmatrix=cmat_ft, rsmat=matrix_s_kp, &
411 ispin=1, xkp=xkp(1:3, ik), &
412 cell_to_index=cell_to_index, sab_nl=sab_nl)
420 para_env_ft, info_ft(indx_ft, 1))
422 para_env_ft, info_ft(indx_ft, 2))
425 para_env_ft, info_ft(indx_ft, 1))
427 para_env_ft, info_ft(indx_ft, 2))
435 CALL cp_cfm_create(snapshot%overlap_cfm_kp(ikp), ao_ao_struct_ft)
437 DO igroup = 1, nkp_grps
438 ik = kp_dist(1, igroup) + ikp - 1
439 my_kpgrp = (ik >= kp_range(1) .AND. ik <= kp_range(2))
440 indx_ft = indx_ft + 1
453 DO indx_ft = 1, kplocal*nkp_grps
466 IF (wf_history%store_frozen_density)
THEN
471 CALL timestop(handle)
473 END SUBROUTINE wfs_update
487 SUBROUTINE wfi_create(wf_history, interpolation_method_nr, extrapolation_order, &
490 INTEGER,
INTENT(in) :: interpolation_method_nr, &
492 LOGICAL,
INTENT(IN) :: has_unit_metric
494 CHARACTER(len=*),
PARAMETER :: routinen =
'wfi_create'
498 CALL timeset(routinen, handle)
500 ALLOCATE (wf_history)
501 wf_history%ref_count = 1
502 wf_history%memory_depth = 0
503 wf_history%snapshot_count = 0
504 wf_history%last_state_index = 1
505 wf_history%store_wf = .false.
506 wf_history%store_rho_r = .false.
507 wf_history%store_rho_g = .false.
508 wf_history%store_rho_ao = .false.
509 wf_history%store_rho_ao_kp = .false.
510 wf_history%store_overlap = .false.
511 wf_history%store_wf_kp = .false.
512 wf_history%store_overlap_kp = .false.
513 wf_history%store_frozen_density = .false.
514 NULLIFY (wf_history%past_states)
516 wf_history%interpolation_method_nr = interpolation_method_nr
518 SELECT CASE (wf_history%interpolation_method_nr)
520 wf_history%memory_depth = 0
522 wf_history%memory_depth = 0
524 wf_history%memory_depth = 1
525 wf_history%store_rho_ao = .true.
527 wf_history%memory_depth = 1
528 wf_history%store_rho_ao = .true.
530 wf_history%memory_depth = 2
531 wf_history%store_wf = .true.
533 wf_history%memory_depth = 2
534 wf_history%store_rho_ao = .true.
536 wf_history%memory_depth = 2
537 wf_history%store_wf = .true.
538 IF (.NOT. has_unit_metric) wf_history%store_overlap = .true.
541 wf_history%memory_depth = extrapolation_order + 1
542 wf_history%store_wf = .true.
543 wf_history%store_wf_kp = .true.
544 IF (.NOT. has_unit_metric)
THEN
545 wf_history%store_overlap = .true.
546 wf_history%store_overlap_kp = .true.
549 wf_history%memory_depth = 1
550 wf_history%store_frozen_density = .true.
552 wf_history%memory_depth = extrapolation_order + 2
553 wf_history%store_wf = .true.
554 wf_history%store_wf_kp = .true.
555 IF (.NOT. has_unit_metric)
THEN
556 wf_history%store_overlap = .true.
557 wf_history%store_overlap_kp = .true.
560 CALL cp_abort(__location__, &
561 "Unknown interpolation method: "// &
565 ALLOCATE (wf_history%past_states(wf_history%memory_depth))
567 DO i = 1,
SIZE(wf_history%past_states)
568 NULLIFY (wf_history%past_states(i)%snapshot)
571 CALL timestop(handle)
586 cpassert(
ASSOCIATED(wf_history))
587 IF (wf_history%store_rho_ao)
THEN
588 wf_history%store_rho_ao_kp = .true.
589 wf_history%store_rho_ao = .false.
592 IF (wf_history%store_wf_kp)
THEN
593 wf_history%store_wf = .false.
594 wf_history%store_overlap = .false.
598 IF (wf_history%store_wf .OR. wf_history%store_overlap)
THEN
599 cpabort(
"Linear WFN-based extrapolation methods not implemented for k-points.")
602 IF (wf_history%store_frozen_density)
THEN
603 cpabort(
"Frozen density initialization method not possible for kpoints.")
617 INTEGER,
INTENT(in) :: method_nr
618 CHARACTER(len=30) :: res
621 SELECT CASE (method_nr)
627 res =
"previous_rho_r"
629 res =
"initial_guess"
639 res =
"frozen density approximation"
643 CALL cp_abort(__location__, &
644 "Unknown interpolation method: "// &
667 REAL(kind=
dp),
INTENT(IN) :: dt
668 INTEGER,
INTENT(OUT),
OPTIONAL :: extrapolation_method_nr
669 LOGICAL,
INTENT(OUT),
OPTIONAL :: orthogonal_wf
671 CHARACTER(len=*),
PARAMETER :: routinen =
'wfi_extrapolate'
673 INTEGER :: actual_extrapolation_method_nr, handle, &
674 i, img, io_unit, ispin, k, n, nmo, &
676 LOGICAL :: do_kpoints, my_orthogonal_wf, use_overlap
677 REAL(kind=
dp) :: alpha, t0, t1, t2
683 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: rho_ao, rho_frozen_ao
684 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: rho_ao_kp
689 NULLIFY (mos, ao_mo_fm_pools, t0_state, t1_state, mo_coeff, &
690 rho, rho_ao, rho_frozen_ao)
692 use_overlap = wf_history%store_overlap
694 CALL timeset(routinen, handle)
696 print_level = logger%iter_info%print_level
700 cpassert(
ASSOCIATED(wf_history))
701 cpassert(wf_history%ref_count > 0)
702 cpassert(
ASSOCIATED(qs_env))
703 CALL get_qs_env(qs_env, mos=mos, rho=rho, do_kpoints=do_kpoints)
704 CALL mpools_get(qs_env%mpools, ao_mo_fm_pools=ao_mo_fm_pools)
706 IF (wf_history%snapshot_count < 1)
THEN
709 actual_extrapolation_method_nr = wf_history%interpolation_method_nr
712 SELECT CASE (actual_extrapolation_method_nr)
714 IF (wf_history%snapshot_count < 2)
THEN
718 IF (wf_history%snapshot_count < 2)
THEN
722 IF (wf_history%snapshot_count < 2)
THEN
727 IF (
PRESENT(extrapolation_method_nr)) &
728 extrapolation_method_nr = actual_extrapolation_method_nr
729 my_orthogonal_wf = .false.
731 SELECT CASE (actual_extrapolation_method_nr)
733 cpassert(.NOT. do_kpoints)
735 cpassert(
ASSOCIATED(t0_state%rho_frozen))
737 nvec = min(wf_history%memory_depth, wf_history%snapshot_count)
738 CALL wfi_set_history_variables(qs_env=qs_env, nvec=nvec)
740 CALL qs_rho_get(t0_state%rho_frozen, rho_ao=rho_frozen_ao)
742 DO ispin = 1,
SIZE(rho_frozen_ao)
744 rho_frozen_ao(ispin)%matrix, &
745 keep_sparsity=.true.)
752 my_orthogonal_wf = .false.
755 nvec = min(wf_history%memory_depth, wf_history%snapshot_count)
756 CALL wfi_set_history_variables(qs_env=qs_env, nvec=nvec)
758 cpassert(
ASSOCIATED(t0_state%rho_ao_kp))
760 DO ispin = 1,
SIZE(t0_state%rho_ao_kp, 1)
761 DO img = 1,
SIZE(t0_state%rho_ao_kp, 2)
762 IF (img >
SIZE(rho_ao_kp, 2))
THEN
763 cpwarn(
"Change in cell neighborlist: might affect quality of initial guess")
765 CALL dbcsr_copy(rho_ao_kp(ispin, img)%matrix, &
766 t0_state%rho_ao_kp(ispin, img)%matrix, &
767 keep_sparsity=.true.)
772 cpassert(
ASSOCIATED(t0_state%rho_ao))
774 DO ispin = 1,
SIZE(t0_state%rho_ao)
776 t0_state%rho_ao(ispin)%matrix, &
777 keep_sparsity=.true.)
785 my_orthogonal_wf = .true.
786 nvec = min(wf_history%memory_depth, wf_history%snapshot_count)
787 CALL wfi_set_history_variables(qs_env=qs_env, nvec=nvec)
790 CALL wfi_use_prev_wf_kp(qs_env, io_unit, print_level)
793 DO ispin = 1,
SIZE(mos)
794 CALL get_mo_set(mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
804 nvec = min(wf_history%memory_depth, wf_history%snapshot_count)
805 CALL wfi_set_history_variables(qs_env=qs_env, nvec=nvec)
807 cpassert(
ASSOCIATED(t0_state%rho_ao_kp))
809 DO ispin = 1,
SIZE(t0_state%rho_ao_kp, 1)
810 DO img = 1,
SIZE(t0_state%rho_ao_kp, 2)
811 IF (img >
SIZE(rho_ao_kp, 2))
THEN
812 cpwarn(
"Change in cell neighborlist: might affect quality of initial guess")
814 CALL dbcsr_copy(rho_ao_kp(ispin, img)%matrix, &
815 t0_state%rho_ao_kp(ispin, img)%matrix, &
816 keep_sparsity=.true.)
821 cpassert(
ASSOCIATED(t0_state%rho_ao))
823 DO ispin = 1,
SIZE(t0_state%rho_ao)
825 t0_state%rho_ao(ispin)%matrix, &
826 keep_sparsity=.true.)
840 nvec = min(wf_history%memory_depth, wf_history%snapshot_count)
841 CALL wfi_set_history_variables(qs_env=qs_env, nvec=nvec)
843 cpassert(.NOT. do_kpoints)
846 cpassert(
ASSOCIATED(t0_state))
847 cpassert(
ASSOCIATED(t1_state))
848 cpassert(
ASSOCIATED(t0_state%wf))
849 cpassert(
ASSOCIATED(t1_state%wf))
850 nvec = min(wf_history%memory_depth, wf_history%snapshot_count)
851 CALL wfi_set_history_variables(qs_env=qs_env, nvec=nvec)
853 my_orthogonal_wf = .true.
858 DO ispin = 1,
SIZE(mos)
859 CALL get_mo_set(mos(ispin), mo_coeff=mo_coeff, &
863 matrix_b=t1_state%wf(ispin), &
864 beta=(t2 - t0)/(t1 - t0))
868 beta=(t1 - t2)/(t1 - t0), matrix_b=t0_state%wf(ispin))
873 density_matrix=rho_ao(ispin)%matrix)
882 cpassert(
ASSOCIATED(t0_state))
883 cpassert(
ASSOCIATED(t1_state))
885 cpassert(
ASSOCIATED(t0_state%rho_ao_kp))
886 cpassert(
ASSOCIATED(t1_state%rho_ao_kp))
888 cpassert(
ASSOCIATED(t0_state%rho_ao))
889 cpassert(
ASSOCIATED(t1_state%rho_ao))
891 nvec = min(wf_history%memory_depth, wf_history%snapshot_count)
892 CALL wfi_set_history_variables(qs_env=qs_env, nvec=nvec)
899 DO ispin = 1,
SIZE(rho_ao_kp, 1)
900 DO img = 1,
SIZE(rho_ao_kp, 2)
901 IF (img >
SIZE(t0_state%rho_ao_kp, 2) .OR. &
902 img >
SIZE(t1_state%rho_ao_kp, 2))
THEN
903 cpwarn(
"Change in cell neighborlist: might affect quality of initial guess")
905 CALL dbcsr_add(rho_ao_kp(ispin, img)%matrix, t1_state%rho_ao_kp(ispin, img)%matrix, &
906 alpha_scalar=0.0_dp, beta_scalar=(t2 - t0)/(t1 - t0))
907 CALL dbcsr_add(rho_ao_kp(ispin, img)%matrix, t0_state%rho_ao_kp(ispin, img)%matrix, &
908 alpha_scalar=1.0_dp, beta_scalar=(t1 - t2)/(t1 - t0))
914 DO ispin = 1,
SIZE(rho_ao)
915 CALL dbcsr_add(rho_ao(ispin)%matrix, t1_state%rho_ao(ispin)%matrix, &
916 alpha_scalar=0.0_dp, beta_scalar=(t2 - t0)/(t1 - t0))
917 CALL dbcsr_add(rho_ao(ispin)%matrix, t0_state%rho_ao(ispin)%matrix, &
918 alpha_scalar=1.0_dp, beta_scalar=(t1 - t2)/(t1 - t0))
926 cpassert(.NOT. do_kpoints)
929 cpassert(
ASSOCIATED(t0_state))
930 cpassert(
ASSOCIATED(t1_state))
931 cpassert(
ASSOCIATED(t0_state%wf))
932 cpassert(
ASSOCIATED(t1_state%wf))
933 IF (wf_history%store_overlap)
THEN
934 cpassert(
ASSOCIATED(t0_state%overlap))
935 cpassert(
ASSOCIATED(t1_state%overlap))
937 nvec = min(wf_history%memory_depth, wf_history%snapshot_count)
938 IF (nvec >= wf_history%memory_depth)
THEN
939 IF ((qs_env%scf_control%max_scf_hist /= 0) .AND. (qs_env%scf_control%eps_scf_hist /= 0))
THEN
940 qs_env%scf_control%max_scf = qs_env%scf_control%max_scf_hist
941 qs_env%scf_control%eps_scf = qs_env%scf_control%eps_scf_hist
942 qs_env%scf_control%outer_scf%have_scf = .false.
943 ELSE IF (qs_env%scf_control%max_scf_hist /= 0)
THEN
944 qs_env%scf_control%max_scf = qs_env%scf_control%max_scf_hist
945 qs_env%scf_control%outer_scf%have_scf = .false.
946 ELSE IF (qs_env%scf_control%eps_scf_hist /= 0)
THEN
947 qs_env%scf_control%eps_scf = qs_env%scf_control%eps_scf_hist
951 my_orthogonal_wf = .true.
955 DO ispin = 1,
SIZE(mos)
956 NULLIFY (mo_coeff, matrix_struct, matrix_struct_new)
957 CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff)
959 matrix_struct=matrix_struct)
961 nrow_global=k, ncol_global=k)
965 IF (use_overlap)
THEN
967 CALL parallel_gemm(
'T',
'N', k, k, n, 1.0_dp, t0_state%wf(ispin), mo_coeff, 0.0_dp, csc)
969 CALL parallel_gemm(
'T',
'N', k, k, n, 1.0_dp, t0_state%wf(ispin), &
970 t1_state%wf(ispin), 0.0_dp, csc)
972 CALL parallel_gemm(
'N',
'N', n, k, k, 1.0_dp, t0_state%wf(ispin), csc, 0.0_dp, mo_coeff)
979 density_matrix=rho_ao(ispin)%matrix)
986 nvec = min(wf_history%memory_depth, wf_history%snapshot_count)
988 IF (nvec >= wf_history%memory_depth)
THEN
989 IF ((qs_env%scf_control%max_scf_hist /= 0) .AND. (qs_env%scf_control%eps_scf_hist /= 0))
THEN
990 qs_env%scf_control%max_scf = qs_env%scf_control%max_scf_hist
991 qs_env%scf_control%eps_scf = qs_env%scf_control%eps_scf_hist
992 qs_env%scf_control%outer_scf%have_scf = .false.
993 ELSE IF (qs_env%scf_control%max_scf_hist /= 0)
THEN
994 qs_env%scf_control%max_scf = qs_env%scf_control%max_scf_hist
995 qs_env%scf_control%outer_scf%have_scf = .false.
996 ELSE IF (qs_env%scf_control%eps_scf_hist /= 0)
THEN
997 qs_env%scf_control%eps_scf = qs_env%scf_control%eps_scf_hist
1001 IF (do_kpoints)
THEN
1002 CALL wfi_extrapolate_ps_aspc_kp(wf_history, qs_env, nvec, io_unit, print_level)
1003 my_orthogonal_wf = .true.
1005 my_orthogonal_wf = .true.
1006 DO ispin = 1,
SIZE(mos)
1007 NULLIFY (mo_coeff, matrix_struct, matrix_struct_new)
1008 CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff)
1010 matrix_struct=matrix_struct)
1013 nrow_global=k, ncol_global=k)
1024 IF (use_overlap)
THEN
1026 CALL parallel_gemm(
'T',
'N', k, k, n, 1.0_dp, t0_state%wf(ispin), fm_tmp, 0.0_dp, csc)
1028 CALL parallel_gemm(
'T',
'N', k, k, n, 1.0_dp, t0_state%wf(ispin), &
1029 t1_state%wf(ispin), 0.0_dp, csc)
1031 CALL parallel_gemm(
'N',
'N', n, k, k, 1.0_dp, t0_state%wf(ispin), csc, 0.0_dp, fm_tmp)
1032 alpha = -1.0_dp*alpha*real(nvec - i + 1,
dp)/real(i,
dp)
1039 v_matrix=mo_coeff, &
1042 density_matrix=rho_ao(ispin)%matrix)
1052 nvec = min(wf_history%memory_depth, wf_history%snapshot_count)
1054 IF (nvec >= wf_history%memory_depth)
THEN
1055 IF ((qs_env%scf_control%max_scf_hist /= 0) .AND. &
1056 (qs_env%scf_control%eps_scf_hist /= 0))
THEN
1057 qs_env%scf_control%max_scf = qs_env%scf_control%max_scf_hist
1058 qs_env%scf_control%eps_scf = qs_env%scf_control%eps_scf_hist
1059 qs_env%scf_control%outer_scf%have_scf = .false.
1060 ELSE IF (qs_env%scf_control%max_scf_hist /= 0)
THEN
1061 qs_env%scf_control%max_scf = qs_env%scf_control%max_scf_hist
1062 qs_env%scf_control%outer_scf%have_scf = .false.
1063 ELSE IF (qs_env%scf_control%eps_scf_hist /= 0)
THEN
1064 qs_env%scf_control%eps_scf = qs_env%scf_control%eps_scf_hist
1068 IF (do_kpoints)
THEN
1069 CALL wfi_extrapolate_ps_aspc_kp(wf_history, qs_env, nvec, io_unit, print_level)
1070 my_orthogonal_wf = .true.
1072 my_orthogonal_wf = .true.
1074 DO ispin = 1,
SIZE(mos)
1075 NULLIFY (mo_coeff, matrix_struct, matrix_struct_new)
1076 CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff)
1080 matrix_struct=matrix_struct)
1081 CALL cp_fm_create(fm_tmp, matrix_struct, set_zero=.true.)
1083 template_fmstruct=matrix_struct, &
1086 CALL cp_fm_create(csc, matrix_struct_new, set_zero=.true.)
1092 alpha = real(4*nvec - 2, kind=
dp)/real(nvec + 1, kind=
dp)
1094 WRITE (unit=io_unit, fmt=
"(/,T2,A,/,/,T3,A,I0,/,/,T3,A2,I0,A4,F10.6)") &
1095 "Parameters for the always stable predictor-corrector (ASPC) method:", &
1096 "ASPC order: ", max(nvec - 2, 0), &
1097 "B(", 1,
") = ", alpha
1103 IF (use_overlap)
THEN
1105 CALL parallel_gemm(
'T',
'N', k, k, n, 1.0_dp, t0_state%wf(ispin), fm_tmp, 0.0_dp, csc)
1107 CALL parallel_gemm(
'T',
'N', k, k, n, 1.0_dp, t0_state%wf(ispin), &
1108 t1_state%wf(ispin), 0.0_dp, csc)
1110 CALL parallel_gemm(
'N',
'N', n, k, k, 1.0_dp, t0_state%wf(ispin), csc, 0.0_dp, fm_tmp)
1111 alpha = (-1.0_dp)**(i + 1)*real(i, kind=
dp)* &
1114 WRITE (unit=io_unit, fmt=
"(T3,A2,I0,A4,F10.6)") &
1115 "B(", i,
") = ", alpha
1122 v_matrix=mo_coeff, &
1125 density_matrix=rho_ao(ispin)%matrix)
1132 CALL cp_abort(__location__, &
1133 "Unknown interpolation method: "// &
1134 trim(adjustl(
cp_to_string(wf_history%interpolation_method_nr))))
1136 IF (
PRESENT(orthogonal_wf)) orthogonal_wf = my_orthogonal_wf
1138 "DFT%SCF%PRINT%PROGRAM_RUN_INFO")
1139 CALL timestop(handle)
1149 SUBROUTINE wfi_use_prev_wf_kp(qs_env, io_unit, print_level)
1151 INTEGER,
INTENT(IN) :: io_unit, print_level
1153 CHARACTER(len=*),
PARAMETER :: routinen =
'wfi_use_prev_wf_kp'
1155 COMPLEX(KIND=dp),
ALLOCATABLE,
DIMENSION(:) :: col_scaling
1156 INTEGER :: chol_info, handle, igroup, ik, ikp, &
1157 indx, ispin, j, kplocal, nao, nkp, &
1158 nkp_groups, nmo, nspin
1159 INTEGER,
DIMENSION(2) :: kp_range
1160 INTEGER,
DIMENSION(:, :),
POINTER :: kp_dist
1161 INTEGER,
DIMENSION(:, :, :),
POINTER :: cell_to_index
1162 LOGICAL :: my_kpgrp, use_real_wfn
1163 REAL(kind=
dp) :: eval_thresh
1164 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: eigenvalues
1165 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: xkp
1167 TYPE(
cp_cfm_type) :: cfm_evecs, cfm_mhalf, cfm_nao_nmo_work, &
1169 TYPE(
cp_cfm_type),
ALLOCATABLE,
DIMENSION(:) :: csmat_cur
1173 TYPE(
cp_fm_type),
POINTER :: imos, mo_coeff, rmos
1174 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrix_s_kp, rho_ao_kp
1175 TYPE(
dbcsr_type),
POINTER :: cmatrix_db, rmatrix, tmpmat
1187 CALL timeset(routinen, handle)
1189 NULLIFY (kpoints, dft_control, matrix_s_kp, scf_env, scf_control, rho, sab_nl, kp, mo_coeff, rmos, imos)
1191 CALL get_qs_env(qs_env, kpoints=kpoints, dft_control=dft_control, &
1192 matrix_s_kp=matrix_s_kp, scf_env=scf_env, &
1193 scf_control=scf_control, rho=rho)
1194 CALL get_kpoint_info(kpoints, nkp=nkp, xkp=xkp, use_real_wfn=use_real_wfn, &
1195 kp_range=kp_range, nkp_groups=nkp_groups, kp_dist=kp_dist, &
1196 sab_nl=sab_nl, cell_to_index=cell_to_index)
1197 kplocal = kp_range(2) - kp_range(1) + 1
1199 IF (use_real_wfn)
THEN
1200 CALL timestop(handle)
1204 kp => kpoints%kp_env(1)%kpoint_env
1205 nspin =
SIZE(kp%mos, 2)
1206 CALL get_mo_set(kp%mos(1, 1), nao=nao, nmo=nmo, mo_coeff=mo_coeff)
1209 WRITE (unit=io_unit, fmt=
"(/,T2,A)") &
1210 "Using previous wavefunctions as initial guess for k-points (with reorthogonalization)"
1214 ALLOCATE (rmatrix, cmatrix_db, tmpmat)
1215 CALL dbcsr_create(rmatrix, template=matrix_s_kp(1, 1)%matrix, matrix_type=dbcsr_type_symmetric)
1216 CALL dbcsr_create(cmatrix_db, template=matrix_s_kp(1, 1)%matrix, matrix_type=dbcsr_type_antisymmetric)
1217 CALL dbcsr_create(tmpmat, template=matrix_s_kp(1, 1)%matrix, matrix_type=dbcsr_type_no_symmetry)
1222 CALL mpools_get(mpools_kp, ao_ao_fm_pools=ao_ao_fm_pools_kp)
1227 CALL cp_cfm_create(cfm_nao_nmo_work, mo_coeff%matrix_struct)
1229 NULLIFY (nmo_nmo_struct)
1231 nrow_global=nmo, ncol_global=nmo)
1235 para_env => kpoints%blacs_env_all%para_env
1236 ALLOCATE (info(kplocal*nkp_groups, 2))
1238 ALLOCATE (csmat_cur(kplocal))
1246 DO igroup = 1, nkp_groups
1247 ik = kp_dist(1, igroup) + ikp - 1
1248 my_kpgrp = (ik >= kp_range(1) .AND. ik <= kp_range(2))
1253 CALL rskp_transform(rmatrix=rmatrix, cmatrix=cmatrix_db, rsmat=matrix_s_kp, &
1254 ispin=1, xkp=xkp(1:3, ik), cell_to_index=cell_to_index, sab_nl=sab_nl)
1273 DO igroup = 1, nkp_groups
1274 ik = kp_dist(1, igroup) + ikp - 1
1275 my_kpgrp = (ik >= kp_range(1) .AND. ik <= kp_range(2))
1286 DO indx = 1, kplocal*nkp_groups
1292 ALLOCATE (eigenvalues(nmo))
1293 eval_thresh = 1.0e-12_dp
1296 kp => kpoints%kp_env(ikp)%kpoint_env
1298 CALL get_mo_set(kp%mos(1, ispin), mo_coeff=rmos)
1299 CALL get_mo_set(kp%mos(2, ispin), mo_coeff=imos)
1303 csmat_cur(ikp), cmos_new,
z_zero, cfm_nao_nmo_work)
1305 cmos_new, cfm_nao_nmo_work,
z_zero, csc_cfm)
1308 IF (chol_info == 0)
THEN
1316 ALLOCATE (col_scaling(nmo))
1318 IF (eigenvalues(j) > eval_thresh)
THEN
1319 col_scaling(j) = cmplx(1.0_dp/sqrt(eigenvalues(j)), 0.0_dp, kind=
dp)
1325 DEALLOCATE (col_scaling)
1335 DEALLOCATE (eigenvalues)
1342 matrix_s_kp(1, 1)%matrix, sab_nl, scf_env%scf_work1)
1350 DEALLOCATE (csmat_cur)
1360 CALL timestop(handle)
1361 END SUBROUTINE wfi_use_prev_wf_kp
1374 SUBROUTINE wfi_extrapolate_ps_aspc_kp(wf_history, qs_env, nvec, io_unit, print_level)
1377 INTEGER,
INTENT(IN) :: nvec, io_unit, print_level
1379 CHARACTER(len=*),
PARAMETER :: routinen =
'wfi_extrapolate_ps_aspc_kp'
1381 INTEGER :: handle, i, ikp, ispin, kplocal, &
1382 method_nr, nao, nmo, nspin
1383 INTEGER,
DIMENSION(2) :: kp_range
1384 LOGICAL :: use_real_wfn
1385 REAL(kind=
dp) :: alpha_coeff
1386 TYPE(
cp_cfm_type) :: cfm_nao_nmo_work, cmos_1, cmos_i, &
1389 TYPE(
cp_fm_type),
POINTER :: imos, mo_coeff, rmos
1394 method_nr = wf_history%interpolation_method_nr
1396 CALL timeset(routinen, handle)
1397 NULLIFY (kpoints, kp, mo_coeff, rmos, imos, t0_state, t1_state, nmo_nmo_struct)
1400 CALL get_kpoint_info(kpoints, use_real_wfn=use_real_wfn, kp_range=kp_range)
1401 kplocal = kp_range(2) - kp_range(1) + 1
1403 IF (use_real_wfn)
THEN
1405 CALL cp_warn(__location__,
"ASPC with k-points requires complex wavefunctions; "// &
1406 "falling back to USE_PREV_WF.")
1408 CALL cp_warn(__location__,
"PS with k-points requires complex wavefunctions; "// &
1409 "falling back to USE_PREV_WF.")
1411 CALL wfi_use_prev_wf_kp(qs_env, io_unit, print_level)
1412 CALL timestop(handle)
1416 kp => kpoints%kp_env(1)%kpoint_env
1417 nspin =
SIZE(kp%mos, 2)
1418 CALL get_mo_set(kp%mos(1, 1), nao=nao, nmo=nmo, mo_coeff=mo_coeff)
1422 WRITE (unit=io_unit, fmt=
"(/,T2,A,/,T3,A,I0)") &
1423 "Parameters for the always stable predictor-corrector (ASPC) method:", &
1424 "ASPC order: ", max(nvec - 2, 0)
1429 CALL cp_cfm_create(cmos_new, mo_coeff%matrix_struct, set_zero=.true.)
1430 CALL cp_cfm_create(cmos_1, mo_coeff%matrix_struct, set_zero=.true.)
1431 CALL cp_cfm_create(cmos_i, mo_coeff%matrix_struct, set_zero=.true.)
1432 CALL cp_cfm_create(cfm_nao_nmo_work, mo_coeff%matrix_struct, set_zero=.true.)
1437 CALL cp_cfm_create(cfm_nao_nmo_work, mo_coeff%matrix_struct)
1441 nrow_global=nmo, ncol_global=nmo)
1443 CALL cp_cfm_create(csc_cfm, nmo_nmo_struct, set_zero=.true.)
1452 alpha_coeff = real(4*nvec - 2, kind=
dp)/real(nvec + 1, kind=
dp)
1454 WRITE (unit=io_unit, fmt=
"(T3,A2,I0,A4,F10.6)")
"B(", 1,
") = ", alpha_coeff
1461 kp => kpoints%kp_env(ikp)%kpoint_env
1463 CALL get_mo_set(kp%mos(1, ispin), mo_coeff=rmos)
1464 CALL get_mo_set(kp%mos(2, ispin), mo_coeff=imos)
1465 CALL cp_fm_to_fm(t1_state%wf_kp(ikp, 1, ispin), rmos)
1466 CALL cp_fm_to_fm(t1_state%wf_kp(ikp, 2, ispin), imos)
1476 alpha_coeff = (-1.0_dp)**(i + 1)*real(i, kind=
dp)* &
1479 WRITE (unit=io_unit, fmt=
"(T3,A2,I0,A4,F10.6)")
"B(", i,
") = ", alpha_coeff
1482 alpha_coeff = -1.0_dp*alpha_coeff*real(nvec - i + 1,
dp)/real(i,
dp)
1486 kp => kpoints%kp_env(ikp)%kpoint_env
1488 CALL cp_fm_to_cfm(t1_state%wf_kp(ikp, 1, ispin), t1_state%wf_kp(ikp, 2, ispin), cmos_1)
1489 CALL cp_fm_to_cfm(t0_state%wf_kp(ikp, 1, ispin), t0_state%wf_kp(ikp, 2, ispin), cmos_i)
1493 t0_state%overlap_cfm_kp(ikp), cmos_1,
z_zero, cfm_nao_nmo_work)
1495 cmos_i, cfm_nao_nmo_work,
z_zero, csc_cfm)
1497 cmos_i, csc_cfm,
z_zero, cfm_nao_nmo_work)
1499 CALL get_mo_set(kp%mos(1, ispin), mo_coeff=rmos)
1500 CALL get_mo_set(kp%mos(2, ispin), mo_coeff=imos)
1517 CALL wfi_use_prev_wf_kp(qs_env, 0, print_level)
1519 CALL timestop(handle)
1521 END SUBROUTINE wfi_extrapolate_ps_aspc_kp
1532 ELEMENTAL SUBROUTINE wfi_set_history_variables(qs_env, nvec)
1534 INTEGER,
INTENT(IN) :: nvec
1536 IF (nvec >= qs_env%wf_history%memory_depth)
THEN
1537 IF ((qs_env%scf_control%max_scf_hist /= 0) .AND. (qs_env%scf_control%eps_scf_hist /= 0))
THEN
1538 qs_env%scf_control%max_scf = qs_env%scf_control%max_scf_hist
1539 qs_env%scf_control%eps_scf = qs_env%scf_control%eps_scf_hist
1540 qs_env%scf_control%outer_scf%have_scf = .false.
1541 ELSE IF (qs_env%scf_control%max_scf_hist /= 0)
THEN
1542 qs_env%scf_control%max_scf = qs_env%scf_control%max_scf_hist
1543 qs_env%scf_control%outer_scf%have_scf = .false.
1544 ELSE IF (qs_env%scf_control%eps_scf_hist /= 0)
THEN
1545 qs_env%scf_control%eps_scf = qs_env%scf_control%eps_scf_hist
1546 qs_env%scf_control%outer_scf%eps_scf = qs_env%scf_control%eps_scf_hist
1550 END SUBROUTINE wfi_set_history_variables
1564 REAL(kind=
dp),
INTENT(in) :: dt
1566 cpassert(
ASSOCIATED(wf_history))
1567 cpassert(wf_history%ref_count > 0)
1568 cpassert(
ASSOCIATED(qs_env))
1570 wf_history%snapshot_count = wf_history%snapshot_count + 1
1571 IF (wf_history%memory_depth > 0)
THEN
1572 wf_history%last_state_index =
modulo(wf_history%snapshot_count, &
1573 wf_history%memory_depth) + 1
1574 CALL wfs_update(snapshot=wf_history%past_states &
1575 (wf_history%last_state_index)%snapshot, wf_history=wf_history, &
1576 qs_env=qs_env, dt=dt)
1592 INTEGER,
INTENT(in),
OPTIONAL :: n_col
1594 CHARACTER(len=*),
PARAMETER :: routinen =
'reorthogonalize_vectors'
1596 INTEGER :: handle, my_n_col
1597 LOGICAL :: has_unit_metric, &
1598 ortho_contains_cholesky, &
1607 NULLIFY (scf_env, scf_control, maxao_maxmo_fm_pool, matrix_s, mpools, dft_control)
1608 CALL timeset(routinen, handle)
1610 cpassert(
ASSOCIATED(qs_env))
1613 IF (
PRESENT(n_col)) my_n_col = n_col
1616 scf_control=scf_control, &
1617 matrix_s=matrix_s, &
1618 dft_control=dft_control)
1619 CALL mpools_get(mpools, maxao_maxmo_fm_pool=maxao_maxmo_fm_pool)
1620 IF (
ASSOCIATED(scf_env))
THEN
1621 ortho_contains_cholesky = (scf_env%method /=
ot_method_nr) .AND. &
1622 (scf_env%cholesky_method > 0) .AND. &
1623 ASSOCIATED(scf_env%ortho)
1625 ortho_contains_cholesky = .false.
1628 CALL get_qs_env(qs_env, has_unit_metric=has_unit_metric)
1629 smearing_is_used = .false.
1630 IF (dft_control%smear)
THEN
1631 smearing_is_used = .true.
1634 IF (has_unit_metric)
THEN
1636 ELSE IF (smearing_is_used)
THEN
1638 matrix_s=matrix_s(1)%matrix)
1639 ELSE IF (ortho_contains_cholesky)
THEN
1641 ortho=scf_env%ortho)
1645 CALL timestop(handle)
1659 CHARACTER(len=*),
PARAMETER :: routinen =
'wfi_purge_history'
1661 INTEGER :: handle, io_unit, print_level
1666 NULLIFY (dft_control, wf_history)
1668 CALL timeset(routinen, handle)
1670 print_level = logger%iter_info%print_level
1672 extension=
".scfLog")
1674 cpassert(
ASSOCIATED(qs_env))
1675 cpassert(
ASSOCIATED(qs_env%wf_history))
1676 cpassert(qs_env%wf_history%ref_count > 0)
1677 CALL get_qs_env(qs_env, dft_control=dft_control)
1679 SELECT CASE (qs_env%wf_history%interpolation_method_nr)
1687 IF (qs_env%wf_history%snapshot_count >= 2)
THEN
1688 IF (debug_this_module .AND. io_unit > 0) &
1689 WRITE (io_unit, fmt=
"(T2,A)")
"QS| Purging WFN history"
1690 CALL wfi_create(wf_history, interpolation_method_nr= &
1691 dft_control%qs_control%wf_interpolation_method_nr, &
1692 extrapolation_order=dft_control%qs_control%wf_extrapolation_order, &
1693 has_unit_metric=qs_env%has_unit_metric)
1695 wf_history=wf_history)
1697 CALL wfi_update(qs_env%wf_history, qs_env=qs_env, dt=1.0_dp)
1700 cpabort(
"Unknown extrapolation method.")
1702 CALL timestop(handle)
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public vandevondele2005a
integer, save, public kuhne2007
integer, save, public kolafa2004
Basic linear algebra operations for complex full matrices.
subroutine, public cp_cfm_scale_and_add(alpha, matrix_a, beta, matrix_b)
Scale and add two BLACS matrices (a = alpha*a + beta*b).
subroutine, public cp_cfm_gemm(transa, transb, m, n, k, alpha, matrix_a, matrix_b, beta, matrix_c, a_first_col, a_first_row, b_first_col, b_first_row, c_first_col, c_first_row)
Performs one of the matrix-matrix operations: matrix_c = alpha * op1( matrix_a ) * op2( matrix_b ) + ...
subroutine, public cp_cfm_scale_and_add_fm(alpha, matrix_a, beta, matrix_b)
Scale and add two BLACS matrices (a = alpha*a + beta*b). where b is a real matrix (adapted from cp_cf...
subroutine, public cp_cfm_triangular_multiply(triangular_matrix, matrix_b, side, transa_tr, invert_tr, uplo_tr, unit_diag_tr, n_rows, n_cols, alpha)
Multiplies in place by a triangular matrix: matrix_b = alpha op(triangular_matrix) matrix_b or (if si...
subroutine, public cp_cfm_column_scale(matrix_a, scaling)
Scales columns of the full matrix by corresponding factors.
various cholesky decomposition related routines
subroutine, public cp_cfm_cholesky_decompose(matrix, n, info_out)
Used to replace a symmetric positive definite matrix M with its Cholesky decomposition U: M = U^T * U...
used for collecting diagonalization schemes available for cp_cfm_type
subroutine, public cp_cfm_heevd(matrix, eigenvectors, eigenvalues)
Perform a diagonalisation of a complex matrix.
Represents a complex full matrix distributed on many processors.
subroutine, public cp_cfm_release(matrix)
Releases a full matrix.
subroutine, public cp_fm_to_cfm(msourcer, msourcei, mtarget)
Construct a complex full matrix by taking its real and imaginary parts from two separate real-value f...
subroutine, public cp_cfm_create(matrix, matrix_struct, name, nrow, ncol, set_zero)
Creates a new full matrix with the given structure.
subroutine, public cp_cfm_set_all(matrix, alpha, beta)
Set all elements of the full matrix to alpha. Besides, set all diagonal matrix elements to beta (if g...
subroutine, public cp_cfm_to_fm(msource, mtargetr, mtargeti)
Copy real and imaginary parts of a complex full matrix into separate real-value full matrices.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_deallocate_matrix(matrix)
...
subroutine, public dbcsr_desymmetrize(matrix_a, matrix_b)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_add(matrix_a, matrix_b, alpha_scalar, beta_scalar)
...
Routines that link DBCSR and CP2K concepts together.
subroutine, public cp_dbcsr_alloc_block_from_nbl(matrix, sab_orb, desymmetrize)
allocate the blocks of a dbcsr based on the neighbor list
DBCSR operations in CP2K.
subroutine, public cp_dbcsr_sm_fm_multiply(matrix, fm_in, fm_out, ncol, alpha, beta)
multiply a dbcsr with a fm matrix
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_scale_and_add(alpha, matrix_a, beta, matrix_b)
calc A <- alpha*A + beta*B optimized for alpha == 1.0 (just add beta*B) and beta == 0....
subroutine, public cp_fm_scale(alpha, matrix_a)
scales a matrix matrix_a = alpha * matrix_b
pool for for elements that are retained and released
subroutine, public fm_pool_create_fm(pool, element, name)
returns an element, allocating it if none is in the pool
subroutine, public fm_pool_give_back_fm(pool, element)
returns the element to the pool
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
subroutine, public cp_fm_start_copy_general(source, destination, para_env, info)
Initiates the copy operation: get distribution data, post MPI isend and irecvs.
subroutine, public cp_fm_cleanup_copy_general(info)
Completes the copy operation: wait for comms clean up MPI state.
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
subroutine, public cp_fm_finish_copy_general(destination, info)
Completes the copy operation: wait for comms, unpack, clean up MPI state.
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)
...
integer, parameter, public low_print_level
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,...
Defines the basic variable types.
integer, parameter, public dp
Routines needed for kpoint calculation.
subroutine, public rskp_transform(rmatrix, cmatrix, rsmat, ispin, xkp, cell_to_index, sab_nl, is_complex, rs_sign)
Transformation of real space matrices to a kpoint.
subroutine, public kpoint_density_transform(kpoint, denmat, wtype, tempmat, sab_nl, fmwork, for_aux_fit, pmat_ext)
generate real space density matrices in DBCSR format
subroutine, public kpoint_density_matrices(kpoint, energy_weighted, for_aux_fit)
Calculate kpoint density matrices (rho(k), owned by kpoint groups)
subroutine, public kpoint_set_mo_occupation(kpoint, smear, probe)
Given the eigenvalues of all kpoints, calculates the occupation numbers.
Types and basic routines needed for a kpoint calculation.
subroutine, public get_kpoint_info(kpoint, kp_scheme, nkp_grid, kp_shift, symmetry, verbose, full_grid, use_real_wfn, eps_geo, parallel_group_size, kp_range, nkp, xkp, wkp, para_env, blacs_env_all, para_env_kp, para_env_inter_kp, blacs_env, kp_env, kp_aux_env, mpools, iogrp, nkp_groups, kp_dist, cell_to_index, index_to_cell, sab_nl, sab_nl_nosym)
Retrieve information from a kpoint environment.
Definition of mathematical constants and functions.
complex(kind=dp), parameter, public z_one
complex(kind=dp), parameter, public gaussi
complex(kind=dp), parameter, public z_zero
Collection of simple mathematical functions and subroutines.
elemental real(kind=dp) function, public binomial(n, k)
The binomial coefficient n over k for 0 <= k <= n is calculated, otherwise zero is returned.
Interface to the message passing library MPI.
basic linear algebra operations for full matrixes
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 ...
collects routines that calculate density matrices
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.
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...
wrapper for the pools of matrixes
subroutine, public mpools_get(mpools, ao_mo_fm_pools, ao_ao_fm_pools, mo_mo_fm_pools, ao_mosub_fm_pools, mosub_mosub_fm_pools, maxao_maxmo_fm_pool, maxao_maxao_fm_pool, maxmo_maxmo_fm_pool)
returns various attributes of the mpools (notably the pools contained in it)
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_lowdin(vmatrix, ncol, matrix_s)
return a set of S orthonormal vectors (C^T S C == 1) where a Loedwin transformation is applied to kee...
subroutine, public make_basis_cholesky(vmatrix, ncol, ortho)
return a set of S orthonormal vectors (C^T S C == 1) where the cholesky decomposed form of S is passe...
subroutine, public make_basis_sm(vmatrix, ncol, matrix_s)
returns an S-orthonormal basis v (v^T S v ==1)
Definition and initialisation of the mo data type.
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.
Define the neighbor list data types and the corresponding functionality.
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_set(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)
...
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...
module that contains the definitions of the scf types
integer, parameter, public ot_method_nr
Storage of past states of the qs_env. Methods to interpolate (or actually normally extrapolate) the n...
subroutine, public wfi_create_for_kp(wf_history)
Adapts wf_history storage flags for k-point calculations. For ASPC, switches from Gamma WFN storage t...
subroutine, public wfi_purge_history(qs_env)
purges wf_history retaining only the latest snapshot
character(len=30) function, public wfi_get_method_label(method_nr)
returns a string describing the interpolation method
subroutine, public reorthogonalize_vectors(qs_env, v_matrix, n_col)
reorthogonalizes the mos
subroutine, public wfi_update(wf_history, qs_env, dt)
updates the snapshot buffer, taking a new snapshot
subroutine, public wfi_extrapolate(wf_history, qs_env, dt, extrapolation_method_nr, orthogonal_wf)
calculates the new starting state for the scf for the next wf optimization
subroutine, public wfi_create(wf_history, interpolation_method_nr, extrapolation_order, has_unit_metric)
...
interpolate the wavefunctions to speed up the convergence when doing MD
type(qs_wf_snapshot_type) function, pointer, public wfi_get_snapshot(wf_history, wf_index)
returns a snapshot, the first being the latest snapshot
subroutine, public wfi_release(wf_history)
releases a wf_history of a wavefunction (see doc/ReferenceCounting.html)
parameters that control an scf iteration
Represent a complex full matrix.
to create arrays of pools
represent a pool of elements with the same structure
keeps the information about the structure of a full matrix
Stores the state of a copy between cp_fm_start_copy_general and cp_fm_finish_copy_general.
type of a logger, at the moment it contains just a print level starting at which level it should be l...
Keeps information about a specific k-point.
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 the pools of matrixes used by qs
keeps the density in various representations, keeping track of which ones are valid.
keeps track of the previous wavefunctions and can extrapolate them for the next step of md
represent a past snapshot of the wavefunction. some elements might not be associated (to spare memory...