72 USE iso_c_binding,
ONLY: c_null_char
174#include "./base/base_uses.f90"
179 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_active_space_methods'
184 INTEGER :: unit_nr = -1, bra_start = -1, ket_start = -1
186 PROCEDURE :: func => eri_fcidump_print_func
187 END TYPE eri_fcidump_print
190 INTEGER :: bra_start = 0, ket_start = 0
191 REAL(KIND=
dp) :: checksum = 0.0_dp
194 PROCEDURE :: func => eri_fcidump_checksum_func
195 END TYPE eri_fcidump_checksum
206 CLASS(eri_fcidump_checksum) :: this
207 INTEGER,
INTENT(IN) :: bra_start, ket_start
208 this%bra_start = bra_start
209 this%ket_start = ket_start
219 CHARACTER(len=*),
PARAMETER :: routinen =
'active_space_main'
221 CHARACTER(len=10) :: cshell, lnam(5)
222 CHARACTER(len=default_path_length) :: qcschema_filename
223 CHARACTER(LEN=default_string_length) :: basis_type, kp_scheme
224 INTEGER :: as_solver, eri_method, eri_operator, eri_print, group_size, handle, i, iatom, &
225 ishell, isp, ispin, iw, j, jm, m, max_orb_ind, mselect, n1, n2, nao, natom, nel, &
226 nelec_active, nelec_inactive, nelec_total, nkp, nmo, nmo_active, nmo_available, &
227 nmo_inactive, nmo_inactive_remaining, nmo_occ, nmo_virtual, nn1, nn2, nrow_global, nspins
228 INTEGER,
DIMENSION(5) :: nshell
229 INTEGER,
DIMENSION(:),
POINTER :: invals
230 LOGICAL :: do_ddapc, do_kpoints, ex_omega, &
231 ex_operator, ex_perd, ex_rcut, &
232 explicit, stop_after_print, store_wfn, &
234 REAL(kind=
dp) :: eri_eps_filter, eri_eps_grid, eri_eps_int, eri_gpw_cutoff, eri_op_omega, &
235 eri_rcut, eri_rel_cutoff, fel, focc, maxocc, nze_percentage
236 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: eigenvalues
237 REAL(kind=
dp),
DIMENSION(:),
POINTER :: evals_virtual
244 TYPE(
cp_fm_type),
POINTER :: fm_target_active, fm_target_inactive, &
245 fmat, mo_coeff, mo_ref, mo_target
247 TYPE(dbcsr_csr_type),
POINTER :: eri_mat
248 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: ks_matrix, rho_ao, s_matrix
249 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: ks_matrix_kp, rho_ao_kp, s_matrix_kp
254 TYPE(
mo_set_type),
POINTER :: mo_set, mo_set_active, mo_set_inactive
260 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
265 hfx_section, input, loc_print, &
266 loc_section, print_orb, xc_section
273 IF (.NOT. explicit)
RETURN
274 CALL timeset(routinen, handle)
280 WRITE (iw,
'(/,T2,A)') &
281 '!-----------------------------------------------------------------------------!'
282 WRITE (iw,
'(T26,A)')
"Active Space Embedding Module"
283 WRITE (iw,
'(T2,A)') &
284 '!-----------------------------------------------------------------------------!'
289 CALL get_qs_env(qs_env, do_kpoints=do_kpoints, dft_control=dft_control, kpoints=kpoints)
291 IF (.NOT.
ASSOCIATED(kpoints)) &
292 CALL cp_abort(__location__,
"Missing Gamma-point environment for active space module")
293 CALL get_kpoint_info(kpoints, kp_scheme=kp_scheme, nkp=nkp, use_real_wfn=use_real_wfn)
294 IF (trim(kp_scheme) /=
"GAMMA" .OR. nkp /= 1 .OR. .NOT. use_real_wfn)
THEN
295 CALL cp_abort(__location__, &
296 "Only Gamma-point DFT%KPOINTS are supported in the active space module")
298 IF (.NOT.
ASSOCIATED(kpoints%kp_env)) &
299 CALL cp_abort(__location__,
"Missing Gamma-point environment for active space module")
300 IF (.NOT.
ASSOCIATED(kpoints%kp_env(1)%kpoint_env)) &
301 CALL cp_abort(__location__,
"Missing Gamma-point environment for active space module")
302 IF (.NOT.
ASSOCIATED(kpoints%kp_env(1)%kpoint_env%mos)) &
303 CALL cp_abort(__location__,
"Missing Gamma-point MOs for active space module")
310 CALL cp_abort(__location__,
"Adiabatic rescaling not supported in active space module")
314 do_ddapc = dft_control%qs_control%ddapc_restraint .OR. &
315 qs_env%cp_ddapc_ewald%do_decoupling .OR. &
316 qs_env%cp_ddapc_ewald%do_qmmm_periodic_decpl .OR. &
317 qs_env%cp_ddapc_ewald%do_solvation
319 CALL cp_abort(__location__,
"DDAPC charges are not supported in the active space module")
321 IF (dft_control%do_sccs)
THEN
322 CALL cp_abort(__location__,
"SCCS is not supported in the active space module")
324 IF (dft_control%correct_surf_dip)
THEN
325 IF (dft_control%surf_dip_correct_switch)
THEN
326 CALL cp_abort(__location__,
"Surface dipole correction not supported in the AS module")
329 IF (dft_control%smeagol_control%smeagol_enabled)
THEN
330 CALL cp_abort(__location__,
"SMEAGOL is not supported in the active space module")
332 IF (dft_control%qs_control%do_kg)
THEN
333 CALL cp_abort(__location__,
"KG correction not supported in the active space module")
336 NULLIFY (active_space_env)
338 active_space_env%energy_total = 0.0_dp
339 active_space_env%energy_ref = 0.0_dp
340 active_space_env%energy_inactive = 0.0_dp
341 active_space_env%energy_active = 0.0_dp
347 active_space_env%fcidump = .true.
352 active_space_env%qcschema = .true.
353 active_space_env%qcschema_filename = qcschema_filename
357 CALL get_qs_env(qs_env, nelectron_total=nelec_total)
359 IF (nelec_active <= 0) cpabort(
"Specify a positive number of active electrons.")
360 IF (nelec_active > nelec_total) cpabort(
"More active electrons than total electrons.")
362 nelec_inactive = nelec_total - nelec_active
363 IF (mod(nelec_inactive, 2) /= 0)
THEN
364 cpabort(
"The remaining number of inactive electrons has to be even.")
368 WRITE (iw,
'(T3,A,T70,I10)')
"Total number of electrons", nelec_total
369 WRITE (iw,
'(T3,A,T70,I10)')
"Number of inactive electrons", nelec_inactive
370 WRITE (iw,
'(T3,A,T70,I10)')
"Number of active electrons", nelec_active
373 CALL get_qs_env(qs_env, dft_control=dft_control)
374 nspins = dft_control%nspins
376 active_space_env%nelec_active = nelec_active
377 active_space_env%nelec_inactive = nelec_inactive
378 active_space_env%nelec_total = nelec_total
379 active_space_env%nspins = nspins
380 active_space_env%multiplicity = dft_control%multiplicity
381 active_space_env%restricted_orbitals = dft_control%roks
385 IF (.NOT. explicit)
THEN
386 CALL cp_abort(__location__,
"Number of Active Orbitals has to be specified.")
388 active_space_env%nmo_active = nmo_active
390 nmo_inactive = nelec_inactive/2
391 active_space_env%nmo_inactive = nmo_inactive
397 SELECT CASE (mselect)
399 cpabort(
"Unknown orbital selection method")
401 WRITE (iw,
'(/,T3,A)') &
402 "Active space orbitals selected using energy ordered canonical orbitals"
404 WRITE (iw,
'(/,T3,A)') &
405 "Active space orbitals selected using projected Wannier orbitals"
407 WRITE (iw,
'(/,T3,A)') &
408 "Active space orbitals selected using modified atomic orbitals (MAO)"
410 WRITE (iw,
'(/,T3,A)') &
411 "Active space orbitals selected manually"
414 WRITE (iw,
'(T3,A,T70,I10)')
"Number of inactive orbitals", nmo_inactive
415 WRITE (iw,
'(T3,A,T70,I10)')
"Number of active orbitals", nmo_active
422 IF (iatom <= 0 .OR. iatom > natom)
THEN
424 WRITE (iw,
'(/,T3,A,I3)')
"ERROR: SUBSPACE_ATOM number is not valid", iatom
426 cpabort(
"Select a valid SUBSPACE_ATOM")
433 cshell = adjustl(cshell)
437 IF (cshell(n1:n1) ==
" ")
THEN
441 READ (cshell(n1:),
"(I1,A1)") nshell(i), lnam(i)
447 SELECT CASE (mselect)
449 cpabort(
"Unknown orbital selection method")
452 mos => kpoints%kp_env(1)%kpoint_env%mos(1, :)
458 nmo_occ = nmo_inactive + nmo_active
461 ALLOCATE (active_space_env%inactive_orbitals(nmo_inactive, nspins))
463 DO i = 1, nmo_inactive
464 active_space_env%inactive_orbitals(i, ispin) = i
469 ALLOCATE (active_space_env%active_orbitals(nmo_active, nspins))
472 active_space_env%active_orbitals(i, ispin) = nmo_inactive + i
481 IF (nspins > 1) maxocc = 1.0_dp
482 ALLOCATE (active_space_env%mos_active(nspins))
483 ALLOCATE (active_space_env%mos_inactive(nspins))
485 CALL get_mo_set(mos(ispin), mo_coeff=mo_ref, nao=nao)
486 CALL cp_fm_get_info(mo_ref, context=context, para_env=para_env, nrow_global=nrow_global)
488 CALL allocate_mo_set(active_space_env%mos_active(ispin), nao, nmo_occ, 0, 0.0_dp, maxocc, 0.0_dp)
490 nrow_global=nrow_global, ncol_global=nmo_occ)
491 CALL init_mo_set(active_space_env%mos_active(ispin), fm_struct=fm_struct_tmp, name=
"Active Space MO")
493 IF (nspins == 2)
THEN
494 nel = nelec_inactive/2
498 CALL allocate_mo_set(active_space_env%mos_inactive(ispin), nao, nmo_occ, nel, &
499 REAL(nel, kind=
dp), maxocc, 0.0_dp)
501 nrow_global=nrow_global, ncol_global=nmo_occ)
502 CALL init_mo_set(active_space_env%mos_inactive(ispin), fm_struct=fm_struct_tmp, name=
"Inactive Space MO")
507 CALL get_qs_env(qs_env, scf_control=scf_control)
508 IF (dft_control%roks .AND. scf_control%roks_scheme /=
high_spin_roks)
THEN
509 CALL cp_abort(__location__, &
510 "Only high-spin ROKS is supported for ACTIVE_SPACE FCI; "// &
511 "general ROKS MO definitions are not implemented.")
513 IF (dft_control%do_admm)
THEN
514 IF (dft_control%do_admm_mo)
THEN
515 cpabort(
"ADMM currently possible only with purification none_dm")
519 ALLOCATE (eigenvalues(nmo_occ, nspins))
522 CALL get_qs_env(qs_env, matrix_ks_kp=ks_matrix_kp, matrix_s_kp=s_matrix_kp, &
523 scf_control=scf_control)
524 ks_matrix => ks_matrix_kp(:, 1)
525 s_matrix => s_matrix_kp(:, 1)
527 CALL get_qs_env(qs_env, matrix_ks=ks_matrix, matrix_s=s_matrix, scf_control=scf_control)
532 WRITE (iw,
'(/,T3,A)')
"Calculating virtual MOs..."
538 CALL get_mo_set(mos(ispin), mo_coeff=mo_ref, nmo=nmo_available)
541 nmo_virtual = nmo_occ - nmo_available
542 nmo_virtual = max(nmo_virtual, 0)
544 NULLIFY (evals_virtual)
545 ALLOCATE (evals_virtual(nmo_virtual))
548 nrow_global=nrow_global)
551 nrow_global=nrow_global, ncol_global=nmo_virtual)
552 CALL cp_fm_create(mo_virtual, fm_struct_tmp, name=
"virtual")
556 NULLIFY (local_preconditioner)
559 CALL ot_eigensolver(matrix_h=ks_matrix(ispin)%matrix, matrix_s=s_matrix(1)%matrix, &
560 matrix_c_fm=mo_virtual, matrix_orthogonal_space_fm=mo_ref, &
561 eps_gradient=scf_control%eps_lumos, &
563 iter_max=scf_control%max_iter_lumos, &
564 size_ortho_space=nmo_available)
573 evals_arg=eigenvalues(:, ispin), do_rotation=.true.)
576 mo_set => active_space_env%mos_inactive(ispin)
578 DO i = 1,
SIZE(active_space_env%inactive_orbitals, 1)
579 m = active_space_env%inactive_orbitals(i, ispin)
581 mo_set%eigenvalues(m) = eigenvalues(m, ispin)
583 mo_set%occupation_numbers(m) = 1.0
585 mo_set%occupation_numbers(m) = 2.0
590 mo_set => active_space_env%mos_active(ispin)
593 IF (nspins == 2)
THEN
595 nel = (nelec_active + active_space_env%multiplicity - 1)/2
597 nel = (nelec_active - active_space_env%multiplicity + 1)/2
602 mo_set%nelectron = nel
603 mo_set%n_el_f = real(nel, kind=
dp)
605 m = active_space_env%active_orbitals(i, ispin)
606 IF (m > nmo_available)
THEN
607 CALL cp_fm_to_fm(mo_virtual, mo_target, 1, m - nmo_available, m)
608 eigenvalues(m, ispin) = evals_virtual(m - nmo_available)
609 mo_set%occupation_numbers(m) = 0.0
612 mo_set%occupation_numbers(m) = mos(ispin)%occupation_numbers(m)
614 mo_set%eigenvalues(m) = eigenvalues(m, ispin)
617 DEALLOCATE (evals_virtual)
624 WRITE (iw,
'(/,T3,A,I3,T66,A)')
"Canonical Orbital Selection for spin", ispin, &
626 DO i = 1, nmo_inactive, 4
627 jm = min(3, nmo_inactive - i)
628 WRITE (iw,
'(T3,4(F14.6,A5))') (eigenvalues(i + j, ispin),
" [I]", j=0, jm)
630 DO i = nmo_inactive + 1, nmo_inactive + nmo_active, 4
631 jm = min(3, nmo_inactive + nmo_active - i)
632 WRITE (iw,
'(T3,4(F14.6,A5))') (eigenvalues(i + j, ispin),
" [A]", j=0, jm)
634 WRITE (iw,
'(/,T3,A,I3)')
"Active Orbital Indices for spin", ispin
635 DO i = 1,
SIZE(active_space_env%active_orbitals, 1), 4
636 jm = min(3,
SIZE(active_space_env%active_orbitals, 1) - i)
637 WRITE (iw,
'(T3,4(I4))') (active_space_env%active_orbitals(i + j, ispin), j=0, jm)
641 DEALLOCATE (eigenvalues)
646 IF (dft_control%roks)
THEN
647 CALL cp_abort(__location__, &
648 "Manual ACTIVE_SPACE orbital selection is not supported for ROKS; "// &
649 "use canonical high-spin ROKS.")
651 IF (dft_control%do_admm)
THEN
656 IF (dft_control%do_admm_mo)
THEN
657 cpabort(
"ADMM currently possible only with purification none_dm")
662 IF (.NOT. explicit)
THEN
663 CALL cp_abort(__location__,
"Manual orbital selection requires to explicitly "// &
664 "set the active orbital indices via ACTIVE_ORBITAL_INDICES")
667 IF (nspins == 1)
THEN
668 cpassert(
SIZE(invals) == nmo_active)
670 cpassert(
SIZE(invals) == 2*nmo_active)
672 ALLOCATE (active_space_env%inactive_orbitals(nmo_inactive, nspins))
673 ALLOCATE (active_space_env%active_orbitals(nmo_active, nspins))
677 active_space_env%active_orbitals(i, ispin) = invals(i + (ispin - 1)*nmo_active)
682 mos => kpoints%kp_env(1)%kpoint_env%mos(1, :)
688 max_orb_ind = maxval(invals)
690 IF (nspins > 1) maxocc = 1.0_dp
691 ALLOCATE (active_space_env%mos_active(nspins))
692 ALLOCATE (active_space_env%mos_inactive(nspins))
695 CALL get_mo_set(mos(ispin), mo_coeff=mo_ref, nao=nao)
696 CALL cp_fm_get_info(mo_ref, context=context, para_env=para_env, nrow_global=nrow_global)
697 CALL allocate_mo_set(active_space_env%mos_active(ispin), nao, max_orb_ind, 0, 0.0_dp, maxocc, 0.0_dp)
699 nrow_global=nrow_global, ncol_global=max_orb_ind)
700 CALL init_mo_set(active_space_env%mos_active(ispin), fm_struct=fm_struct_tmp, name=
"Active Space MO")
704 IF (nspins == 2)
THEN
705 nel = nelec_inactive/2
709 CALL allocate_mo_set(active_space_env%mos_inactive(ispin), nao, max_orb_ind, nel, real(nel, kind=
dp), maxocc, 0.0_dp)
711 nrow_global=nrow_global, ncol_global=max_orb_ind)
712 CALL init_mo_set(active_space_env%mos_inactive(ispin), fm_struct=fm_struct_tmp, name=
"Inactive Space MO")
714 active_space_env%mos_inactive(ispin)%occupation_numbers = 0.0_dp
718 ALLOCATE (eigenvalues(max_orb_ind, nspins))
721 CALL get_qs_env(qs_env, matrix_ks_kp=ks_matrix_kp, matrix_s_kp=s_matrix_kp, &
722 scf_control=scf_control)
723 ks_matrix => ks_matrix_kp(:, 1)
724 s_matrix => s_matrix_kp(:, 1)
726 CALL get_qs_env(qs_env, matrix_ks=ks_matrix, matrix_s=s_matrix, scf_control=scf_control)
731 WRITE (iw,
'(/,T3,A)')
"Calculating virtual MOs..."
734 CALL get_mo_set(mos(ispin), mo_coeff=mo_ref, nmo=nmo_available)
735 nmo_virtual = max_orb_ind - nmo_available
736 nmo_virtual = max(nmo_virtual, 0)
738 NULLIFY (evals_virtual)
739 ALLOCATE (evals_virtual(nmo_virtual))
742 nrow_global=nrow_global)
745 nrow_global=nrow_global, ncol_global=nmo_virtual)
746 CALL cp_fm_create(mo_virtual, fm_struct_tmp, name=
"virtual")
750 NULLIFY (local_preconditioner)
752 CALL ot_eigensolver(matrix_h=ks_matrix(ispin)%matrix, matrix_s=s_matrix(1)%matrix, &
753 matrix_c_fm=mo_virtual, matrix_orthogonal_space_fm=mo_ref, &
754 eps_gradient=scf_control%eps_lumos, &
756 iter_max=scf_control%max_iter_lumos, &
757 size_ortho_space=nmo_available)
767 evals_arg=eigenvalues(:, ispin), do_rotation=.true.)
769 mo_set_active => active_space_env%mos_active(ispin)
770 CALL get_mo_set(mo_set_active, mo_coeff=fm_target_active)
771 mo_set_inactive => active_space_env%mos_inactive(ispin)
772 CALL get_mo_set(mo_set_inactive, mo_coeff=fm_target_inactive)
775 nmo_inactive_remaining = nmo_inactive
776 DO i = 1, max_orb_ind
778 IF (any(active_space_env%active_orbitals(:, ispin) == i))
THEN
779 IF (i > nmo_available)
THEN
780 CALL cp_fm_to_fm(mo_virtual, fm_target_active, 1, i - nmo_available, i)
781 eigenvalues(i, ispin) = evals_virtual(i - nmo_available)
782 mo_set_active%occupation_numbers(i) = 0.0
784 CALL cp_fm_to_fm(fm_dummy, fm_target_active, 1, i, i)
785 mo_set_active%occupation_numbers(i) = mos(ispin)%occupation_numbers(i)
787 mo_set_active%eigenvalues(i) = eigenvalues(i, ispin)
789 ELSEIF (nmo_inactive_remaining > 0)
THEN
790 CALL cp_fm_to_fm(fm_dummy, fm_target_inactive, 1, i, i)
792 active_space_env%inactive_orbitals(nmo_inactive - nmo_inactive_remaining + 1, ispin) = i
793 mo_set_inactive%eigenvalues(i) = eigenvalues(i, ispin)
794 mo_set_inactive%occupation_numbers(i) = mos(ispin)%occupation_numbers(i)
796 IF (nmo_inactive_remaining == 1)
THEN
797 mo_set_inactive%homo = i
798 mo_set_inactive%lfomo = i + 1
800 nmo_inactive_remaining = nmo_inactive_remaining - 1
807 DEALLOCATE (evals_virtual)
814 WRITE (iw,
'(/,T3,A,I3,T66,A)')
"Orbital Energies and Selection for spin", ispin,
"[atomic units]"
816 DO i = 1, max_orb_ind, 4
817 jm = min(3, max_orb_ind - i)
818 WRITE (iw,
'(T4)', advance=
"no")
820 IF (any(active_space_env%active_orbitals(:, ispin) == i + j))
THEN
821 WRITE (iw,
'(T3,F12.6,A5)', advance=
"no") eigenvalues(i + j, ispin),
" [A]"
822 ELSEIF (any(active_space_env%inactive_orbitals(:, ispin) == i + j))
THEN
823 WRITE (iw,
'(T3,F12.6,A5)', advance=
"no") eigenvalues(i + j, ispin),
" [I]"
825 WRITE (iw,
'(T3,F12.6,A5)', advance=
"no") eigenvalues(i + j, ispin),
" [V]"
830 WRITE (iw,
'(/,T3,A,I3)')
"Active Orbital Indices for spin", ispin
831 DO i = 1,
SIZE(active_space_env%active_orbitals, 1), 4
832 jm = min(3,
SIZE(active_space_env%active_orbitals, 1) - i)
833 WRITE (iw,
'(T3,4(I4))') (active_space_env%active_orbitals(i + j, ispin), j=0, jm)
837 DEALLOCATE (eigenvalues)
841 NULLIFY (loc_section, loc_print)
843 cpassert(
ASSOCIATED(loc_section))
846 cpabort(
"not yet available")
850 cpabort(
"not yet available")
860 CALL print_orbital_cubes(print_orb, qs_env, active_space_env%mos_active)
862 IF (stop_after_print)
THEN
865 WRITE (iw,
'(/,T2,A)') &
866 '!----------------- Early End of Active Space Interface -----------------------!'
869 CALL timestop(handle)
879 rho_ao => rho_ao_kp(:, 1)
883 cpassert(
ASSOCIATED(rho_ao))
888 mo_set => active_space_env%mos_inactive(ispin)
890 active_space_env%pmat_inactive(ispin)%matrix => denmat
895 active_space_env%eri%method = eri_method
897 active_space_env%eri%operator = eri_operator
899 active_space_env%eri%omega = eri_op_omega
901 active_space_env%eri%cutoff_radius = eri_rcut
904 active_space_env%eri%eps_integral = eri_eps_int
907 IF (
SIZE(invals) == 1)
THEN
908 active_space_env%eri%periodicity(1:3) = invals(1)
910 active_space_env%eri%periodicity(1:3) = invals(1:3)
914 active_space_env%eri%periodicity(1:3) = cell%perd(1:3)
917 WRITE (iw,
'(/,T3,A)')
"Calculation of Electron Repulsion Integrals"
919 SELECT CASE (eri_method)
921 WRITE (iw,
'(T3,A,T50,A)')
"Integration method",
"GPW Fourier transform over MOs"
923 WRITE (iw,
'(T3,A,T44,A)')
"Integration method",
"Half transformed integrals from GPW"
925 cpabort(
"Unknown ERI method")
928 SELECT CASE (eri_operator)
930 WRITE (iw,
'(T3,A,T73,A)')
"ERI operator",
"Coulomb"
933 WRITE (iw,
'(T3,A,T74,A)')
"ERI operator",
"Yukawa"
934 IF (.NOT. ex_omega)
CALL cp_abort(__location__, &
935 "Yukawa operator requires OMEGA to be explicitly set")
936 WRITE (iw,
'(T3,A,T66,F14.3)')
"ERI operator parameter OMEGA", eri_op_omega
939 WRITE (iw,
'(T3,A,T63,A)')
"ERI operator",
"Longrange Coulomb"
940 IF (.NOT. ex_omega)
CALL cp_abort(__location__, &
941 "Longrange operator requires OMEGA to be explicitly set")
942 WRITE (iw,
'(T3,A,T66,F14.3)')
"ERI operator parameter OMEGA", eri_op_omega
945 WRITE (iw,
'(T3,A,T62,A)')
"ERI operator",
"Shortrange Coulomb"
946 IF (.NOT. ex_omega)
CALL cp_abort(__location__, &
947 "Shortrange operator requires OMEGA to be explicitly set")
948 WRITE (iw,
'(T3,A,T66,F14.3)')
"ERI operator parameter OMEGA", eri_op_omega
951 WRITE (iw,
'(T3,A,T63,A)')
"ERI operator",
"Truncated Coulomb"
952 IF (.NOT. ex_rcut)
CALL cp_abort(__location__, &
953 "Cutoff radius not specified for trunc. Coulomb operator")
954 WRITE (iw,
'(T3,A,T66,F14.3)')
"ERI operator cutoff radius (au)", eri_rcut
957 WRITE (iw,
'(T3,A,T53,A)')
"ERI operator",
"Longrange truncated Coulomb"
958 IF (.NOT. ex_rcut)
CALL cp_abort(__location__, &
959 "Cutoff radius not specified for trunc. longrange operator")
960 WRITE (iw,
'(T3,A,T66,F14.3)')
"ERI operator cutoff radius (au)", eri_rcut
961 IF (.NOT. ex_omega)
CALL cp_abort(__location__, &
962 "LR truncated operator requires OMEGA to be explicitly set")
963 WRITE (iw,
'(T3,A,T66,F14.3)')
"ERI operator parameter OMEGA", eri_op_omega
964 IF (eri_op_omega < 0.01_dp)
THEN
965 cpabort(
"LR truncated operator requires OMEGA >= 0.01 to be stable")
969 cpabort(
"Unknown ERI operator")
973 WRITE (iw,
'(T3,A,T68,E12.4)')
"Accuracy of ERIs", eri_eps_int
974 WRITE (iw,
'(T3,A,T71,3I3)')
"Periodicity", active_space_env%eri%periodicity(1:3)
978 WRITE (iw,
'(T3,A,T68,I12)')
"Total Number of ERI", (nmo_active**4)/8
980 WRITE (iw,
'(T3,A,T68,I12)')
"Total Number of ERI (aa|aa)", (nmo_active**4)/8
981 WRITE (iw,
'(T3,A,T68,I12)')
"Total Number of ERI (bb|bb)", (nmo_active**4)/8
982 WRITE (iw,
'(T3,A,T68,I12)')
"Total Number of ERI (aa|bb)", (nmo_active**4)/4
988 m = (nspins*(nspins + 1))/2
990 IF (dft_control%roks) m = 1
991 ALLOCATE (active_space_env%eri%eri(m))
993 CALL get_mo_set(active_space_env%mos_active(1), nmo=nmo)
994 ALLOCATE (active_space_env%eri%eri(i)%csr_mat)
995 eri_mat => active_space_env%eri%eri(i)%csr_mat
1006 nn1 = (n1*(n1 + 1))/2
1007 nn2 = (n2*(n2 + 1))/2
1008 CALL dbcsr_csr_create(eri_mat, nn1, nn2, 0_int_8, 0, 0, para_env%get_handle())
1009 active_space_env%eri%norb = nmo
1012 SELECT CASE (eri_method)
1015 active_space_env%eri%eri_gpw%eps_grid = eri_eps_grid
1017 active_space_env%eri%eri_gpw%eps_filter = eri_eps_filter
1019 active_space_env%eri%eri_gpw%cutoff = eri_gpw_cutoff
1021 active_space_env%eri%eri_gpw%rel_cutoff = eri_rel_cutoff
1023 active_space_env%eri%eri_gpw%print_level = eri_print
1025 active_space_env%eri%eri_gpw%store_wfn = store_wfn
1027 active_space_env%eri%eri_gpw%group_size = group_size
1029 active_space_env%eri%eri_gpw%redo_poisson = .true.
1032 WRITE (iw,
'(/,T2,A,T71,F10.1)')
"ERI_GPW| Energy cutoff [Ry]", eri_gpw_cutoff
1033 WRITE (iw,
'(T2,A,T71,F10.1)')
"ERI_GPW| Relative energy cutoff [Ry]", eri_rel_cutoff
1036 CALL calculate_eri_gpw(active_space_env%mos_active, active_space_env%active_orbitals, active_space_env%eri, qs_env, iw, &
1040 cpabort(
"Unknown ERI method")
1043 DO isp = 1,
SIZE(active_space_env%eri%eri)
1044 eri_mat => active_space_env%eri%eri(isp)%csr_mat
1045 nze_percentage = 100.0_dp*(real(eri_mat%nze_total, kind=
dp) &
1046 /real(eri_mat%nrows_total, kind=
dp))/real(eri_mat%ncols_total, kind=
dp)
1047 WRITE (iw,
'(/,T2,A,I2,T30,A,T68,I12)')
"ERI_GPW| Spinmatrix:", isp, &
1048 "Number of CSR non-zero elements:", eri_mat%nze_total
1049 WRITE (iw,
'(T2,A,I2,T30,A,T68,F12.4)')
"ERI_GPW| Spinmatrix:", isp, &
1050 "Percentage CSR non-zero elements:", nze_percentage
1051 WRITE (iw,
'(T2,A,I2,T30,A,T68,I12)')
"ERI_GPW| Spinmatrix:", isp, &
1052 "nrows_total", eri_mat%nrows_total
1053 WRITE (iw,
'(T2,A,I2,T30,A,T68,I12)')
"ERI_GPW| Spinmatrix:", isp, &
1054 "ncols_total", eri_mat%ncols_total
1055 WRITE (iw,
'(T2,A,I2,T30,A,T68,I12)')
"ERI_GPW| Spinmatrix:", isp, &
1056 "nrows_local", eri_mat%nrows_local
1060 CALL para_env%sync()
1063 nspins = active_space_env%nspins
1064 ALLOCATE (active_space_env%p_active(nspins))
1066 mo_set => active_space_env%mos_active(isp)
1067 CALL get_mo_set(mo_set, mo_coeff=mo_coeff, nmo=nmo)
1068 CALL create_subspace_matrix(mo_coeff, active_space_env%p_active(isp), nmo)
1070 SELECT CASE (mselect)
1072 cpabort(
"Unknown orbital selection method")
1075 IF (nspins == 2) focc = 1.0_dp
1077 fmat => active_space_env%p_active(isp)
1079 IF (nspins == 2)
THEN
1081 nel = (active_space_env%nelec_active + active_space_env%multiplicity - 1)/2
1083 nel = (active_space_env%nelec_active - active_space_env%multiplicity + 1)/2
1086 nel = active_space_env%nelec_active
1088 DO i = 1, nmo_active
1089 m = active_space_env%active_orbitals(i, isp)
1090 fel = min(focc, real(nel, kind=
dp))
1092 nel = nel - nint(fel)
1097 cpabort(
"NOT IMPLEMENTED")
1099 cpabort(
"NOT IMPLEMENTED")
1103 CALL calculate_spin_pol_overlap(active_space_env%mos_active, qs_env, active_space_env)
1108 IF (
ASSOCIATED(xc_section))
CALL section_vals_get(xc_section, explicit=explicit)
1113 IF (
ASSOCIATED(qs_env%x_data))
CALL hfx_release(qs_env%x_data)
1117 CALL get_qs_env(qs_env, atomic_kind_set=atomic_kind_set, qs_kind_set=qs_kind_set, &
1118 particle_set=particle_set, cell=cell, ks_env=ks_env)
1119 IF (dft_control%do_admm)
THEN
1120 basis_type =
'AUX_FIT'
1125 CALL hfx_create(qs_env%x_data, para_env, hfx_section, atomic_kind_set, &
1126 qs_kind_set, particle_set, dft_control, cell, orb_basis=basis_type, &
1127 nelectron_total=nelec_total)
1129 qs_env%requires_matrix_vxc = .true.
1132 CALL set_ks_env(ks_env, s_mstruct_changed=.true.)
1134 just_energy=.false., &
1135 ext_xc_section=xc_section)
1137 CALL set_ks_env(ks_env, s_mstruct_changed=.false.)
1142 active_space_env%xc_section => xc_section
1146 CALL calculate_operators(active_space_env%mos_active, qs_env, active_space_env)
1148 active_space_env%energy_ref = energy%total
1150 CALL subspace_fock_matrix(active_space_env, dft_control%roks)
1153 CALL set_qs_env(qs_env, active_space=active_space_env)
1157 SELECT CASE (as_solver)
1160 WRITE (iw,
'(/,T3,A)')
"No active space solver specified, skipping embedding calculation"
1163 CALL para_env%sync()
1165 CALL rsdft_embedding(qs_env, active_space_env, as_input)
1168 CALL local_fci_embedding(qs_env, active_space_env, as_input)
1171 cpabort(
"Unknown active space solver")
1175 IF (active_space_env%fcidump)
CALL fcidump(active_space_env, as_input, dft_control%roks)
1178 IF (active_space_env%qcschema)
THEN
1185 WRITE (iw,
'(/,T2,A)') &
1186 '!-------------------- End of Active Space Interface --------------------------!'
1189 CALL para_env%sync()
1191 CALL timestop(handle)
1203 SUBROUTINE calculate_spin_pol_overlap(mos, qs_env, active_space_env)
1205 TYPE(
mo_set_type),
DIMENSION(:),
INTENT(IN) :: mos
1209 CHARACTER(len=*),
PARAMETER :: routinen =
'calculate_spin_pol_overlap'
1211 INTEGER :: handle, nmo, nspins
1212 LOGICAL :: do_kpoints
1213 TYPE(
cp_fm_type),
POINTER :: mo_coeff_a, mo_coeff_b
1215 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: s_matrix_kp
1217 CALL timeset(routinen, handle)
1219 nspins = active_space_env%nspins
1222 IF (nspins > 1)
THEN
1223 CALL get_qs_env(qs_env, do_kpoints=do_kpoints)
1224 IF (do_kpoints)
THEN
1225 CALL get_qs_env(qs_env, matrix_s_kp=s_matrix_kp)
1226 s_matrix => s_matrix_kp(:, 1)
1230 ALLOCATE (active_space_env%sab_sub(1))
1232 CALL get_mo_set(mo_set=mos(1), mo_coeff=mo_coeff_a, nmo=nmo)
1233 CALL get_mo_set(mo_set=mos(2), mo_coeff=mo_coeff_b, nmo=nmo)
1234 CALL subspace_operator(mo_coeff_a, nmo, s_matrix(1)%matrix, active_space_env%sab_sub(1), mo_coeff_b)
1237 CALL timestop(handle)
1239 END SUBROUTINE calculate_spin_pol_overlap
1249 SUBROUTINE calculate_operators(mos, qs_env, active_space_env)
1251 TYPE(
mo_set_type),
DIMENSION(:),
INTENT(IN) :: mos
1255 CHARACTER(len=*),
PARAMETER :: routinen =
'calculate_operators'
1257 INTEGER :: handle, ispin, nmo, nspins
1259 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: h_matrix, ks_matrix
1261 CALL timeset(routinen, handle)
1263 nspins = active_space_env%nspins
1267 CALL get_qs_env(qs_env, matrix_ks_kp=ks_matrix)
1268 ALLOCATE (active_space_env%ks_sub(nspins))
1269 DO ispin = 1, nspins
1270 CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
1271 CALL subspace_operator(mo_coeff, nmo, ks_matrix(ispin, 1)%matrix, active_space_env%ks_sub(ispin))
1278 CALL get_qs_env(qs_env=qs_env, matrix_h_kp=h_matrix)
1279 ALLOCATE (active_space_env%h_sub(nspins))
1280 DO ispin = 1, nspins
1281 CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
1282 CALL subspace_operator(mo_coeff, nmo, h_matrix(1, 1)%matrix, active_space_env%h_sub(ispin))
1285 CALL timestop(handle)
1287 END SUBROUTINE calculate_operators
1299 SUBROUTINE subspace_operator(mo_coeff, nmo, op_matrix, op_sub, mo_coeff_b)
1302 INTEGER,
INTENT(IN) :: nmo
1305 TYPE(
cp_fm_type),
INTENT(IN),
OPTIONAL :: mo_coeff_b
1307 CHARACTER(len=*),
PARAMETER :: routinen =
'subspace_operator'
1309 INTEGER :: handle, ncol, nrow
1312 CALL timeset(routinen, handle)
1314 CALL cp_fm_get_info(matrix=mo_coeff, ncol_global=ncol, nrow_global=nrow)
1315 cpassert(nmo <= ncol)
1318 CALL cp_fm_create(vectors, mo_coeff%matrix_struct,
"vectors")
1319 CALL create_subspace_matrix(mo_coeff, op_sub, nmo)
1321 IF (
PRESENT(mo_coeff_b))
THEN
1329 CALL parallel_gemm(
'T',
'N', nmo, nmo, nrow, 1.0_dp, mo_coeff, vectors, 0.0_dp, op_sub)
1333 CALL timestop(handle)
1335 END SUBROUTINE subspace_operator
1345 SUBROUTINE create_subspace_matrix(orbitals, op_sub, n)
1349 INTEGER,
INTENT(IN) :: n
1357 para_env=orbitals%matrix_struct%para_env, &
1358 context=orbitals%matrix_struct%context)
1359 CALL cp_fm_create(op_sub, fm_struct, name=
"Subspace operator")
1364 END SUBROUTINE create_subspace_matrix
1377 SUBROUTINE calculate_eri_gpw(mos, orbitals, eri_env, qs_env, iw, restricted)
1379 TYPE(
mo_set_type),
DIMENSION(:),
INTENT(IN) :: mos
1380 INTEGER,
DIMENSION(:, :),
POINTER :: orbitals
1383 INTEGER,
INTENT(IN) :: iw
1384 LOGICAL,
INTENT(IN) :: restricted
1386 CHARACTER(len=*),
PARAMETER :: routinen =
'calculate_eri_gpw'
1388 INTEGER :: col_local, color, handle, i1, i2, i3, i4, i_multigrid, icount2, intcount, isp, &
1389 isp1, isp2, ispin, iwa1, iwa12, iwa2, iwb1, iwb12, iwb2, iwbs, iwbt, iwfn, n_multigrid, &
1390 ncol_global, ncol_local, nmm, nmo, nmo1, nmo2, nrow_global, nrow_local, nspins, &
1391 number_of_subgroups, nx, row_local, stored_integrals
1392 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: eri_index
1393 INTEGER,
DIMENSION(:),
POINTER :: col_indices, row_indices
1394 LOGICAL :: print1, print2, &
1395 skip_load_balance_distributed
1396 REAL(kind=
dp) :: dvol, erint, pair_int, &
1397 progression_factor, rc, rsize, t1, t2
1398 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: eri
1403 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:) :: fm_matrix_pq_rnu, fm_matrix_pq_rs, &
1407 TYPE(
dbcsr_type),
ALLOCATABLE,
DIMENSION(:) :: matrix_pq_rnu, mo_coeff_as
1411 POINTER :: sab_orb_sub
1419 DIMENSION(:, :),
TARGET :: wfn_a
1422 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
1426 CALL timeset(routinen, handle)
1431 SELECT CASE (eri_env%eri_gpw%print_level)
1452 CALL get_qs_env(qs_env, para_env=para_env, blacs_env=blacs_env)
1453 IF (eri_env%eri_gpw%group_size < 1) eri_env%eri_gpw%group_size = para_env%num_pe
1454 IF (mod(para_env%num_pe, eri_env%eri_gpw%group_size) /= 0) &
1455 cpabort(
"Group size must be a divisor of the total number of processes!")
1457 IF (eri_env%eri_gpw%group_size == para_env%num_pe)
THEN
1458 eri_env%para_env_sub => para_env
1459 CALL eri_env%para_env_sub%retain()
1460 blacs_env_sub => blacs_env
1461 CALL blacs_env_sub%retain()
1462 number_of_subgroups = 1
1465 number_of_subgroups = para_env%num_pe/eri_env%eri_gpw%group_size
1466 color = para_env%mepos/eri_env%eri_gpw%group_size
1467 ALLOCATE (eri_env%para_env_sub)
1468 CALL eri_env%para_env_sub%from_split(para_env, color)
1469 NULLIFY (blacs_env_sub)
1472 CALL eri_env%comm_exchange%from_split(para_env, eri_env%para_env_sub%mepos)
1475 CALL get_qs_env(qs_env, dft_control=dft_control)
1476 ALLOCATE (qs_control)
1477 qs_control_old => dft_control%qs_control
1478 qs_control = qs_control_old
1479 dft_control%qs_control => qs_control
1480 progression_factor = qs_control%progression_factor
1481 n_multigrid =
SIZE(qs_control%e_cutoff)
1485 IF (restricted) nspins = 1
1487 ALLOCATE (qs_control%e_cutoff(n_multigrid))
1489 qs_control%cutoff = eri_env%eri_gpw%cutoff*0.5_dp
1490 qs_control%e_cutoff(1) = qs_control%cutoff
1491 DO i_multigrid = 2, n_multigrid
1492 qs_control%e_cutoff(i_multigrid) = qs_control%e_cutoff(i_multigrid - 1) &
1495 qs_control%relative_cutoff = eri_env%eri_gpw%rel_cutoff*0.5_dp
1499 CALL create_mat_munu(mat_munu, qs_env, eri_env%eri_gpw%eps_grid, blacs_env_sub, sab_orb_sub=sab_orb_sub, &
1500 do_alloc_blocks_from_nbl=.true., dbcsr_sym_type=dbcsr_type_symmetric)
1504 NULLIFY (pw_env_sub)
1506 CALL pw_env_rebuild(pw_env_sub, qs_env, external_para_env=eri_env%para_env_sub)
1507 CALL pw_env_get(pw_env_sub, auxbas_pw_pool=auxbas_pw_pool, poisson_env=poisson_env)
1510 IF (eri_env%eri_gpw%redo_poisson)
THEN
1512 IF (sum(eri_env%periodicity) /= 0)
THEN
1517 poisson_env%parameters%periodic = eri_env%periodicity
1526 rc = cell%hmat(1, 1)
1529 rc = min(rc, 0.5_dp*cell%hmat(iwa1, iwa1))
1531 poisson_env%green_fft%radius = rc
1534 CALL pw_eri_green_create(poisson_env%green_fft, eri_env)
1538 IF (sum(cell%perd) /= sum(eri_env%periodicity))
THEN
1539 IF (sum(eri_env%periodicity) /= 0)
THEN
1540 WRITE (unit=iw, fmt=
"(/,T2,A,T51,A30)") &
1541 "ERI_GPW| Switching Poisson solver to",
"PERIODIC"
1543 WRITE (unit=iw, fmt=
"(/,T2,A,T51,A30)") &
1544 "ERI_GPW| Switching Poisson solver to",
"ANALYTIC"
1548 SELECT CASE (poisson_env%green_fft%method)
1550 WRITE (unit=iw, fmt=
"(T2,A,T51,A30)") &
1551 "ERI_GPW| Poisson Greens function",
"PERIODIC"
1553 WRITE (unit=iw, fmt=
"(T2,A,T51,A30)") &
1554 "ERI_GPW| Poisson Greens function",
"ANALYTIC"
1555 WRITE (unit=iw, fmt=
"(T2,A,T71,F10.4)")
"ERI_GPW| Poisson cutoff radius", &
1556 poisson_env%green_fft%radius*
angstrom
1558 cpabort(
"Wrong Greens function setup")
1563 ALLOCATE (mo_coeff_as(nspins), fm_mo_coeff_as(nspins))
1564 DO ispin = 1, nspins
1566 REAL(kind=
dp),
DIMENSION(:, :),
ALLOCATABLE :: c, c_active
1570 CALL get_mo_set(mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
1573 ALLOCATE (c_active(
SIZE(c, 1),
SIZE(orbitals, 1)))
1574 DO i1 = 1,
SIZE(orbitals, 1)
1575 c_active(:, i1) = c(:, orbitals(i1, ispin))
1578 c_active, mat_munu%matrix, gd_array, eri_env%eri_gpw%eps_filter)
1580 DEALLOCATE (c, c_active)
1583 CALL dbcsr_get_info(mo_coeff_as(ispin), nfullrows_total=nrow_global, nfullcols_total=ncol_global)
1586 CALL cp_fm_struct_create(fm_struct, context=blacs_env_sub, para_env=eri_env%para_env_sub, &
1587 nrow_global=nrow_global, ncol_global=ncol_global)
1596 NULLIFY (task_list_sub)
1597 skip_load_balance_distributed = dft_control%qs_control%skip_load_balance_distributed
1600 reorder_rs_grid_ranks=.true., &
1601 skip_load_balance_distributed=skip_load_balance_distributed, &
1602 pw_env_external=pw_env_sub, sab_orb_external=sab_orb_sub)
1607 ALLOCATE (matrix_pq_rnu(nspins), fm_matrix_pq_rnu(nspins), fm_matrix_pq_rs(nspins))
1608 DO ispin = 1, nspins
1609 CALL dbcsr_create(matrix_pq_rnu(ispin), template=mo_coeff_as(ispin))
1610 CALL dbcsr_set(matrix_pq_rnu(ispin), 0.0_dp)
1612 CALL dbcsr_get_info(matrix_pq_rnu(ispin), nfullrows_total=nrow_global, nfullcols_total=ncol_global)
1615 CALL cp_fm_struct_create(fm_struct, context=blacs_env_sub, para_env=eri_env%para_env_sub, &
1616 nrow_global=nrow_global, ncol_global=ncol_global)
1621 CALL cp_fm_struct_create(fm_struct, context=blacs_env_sub, para_env=eri_env%para_env_sub, &
1622 nrow_global=ncol_global, ncol_global=ncol_global)
1630 CALL auxbas_pw_pool%create_pw(wfn_r)
1631 CALL auxbas_pw_pool%create_pw(rho_g)
1632 CALL get_qs_env(qs_env, qs_kind_set=qs_kind_set, cell=cell, &
1633 particle_set=particle_set, atomic_kind_set=atomic_kind_set)
1639 IF (restricted) nspins = 1
1640 IF (eri_env%eri_gpw%store_wfn)
THEN
1644 DO ispin = 1, nspins
1647 rsize = real(
SIZE(wfn_r%array), kind=
dp)*nx
1649 IF (print1 .AND. iw > 0)
THEN
1650 rsize = rsize*8._dp/1000000._dp
1651 WRITE (iw,
"(T2,'ERI_GPW|',' Store active orbitals on real space grid ',T66,F12.3,' MB')") rsize
1653 ALLOCATE (wfn_a(nmo, nspins))
1654 DO ispin = 1, nspins
1655 CALL get_mo_set(mo_set=mos(ispin), mo_coeff=mo_coeff, nmo=nmo)
1656 DO i1 = 1,
SIZE(orbitals, 1)
1657 iwfn = orbitals(i1, ispin)
1658 CALL auxbas_pw_pool%create_pw(wfn_a(iwfn, ispin))
1660 qs_kind_set, cell, dft_control, particle_set, pw_env_sub)
1661 IF (print2 .AND. iw > 0)
THEN
1662 WRITE (iw,
"(T2,'ERI_GPW|',' Orbital stored ',I4,' Spin ',i1)") iwfn, ispin
1668 ALLOCATE (wfn1, wfn2)
1669 CALL auxbas_pw_pool%create_pw(wfn1)
1670 CALL auxbas_pw_pool%create_pw(wfn2)
1672 ALLOCATE (wfn3, wfn4)
1673 CALL auxbas_pw_pool%create_pw(wfn3)
1674 CALL auxbas_pw_pool%create_pw(wfn4)
1679 CALL auxbas_pw_pool%create_pw(rho_r)
1680 CALL auxbas_pw_pool%create_pw(pot_g)
1685 dvol = rho_r%pw_grid%dvol
1691 stored_integrals = 0
1694 nmm = (nmo1*(nmo1 + 1))/2
1695 DO i1 = 1,
SIZE(orbitals, 1)
1696 iwa1 = orbitals(i1, isp1)
1697 IF (eri_env%eri_gpw%store_wfn)
THEN
1698 wfn1 => wfn_a(iwa1, isp1)
1701 qs_kind_set, cell, dft_control, particle_set, pw_env_sub)
1703 DO i2 = i1,
SIZE(orbitals, 1)
1704 iwa2 = orbitals(i2, isp1)
1707 IF (mod(iwa12 - 1, eri_env%comm_exchange%num_pe) /= eri_env%comm_exchange%mepos) cycle
1708 iwa12 = (iwa12 - 1)/eri_env%comm_exchange%num_pe + 1
1709 IF (eri_env%eri_gpw%store_wfn)
THEN
1710 wfn2 => wfn_a(iwa2, isp1)
1713 qs_kind_set, cell, dft_control, particle_set, pw_env_sub)
1722 IF (pair_int < eri_env%eps_integral) cycle
1727 DO isp2 = isp1, nspins
1729 nx = (nmo2*(nmo2 + 1))/2
1730 ALLOCATE (eri(nx), eri_index(nx))
1732 CALL integrate_v_rspace(rho_r, hmat=mat_munu, qs_env=qs_env, &
1733 calculate_forces=.false., compute_tau=.false., gapw=.false., &
1734 pw_env_external=pw_env_sub, task_list_external=task_list_sub)
1736 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, mat_munu%matrix, mo_coeff_as(isp2), &
1737 0.0_dp, matrix_pq_rnu(isp2), filter_eps=eri_env%eri_gpw%eps_filter)
1740 CALL cp_fm_get_info(fm_matrix_pq_rnu(isp2), ncol_global=ncol_global, nrow_global=nrow_global)
1742 CALL parallel_gemm(
"T",
"N", ncol_global, ncol_global, nrow_global, 0.5_dp, &
1743 fm_matrix_pq_rnu(isp2), fm_mo_coeff_as(isp2), &
1744 0.0_dp, fm_matrix_pq_rs(isp2))
1745 CALL parallel_gemm(
"T",
"N", ncol_global, ncol_global, nrow_global, 0.5_dp, &
1746 fm_mo_coeff_as(isp2), fm_matrix_pq_rnu(isp2), &
1747 1.0_dp, fm_matrix_pq_rs(isp2))
1749 CALL cp_fm_get_info(fm_matrix_pq_rs(isp2), ncol_local=ncol_local, nrow_local=nrow_local, &
1750 col_indices=col_indices, row_indices=row_indices)
1753 DO col_local = 1, ncol_local
1754 iwb1 = orbitals(col_indices(col_local), isp2)
1755 IF (isp1 == isp2 .AND. iwb1 < iwa1) cycle
1756 DO row_local = 1, nrow_local
1757 iwb2 = orbitals(row_indices(row_local), isp2)
1758 IF (iwb2 < iwb1) cycle
1759 IF (isp1 == isp2 .AND. iwa1 == iwb1 .AND. iwb2 < iwa2) cycle
1762 erint = fm_matrix_pq_rs(isp2)%local_data(row_local, col_local)
1763 IF (abs(erint) > eri_env%eps_integral)
THEN
1764 icount2 = icount2 + 1
1765 eri(icount2) = erint
1766 eri_index(icount2) = iwb12
1770 stored_integrals = stored_integrals + icount2
1772 isp = (isp1 - 1)*isp2 + (isp2 - isp1 + 1)
1773 CALL update_csr_matrix(eri_env%eri(isp)%csr_mat, icount2, eri, eri_index, iwa12)
1775 DEALLOCATE (eri, eri_index)
1778 DO isp2 = isp1, nspins
1780 nx = (nmo2*(nmo2 + 1))/2
1781 ALLOCATE (eri(nx), eri_index(nx))
1784 IF (isp1 == isp2) iwbs = i1
1785 isp = (isp1 - 1)*isp2 + (isp2 - isp1 + 1)
1786 DO i3 = iwbs,
SIZE(orbitals, 1)
1787 iwb1 = orbitals(i3, isp2)
1788 IF (eri_env%eri_gpw%store_wfn)
THEN
1789 wfn3 => wfn_a(iwb1, isp2)
1792 qs_kind_set, cell, dft_control, particle_set, pw_env_sub)
1797 IF (isp1 == isp2 .AND. i1 == i3) iwbt = i2
1798 DO i4 = iwbt,
SIZE(orbitals, 1)
1799 iwb2 = orbitals(i4, isp2)
1800 IF (eri_env%eri_gpw%store_wfn)
THEN
1801 wfn4 => wfn_a(iwb2, isp2)
1804 qs_kind_set, cell, dft_control, particle_set, pw_env_sub)
1808 icount2 = icount2 + 1
1809 eri(icount2) = erint
1814 CALL eri_env%para_env_sub%sum(eri)
1819 IF (isp1 == isp2) iwbs = i1
1820 isp = (isp1 - 1)*isp2 + (isp2 - isp1 + 1)
1821 DO i3 = iwbs,
SIZE(orbitals, 1)
1822 iwb1 = orbitals(i3, isp2)
1824 IF (isp1 == isp2 .AND. i1 == i3) iwbt = i2
1825 DO i4 = iwbt,
SIZE(orbitals, 1)
1826 iwb2 = orbitals(i4, isp2)
1827 intcount = intcount + 1
1828 erint = eri(intcount)
1829 IF (abs(erint) > eri_env%eps_integral)
THEN
1830 IF (mod(intcount, eri_env%para_env_sub%num_pe) == eri_env%para_env_sub%mepos)
THEN
1831 icount2 = icount2 + 1
1832 eri(icount2) = erint
1833 eri_index(icount2) = eri_index(intcount)
1838 stored_integrals = stored_integrals + icount2
1840 CALL update_csr_matrix(eri_env%eri(isp)%csr_mat, icount2, eri, eri_index, iwa12)
1842 DEALLOCATE (eri, eri_index)
1845 cpabort(
"Unknown option")
1851 IF (print1 .AND. iw > 0)
THEN
1852 WRITE (iw,
"(T2,'ERI_GPW|',' Number of Integrals stored locally',T71,I10)") stored_integrals
1855 IF (eri_env%eri_gpw%store_wfn)
THEN
1856 DO ispin = 1, nspins
1857 DO i1 = 1,
SIZE(orbitals, 1)
1858 iwfn = orbitals(i1, ispin)
1859 CALL wfn_a(iwfn, ispin)%release()
1866 DEALLOCATE (wfn1, wfn2)
1870 DEALLOCATE (wfn3, wfn4)
1873 CALL auxbas_pw_pool%give_back_pw(wfn_r)
1874 CALL auxbas_pw_pool%give_back_pw(rho_g)
1875 CALL auxbas_pw_pool%give_back_pw(rho_r)
1876 CALL auxbas_pw_pool%give_back_pw(pot_g)
1879 DO ispin = 1, nspins
1885 DEALLOCATE (matrix_pq_rnu, fm_matrix_pq_rnu, fm_matrix_pq_rs)
1888 DO ispin = 1, nspins
1892 DEALLOCATE (mo_coeff_as, fm_mo_coeff_as)
1896 DEALLOCATE (mat_munu%matrix)
1899 dft_control%qs_control => qs_control_old
1900 DEALLOCATE (qs_control%e_cutoff)
1901 DEALLOCATE (qs_control)
1906 WRITE (iw,
'(/,T2,A,T66,F14.2)')
"ERI_GPW| ERI calculation took (sec)", t2 - t1
1910 CALL timestop(handle)
1912 END SUBROUTINE calculate_eri_gpw
1922 SUBROUTINE pw_eri_green_create(green, eri_env)
1924 TYPE(greens_fn_type),
INTENT(INOUT) :: green
1925 TYPE(eri_type) :: eri_env
1927 COMPLEX(KIND=dp) :: erf_fac_p, z_p
1929 REAL(kind=dp) :: cossin_fac, ea, erfcos_fac, exp_prefac, &
1930 g, g0, g2, g3d, ga, ginf, omega, &
1934 associate(gf => green%influence_fn, grid => green%influence_fn%pw_grid)
1935 SELECT CASE (green%method)
1938 SELECT CASE (eri_env%operator)
1939 CASE (eri_operator_coulomb)
1940 DO ig = grid%first_gne0, grid%ngpts_cut_local
1942 gf%array(ig) = fourpi/g2
1944 IF (grid%have_g0) gf%array(1) = 0.0_dp
1946 CASE (eri_operator_yukawa)
1947 CALL cp_warn(__location__,
"Yukawa operator has not been tested")
1948 omega2 = eri_env%omega**2
1949 DO ig = grid%first_gne0, grid%ngpts_cut_local
1951 gf%array(ig) = fourpi/(omega2 + g2)
1953 IF (grid%have_g0) gf%array(1) = fourpi/omega2
1955 CASE (eri_operator_erf)
1956 omega2 = eri_env%omega**2
1957 DO ig = grid%first_gne0, grid%ngpts_cut_local
1959 gf%array(ig) = fourpi/g2*exp(-0.25_dp*g2/omega2)
1961 IF (grid%have_g0) gf%array(1) = 0.0_dp
1963 CASE (eri_operator_erfc)
1964 omega2 = eri_env%omega**2
1965 DO ig = grid%first_gne0, grid%ngpts_cut_local
1967 gf%array(ig) = fourpi/g2*(1.0_dp - exp(-0.25_dp*g2/omega2))
1969 IF (grid%have_g0) gf%array(1) = pi/omega2
1971 CASE (eri_operator_trunc)
1972 rc = eri_env%cutoff_radius
1973 DO ig = grid%first_gne0, grid%ngpts_cut_local
1977 IF (g*rc >= 0.005_dp)
THEN
1978 gf%array(ig) = fourpi/g2*(1.0_dp - cos(g*rc))
1980 gf%array(ig) = fourpi/g2*(g*rc)**2/2.0_dp*(1.0_dp - (g*rc)**2/12.0_dp)
1983 IF (grid%have_g0) gf%array(1) = twopi*rc**2
1985 CASE (eri_operator_lr_trunc)
1986 omega = eri_env%omega
1988 rc = eri_env%cutoff_radius
1992 DO ig = grid%first_gne0, grid%ngpts_cut_local
1995 IF (g <= 2.0_dp*g0)
THEN
1996 gf%array(ig) = -pi/omega2*erf(omega*rc) &
1997 + twopi*rc2*erf(omega*rc) &
1998 + 2*rootpi*rc*exp(-omega2*rc2)/omega
1999 ELSE IF (g >= 2.0_dp*ginf*omega)
THEN
2001 exp_prefac = exp(-omega2*rc2)/(rootpi*(omega2*rc2 + 0.25_dp*g2/omega2))
2003 cossin_fac = omega*rc*cos(g*rc) - 0.5_dp*g/omega*sin(g*rc)
2005 erfcos_fac = erf(omega*rc)*cos(g*rc)
2007 gf%array(ig) = fourpi/g2*(-exp_prefac*cossin_fac - erfcos_fac)
2010 exp_prefac = twopi/g2*exp(-0.25_dp*g2/omega2)
2012 z_p = cmplx(omega*rc, 0.5_dp*g/omega, kind=dp)
2014 erf_fac_p = 2.0_dp*real(erfz_fast(z_p))
2016 erfcos_fac = fourpi/g2*erf(omega*rc)*cos(g*rc)
2018 gf%array(ig) = exp_prefac*erf_fac_p - erfcos_fac
2021 IF (grid%have_g0)
THEN
2022 gf%array(1) = -pi/omega2*erf(omega*rc) &
2023 + twopi*rc2*erf(omega*rc) &
2024 + 2*rootpi*rc*exp(-omega2*rc2)/omega
2028 cpabort(
"Please specify a valid operator for the periodic Poisson solver")
2036 SELECT CASE (eri_env%operator)
2040 CASE (eri_operator_coulomb, eri_operator_trunc)
2041 IF (eri_env%operator == eri_operator_coulomb)
THEN
2044 rc = eri_env%cutoff_radius
2046 DO ig = grid%first_gne0, grid%ngpts_cut_local
2050 IF (g*rc >= 0.005_dp)
THEN
2051 gf%array(ig) = fourpi/g2*(1.0_dp - cos(g*rc))
2053 gf%array(ig) = fourpi/g2*(g*rc)**2/2.0_dp*(1.0_dp - (g*rc)**2/12.0_dp)
2056 IF (grid%have_g0) gf%array(1) = twopi*rc**2
2059 CASE (eri_operator_yukawa)
2060 CALL cp_warn(__location__,
"Yukawa operator has not been tested")
2062 omega = eri_env%omega
2064 DO ig = grid%first_gne0, grid%ngpts_cut_local
2067 g3d = fourpi/(omega**2 + g2)
2068 gf%array(ig) = g3d*(1.0_dp - ea*(cos(g*rc) + omega/g*sin(g*rc)))
2070 IF (grid%have_g0) gf%array(1) = fourpi/(omega**2)*(1.0_dp - ea*(1.0_dp + omega*rc))
2074 CASE (eri_operator_erf, eri_operator_lr_trunc)
2075 IF (eri_env%operator == eri_operator_erf)
THEN
2078 rc = eri_env%cutoff_radius
2080 omega2 = eri_env%omega**2
2081 DO ig = grid%first_gne0, grid%ngpts_cut_local
2084 ga = -0.25_dp*g2/omega2
2085 gf%array(ig) = fourpi/g2*exp(ga)*(1.0_dp - cos(g*rc))
2087 IF (grid%have_g0) gf%array(1) = twopi*rc**2
2091 CASE (eri_operator_erfc)
2092 CALL cp_warn(__location__, &
2093 "Short-range Coulomb operator may be incorrect with ANALYTIC0D Poisson solver")
2095 omega2 = eri_env%omega**2
2096 DO ig = grid%first_gne0, grid%ngpts_cut_local
2099 ga = -0.25_dp*g2/omega2
2100 gf%array(ig) = fourpi/g2*(1.0_dp - exp(ga))*(1.0_dp - cos(g*rc))
2102 IF (grid%have_g0) gf%array(1) = pi/omega2
2105 cpabort(
"Unsupported operator")
2109 cpabort(
"Unsupported Poisson solver")
2113 END SUBROUTINE pw_eri_green_create
2125 SUBROUTINE update_csr_matrix(csr_mat, nnz, rdat, rind, irow)
2127 TYPE(dbcsr_csr_type),
INTENT(INOUT) :: csr_mat
2128 INTEGER,
INTENT(IN) :: nnz
2129 REAL(kind=dp),
DIMENSION(:),
INTENT(IN) :: rdat
2130 INTEGER,
DIMENSION(:),
INTENT(IN) :: rind
2131 INTEGER,
INTENT(IN) :: irow
2133 INTEGER :: k, nrow, nze, nze_new
2136 nze = csr_mat%nze_local
2139 CALL reallocate(csr_mat%nzval_local%r_dp, 1, nze_new)
2140 csr_mat%nzval_local%r_dp(nze + 1:nze_new) = rdat(1:nnz)
2142 CALL reallocate(csr_mat%colind_local, 1, nze_new)
2143 csr_mat%colind_local(nze + 1:nze_new) = rind(1:nnz)
2145 nrow = csr_mat%nrows_local
2146 CALL reallocate(csr_mat%rowptr_local, 1, irow + 1)
2147 csr_mat%rowptr_local(nrow + 1:irow) = nze + 1
2148 csr_mat%rowptr_local(irow + 1) = nze_new + 1
2150 CALL reallocate(csr_mat%nzerow_local, 1, irow)
2151 DO k = nrow + 1, irow
2152 csr_mat%nzerow_local(k) = csr_mat%rowptr_local(k + 1) - csr_mat%rowptr_local(k)
2154 csr_mat%nrows_local = irow
2155 csr_mat%nze_local = csr_mat%nze_local + nnz
2157 csr_mat%nze_total = csr_mat%nze_total + nnz
2158 csr_mat%has_indices = .true.
2160 END SUBROUTINE update_csr_matrix
2168 SUBROUTINE print_orbital_cubes(input, qs_env, mos)
2169 TYPE(section_vals_type),
POINTER :: input
2170 TYPE(qs_environment_type),
POINTER :: qs_env
2171 TYPE(mo_set_type),
DIMENSION(:),
INTENT(IN) :: mos
2173 CHARACTER(LEN=default_path_length) :: filebody, filename, title
2174 INTEGER :: i, imo, isp, nmo, str(3), unit_nr
2175 INTEGER,
DIMENSION(:),
POINTER :: alist, blist, istride
2176 LOGICAL :: do_mo, explicit_a, explicit_b
2177 TYPE(atomic_kind_type),
DIMENSION(:),
POINTER :: atomic_kind_set
2178 TYPE(cell_type),
POINTER :: cell
2179 TYPE(cp_fm_type),
POINTER :: mo_coeff
2180 TYPE(dft_control_type),
POINTER :: dft_control
2181 TYPE(mp_para_env_type),
POINTER :: para_env
2182 TYPE(particle_list_type),
POINTER :: particles
2183 TYPE(particle_type),
DIMENSION(:),
POINTER :: particle_set
2184 TYPE(pw_c1d_gs_type) :: wf_g
2185 TYPE(pw_env_type),
POINTER :: pw_env
2186 TYPE(pw_pool_type),
POINTER :: auxbas_pw_pool
2187 TYPE(pw_r3d_rs_type) :: wf_r
2188 TYPE(qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
2189 TYPE(qs_subsys_type),
POINTER :: subsys
2190 TYPE(section_vals_type),
POINTER :: dft_section, scf_input
2192 CALL section_vals_val_get(input,
"FILENAME", c_val=filebody)
2193 CALL section_vals_val_get(input,
"STRIDE", i_vals=istride)
2194 IF (
SIZE(istride) == 1)
THEN
2195 str(1:3) = istride(1)
2196 ELSEIF (
SIZE(istride) == 3)
THEN
2197 str(1:3) = istride(1:3)
2199 cpabort(
"STRIDE arguments inconsistent")
2201 CALL section_vals_val_get(input,
"ALIST", i_vals=alist, explicit=explicit_a)
2202 CALL section_vals_val_get(input,
"BLIST", i_vals=blist, explicit=explicit_b)
2204 CALL get_qs_env(qs_env=qs_env, &
2205 dft_control=dft_control, &
2206 para_env=para_env, &
2208 atomic_kind_set=atomic_kind_set, &
2209 qs_kind_set=qs_kind_set, &
2211 particle_set=particle_set, &
2215 CALL qs_subsys_get(subsys, particles=particles)
2217 CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
2218 CALL auxbas_pw_pool%create_pw(wf_r)
2219 CALL auxbas_pw_pool%create_pw(wf_g)
2221 dft_section => section_vals_get_subs_vals(scf_input,
"DFT")
2223 DO isp = 1,
SIZE(mos)
2224 CALL get_mo_set(mo_set=mos(isp), mo_coeff=mo_coeff, nmo=nmo)
2226 IF (
SIZE(mos) > 1)
THEN
2229 CALL write_mo_set_to_output_unit(mos(isp), qs_kind_set, particle_set, &
2230 dft_section, 4, 0, final_mos=.true., spin=
"ALPHA")
2232 CALL write_mo_set_to_output_unit(mos(isp), qs_kind_set, particle_set, &
2233 dft_section, 4, 0, final_mos=.true., spin=
"BETA")
2235 cpabort(
"Invalid spin")
2238 CALL write_mo_set_to_output_unit(mos(isp), qs_kind_set, particle_set, &
2239 dft_section, 4, 0, final_mos=.true.)
2243 IF (isp == 1 .AND. explicit_a)
THEN
2244 IF (alist(1) == -1)
THEN
2248 DO i = 1,
SIZE(alist)
2249 IF (imo == alist(i)) do_mo = .true.
2252 ELSE IF (isp == 2 .AND. explicit_b)
THEN
2253 IF (blist(1) == -1)
THEN
2257 DO i = 1,
SIZE(blist)
2258 IF (imo == blist(i)) do_mo = .true.
2264 IF (.NOT. do_mo) cycle
2265 CALL calculate_wavefunction(mo_coeff, imo, wf_r, wf_g, atomic_kind_set, &
2266 qs_kind_set, cell, dft_control, particle_set, pw_env)
2267 IF (para_env%is_source())
THEN
2268 WRITE (filename,
'(A,A1,I4.4,A1,I1.1,A)') trim(filebody),
"_", imo,
"_", isp,
".cube"
2269 CALL open_file(filename, unit_number=unit_nr, file_status=
"UNKNOWN", file_action=
"WRITE")
2270 WRITE (title, *)
"Active Orbital ", imo,
" spin ", isp
2274 CALL cp_pw_to_cube(wf_r, unit_nr, title, particles=particles, stride=istride)
2275 IF (para_env%is_source())
THEN
2276 CALL close_file(unit_nr)
2281 CALL auxbas_pw_pool%give_back_pw(wf_r)
2282 CALL auxbas_pw_pool%give_back_pw(wf_g)
2284 END SUBROUTINE print_orbital_cubes
2294 SUBROUTINE fcidump(active_space_env, as_input, restricted)
2296 TYPE(active_space_type),
POINTER :: active_space_env
2297 TYPE(section_vals_type),
POINTER :: as_input
2298 LOGICAL,
INTENT(IN) :: restricted
2300 INTEGER :: i, i1, i2, i3, i4, isym, iw, m1, m2, &
2301 ms2, nmo, norb, nspins
2302 REAL(kind=dp) :: checksum, esub
2303 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:, :) :: fmat
2304 TYPE(cp_logger_type),
POINTER :: logger
2305 TYPE(eri_fcidump_checksum) :: eri_checksum
2309 logger => cp_get_default_logger()
2310 iw = cp_print_key_unit_nr(logger, as_input,
"FCIDUMP", &
2311 extension=
".fcidump", file_status=
"REPLACE", file_action=
"WRITE", file_form=
"FORMATTED")
2313 nspins = active_space_env%nspins
2314 norb =
SIZE(active_space_env%active_orbitals, 1)
2315 ms2 = active_space_env%multiplicity - 1
2316 IF (nspins == 1 .OR. restricted)
THEN
2318 associate(nelec => active_space_env%nelec_active)
2321 WRITE (iw,
"(A,A,I4,A,I4,A,I2,A)")
"&FCI",
" NORB=", norb,
",NELEC=", nelec,
",MS2=", ms2,
","
2323 WRITE (iw,
"(A,1000(I1,','))")
" ORBSYM=", (isym, i=1, norb)
2325 WRITE (iw,
"(A,I1,A)")
" ISYM=", isym,
","
2326 IF (restricted)
WRITE (iw,
"(A,I1,A)")
" UHF=", 0,
","
2327 WRITE (iw,
"(A)")
" /"
2331 CALL active_space_env%eri%eri_foreach(1, active_space_env%active_orbitals, &
2332 eri_fcidump_print(iw, 1, 1), 1, 1)
2333 CALL eri_checksum%set(1, 1)
2334 CALL active_space_env%eri%eri_foreach(1, active_space_env%active_orbitals, eri_checksum, 1, 1)
2338 nmo = active_space_env%eri%norb
2339 ALLOCATE (fmat(nmo, nmo))
2340 CALL replicate_and_symmetrize_matrix(nmo, active_space_env%fock_sub(1), fmat)
2343 DO m1 = 1,
SIZE(active_space_env%active_orbitals, 1)
2344 i1 = active_space_env%active_orbitals(m1, 1)
2345 DO m2 = m1,
SIZE(active_space_env%active_orbitals, 1)
2346 i2 = active_space_env%active_orbitals(m2, 1)
2347 checksum = checksum + abs(fmat(i1, i2))
2348 WRITE (iw,
"(ES23.16,4I4)") fmat(i1, i2), m1, m2, i3, i4
2354 esub = active_space_env%energy_inactive
2355 i1 = 0; i2 = 0; i3 = 0; i4 = 0
2356 checksum = checksum + abs(esub)
2357 IF (iw > 0)
WRITE (iw,
"(ES23.16,4I4)") esub, i1, i2, i3, i4
2361 associate(nelec => active_space_env%nelec_active)
2364 WRITE (iw,
"(A,A,I4,A,I4,A,I2,A)")
"&FCI",
" NORB=", norb,
",NELEC=", nelec,
",MS2=", ms2,
","
2366 WRITE (iw,
"(A,1000(I1,','))")
" ORBSYM=", (isym, i=1, norb)
2368 WRITE (iw,
"(A,I1,A)")
" ISYM=", isym,
","
2369 WRITE (iw,
"(A,I1,A)")
" UHF=", 1,
","
2370 WRITE (iw,
"(A)")
" /"
2375 CALL active_space_env%eri%eri_foreach(1, active_space_env%active_orbitals, &
2376 eri_fcidump_print(iw, 1, 1), 1, 1)
2377 CALL eri_checksum%set(1, 1)
2378 CALL active_space_env%eri%eri_foreach(1, active_space_env%active_orbitals, eri_checksum, 1, 1)
2380 CALL active_space_env%eri%eri_foreach(2, active_space_env%active_orbitals, &
2381 eri_fcidump_print(iw, 1, norb + 1), 1, 2)
2382 CALL eri_checksum%set(1, norb + 1)
2383 CALL active_space_env%eri%eri_foreach(2, active_space_env%active_orbitals, eri_checksum, 1, 2)
2385 CALL active_space_env%eri%eri_foreach(3, active_space_env%active_orbitals, &
2386 eri_fcidump_print(iw, norb + 1, norb + 1), 2, 2)
2387 CALL eri_checksum%set(norb + 1, norb + 1)
2388 CALL active_space_env%eri%eri_foreach(3, active_space_env%active_orbitals, eri_checksum, 2, 2)
2391 nmo = active_space_env%eri%norb
2392 ALLOCATE (fmat(nmo, nmo))
2393 CALL replicate_and_symmetrize_matrix(nmo, active_space_env%fock_sub(1), fmat)
2397 i1 = active_space_env%active_orbitals(m1, 1)
2399 i2 = active_space_env%active_orbitals(m2, 1)
2400 checksum = checksum + abs(fmat(i1, i2))
2401 WRITE (iw,
"(ES23.16,4I4)") fmat(i1, i2), m1, m2, i3, i4
2407 ALLOCATE (fmat(nmo, nmo))
2408 CALL replicate_and_symmetrize_matrix(nmo, active_space_env%fock_sub(2), fmat)
2411 DO m1 = 1,
SIZE(active_space_env%active_orbitals, 1)
2412 i1 = active_space_env%active_orbitals(m1, 2)
2413 DO m2 = m1,
SIZE(active_space_env%active_orbitals, 1)
2414 i2 = active_space_env%active_orbitals(m2, 2)
2415 checksum = checksum + abs(fmat(i1, i2))
2416 WRITE (iw,
"(ES23.16,4I4)") fmat(i1, i2), m1 + norb, m2 + norb, i3, i4
2422 esub = active_space_env%energy_inactive
2423 i1 = 0; i2 = 0; i3 = 0; i4 = 0
2424 checksum = checksum + abs(esub)
2425 IF (iw > 0)
WRITE (iw,
"(ES23.16,4I4)") esub, i1, i2, i3, i4
2429 CALL cp_print_key_finished_output(iw, logger, as_input,
"FCIDUMP")
2432 iw = cp_logger_get_default_io_unit(logger)
2433 IF (iw > 0)
WRITE (iw,
'(T4,A,T66,F12.8)')
"FCIDUMP| Checksum:", eri_checksum%checksum + checksum
2436 END SUBROUTINE fcidump
2444 SUBROUTINE replicate_and_symmetrize_matrix(norb, distributed_matrix, replicated_matrix)
2445 INTEGER,
INTENT(IN) :: norb
2446 TYPE(cp_fm_type),
INTENT(IN) :: distributed_matrix
2447 REAL(dp),
DIMENSION(:, :),
INTENT(INOUT) :: replicated_matrix
2452 replicated_matrix(:, :) = 0.0_dp
2455 CALL cp_fm_get_element(distributed_matrix, i1, i2, mval)
2456 replicated_matrix(i1, i2) = mval
2457 replicated_matrix(i2, i1) = mval
2460 END SUBROUTINE replicate_and_symmetrize_matrix
2469 SUBROUTINE subspace_fock_matrix(active_space_env, restricted)
2471 TYPE(active_space_type),
POINTER :: active_space_env
2472 LOGICAL,
INTENT(IN) :: restricted
2474 INTEGER :: i1, i2, is, norb, nspins
2475 REAL(kind=dp) :: eeri, eref, esub, mval
2476 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:, :) :: ks_a_mat, ks_a_ref, ks_b_mat, ks_b_ref, &
2477 ks_mat, ks_ref, p_a_mat, p_b_mat, p_mat
2478 TYPE(cp_fm_type),
POINTER :: matrix, mo_coef
2479 TYPE(dbcsr_csr_type),
POINTER :: eri, eri_aa, eri_ab, eri_bb
2481 eref = active_space_env%energy_ref
2482 nspins = active_space_env%nspins
2484 IF (nspins == 1)
THEN
2485 CALL get_mo_set(active_space_env%mos_active(1), nmo=norb, mo_coeff=mo_coef)
2490 ALLOCATE (ks_mat(norb, norb), ks_ref(norb, norb), p_mat(norb, norb))
2494 CALL replicate_and_symmetrize_matrix(norb, active_space_env%ks_sub(1), ks_mat)
2495 CALL replicate_and_symmetrize_matrix(norb, active_space_env%p_active(1), p_mat)
2498 eri => active_space_env%eri%eri(1)%csr_mat
2499 CALL build_subspace_fock_matrix(active_space_env%active_orbitals, eri, p_mat, ks_ref, &
2500 active_space_env%eri%comm_exchange)
2504 eeri = 0.5_dp*sum(ks_ref*p_mat)
2510 esub = eref - sum(ks_mat(1:norb, 1:norb)*p_mat(1:norb, 1:norb)) + eeri
2513 ks_mat(1:norb, 1:norb) = ks_mat(1:norb, 1:norb) - ks_ref(1:norb, 1:norb)
2516 active_space_env%energy_inactive = esub
2518 CALL cp_fm_release(active_space_env%fock_sub)
2519 ALLOCATE (active_space_env%fock_sub(nspins))
2521 matrix => active_space_env%ks_sub(is)
2522 CALL cp_fm_create(active_space_env%fock_sub(is), matrix%matrix_struct, &
2523 name=
"Active Fock operator")
2525 matrix => active_space_env%fock_sub(1)
2528 mval = ks_mat(i1, i2)
2529 CALL cp_fm_set_element(matrix, i1, i2, mval)
2534 CALL get_mo_set(active_space_env%mos_active(1), nmo=norb)
2539 ALLOCATE (ks_a_mat(norb, norb), ks_b_mat(norb, norb), &
2540 & ks_a_ref(norb, norb), ks_b_ref(norb, norb), &
2541 & p_a_mat(norb, norb), p_b_mat(norb, norb))
2542 ks_a_ref(:, :) = 0.0_dp; ks_b_ref(:, :) = 0.0_dp
2544 CALL replicate_and_symmetrize_matrix(norb, active_space_env%p_active(1), p_a_mat)
2545 CALL replicate_and_symmetrize_matrix(norb, active_space_env%p_active(2), p_b_mat)
2546 CALL replicate_and_symmetrize_matrix(norb, active_space_env%ks_sub(1), ks_a_mat)
2547 CALL replicate_and_symmetrize_matrix(norb, active_space_env%ks_sub(2), ks_b_mat)
2550 IF (restricted)
THEN
2552 eri_aa => active_space_env%eri%eri(1)%csr_mat
2553 CALL build_subspace_spin_fock_matrix(active_space_env%active_orbitals, eri_aa, eri_aa, p_a_mat, p_b_mat, ks_a_ref, &
2554 tr_mixed_eri=.false., comm_exchange=active_space_env%eri%comm_exchange)
2555 CALL build_subspace_spin_fock_matrix(active_space_env%active_orbitals, eri_aa, eri_aa, p_b_mat, p_a_mat, ks_b_ref, &
2556 tr_mixed_eri=.true., comm_exchange=active_space_env%eri%comm_exchange)
2558 eri_aa => active_space_env%eri%eri(1)%csr_mat
2559 eri_ab => active_space_env%eri%eri(2)%csr_mat
2560 eri_bb => active_space_env%eri%eri(3)%csr_mat
2561 CALL build_subspace_spin_fock_matrix(active_space_env%active_orbitals, eri_aa, eri_ab, p_a_mat, p_b_mat, ks_a_ref, &
2562 tr_mixed_eri=.false., comm_exchange=active_space_env%eri%comm_exchange)
2563 CALL build_subspace_spin_fock_matrix(active_space_env%active_orbitals, eri_bb, eri_ab, p_b_mat, p_a_mat, ks_b_ref, &
2564 tr_mixed_eri=.true., comm_exchange=active_space_env%eri%comm_exchange)
2569 eeri = 0.5_dp*(sum(ks_a_ref*p_a_mat) + sum(ks_b_ref*p_b_mat))
2570 esub = eref - sum(ks_a_mat*p_a_mat) - sum(ks_b_mat*p_b_mat) + eeri
2571 ks_a_mat(:, :) = ks_a_mat(:, :) - ks_a_ref(:, :)
2572 ks_b_mat(:, :) = ks_b_mat(:, :) - ks_b_ref(:, :)
2574 active_space_env%energy_inactive = esub
2576 CALL cp_fm_release(active_space_env%fock_sub)
2577 ALLOCATE (active_space_env%fock_sub(nspins))
2579 matrix => active_space_env%ks_sub(is)
2580 CALL cp_fm_create(active_space_env%fock_sub(is), matrix%matrix_struct, &
2581 name=
"Active Fock operator")
2584 matrix => active_space_env%fock_sub(1)
2587 mval = ks_a_mat(i1, i2)
2588 CALL cp_fm_set_element(matrix, i1, i2, mval)
2591 matrix => active_space_env%fock_sub(2)
2594 mval = ks_b_mat(i1, i2)
2595 CALL cp_fm_set_element(matrix, i1, i2, mval)
2601 END SUBROUTINE subspace_fock_matrix
2611 SUBROUTINE build_subspace_fock_matrix(active_orbitals, eri, p_mat, ks_ref, comm_exchange)
2612 INTEGER,
DIMENSION(:, :),
INTENT(IN) :: active_orbitals
2613 TYPE(dbcsr_csr_type),
INTENT(IN) :: eri
2614 REAL(dp),
DIMENSION(:, :),
INTENT(IN) :: p_mat
2615 REAL(dp),
DIMENSION(:, :),
INTENT(INOUT) :: ks_ref
2616 TYPE(mp_comm_type),
INTENT(IN) :: comm_exchange
2618 CHARACTER(LEN=*),
PARAMETER :: routinen =
'build_subspace_fock_matrix'
2620 INTEGER :: handle, i1, i12, i12l, i2, i3, i34, &
2621 i34l, i4, irptr, m1, m2, nindex, &
2624 TYPE(mp_comm_type) :: mp_group
2626 CALL timeset(routinen, handle)
2629 norb =
SIZE(active_orbitals, 1)
2630 nmo_total =
SIZE(p_mat, 1)
2631 nindex = (nmo_total*(nmo_total + 1))/2
2632 CALL mp_group%set_handle(eri%mp_group%get_handle())
2634 i1 = active_orbitals(m1, 1)
2636 i2 = active_orbitals(m2, 1)
2637 i12 = csr_idx_to_combined(i1, i2, nmo_total)
2638 IF (mod(i12 - 1, comm_exchange%num_pe) == comm_exchange%mepos)
THEN
2639 i12l = (i12 - 1)/comm_exchange%num_pe + 1
2640 irptr = eri%rowptr_local(i12l) - 1
2641 DO i34l = 1, eri%nzerow_local(i12l)
2642 i34 = eri%colind_local(irptr + i34l)
2643 CALL csr_idx_from_combined(i34, nmo_total, i3, i4)
2644 erint = eri%nzval_local%r_dp(irptr + i34l)
2646 ks_ref(i1, i2) = ks_ref(i1, i2) + erint*p_mat(i3, i4)
2648 ks_ref(i1, i2) = ks_ref(i1, i2) + erint*p_mat(i3, i4)
2650 IF (i12 /= i34)
THEN
2651 ks_ref(i3, i4) = ks_ref(i3, i4) + erint*p_mat(i1, i2)
2653 ks_ref(i3, i4) = ks_ref(i3, i4) + erint*p_mat(i1, i2)
2657 erint = -0.5_dp*erint
2658 ks_ref(i1, i3) = ks_ref(i1, i3) + erint*p_mat(i2, i4)
2660 ks_ref(i2, i3) = ks_ref(i2, i3) + erint*p_mat(i1, i4)
2663 ks_ref(i1, i4) = ks_ref(i1, i4) + erint*p_mat(i2, i3)
2665 IF (i1 /= i2 .AND. i3 /= i4)
THEN
2666 ks_ref(i2, i4) = ks_ref(i2, i4) + erint*p_mat(i1, i3)
2674 i1 = active_orbitals(m1, 1)
2676 i2 = active_orbitals(m2, 1)
2677 ks_ref(i2, i1) = ks_ref(i1, i2)
2680 CALL mp_group%sum(ks_ref)
2682 CALL timestop(handle)
2684 END SUBROUTINE build_subspace_fock_matrix
2697 SUBROUTINE build_subspace_spin_fock_matrix(active_orbitals, eri_aa, eri_ab, p_a_mat, p_b_mat, ks_a_ref, tr_mixed_eri, &
2699 INTEGER,
DIMENSION(:, :),
INTENT(IN) :: active_orbitals
2700 TYPE(dbcsr_csr_type),
INTENT(IN) :: eri_aa, eri_ab
2701 REAL(dp),
DIMENSION(:, :),
INTENT(IN) :: p_a_mat, p_b_mat
2702 REAL(dp),
DIMENSION(:, :),
INTENT(INOUT) :: ks_a_ref
2703 LOGICAL,
INTENT(IN) :: tr_mixed_eri
2704 TYPE(mp_comm_type),
INTENT(IN) :: comm_exchange
2706 CHARACTER(LEN=*),
PARAMETER :: routinen =
'build_subspace_spin_fock_matrix'
2708 INTEGER :: handle, i1, i12, i12l, i2, i3, i34, &
2709 i34l, i4, irptr, m1, m2, nindex, &
2710 nmo_total, norb, spin1, spin2
2712 TYPE(mp_comm_type) :: mp_group
2714 CALL timeset(routinen, handle)
2716 norb =
SIZE(active_orbitals, 1)
2717 nmo_total =
SIZE(p_a_mat, 1)
2718 nindex = (nmo_total*(nmo_total + 1))/2
2719 IF (tr_mixed_eri)
THEN
2727 i1 = active_orbitals(m1, spin1)
2729 i2 = active_orbitals(m2, spin1)
2730 i12 = csr_idx_to_combined(i1, i2, nmo_total)
2731 IF (mod(i12 - 1, comm_exchange%num_pe) == comm_exchange%mepos)
THEN
2732 i12l = (i12 - 1)/comm_exchange%num_pe + 1
2733 irptr = eri_aa%rowptr_local(i12l) - 1
2734 DO i34l = 1, eri_aa%nzerow_local(i12l)
2735 i34 = eri_aa%colind_local(irptr + i34l)
2736 CALL csr_idx_from_combined(i34, nmo_total, i3, i4)
2737 erint = eri_aa%nzval_local%r_dp(irptr + i34l)
2740 ks_a_ref(i1, i2) = ks_a_ref(i1, i2) + erint*p_a_mat(i3, i4)
2741 IF (i12 /= i34)
THEN
2743 ks_a_ref(i3, i4) = ks_a_ref(i3, i4) + erint*p_a_mat(i1, i2)
2746 erint = -1.0_dp*erint
2748 ks_a_ref(i1, i3) = ks_a_ref(i1, i3) + erint*p_a_mat(i2, i4)
2751 ks_a_ref(i2, i3) = ks_a_ref(i2, i3) + erint*p_a_mat(i1, i4)
2755 ks_a_ref(i1, i4) = ks_a_ref(i1, i4) + erint*p_a_mat(i2, i3)
2757 IF (i1 /= i2 .AND. i3 /= i4)
THEN
2759 ks_a_ref(i2, i4) = ks_a_ref(i2, i4) + erint*p_a_mat(i1, i3)
2768 i1 = active_orbitals(m1, 1)
2770 i2 = active_orbitals(m2, 1)
2771 i12 = csr_idx_to_combined(i1, i2, nmo_total)
2772 IF (mod(i12 - 1, comm_exchange%num_pe) == comm_exchange%mepos)
THEN
2773 i12l = (i12 - 1)/comm_exchange%num_pe + 1
2774 irptr = eri_ab%rowptr_local(i12l) - 1
2775 DO i34l = 1, eri_ab%nzerow_local(i12l)
2776 i34 = eri_ab%colind_local(irptr + i34l)
2777 CALL csr_idx_from_combined(i34, nmo_total, i3, i4)
2778 erint = eri_ab%nzval_local%r_dp(irptr + i34l)
2780 IF (tr_mixed_eri)
THEN
2782 ks_a_ref(i3, i4) = ks_a_ref(i3, i4) + erint*p_b_mat(i1, i2)
2785 ks_a_ref(i1, i2) = ks_a_ref(i1, i2) + erint*p_b_mat(i3, i4)
2793 i1 = active_orbitals(m1, spin1)
2795 i2 = active_orbitals(m2, spin1)
2796 ks_a_ref(i2, i1) = ks_a_ref(i1, i2)
2799 CALL mp_group%set_handle(eri_aa%mp_group%get_handle())
2800 CALL mp_group%sum(ks_a_ref)
2802 CALL timestop(handle)
2804 END SUBROUTINE build_subspace_spin_fock_matrix
2816 SUBROUTINE create_pro_basis(pro_basis_set, zval, ishell, nshell, lnam)
2817 TYPE(gto_basis_set_type),
POINTER :: pro_basis_set
2818 INTEGER,
INTENT(IN) :: zval, ishell
2819 INTEGER,
DIMENSION(:),
INTENT(IN) :: nshell
2820 CHARACTER(len=*),
DIMENSION(:),
INTENT(IN) :: lnam
2822 CHARACTER(len=6),
DIMENSION(:),
POINTER :: sym
2824 INTEGER,
DIMENSION(4, 7) :: ne
2825 INTEGER,
DIMENSION(:),
POINTER :: lq, nq
2826 REAL(kind=dp),
DIMENSION(:),
POINTER :: zet
2827 TYPE(sto_basis_set_type),
POINTER :: sto_basis_set
2829 cpassert(.NOT.
ASSOCIATED(pro_basis_set))
2830 NULLIFY (sto_basis_set)
2837 ne(l, i) = ptable(zval)%e_conv(l - 1) - 2*nj*(i - l)
2838 ne(l, i) = max(ne(l, i), 0)
2839 ne(l, i) = min(ne(l, i), 2*nj)
2842 ALLOCATE (nq(ishell), lq(ishell), zet(ishell), sym(ishell))
2845 SELECT CASE (lnam(i))
2855 cpabort(
"Wrong l QN")
2858 zet(i) = srules(zval, ne, nq(1), lq(1))
2860 CALL allocate_sto_basis_set(sto_basis_set)
2861 CALL set_sto_basis_set(sto_basis_set, nshell=1, nq=nq, lq=lq, zet=zet, symbol=sym)
2862 CALL create_gto_from_sto_basis(sto_basis_set, pro_basis_set, 6)
2863 pro_basis_set%norm_type = 2
2864 CALL init_orb_basis_set(pro_basis_set)
2865 CALL deallocate_sto_basis_set(sto_basis_set)
2867 END SUBROUTINE create_pro_basis
2874 SUBROUTINE update_density_ao(active_space_env, rho_ao)
2875 TYPE(active_space_type),
POINTER :: active_space_env
2876 TYPE(dbcsr_p_type),
DIMENSION(:),
POINTER :: rho_ao
2878 INTEGER :: ispin, nao, nmo, nspins
2879 TYPE(cp_fm_type) :: r, u
2880 TYPE(cp_fm_type),
POINTER :: c_active, p_active_mo
2881 TYPE(dbcsr_type),
POINTER :: p_inactive_ao
2882 TYPE(mo_set_type),
DIMENSION(:),
POINTER :: mos_active
2886 nspins = active_space_env%nspins
2887 mos_active => active_space_env%mos_active
2888 DO ispin = 1, nspins
2890 p_inactive_ao => active_space_env%pmat_inactive(ispin)%matrix
2893 CALL dbcsr_copy(rho_ao(ispin)%matrix, p_inactive_ao)
2896 p_active_mo => active_space_env%p_active(ispin)
2899 CALL cp_fm_create(r, p_active_mo%matrix_struct)
2900 CALL cp_fm_to_fm(p_active_mo, r)
2903 CALL get_mo_set(mos_active(ispin), mo_coeff=c_active, nao=nao, nmo=nmo)
2904 CALL cp_fm_create(u, c_active%matrix_struct)
2905 CALL parallel_gemm(
"N",
"N", nao, nmo, nmo, 1.0_dp, c_active, r, 0.0_dp, u)
2907 CALL cp_dbcsr_plus_fm_fm_t(sparse_matrix=rho_ao(ispin)%matrix, &
2908 matrix_v=u, matrix_g=c_active, ncol=nmo, alpha=1.0_dp)
2910 CALL cp_fm_release(r)
2911 CALL cp_fm_release(u)
2914 END SUBROUTINE update_density_ao
2926 LOGICAL FUNCTION eri_fcidump_print_func(this, i, j, k, l, val)
RESULT(cont)
2927 CLASS(eri_fcidump_print),
INTENT(inout) :: this
2928 INTEGER,
INTENT(in) :: i, j, k, l
2929 REAL(kind=dp),
INTENT(in) :: val
2932 IF (this%unit_nr > 0)
THEN
2933 WRITE (this%unit_nr,
"(ES23.16,4I4)") val, i + this%bra_start - 1, j + this%bra_start - 1, &
2934 & k + this%ket_start - 1, l + this%ket_start - 1
2938 END FUNCTION eri_fcidump_print_func
2950 LOGICAL FUNCTION eri_fcidump_checksum_func(this, i, j, k, l, val)
RESULT(cont)
2951 CLASS(eri_fcidump_checksum),
INTENT(inout) :: this
2952 INTEGER,
INTENT(in) :: i, j, k, l
2953 REAL(kind=dp),
INTENT(in) :: val
2959 this%checksum = this%checksum + abs(val)
2962 END FUNCTION eri_fcidump_checksum_func
2970 SUBROUTINE print_pmat_noon(active_space_env, iw)
2971 TYPE(active_space_type),
POINTER :: active_space_env
2974 INTEGER :: i1, i2, ii, ispin, jm, m1, m2, &
2976 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:, :) :: noon, pmat
2977 TYPE(cp_fm_type),
POINTER :: p_active
2979 nspins = active_space_env%nspins
2980 nmo_active = active_space_env%nmo_active
2982 ALLOCATE (noon(nmo_active, nspins))
2983 ALLOCATE (pmat(nmo_active, nmo_active))
2985 DO ispin = 1, nspins
2986 p_active => active_space_env%p_active(ispin)
2987 noon(:, ispin) = 0.0_dp
2990 DO i1 = 1, nmo_active
2991 m1 = active_space_env%active_orbitals(i1, ispin)
2992 DO i2 = 1, nmo_active
2993 m2 = active_space_env%active_orbitals(i2, ispin)
2994 CALL cp_fm_get_element(p_active, m1, m2, pmat(i1, i2))
2999 WRITE (iw,
'(/,T3,A,I2,A)')
"Active space density matrix for spin ", ispin
3000 DO i1 = 1, nmo_active
3001 DO ii = 1, nmo_active, 8
3002 jm = min(7, nmo_active - ii)
3003 WRITE (iw,
'(T3,6(F9.4))') (pmat(i1, ii + i2), i2=0, jm)
3009 CALL diamat_all(pmat, noon(:, ispin))
3012 WRITE (iw,
'(/,T3,A,I2,A)')
"Natural orbitals occupation numbers for spin ", ispin
3013 DO i1 = 1, nmo_active, 8
3014 jm = min(7, nmo_active - i1)
3016 WRITE (iw,
'(T3,6(F9.4))') (noon(nmo_active - i1 - i2 + 1, ispin), i2=0, jm)
3025 END SUBROUTINE print_pmat_noon
3033 SUBROUTINE local_fci_embedding(qs_env, active_space_env, as_input)
3034 TYPE(qs_environment_type),
POINTER :: qs_env
3035 TYPE(active_space_type),
POINTER :: active_space_env
3036 TYPE(section_vals_type),
POINTER :: as_input
3038 CHARACTER(len=*),
PARAMETER :: routinen =
'local_fci_embedding'
3040 INTEGER :: handle, iter, iw, max_iter
3041 LOGICAL :: converged, do_scf_embedding
3042 REAL(kind=dp) :: delta_e, energy_corr, energy_new, &
3043 energy_old, energy_scf, eps_iter, t1, &
3045 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:) :: p_act_mo_a, p_act_mo_b
3046 TYPE(cp_logger_type),
POINTER :: logger
3047 TYPE(dbcsr_p_type),
DIMENSION(:),
POINTER :: rho_ao
3048 TYPE(dft_control_type),
POINTER :: dft_control
3049 TYPE(mo_set_type),
DIMENSION(:),
POINTER :: mos_active
3050 TYPE(mp_para_env_type),
POINTER :: para_env
3051 TYPE(qs_energy_type),
POINTER :: energy
3052 TYPE(qs_ks_env_type),
POINTER :: ks_env
3053 TYPE(qs_rho_type),
POINTER :: rho
3055 CALL timeset(routinen, handle)
3058 logger => cp_get_default_logger()
3059 iw = cp_logger_get_default_io_unit(logger)
3061 CALL get_qs_env(qs_env, para_env=para_env, dft_control=dft_control)
3063 CALL section_vals_val_get(as_input,
"SCF_EMBEDDING", l_val=do_scf_embedding)
3064 active_space_env%do_scf_embedding = do_scf_embedding
3065 CALL section_vals_val_get(as_input,
"MAX_ITER", i_val=max_iter)
3066 IF (max_iter < 0) cpabort(
"Specify a non-negative number of max iterations.")
3067 CALL section_vals_val_get(as_input,
"EPS_ITER", r_val=eps_iter)
3068 IF (eps_iter < 0.0) cpabort(
"Specify a non-negative convergence threshold.")
3070 CALL get_qs_env(qs_env, rho=rho, energy=energy, ks_env=ks_env)
3071 CALL qs_rho_get(rho, rho_ao=rho_ao)
3074 WRITE (unit=iw, fmt=
"(/,T2,A,/)") &
3075 "RANGE-SEPARATED DFT EMBEDDING WITH LIBFCI SOLVER"
3076 WRITE (iw,
'(T3,A,T68,I12)')
"Max. iterations", max_iter
3077 WRITE (iw,
'(T3,A,T68,E12.4)')
"Conv. threshold", eps_iter
3078 WRITE (iw,
'(T3,A,T68,A)')
"Density mixer", trim(active_space_mixing_label(active_space_env))
3079 WRITE (iw,
'(T3,A,T66,F14.2)')
"Mixing alpha", active_space_env%alpha
3080 WRITE (unit=iw, fmt=
"(/,T3,A,T11,A,T21,A,T34,A,T55,A,T75,A,/,T3,A)") &
3081 "Iter",
"Update",
"Time",
"Corr. energy",
"Total energy",
"Change", repeat(
"-", 78)
3086 energy_scf = active_space_env%energy_ref
3087 energy_new = energy_scf
3088 mos_active => active_space_env%mos_active
3090 DO WHILE (iter < max_iter)
3093 IF (active_space_env%nspins == 2)
THEN
3094 CALL solve_active_space_fci(active_space_env, para_env, p_act_mo_a, p_act_mo_b)
3095 active_space_env%energy_total = active_space_env%energy_inactive + active_space_env%energy_active
3096 CALL update_active_density(p_act_mo_a, active_space_env, p_act_mo_b)
3097 DEALLOCATE (p_act_mo_a, p_act_mo_b)
3099 CALL solve_active_space_fci(active_space_env, para_env, p_act_mo_a)
3100 active_space_env%energy_total = active_space_env%energy_inactive + active_space_env%energy_active
3101 CALL update_active_density(p_act_mo_a, active_space_env)
3102 DEALLOCATE (p_act_mo_a)
3105 energy_old = energy_new
3106 energy_new = active_space_env%energy_total
3107 energy_corr = energy_new - energy_scf
3108 delta_e = energy_new - energy_old
3114 fmt=
"(T3,I4,T11,A,T19,F6.1,T28,F18.10,T49,F18.10,T70,ES11.2)") &
3115 iter, trim(active_space_mixing_label(active_space_env)), &
3116 t1 - t2, energy_corr, energy_new, delta_e
3120 CALL update_density_ao(active_space_env, rho_ao)
3121 CALL qs_rho_update_rho(rho, qs_env=qs_env)
3122 CALL qs_ks_did_change(qs_env%ks_env, rho_changed=.true.)
3123 CALL evaluate_core_matrix_traces(qs_env)
3124 CALL qs_ks_build_kohn_sham_matrix(qs_env, calculate_forces=.false., &
3125 just_energy=.false., &
3126 ext_xc_section=active_space_env%xc_section)
3128 active_space_env%energy_ref = energy%total
3129 CALL calculate_operators(mos_active, qs_env, active_space_env)
3130 CALL subspace_fock_matrix(active_space_env, dft_control%roks)
3132 IF (.NOT. active_space_env%do_scf_embedding)
THEN
3134 WRITE (unit=iw, fmt=
"(/,T3,A,I5,A)") &
3135 "*** one-shot embedding correction finished ***"
3139 ELSEIF (abs(delta_e) <= eps_iter)
THEN
3141 WRITE (unit=iw, fmt=
"(/,T3,A,I5,A)") &
3142 "*** rs-DFT embedding run converged in ", iter,
" iteration(s) ***"
3149 IF (.NOT. converged)
THEN
3151 WRITE (unit=iw, fmt=
"(/,T3,A,I5,A)") &
3152 "*** rs-DFT embedding did not converged after ", iter,
" iteration(s) ***"
3156 energy%total = active_space_env%energy_total
3159 WRITE (unit=iw, fmt=
"(/,T3,A)") &
3160 "Final energy contributions:"
3161 WRITE (unit=iw, fmt=
"(T6,A,T56,F20.10)") &
3162 "Inactive energy:", active_space_env%energy_inactive
3163 WRITE (unit=iw, fmt=
"(T6,A,T56,F20.10)") &
3164 "Active energy:", active_space_env%energy_active
3165 WRITE (unit=iw, fmt=
"(T6,A,T56,F20.10)") &
3166 "Correlation energy:", energy_corr
3167 WRITE (unit=iw, fmt=
"(T6,A,T56,F20.10)") &
3168 "Total rs-DFT energy:", active_space_env%energy_total
3171 CALL print_pmat_noon(active_space_env, iw)
3172 CALL para_env%sync()
3173 CALL timestop(handle)
3175 END SUBROUTINE local_fci_embedding
3183 SUBROUTINE rsdft_embedding(qs_env, active_space_env, as_input)
3184 TYPE(qs_environment_type),
POINTER :: qs_env
3185 TYPE(active_space_type),
POINTER :: active_space_env
3186 TYPE(section_vals_type),
POINTER :: as_input
3188 CHARACTER(len=*),
PARAMETER :: routinen =
'rsdft_embedding'
3192 CALL timeset(routinen, handle)
3193 cpabort(
"CP2K was compiled with the __NO_SOCKETS option!")
3195 mark_used(active_space_env)
3199 INTEGER :: iw, client_fd, socket_fd, iter, max_iter
3200 LOGICAL :: converged, do_scf_embedding, ionode
3201 REAL(kind=dp) :: delta_e, energy_corr, energy_new, &
3202 energy_old, energy_scf, eps_iter, t1, t2
3203 TYPE(cp_logger_type),
POINTER :: logger
3204 TYPE(dbcsr_p_type),
DIMENSION(:),
POINTER :: rho_ao
3205 TYPE(mo_set_type),
DIMENSION(:),
POINTER :: mos_active
3206 TYPE(mp_para_env_type),
POINTER :: para_env
3207 TYPE(qs_energy_type),
POINTER :: energy
3208 TYPE(qs_ks_env_type),
POINTER :: ks_env
3209 TYPE(qs_rho_type),
POINTER :: rho
3210 TYPE(dft_control_type),
POINTER :: dft_control
3212 CALL timeset(routinen, handle)
3216 logger => cp_get_default_logger()
3217 iw = cp_logger_get_default_io_unit(logger)
3219 CALL get_qs_env(qs_env, para_env=para_env, dft_control=dft_control)
3220 ionode = para_env%is_source()
3223 CALL section_vals_val_get(as_input,
"SCF_EMBEDDING", l_val=do_scf_embedding)
3224 active_space_env%do_scf_embedding = do_scf_embedding
3225 CALL section_vals_val_get(as_input,
"MAX_ITER", i_val=max_iter)
3226 IF (max_iter < 0) cpabort(
"Specify a non-negative number of max iterations.")
3227 CALL section_vals_val_get(as_input,
"EPS_ITER", r_val=eps_iter)
3228 IF (eps_iter < 0.0) cpabort(
"Specify a non-negative convergence threshold.")
3231 CALL initialize_socket(socket_fd, client_fd, as_input, ionode)
3232 CALL para_env%sync()
3235 CALL send_eri_to_client(client_fd, active_space_env, para_env)
3238 CALL get_qs_env(qs_env, rho=rho, energy=energy, ks_env=ks_env)
3239 CALL qs_rho_get(rho, rho_ao=rho_ao)
3242 WRITE (unit=iw, fmt=
"(/,T2,A,/)") &
3243 "RANGE-SEPARATED DFT EMBEDDING SELF-CONSISTENT OPTIMIZATION"
3245 WRITE (iw,
'(T3,A,T68,I12)')
"Max. iterations", max_iter
3246 WRITE (iw,
'(T3,A,T68,E12.4)')
"Conv. threshold", eps_iter
3247 WRITE (iw,
'(T3,A,T68,A)')
"Density mixer", trim(active_space_mixing_label(active_space_env))
3248 WRITE (iw,
'(T3,A,T66,F14.2)')
"Mixing alpha", active_space_env%alpha
3250 WRITE (unit=iw, fmt=
"(/,T3,A,T11,A,T21,A,T34,A,T55,A,T75,A,/,T3,A)") &
3251 "Iter",
"Update",
"Time",
"Corr. energy",
"Total energy",
"Change", repeat(
"-", 78)
3258 energy_scf = active_space_env%energy_ref
3259 energy_new = energy_scf
3260 mos_active => active_space_env%mos_active
3264 DO WHILE (iter < max_iter)
3269 CALL send_fock_to_client(client_fd, active_space_env, para_env)
3272 energy_old = energy_new
3273 energy_new = active_space_env%energy_total
3274 energy_corr = energy_new - energy_scf
3275 delta_e = energy_new - energy_old
3283 fmt=
"(T3,I4,T11,A,T19,F6.1,T28,F18.10,T49,F18.10,T70,ES11.2)") &
3284 iter, trim(active_space_mixing_label(active_space_env)), &
3285 t1 - t2, energy_corr, energy_new, delta_e
3290 CALL update_density_ao(active_space_env, rho_ao)
3293 CALL qs_rho_update_rho(rho, qs_env=qs_env)
3294 CALL qs_ks_did_change(qs_env%ks_env, rho_changed=.true.)
3296 CALL evaluate_core_matrix_traces(qs_env)
3299 CALL qs_ks_build_kohn_sham_matrix(qs_env, calculate_forces=.false., &
3300 just_energy=.false., &
3301 ext_xc_section=active_space_env%xc_section)
3304 active_space_env%energy_ref = energy%total
3307 CALL calculate_operators(mos_active, qs_env, active_space_env)
3310 CALL subspace_fock_matrix(active_space_env, dft_control%roks)
3313 IF (.NOT. active_space_env%do_scf_embedding)
THEN
3315 WRITE (unit=iw, fmt=
"(/,T3,A,I5,A)") &
3316 "*** one-shot embedding correction finished ***"
3321 ELSEIF (abs(delta_e) <= eps_iter)
THEN
3323 WRITE (unit=iw, fmt=
"(/,T3,A,I5,A)") &
3324 "*** rs-DFT embedding run converged in ", iter,
" iteration(s) ***"
3331 IF (.NOT. converged)
THEN
3333 WRITE (unit=iw, fmt=
"(/,T3,A,I5,A)") &
3334 "*** rs-DFT embedding did not converged after ", iter,
" iteration(s) ***"
3339 energy%total = active_space_env%energy_total
3343 WRITE (unit=iw, fmt=
"(/,T3,A)") &
3344 "Final energy contributions:"
3345 WRITE (unit=iw, fmt=
"(T6,A,T56,F20.10)") &
3346 "Inactive energy:", active_space_env%energy_inactive
3347 WRITE (unit=iw, fmt=
"(T6,A,T56,F20.10)") &
3348 "Active energy:", active_space_env%energy_active
3349 WRITE (unit=iw, fmt=
"(T6,A,T56,F20.10)") &
3350 "Correlation energy:", energy_corr
3351 WRITE (unit=iw, fmt=
"(T6,A,T56,F20.10)") &
3352 "Total rs-DFT energy:", active_space_env%energy_total
3356 CALL print_pmat_noon(active_space_env, iw)
3358 CALL finalize_socket(socket_fd, client_fd, as_input, ionode)
3359 CALL para_env%sync()
3362 CALL timestop(handle)
3364 END SUBROUTINE rsdft_embedding
3374 SUBROUTINE initialize_socket(socket_fd, client_fd, as_input, ionode)
3375 INTEGER,
INTENT(OUT) :: socket_fd, client_fd
3376 TYPE(section_vals_type),
INTENT(IN),
POINTER :: as_input
3377 LOGICAL,
INTENT(IN) :: ionode
3379 CHARACTER(len=*),
PARAMETER :: routinen =
'initialize_socket'
3380 INTEGER,
PARAMETER :: backlog = 10
3382 CHARACTER(len=default_path_length) :: hostname
3383 INTEGER :: handle, iw, port, protocol
3385 TYPE(cp_logger_type),
POINTER :: logger
3387 CALL timeset(routinen, handle)
3389 logger => cp_get_default_logger()
3390 iw = cp_logger_get_default_io_unit(logger)
3393 CALL section_vals_val_get(as_input,
"SOCKET%INET", l_val=inet)
3399 CALL section_vals_val_get(as_input,
"SOCKET%HOST", c_val=hostname)
3400 CALL section_vals_val_get(as_input,
"SOCKET%PORT", i_val=port)
3403 CALL open_bind_socket(socket_fd, protocol, port, trim(hostname)//c_null_char)
3404 WRITE (iw,
'(/,T2,A,A)')
"@SERVER: Created socket with address ", trim(hostname)
3408 WRITE (iw,
'(T2,A)')
"@SERVER: Waiting for requests..."
3410 WRITE (iw,
'(T2,A,I2)')
"@SERVER: Accepted socket with fd ", client_fd
3413 CALL timestop(handle)
3415 END SUBROUTINE initialize_socket
3424 SUBROUTINE finalize_socket(socket_fd, client_fd, as_input, ionode)
3425 INTEGER,
INTENT(IN) :: socket_fd, client_fd
3426 TYPE(section_vals_type),
INTENT(IN),
POINTER :: as_input
3427 LOGICAL,
INTENT(IN) :: ionode
3429 CHARACTER(len=*),
PARAMETER :: routinen =
'finalize_socket'
3430 INTEGER,
PARAMETER :: header_len = 12
3432 CHARACTER(len=default_path_length) :: hostname
3435 CALL timeset(routinen, handle)
3437 CALL section_vals_val_get(as_input,
"SOCKET%HOST", c_val=hostname)
3447 IF (file_exists(trim(hostname)))
THEN
3452 CALL timestop(handle)
3454 END SUBROUTINE finalize_socket
3462 SUBROUTINE send_eri_to_client(client_fd, active_space_env, para_env)
3463 INTEGER,
INTENT(IN) :: client_fd
3464 TYPE(active_space_type),
INTENT(IN),
POINTER :: active_space_env
3465 TYPE(mp_para_env_type),
INTENT(IN),
POINTER :: para_env
3467 CHARACTER(len=*),
PARAMETER :: routinen =
'send_eri_to_client'
3468 INTEGER,
PARAMETER :: header_len = 12
3470 CHARACTER(len=default_string_length) ::
header
3471 INTEGER :: handle, iw
3472 LOGICAL :: ionode, restricted_orbitals
3473 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:) :: eri_aa, eri_ab, eri_bb, s_ab
3474 TYPE(cp_logger_type),
POINTER :: logger
3476 CALL timeset(routinen, handle)
3478 logger => cp_get_default_logger()
3479 iw = cp_logger_get_default_io_unit(logger)
3480 ionode = para_env%is_source()
3481 restricted_orbitals = active_space_env%restricted_orbitals
3483 ALLOCATE (eri_aa(active_space_env%nmo_active**4))
3484 CALL eri_to_array(active_space_env%eri, eri_aa, active_space_env%active_orbitals, 1, 1)
3485 IF (active_space_env%nspins == 2)
THEN
3486 ALLOCATE (eri_ab(active_space_env%nmo_active**4))
3487 ALLOCATE (eri_bb(active_space_env%nmo_active**4))
3488 IF (restricted_orbitals)
THEN
3492 CALL eri_to_array(active_space_env%eri, eri_ab, active_space_env%active_orbitals, 1, 2)
3493 CALL eri_to_array(active_space_env%eri, eri_bb, active_space_env%active_orbitals, 2, 2)
3496 ALLOCATE (s_ab(active_space_env%nmo_active**2))
3497 associate(act_indices_a => active_space_env%active_orbitals(:, 1), &
3498 act_indices_b => active_space_env%active_orbitals(:, 2))
3499 CALL subspace_matrix_to_array(active_space_env%sab_sub(1), s_ab, act_indices_a, act_indices_b)
3504 IF (ionode)
CALL writebuffer(client_fd,
"STATUS ", header_len)
3507 CALL para_env%sync()
3512 CALL para_env%bcast(
header, para_env%source)
3516 IF (trim(
header) ==
"READY")
THEN
3518 CALL para_env%sync()
3520 CALL writebuffer(client_fd,
"TWOBODY ", header_len)
3521 CALL writebuffer(client_fd, active_space_env%nspins)
3522 CALL writebuffer(client_fd, active_space_env%nmo_active)
3523 CALL writebuffer(client_fd, active_space_env%nelec_active)
3524 CALL writebuffer(client_fd, active_space_env%multiplicity)
3528 IF (active_space_env%nspins == 2)
THEN
3534 ELSE IF (trim(
header) ==
"RECEIVED")
THEN
3540 IF (active_space_env%nspins == 2)
THEN
3546 CALL para_env%sync()
3548 CALL timestop(handle)
3550 END SUBROUTINE send_eri_to_client
3558 SUBROUTINE send_fock_to_client(client_fd, active_space_env, para_env)
3559 INTEGER,
INTENT(IN) :: client_fd
3560 TYPE(active_space_type),
INTENT(INOUT),
POINTER :: active_space_env
3561 TYPE(mp_para_env_type),
INTENT(IN),
POINTER :: para_env
3563 CHARACTER(len=*),
PARAMETER :: routinen =
'send_fock_to_client'
3564 INTEGER,
PARAMETER :: header_len = 12
3566 CHARACTER(len=default_string_length) ::
header
3567 INTEGER :: handle, iw
3568 LOGICAL :: debug, ionode
3569 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:) :: fock_a, fock_b, p_act_mo_a, p_act_mo_b
3570 TYPE(cp_logger_type),
POINTER :: logger
3572 CALL timeset(routinen, handle)
3577 logger => cp_get_default_logger()
3578 iw = cp_logger_get_default_io_unit(logger)
3579 ionode = para_env%is_source()
3581 ALLOCATE (p_act_mo_a(active_space_env%nmo_active**2))
3582 ALLOCATE (fock_a(active_space_env%nmo_active**2))
3583 IF (active_space_env%nspins == 2)
THEN
3584 ALLOCATE (p_act_mo_b(active_space_env%nmo_active**2))
3585 ALLOCATE (fock_b(active_space_env%nmo_active**2))
3589 associate(act_indices => active_space_env%active_orbitals(:, 1))
3590 CALL subspace_matrix_to_array(active_space_env%fock_sub(1), fock_a, act_indices, act_indices)
3593 IF (active_space_env%nspins == 2)
THEN
3594 associate(act_indices => active_space_env%active_orbitals(:, 2))
3595 CALL subspace_matrix_to_array(active_space_env%fock_sub(2), fock_b, act_indices, act_indices)
3600 IF (ionode)
CALL writebuffer(client_fd,
"STATUS ", header_len)
3604 CALL para_env%sync()
3606 IF (debug .AND. iw > 0)
WRITE (iw, *)
"@SERVER: Waiting for messages..."
3609 CALL para_env%bcast(
header, para_env%source)
3611 IF (debug .AND. iw > 0)
WRITE (iw, *)
"@SERVER: Message from client: ", trim(
header)
3613 IF (trim(
header) ==
"READY")
THEN
3615 CALL para_env%sync()
3617 CALL writebuffer(client_fd,
"ONEBODY ", header_len)
3618 CALL writebuffer(client_fd, active_space_env%energy_inactive)
3622 IF (active_space_env%nspins == 2)
THEN
3627 ELSE IF (trim(
header) ==
"HAVEDATA")
THEN
3629 CALL para_env%sync()
3631 IF (debug .AND. iw > 0)
WRITE (iw, *)
"@SERVER: Qiskit has data to transfer"
3632 CALL writebuffer(client_fd,
"GETDENSITY ", header_len)
3635 CALL readbuffer(client_fd, active_space_env%energy_active)
3636 CALL readbuffer(client_fd, p_act_mo_a,
SIZE(p_act_mo_a))
3637 IF (active_space_env%nspins == 2)
THEN
3638 CALL readbuffer(client_fd, p_act_mo_b,
SIZE(p_act_mo_b))
3643 CALL para_env%bcast(active_space_env%energy_active, para_env%source)
3644 CALL para_env%bcast(p_act_mo_a, para_env%source)
3645 IF (active_space_env%nspins == 2)
THEN
3646 CALL para_env%bcast(p_act_mo_b, para_env%source)
3650 active_space_env%energy_total = active_space_env%energy_inactive + active_space_env%energy_active
3653 IF (active_space_env%nspins == 2)
THEN
3654 CALL update_active_density(p_act_mo_a, active_space_env, p_act_mo_b)
3656 CALL update_active_density(p_act_mo_a, active_space_env)
3665 DEALLOCATE (p_act_mo_a)
3667 IF (active_space_env%nspins == 2)
THEN
3668 DEALLOCATE (p_act_mo_b)
3672 CALL para_env%sync()
3674 CALL timestop(handle)
3676 END SUBROUTINE send_fock_to_client
Types and set/get functions for auxiliary density matrix methods.
subroutine, public get_admm_env(admm_env, mo_derivs_aux_fit, mos_aux_fit, sab_aux_fit, sab_aux_fit_asymm, sab_aux_fit_vs_orb, matrix_s_aux_fit, matrix_s_aux_fit_kp, matrix_s_aux_fit_vs_orb, matrix_s_aux_fit_vs_orb_kp, task_list_aux_fit, matrix_ks_aux_fit, matrix_ks_aux_fit_kp, matrix_ks_aux_fit_im, matrix_ks_aux_fit_dft, matrix_ks_aux_fit_hfx, matrix_ks_aux_fit_dft_kp, matrix_ks_aux_fit_hfx_kp, rho_aux_fit, rho_aux_fit_buffer, admm_dm)
Get routine for the ADMM env.
subroutine, public admm_env_release(admm_env)
releases the ADMM environment, cleans up all types
Define the atomic kind types and their sub types.
pure real(dp) function, public srules(z, ne, n, l)
...
subroutine, public deallocate_sto_basis_set(sto_basis_set)
...
subroutine, public allocate_sto_basis_set(sto_basis_set)
...
subroutine, public create_gto_from_sto_basis(sto_basis_set, gto_basis_set, ngauss, ortho)
...
subroutine, public set_sto_basis_set(sto_basis_set, name, nshell, symbol, nq, lq, zet)
...
subroutine, public init_orb_basis_set(gto_basis_set)
Initialise a Gaussian-type orbital (GTO) basis set data set.
Handles all functions related to the CELL.
subroutine, public write_cell_low(cell, unit_str, output_unit, label)
Write the cell parameters to the output unit.
subroutine, public set_cell_param(cell, cell_length, cell_angle, periodic, do_init_cell)
Sets the cell using the internal parameters (a,b,c) (alpha,beta,gamma) using the convention: a parall...
subroutine, public init_cell(cell, hmat, periodic)
Initialise/readjust a simulation cell after hmat has been changed.
Handles all functions related to the CELL.
integer, parameter, public use_perd_xyz
integer, parameter, public use_perd_none
methods related to the blacs parallel environment
integer, parameter, public blacs_grid_square
subroutine, public cp_blacs_env_release(blacs_env)
releases the given blacs_env
subroutine, public cp_blacs_env_create(blacs_env, para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
allocates and initializes a type that represent a blacs context
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_multiply(transa, transb, alpha, matrix_a, matrix_b, beta, matrix_c, first_row, last_row, first_column, last_column, first_k, last_k, retain_sparsity, filter_eps, flop)
...
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_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset, transposed)
...
real(kind=dp) function, public dbcsr_get_occupation(matrix)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_release(matrix)
...
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.
subroutine, public cp_dbcsr_plus_fm_fm_t(sparse_matrix, matrix_v, matrix_g, ncol, alpha, keep_sparsity, symmetry_mode)
performs the multiplication sparse_matrix+dense_mat*dens_mat^T if matrix_g is not explicitly given,...
subroutine, public cp_dbcsr_m_by_n_from_template(matrix, template, m, n, sym)
Utility function to create an arbitrary shaped dbcsr matrix with the same processor grid as the templ...
subroutine, public cp_dbcsr_write_sparse_matrix(sparse_matrix, before, after, qs_env, para_env, first_row, last_row, first_col, last_col, scale, output_unit, omit_headers, cartesian_basis)
...
Utility routines to open and close files. Tracking of preconnections.
subroutine, public open_file(file_name, file_status, file_form, file_action, file_position, file_pad, unit_number, debug, skip_get_unit_number, file_access)
Opens the requested file using a free unit number.
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.
logical function, public file_exists(file_name)
Checks if file exists, considering also the file discovery mechanism.
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_column_scale(matrixa, scaling)
scales column i of matrix a with scaling(i)
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_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_get_element(matrix, irow_global, icol_global, alpha, local)
returns an element of a fm this value is valid on every cpu using this call is expensive
subroutine, public cp_fm_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
subroutine, public cp_fm_set_element(matrix, irow_global, icol_global, alpha)
sets an element of a matrix
subroutine, public cp_fm_init_random(matrix, ncol, start_col)
fills a matrix with random numbers
subroutine, public cp_fm_write_formatted(fm, unit, header, value_format)
Write out a full matrix in plain text.
various routines to log and control the output. The idea is that decisions about where to log should ...
integer function, public cp_logger_get_default_io_unit(logger)
returns the unit nr for the ionode (-1 on all other processors) skips as well checks if the procs cal...
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, parameter, public debug_print_level
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
integer, parameter, public medium_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,...
integer, parameter, public cp_p_file
integer, parameter, public high_print_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...
integer, parameter, public silent_print_level
A wrapper around pw_to_cube() which accepts particle_list_type.
subroutine, public cp_pw_to_cube(pw, unit_nr, title, particles, zeff, stride, max_file_size_mb, zero_tails, silent, mpi_io)
...
Module to compute the error function of a complex argument.
elemental complex(kind=dp) function, public erfz_fast(z)
Computes the error function of a complex argument using the Poppe and Wijers algorithm.
Types to describe group distributions.
Types and set/get functions for HFX.
subroutine, public hfx_create(x_data, para_env, hfx_section, atomic_kind_set, qs_kind_set, particle_set, dft_control, cell, orb_basis, ri_basis, nelectron_total, nkp_grid)
This routine allocates and initializes all types in hfx_data
subroutine, public hfx_release(x_data)
This routine deallocates all data structures
Defines the basic variable types.
integer, parameter, public int_8
integer, parameter, public dp
integer, parameter, public default_string_length
integer, parameter, public default_path_length
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, inversion_symmetry_only, symmetry_backend, symmetry_reduction_method, gamma_centered)
Retrieve information from a kpoint environment.
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
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
real(kind=dp), parameter, public rootpi
real(kind=dp), parameter, public fourpi
real(kind=dp), parameter, public twopi
Collection of simple mathematical functions and subroutines.
subroutine, public diamat_all(a, eigval, dac)
Diagonalize the symmetric n by n matrix a using the LAPACK library. Only the upper triangle of matrix...
Utility routines for the memory handling.
Interface to the message passing library MPI.
subroutine, public mp_para_env_release(para_env)
releases the para object (to be called when you don't want anymore the shared copy of this object)
Calls routines to get RI integrals and calculate total energies.
subroutine, public grep_rows_in_subgroups(para_env, para_env_sub, mo_coeff, gd_array, c)
...
subroutine, public build_dbcsr_from_rows(para_env_sub, mo_coeff_to_build, cread, mat_munu, gd_array, eps_filter)
Encapsulate the building of dbcsr_matrices mo_coeff_(v,o,all)
subroutine, public create_mat_munu(mat_munu, qs_env, eps_grid, blacs_env_sub, do_ri_aux_basis, do_mixed_basis, group_size_prim, do_alloc_blocks_from_nbl, do_kpoints, sab_orb_sub, dbcsr_sym_type)
Encapsulate the building of dbcsr_matrix mat_munu.
integer, parameter, public mt0d
basic linear algebra operations for full matrixes
represent a simple array based list of the given type
Define the data structure for the particle information.
Periodic Table related data definitions.
type(atom), dimension(0:nelem), public ptable
Definition of physical constants:
real(kind=dp), parameter, public angstrom
real(kind=dp), parameter, public bohr
computes preconditioners, and implements methods to apply them currently used in qs_ot
methods of pw_env that have dependence on qs_env
subroutine, public pw_env_rebuild(pw_env, qs_env, external_para_env)
rebuilds the pw_env data (necessary if cell or cutoffs change)
subroutine, public pw_env_create(pw_env)
creates a pw_env, if qs_env is given calls pw_env_rebuild
container for various plainwaves related things
subroutine, public pw_env_release(pw_env, para_env)
releases the given pw_env (see doc/ReferenceCounting.html)
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
subroutine, public pw_gauss_damp(pw, omega)
Multiply all data points with a Gaussian damping factor Needed for longrange Coulomb potential V(\vec...
subroutine, public pw_compl_gauss_damp(pw, omega)
Multiply all data points with a Gaussian damping factor Needed for longrange Coulomb potential V(\vec...
functions related to the poisson solver on regular grids
integer, parameter, public pw_poisson_periodic
integer, parameter, public periodic3d
integer, parameter, public analytic0d
integer, parameter, public pw_poisson_analytic
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
The module to read/write QCSchema HDF5 files for interfacing CP2K with other programs.
subroutine, public qcschema_env_release(qcschema_env)
Releases the allocated memory of a qcschema environment.
subroutine, public qcschema_env_create(qcschema_env, qs_env)
Create and initialize a qcschema object from a quickstep environment.
subroutine, public qcschema_to_hdf5(qcschema_env, filename)
Writes a qcschema object to an hdf5 file.
Local FCI solver interface for active-space embedding.
subroutine, public solve_active_space_fci(active_space_env, para_env, p_act_mo_a, p_act_mo_b)
Solve the current active-space Hamiltonian with the local FCI kernel.
Determine active space Hamiltonian.
subroutine eri_fcidump_set(this, bra_start, ket_start)
Sets the starting indices of the bra and ket.
subroutine, public active_space_main(qs_env)
Main method for determining the active space Hamiltonian.
Dense density mixing for active-space embedding.
character(len=15) function, public active_space_mixing_label(active_space_env)
Return the current active-space mixer label for iteration output.
subroutine, public update_active_density(p_act_mo_a, active_space_env, p_act_mo_b)
Update active space density matrix from Fortran arrays.
subroutine, public initialize_active_space_mixing(active_space_env, as_input)
Initialize the density mixer used in the self-consistent active-space embedding loop.
The types needed for the calculation of active space Hamiltonians.
subroutine, public csr_idx_from_combined(ij, n, i, j)
extracts indices i and j from combined index ij
integer function, public csr_idx_to_combined(i, j, n)
calculates combined index (ij)
subroutine, public create_active_space_type(active_space_env)
Creates an active space environment type, nullifying all quantities.
Contains utility routines for the active space module.
subroutine, public eri_to_array(eri_env, array, active_orbitals, spin1, spin2)
Copy the eri tensor for spins isp1 and isp2 to a standard 1D Fortran array.
subroutine, public subspace_matrix_to_array(source_matrix, target_array, row_index, col_index)
Copy a (square portion) of a cp_fm_type matrix to a standard 1D Fortran array.
Calculate the plane wave density by collocating the primitive Gaussian functions (pgf).
subroutine, public calculate_wavefunction(mo_vectors, ivector, rho, rho_gspace, atomic_kind_set, qs_kind_set, cell, dft_control, particle_set, pw_env, basis_type)
maps a given wavefunction on the grid
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.
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_build_kohn_sham_matrix(qs_env, calculate_forces, just_energy, print_active, ext_ks_matrix, ext_xc_section)
routine where the real calculations are made: the KS matrix is calculated
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 set_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, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, vee, neighbor_list_id, kpoints, 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, subsys, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env)
...
Definition and initialisation of the mo data type.
subroutine, public write_mo_set_to_output_unit(mo_set, qs_kind_set, particle_set, dft_section, before, kpoint, final_mos, spin, solver_method, rtp, cpart, sim_step, umo_set, qs_env)
Write MO information to output file (eigenvalues, occupation numbers, coefficients)
collects routines that perform operations directly related to MOs
Definition and initialisation of the mo data type.
subroutine, public init_mo_set(mo_set, fm_pool, fm_ref, fm_struct, name, counter)
initializes an allocated mo_set. eigenvalues, mo_coeff, occupation_numbers are valid only after this ...
subroutine, public allocate_mo_set(mo_set, nao, nmo, nelectron, n_el_f, maxocc, flexible_electron_count)
Allocates a mo set and partially initializes it (nao,nmo,nelectron, and flexible_electron_count are v...
subroutine, public 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.
subroutine, public release_neighbor_list_sets(nlists)
releases an array of neighbor_list_sets
an eigen-space solver for the generalised symmetric eigenvalue problem for sparse matrices,...
subroutine, public ot_eigensolver(matrix_h, matrix_s, matrix_orthogonal_space_fm, matrix_c_fm, preconditioner, eps_gradient, iter_max, size_ortho_space, silent, ot_settings)
...
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_compute_properties(qs_env, wf_type, do_mp2)
computes properties for a given hamilonian using the current wfn
types that represent a quickstep subsys
subroutine, public qs_subsys_get(subsys, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell, cell_ref, use_ref_cell, energy, force, qs_kind_set, cp_subsys, nelectron_total, nelectron_spin)
...
parameters that control an scf iteration
Implements UNIX and INET sockets.
generate the tasks lists used by collocate and integrate routines
subroutine, public generate_qs_task_list(ks_env, task_list, basis_type, reorder_rs_grid_ranks, skip_load_balance_distributed, pw_env_external, sab_orb_external, ext_kpoints)
...
subroutine, public deallocate_task_list(task_list)
deallocates the components and the object itself
subroutine, public allocate_task_list(task_list)
allocates and initialised the components of the task_list_type
All kind of helpful little routines.
pure integer function, dimension(2), public get_limit(m, n, me)
divide m entries into n parts, return size of part me
void open_bind_socket(int *psockfd, int *inet, int *port, char *host)
Opens and binds a socket.
void writebuffer(int *psockfd, char *data, int *plen)
Writes to a socket.
void accept_socket(int *psockfd, int *pclientfd)
Listens to a socket.
void close_socket(int *psockfd)
Closes a socket.
void listen_socket(int *psockfd, int *backlog)
Listens to a socket.
void readbuffer(int *psockfd, char *data, int *plen)
Reads from a socket.
void remove_socket_file(char *host)
Removes a socket file.
stores some data used in wavefunction fitting
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
keeps the information about the structure of a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
Contains information about kpoints.
stores all the informations relevant to an mpi environment
represent a list of objects
contained for different pw related things
contains all the informations needed by the fft based poisson solvers
environment for the poisson solver
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
The full QCSchema output type. For more information refer to: https://molssi-qc-schema....
Abstract function object for the eri_type_eri_foreach method.
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.