90#include "./base/base_uses.f90"
96 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'eeq_method'
98 INTEGER,
PARAMETER :: maxElem = 86
101 REAL(KIND=
dp),
PARAMETER :: rcov(1:maxelem) = [&
102 & 0.32_dp, 0.46_dp, 1.20_dp, 0.94_dp, 0.77_dp, 0.75_dp, 0.71_dp, 0.63_dp, &
103 & 0.64_dp, 0.67_dp, 1.40_dp, 1.25_dp, 1.13_dp, 1.04_dp, 1.10_dp, 1.02_dp, &
104 & 0.99_dp, 0.96_dp, 1.76_dp, 1.54_dp, 1.33_dp, 1.22_dp, 1.21_dp, 1.10_dp, &
105 & 1.07_dp, 1.04_dp, 1.00_dp, 0.99_dp, 1.01_dp, 1.09_dp, 1.12_dp, 1.09_dp, &
106 & 1.15_dp, 1.10_dp, 1.14_dp, 1.17_dp, 1.89_dp, 1.67_dp, 1.47_dp, 1.39_dp, &
107 & 1.32_dp, 1.24_dp, 1.15_dp, 1.13_dp, 1.13_dp, 1.08_dp, 1.15_dp, 1.23_dp, &
108 & 1.28_dp, 1.26_dp, 1.26_dp, 1.23_dp, 1.32_dp, 1.31_dp, 2.09_dp, 1.76_dp, &
109 & 1.62_dp, 1.47_dp, 1.58_dp, 1.57_dp, 1.56_dp, 1.55_dp, 1.51_dp, 1.52_dp, &
110 & 1.51_dp, 1.50_dp, 1.49_dp, 1.49_dp, 1.48_dp, 1.53_dp, 1.46_dp, 1.37_dp, &
111 & 1.31_dp, 1.23_dp, 1.18_dp, 1.16_dp, 1.11_dp, 1.12_dp, 1.13_dp, 1.32_dp, &
112 & 1.30_dp, 1.30_dp, 1.36_dp, 1.31_dp, 1.38_dp, 1.42_dp]
129 INTEGER,
INTENT(IN) :: iounit, print_level
130 LOGICAL,
INTENT(IN) :: ext
132 CHARACTER(LEN=2) :: element_symbol
133 INTEGER :: enshift_type, iatom, ikind, natom
134 REAL(kind=
dp),
DIMENSION(:),
POINTER :: charges
139 mark_used(print_level)
141 CALL get_qs_env(qs_env, natom=natom, particle_set=particle_set, cell=cell)
145 cpassert(
ASSOCIATED(charges))
148 ALLOCATE (charges(natom))
153 CALL eeq_charges(qs_env, charges, eeq_sparam, 2, enshift_type)
158 IF (enshift_type == 0)
THEN
159 WRITE (unit=iounit, fmt=
"(/,T2,A)")
"EEQ Charges (External)"
160 ELSE IF (enshift_type == 1)
THEN
161 WRITE (unit=iounit, fmt=
"(/,T2,A)")
"EEQ Charges (Parametrization 2019 (Molecules))"
162 ELSE IF (enshift_type == 2)
THEN
163 WRITE (unit=iounit, fmt=
"(/,T2,A)")
"EEQ Charges (Parametrization 2019 (Crystals))"
165 cpabort(
"Unknown enshift_type")
167 WRITE (unit=iounit, fmt=
"(/,T2,A)") &
168 "# Atom Element Kind Atomic Charge"
172 element_symbol=element_symbol, &
174 WRITE (unit=iounit, fmt=
"(T4,I8,T18,A2,I10,T43,F12.4)") &
175 iatom, element_symbol, ikind, charges(iatom)
180 IF (.NOT. ext)
DEALLOCATE (charges)
194 SUBROUTINE eeq_charges(qs_env, charges, eeq_sparam, eeq_model, enshift_type, exclude, cn_max)
197 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: charges
199 INTEGER,
INTENT(IN) :: eeq_model, enshift_type
200 LOGICAL,
DIMENSION(:),
INTENT(IN),
OPTIONAL :: exclude
201 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: cn_max
203 CHARACTER(len=*),
PARAMETER :: routinen =
'eeq_charges'
205 INTEGER :: handle, iatom, ikind, iunit, jkind, &
207 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: kind_of
208 INTEGER,
DIMENSION(3) :: periodic
210 REAL(kind=
dp) :: ala, alb, eeq_energy, esg, kappa, &
211 lambda, scn, sgamma, totalcharge, xi
212 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: chia, cnumbers, efr, gam
213 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: gab
215 TYPE(
cell_type),
POINTER :: cell, cell_ref
218 TYPE(
dcnum_type),
ALLOCATABLE,
DIMENSION(:) :: dcnum
224 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
228 CALL timeset(routinen, handle)
231 qs_kind_set=qs_kind_set, &
232 atomic_kind_set=atomic_kind_set, &
233 particle_set=particle_set, &
235 blacs_env=blacs_env, &
237 dft_control=dft_control)
238 CALL get_qs_env(qs_env, nkind=nkind, natom=natom)
241 IF (para_env%is_source() .AND. logger%iter_info%print_level >=
medium_print_level)
THEN
247 totalcharge = dft_control%charge
249 CALL get_cnumbers(qs_env, cnumbers, dcnum, .false.)
252 IF (
PRESENT(cn_max))
THEN
254 cnumbers(iatom) = log(1.0_dp + exp(cn_max)) - log(1.0_dp + exp(cn_max - cnumbers(iatom)))
259 ALLOCATE (gab(nkind, nkind), gam(nkind))
264 CALL get_eeq_data(za, eeq_model, eta=gam(ikind), rad=ala)
269 gab(ikind, jkind) = sqrt(1._dp/(ala*ala + alb*alb))
276 IF (
PRESENT(exclude))
THEN
278 IF (exclude(ikind))
THEN
279 gam(ikind) = 1.0e30_dp
280 gab(ikind, :) = 0.0_dp
281 gab(:, ikind) = 0.0_dp
288 esg = 1.0_dp + exp(sgamma)
289 ALLOCATE (chia(natom))
292 ikind = kind_of(iatom)
296 IF (enshift_type == 1)
THEN
297 scn = cnumbers(iatom)/sqrt(cnumbers(iatom) + 1.0e-14_dp)
298 ELSE IF (enshift_type == 2)
THEN
299 scn = log(esg/(esg - cnumbers(iatom)))
301 cpabort(
"Unknown enshift_type")
303 chia(iatom) = xi - kappa*scn
308 IF (
PRESENT(exclude))
THEN
310 ikind = kind_of(iatom)
311 IF (exclude(ikind)) chia(iatom) = 0.0_dp
316 IF (dft_control%apply_period_efield .OR. dft_control%apply_efield .OR. &
317 dft_control%apply_efield_field)
THEN
318 ALLOCATE (efr(natom))
319 efr(1:natom) = 0.0_dp
321 chia(1:natom) = chia(1:natom) + efr(1:natom)
327 CALL get_cell(cell, periodic=periodic)
328 do_ewald = .NOT. all(periodic == 0)
333 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
338 silent=.true., pset=
"EEQ")
340 CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section=print_section)
342 CALL eeq_solver(charges, lambda, eeq_energy, &
343 particle_set, kind_of, cell, chia, gam, gab, &
344 para_env, blacs_env, dft_control, eeq_sparam, &
345 totalcharge=totalcharge, ewald=do_ewald, &
346 ewald_env=ewald_env, ewald_pw=ewald_pw, iounit=iunit)
350 DEALLOCATE (ewald_env, ewald_pw)
352 CALL eeq_solver(charges, lambda, eeq_energy, &
353 particle_set, kind_of, cell, chia, gam, gab, &
354 para_env, blacs_env, dft_control, eeq_sparam, &
355 totalcharge=totalcharge, iounit=iunit)
358 DEALLOCATE (gab, gam, chia)
360 CALL timestop(handle)
378 SUBROUTINE eeq_forces(qs_env, charges, dcharges, gradient, stress, &
379 eeq_sparam, eeq_model, enshift_type, response_only, exclude, cn_max)
382 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: charges, dcharges
383 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(INOUT) :: gradient
384 REAL(kind=
dp),
DIMENSION(3, 3),
INTENT(INOUT) :: stress
386 INTEGER,
INTENT(IN) :: eeq_model, enshift_type
387 LOGICAL,
INTENT(IN) :: response_only
388 LOGICAL,
DIMENSION(:),
INTENT(IN),
OPTIONAL :: exclude
389 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: cn_max
391 CHARACTER(len=*),
PARAMETER :: routinen =
'eeq_forces'
393 INTEGER :: handle, i, ia, iatom, ikind, iunit, &
394 jatom, jkind, katom, natom, nkind, za, &
396 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: kind_of
397 INTEGER,
DIMENSION(3) :: periodic
398 LOGICAL :: do_ewald, use_virial
399 LOGICAL,
ALLOCATABLE,
DIMENSION(:) :: default_present
400 REAL(kind=
dp) :: ala, alb, alpha, cn, ctot, dcnpdcn, dr, dr2, drk, elag, esg, fe, gam2, &
401 gama, grc, kappa, qlam, qq, qq1, qq2, rcut, scn, sgamma, subcells, totalcharge, xi
402 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: c_radius, cnumbers, gam, qlag
403 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: epforce, gab, pair_radius
404 REAL(kind=
dp),
DIMENSION(3) :: fdik, ri, rij, rik, rj
405 REAL(kind=
dp),
DIMENSION(3, 3) :: pvir
406 REAL(kind=
dp),
DIMENSION(:),
POINTER :: chrgx, dchia
409 TYPE(
cell_type),
POINTER :: cell, cell_ref
412 TYPE(
dcnum_type),
ALLOCATABLE,
DIMENSION(:) :: dcnum
422 DIMENSION(:),
POINTER :: nl_iterator
427 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
432 CALL timeset(routinen, handle)
435 qs_kind_set=qs_kind_set, &
436 atomic_kind_set=atomic_kind_set, &
437 particle_set=particle_set, &
439 blacs_env=blacs_env, &
444 dft_control=dft_control)
447 IF (para_env%is_source() .AND. logger%iter_info%print_level >=
medium_print_level)
THEN
453 CALL get_qs_env(qs_env, nkind=nkind, natom=natom)
454 use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
458 totalcharge = dft_control%charge
460 CALL get_cnumbers(qs_env, cnumbers, dcnum, .true.)
464 IF (
PRESENT(cn_max))
THEN
466 dcnpdcn = exp(cn_max)/(exp(cn_max) + exp(cnumbers(iatom)))
467 cnumbers(iatom) = log(1.0_dp + exp(cn_max)) - log(1.0_dp + exp(cn_max - cnumbers(iatom)))
468 DO i = 1, dcnum(iatom)%neighbors
469 dcnum(iatom)%dvals(i) = dcnum(iatom)%dvals(i)*dcnpdcn
475 ALLOCATE (gab(nkind, nkind), gam(nkind))
480 CALL get_eeq_data(za, eeq_model, eta=gam(ikind), rad=ala)
485 gab(ikind, jkind) = sqrt(1._dp/(ala*ala + alb*alb))
491 IF (
PRESENT(exclude))
THEN
493 IF (exclude(ikind))
THEN
494 gam(ikind) = 1.0e30_dp
495 gab(ikind, :) = 0.0_dp
496 gab(:, ikind) = 0.0_dp
501 ALLOCATE (qlag(natom))
503 CALL get_cell(cell, periodic=periodic)
504 do_ewald = .NOT. all(periodic == 0)
509 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
514 silent=.true., pset=
"EEQ")
516 CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section=print_section)
519 particle_set, kind_of, cell, -dcharges, gam, gab, &
520 para_env, blacs_env, dft_control, eeq_sparam, &
521 ewald=do_ewald, ewald_env=ewald_env, ewald_pw=ewald_pw, iounit=iunit)
524 particle_set, kind_of, cell, -dcharges, gam, gab, &
525 para_env, blacs_env, dft_control, eeq_sparam, iounit=iunit)
529 esg = 1.0_dp + exp(sgamma)
530 ALLOCATE (chrgx(natom), dchia(natom))
532 ikind = kind_of(iatom)
536 IF (response_only)
THEN
537 ctot = -0.5_dp*qlag(iatom)
539 ctot = 0.5_dp*(charges(iatom) - qlag(iatom))
541 IF (enshift_type == 1)
THEN
542 scn = sqrt(cnumbers(iatom)) + 1.0e-14_dp
543 dchia(iatom) = -ctot*kappa/scn
544 ELSE IF (enshift_type == 2)
THEN
546 scn = 1.0_dp/(esg - cn)
547 dchia(iatom) = -ctot*kappa*scn
549 cpabort(
"Unknown enshift_type")
554 IF (dft_control%apply_period_efield)
THEN
556 ELSE IF (dft_control%apply_efield)
THEN
558 ELSE IF (dft_control%apply_efield_field)
THEN
559 cpabort(
"apply field")
563 CALL get_qs_env(qs_env=qs_env, local_particles=local_particles)
565 DO ia = 1, local_particles%n_el(ikind)
566 iatom = local_particles%list(ikind)%array(ia)
567 DO i = 1, dcnum(iatom)%neighbors
568 katom = dcnum(iatom)%nlist(i)
569 rik = dcnum(iatom)%rik(:, i)
570 drk = sqrt(sum(rik(:)**2))
571 IF (drk > 1.e-3_dp)
THEN
572 fdik(:) = dchia(iatom)*dcnum(iatom)%dvals(i)*rik(:)/drk
573 gradient(:, iatom) = gradient(:, iatom) - fdik(:)
574 gradient(:, katom) = gradient(:, katom) + fdik(:)
588 distribution_2d=distribution_2d, &
589 local_particles=distribution_1d, &
590 molecule_set=molecule_set)
595 ALLOCATE (c_radius(nkind), default_present(nkind), pair_radius(nkind, nkind))
597 default_present = .true.
598 ALLOCATE (atom2d(nkind))
599 CALL atom2d_build(atom2d, distribution_1d, distribution_2d, atomic_kind_set, &
600 molecule_set, .false., particle_set=particle_set)
601 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
603 subcells=subcells, operator_type=
"PP", nlname=
"sab_ew")
604 DEALLOCATE (c_radius, pair_radius, default_present)
610 iatom=iatom, jatom=jatom, r=rij)
614 IF (dr > rcut .OR. dr < 1.e-6_dp) cycle
616 IF (iatom == jatom) fe = 0.5_dp
617 IF (response_only)
THEN
618 qq = -qlag(iatom)*charges(jatom)
620 qq = (0.5_dp*charges(iatom) - qlag(iatom))*charges(jatom)
622 gama = gab(ikind, jkind)
624 grc = 2._dp*gama*exp(-gam2*dr2)*
oorootpi/dr - erf(gama*dr)/dr2 &
625 - 2._dp*alpha*exp(-alpha**2*dr2)*
oorootpi/dr + erf(alpha*dr)/dr2
626 IF (response_only)
THEN
627 qq1 = -qlag(iatom)*charges(jatom)
628 qq2 = -qlag(jatom)*charges(iatom)
630 qq1 = (0.5_dp*charges(iatom) - qlag(iatom))*charges(jatom)
631 qq2 = (0.5_dp*charges(jatom) - qlag(jatom))*charges(iatom)
633 fdik(:) = -qq1*grc*rij(:)/dr
634 gradient(:, iatom) = gradient(:, iatom) + fdik(:)
635 gradient(:, jatom) = gradient(:, jatom) - fdik(:)
639 fdik(:) = qq2*grc*rij(:)/dr
640 gradient(:, iatom) = gradient(:, iatom) - fdik(:)
641 gradient(:, jatom) = gradient(:, jatom) + fdik(:)
651 DO ia = 1, local_particles%n_el(ikind)
652 iatom = local_particles%list(ikind)%array(ia)
653 ri(1:3) = particle_set(iatom)%r(1:3)
655 IF (iatom == jatom) cycle
656 jkind = kind_of(jatom)
657 IF (response_only)
THEN
658 qq = -qlag(iatom)*charges(jatom)
660 qq = (0.5_dp*charges(iatom) - qlag(iatom))*charges(jatom)
662 rj(1:3) = particle_set(jatom)%r(1:3)
663 rij(1:3) = ri(1:3) - rj(1:3)
667 gama = gab(ikind, jkind)
669 grc = 2._dp*gama*exp(-gam2*dr2)*
oorootpi/dr - erf(gama*dr)/dr2
670 fdik(:) = qq*grc*rij(:)/dr
671 gradient(:, iatom) = gradient(:, iatom) + fdik(:)
672 gradient(:, jatom) = gradient(:, jatom) - fdik(:)
680 ALLOCATE (epforce(3, natom))
682 IF (response_only)
THEN
683 dchia(1:natom) = qlag(1:natom)
685 dchia(1:natom) = -charges(1:natom) + qlag(1:natom)
687 chrgx(1:natom) = charges(1:natom)
688 CALL spme_forces(ewald_env, ewald_pw, cell, particle_set, chrgx, &
689 particle_set, dchia, epforce)
690 dchia(1:natom) = charges(1:natom)
691 chrgx(1:natom) = qlag(1:natom)
692 CALL spme_forces(ewald_env, ewald_pw, cell, particle_set, chrgx, &
693 particle_set, dchia, epforce)
694 gradient(1:3, 1:natom) = gradient(1:3, 1:natom) + epforce(1:3, 1:natom)
699 chrgx(1:natom) = charges(1:natom) - qlag(1:natom)
700 CALL spme_virial(ewald_env, ewald_pw, particle_set, cell, chrgx, pvir)
701 stress = stress - pvir
702 chrgx(1:natom) = qlag(1:natom)
703 CALL spme_virial(ewald_env, ewald_pw, particle_set, cell, chrgx, pvir)
704 stress = stress + pvir
705 IF (response_only)
THEN
706 chrgx(1:natom) = charges(1:natom)
707 CALL spme_virial(ewald_env, ewald_pw, particle_set, cell, chrgx, pvir)
708 stress = stress + pvir
714 DEALLOCATE (ewald_env, ewald_pw)
719 DEALLOCATE (gab, gam, qlag, chrgx, dchia)
721 CALL timestop(handle)
732 SUBROUTINE get_cnumbers(qs_env, cnumbers, dcnum, calculate_forces)
735 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: cnumbers
736 TYPE(
dcnum_type),
ALLOCATABLE,
DIMENSION(:) :: dcnum
737 LOGICAL,
INTENT(IN) :: calculate_forces
739 INTEGER :: ikind, natom, nkind, za
740 LOGICAL,
ALLOCATABLE,
DIMENSION(:) :: default_present
741 REAL(kind=
dp) :: subcells
742 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: c_radius
743 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: pair_radius
754 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
757 qs_kind_set=qs_kind_set, &
758 atomic_kind_set=atomic_kind_set, &
759 particle_set=particle_set, &
761 distribution_2d=distribution_2d, &
762 local_particles=distribution_1d, &
763 molecule_set=molecule_set)
764 CALL get_qs_env(qs_env, nkind=nkind, natom=natom)
769 disp%k2 = 4._dp/3._dp
770 disp%eps_cn = 1.e-6_dp
771 disp%max_elem = maxelem
772 ALLOCATE (disp%rcov(maxelem))
773 disp%rcov(1:maxelem) =
bohr*disp%k2*rcov(1:maxelem)
777 ALLOCATE (c_radius(nkind), default_present(nkind), pair_radius(nkind, nkind))
779 default_present = .true.
782 c_radius(ikind) = 4._dp*rcov(za)*
bohr
784 ALLOCATE (atom2d(nkind))
785 CALL atom2d_build(atom2d, distribution_1d, distribution_2d, atomic_kind_set, &
786 molecule_set, .false., particle_set=particle_set)
787 CALL pair_radius_setup(default_present, default_present, c_radius, c_radius, pair_radius)
789 subcells=subcells, operator_type=
"PP", nlname=
"sab_cn")
790 disp%sab_cn => sab_cn
791 DEALLOCATE (c_radius, pair_radius, default_present)
795 CALL cnumber_init(qs_env, cnumbers, dcnum, 2, calculate_forces, disp_env=disp)
799 END SUBROUTINE get_cnumbers
822 SUBROUTINE eeq_solver(charges, lambda, eeq_energy, particle_set, kind_of, cell, &
823 chia, gam, gab, para_env, blacs_env, dft_control, eeq_sparam, &
824 totalcharge, ewald, ewald_env, ewald_pw, iounit)
826 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: charges
827 REAL(kind=
dp),
INTENT(INOUT) :: lambda, eeq_energy
828 TYPE(
particle_type),
DIMENSION(:),
INTENT(IN) :: particle_set
829 INTEGER,
DIMENSION(:),
INTENT(IN) :: kind_of
831 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: chia, gam
832 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: gab
837 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: totalcharge
838 LOGICAL,
INTENT(IN),
OPTIONAL :: ewald
841 INTEGER,
INTENT(IN),
OPTIONAL :: iounit
843 CHARACTER(len=*),
PARAMETER :: routinen =
'eeq_solver'
845 INTEGER :: handle, ierror, iunit, natom, nkind, ns
846 LOGICAL :: do_direct, do_displ, do_ewald, do_sparse
847 REAL(kind=
dp) :: alpha, deth, ftime, qtot
851 CALL timeset(routinen, handle)
854 IF (
PRESENT(ewald)) do_ewald = ewald
857 IF (
PRESENT(totalcharge)) qtot = totalcharge
860 IF (
PRESENT(iounit)) iunit = iounit
863 do_direct = eeq_sparam%direct
864 do_sparse = eeq_sparam%sparse
867 IF (dft_control%apply_period_efield .AND.
ASSOCIATED(dft_control%period_efield))
THEN
868 do_displ = dft_control%period_efield%displacement_field
873 CALL cp_abort(__location__,
"EEQ: Sparse option not yet available")
876 natom =
SIZE(particle_set)
880 nrow_global=ns, ncol_global=ns)
885 cpassert(
PRESENT(ewald_env))
886 cpassert(
PRESENT(ewald_pw))
891 CALL fpbc_solver(charges, lambda, eeq_energy, eeq_mat, particle_set, &
892 kind_of, cell, chia, gam, gab, qtot, &
893 ewald_env, ewald_pw, iounit)
900 CALL pbc_solver(charges, lambda, eeq_energy, eeq_mat, particle_set, &
901 kind_of, cell, chia, gam, gab, qtot, &
902 ewald_env, ewald_pw, eeq_sparam, ierror, iounit)
903 IF (ierror /= 0)
THEN
905 CALL fpbc_solver(charges, lambda, eeq_energy, eeq_mat, particle_set, &
906 kind_of, cell, chia, gam, gab, qtot, &
907 ewald_env, ewald_pw, iounit)
911 IF (qtot /= 0._dp)
THEN
914 eeq_energy = eeq_energy - 0.5_dp*qtot**2/alpha**2/deth
917 CALL mi_solver(charges, lambda, eeq_energy, eeq_mat, particle_set, kind_of, &
918 cell, chia, gam, gab, qtot, ftime)
920 WRITE (iunit,
'(A,T67,F14.3)')
" EEQ| Molecular solver time[s]", ftime
926 CALL timestop(handle)
945 SUBROUTINE mi_solver(charges, lambda, eeq_energy, eeq_mat, particle_set, kind_of, cell, &
946 chia, gam, gab, qtot, ftime)
948 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: charges
949 REAL(kind=
dp),
INTENT(INOUT) :: lambda, eeq_energy
951 TYPE(
particle_type),
DIMENSION(:),
INTENT(IN) :: particle_set
952 INTEGER,
DIMENSION(:),
INTENT(IN) :: kind_of
954 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: chia, gam
955 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: gab
956 REAL(kind=
dp),
INTENT(IN) :: qtot
957 REAL(kind=
dp),
INTENT(OUT) :: ftime
959 CHARACTER(len=*),
PARAMETER :: routinen =
'mi_solver'
961 INTEGER :: handle, ia, iac, iar, ic, ikind, ir, &
962 jkind, natom, ncloc, ncvloc, nkind, &
964 INTEGER,
DIMENSION(:),
POINTER :: cind, cvind, rind, rvind
965 REAL(kind=
dp) :: dr, grc, te, ti, xr
966 REAL(kind=
dp),
DIMENSION(3) :: ri, rij, rj
971 CALL timeset(routinen, handle)
974 natom =
SIZE(particle_set)
978 CALL cp_fm_get_info(eeq_mat, matrix_struct=mat_struct, para_env=para_env)
979 CALL cp_fm_get_info(eeq_mat, nrow_local=nrloc, ncol_local=ncloc, &
980 row_indices=rind, col_indices=cind)
982 nrow_global=ns, ncol_global=1)
984 CALL cp_fm_get_info(rhs_vec, nrow_local=nrvloc, ncol_local=ncvloc, &
985 row_indices=rvind, col_indices=cvind)
992 IF (iar > natom) cycle
994 ri(1:3) = particle_set(iar)%r(1:3)
997 IF (iac > natom) cycle
999 rj(1:3) = particle_set(iac)%r(1:3)
1000 IF (iar == iac)
THEN
1001 grc = gam(ikind) + 2.0_dp*gab(ikind, ikind)*
oorootpi
1003 rij(1:3) = ri(1:3) - rj(1:3)
1004 rij =
pbc(rij, cell)
1005 dr = sqrt(sum(rij**2))
1006 grc = erf(gab(ikind, jkind)*dr)/dr
1008 eeq_mat%local_data(ir, ic) = grc
1017 IF (ia > natom)
THEN
1022 rhs_vec%local_data(ir, ic) = xr
1035 IF (ia <= natom)
THEN
1036 xr = rhs_vec%local_data(ir, ic)
1039 lambda = rhs_vec%local_data(ir, ic)
1043 CALL para_env%sum(lambda)
1044 CALL para_env%sum(charges)
1047 eeq_energy = 0.5*sum(charges(1:natom)*chia(1:natom)) - 0.5_dp*lambda*qtot
1054 CALL timestop(handle)
1056 END SUBROUTINE mi_solver
1077 SUBROUTINE pbc_solver(charges, lambda, eeq_energy, eeq_mat, particle_set, &
1078 kind_of, cell, chia, gam, gab, qtot, &
1079 ewald_env, ewald_pw, eeq_sparam, ierror, iounit)
1081 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: charges
1082 REAL(kind=
dp),
INTENT(INOUT) :: lambda, eeq_energy
1084 TYPE(
particle_type),
DIMENSION(:),
INTENT(IN) :: particle_set
1085 INTEGER,
DIMENSION(:),
INTENT(IN) :: kind_of
1087 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: chia, gam
1088 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: gab
1089 REAL(kind=
dp),
INTENT(IN) :: qtot
1093 INTEGER,
INTENT(OUT) :: ierror
1094 INTEGER,
OPTIONAL :: iounit
1096 CHARACTER(len=*),
PARAMETER :: routinen =
'pbc_solver'
1098 INTEGER :: ewald_type, handle, i, iac, iar, ic, ikind, info, ir, iunit, iv, ix, iy, iz, &
1099 jkind, max_diis, mdiis, natom, ncloc, ndiis, nkind, now, nrloc, ns, sdiis
1100 INTEGER,
DIMENSION(3) :: cvec, ncell, periodic
1101 INTEGER,
DIMENSION(:),
POINTER :: cind, rind
1102 REAL(kind=
dp) :: ad, alpha, astep, deth, dr, eeqn, &
1103 eps_diis, ftime, grc1, grc2, rcut, &
1104 res, resin, rmax, te, ti
1105 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: bvec, dvec
1106 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: dmat, fvec, vmat, xvec
1107 REAL(kind=
dp),
DIMENSION(3) :: ri, rij, rijl, rj
1108 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
1109 REAL(kind=
dp),
DIMENSION(:),
POINTER :: rhs, rv0, xv0
1114 CALL timeset(routinen, handle)
1120 IF (
PRESENT(iounit)) iunit = iounit
1122 natom =
SIZE(particle_set)
1125 CALL get_cell(cell=cell, deth=deth)
1126 CALL ewald_env_get(ewald_env, alpha=alpha, rcut=rcut, ewald_type=ewald_type)
1129 CALL cp_abort(__location__,
"Only SPME Ewald method available with EEQ.")
1134 CALL get_cell(cell, h=hmat, periodic=periodic)
1138 IF (periodic(1) == 0) ncell(1) = 0
1139 IF (periodic(2) == 0) ncell(2) = 0
1140 IF (periodic(3) == 0) ncell(3) = 0
1142 CALL mi_solver(charges, lambda, eeqn, eeq_mat, particle_set, kind_of, cell, &
1143 chia, gam, gab, qtot, ftime)
1145 WRITE (iunit,
'(A,T67,F14.3)')
" EEQ| Iterative PBC guess time[s]", ftime
1147 CALL cp_fm_get_info(eeq_mat, matrix_struct=mat_struct, para_env=para_env)
1152 CALL cp_fm_get_info(eeq_mat, nrow_local=nrloc, ncol_local=ncloc, &
1153 row_indices=rind, col_indices=cind)
1158 IF (iar <= natom)
THEN
1159 ikind = kind_of(iar)
1160 ri(1:3) = particle_set(iar)%r(1:3)
1164 IF (iac > natom .AND. iar > natom)
THEN
1165 eeq_mat%local_data(ir, ic) = 0.0_dp
1167 ELSE IF ((iac > natom) .OR. (iar > natom))
THEN
1168 eeq_mat%local_data(ir, ic) = 1.0_dp
1171 jkind = kind_of(iac)
1172 rj(1:3) = particle_set(iac)%r(1:3)
1173 rij(1:3) = ri(1:3) - rj(1:3)
1174 rij =
pbc(rij, cell)
1175 DO ix = -ncell(1), ncell(1)
1176 DO iy = -ncell(2), ncell(2)
1177 DO iz = -ncell(3), ncell(3)
1179 rijl = rij + matmul(hmat, cvec)
1180 dr = sqrt(sum(rijl**2))
1181 IF (dr > rmax) cycle
1182 IF (iar == iac .AND. dr < 0.00001_dp)
THEN
1183 grc1 = gam(ikind) + 2.0_dp*gab(ikind, ikind)*
oorootpi - ad
1185 grc1 = erf(gab(ikind, jkind)*dr)/dr - erf(alpha*dr)/dr
1187 eeq_mat%local_data(ir, ic) = eeq_mat%local_data(ir, ic) + grc1
1196 row_indices=rind, col_indices=cind)
1201 IF (iar <= natom)
THEN
1202 ikind = kind_of(iar)
1203 ri(1:3) = particle_set(iar)%r(1:3)
1207 IF (iac > natom .AND. iar > natom)
THEN
1208 pmat%local_data(ir, ic) = 0.0_dp
1210 ELSE IF ((iac > natom) .OR. (iar > natom))
THEN
1211 pmat%local_data(ir, ic) = 1.0_dp
1214 jkind = kind_of(iac)
1215 rj(1:3) = particle_set(iac)%r(1:3)
1216 rij(1:3) = ri(1:3) - rj(1:3)
1217 rij =
pbc(rij, cell)
1218 IF (iar == iac)
THEN
1219 grc2 = gam(ikind) + 2.0_dp*gab(ikind, ikind)*
oorootpi
1221 grc2 = erf(gab(ikind, jkind)*dr)/dr
1223 pmat%local_data(ir, ic) = grc2
1233 rhs(1:natom) = chia(1:natom)
1236 ALLOCATE (xv0(ns), rv0(ns))
1238 xv0(1:natom) = charges(1:natom)
1241 max_diis = eeq_sparam%max_diis
1242 mdiis = eeq_sparam%mdiis
1243 sdiis = eeq_sparam%sdiis
1244 eps_diis = eeq_sparam%eps_diis
1245 astep = eeq_sparam%alpha
1246 ALLOCATE (xvec(ns, mdiis), fvec(ns, mdiis), bvec(ns))
1247 xvec = 0.0_dp; fvec = 0.0_dp
1248 ALLOCATE (vmat(mdiis, mdiis), dmat(mdiis + 1, mdiis + 1), dvec(mdiis + 1))
1249 dmat = 0.0_dp; dvec = 0.0_dp
1252 CALL get_energy_gradient(eeqn, eeq_mat, mmat, ewald_env, ewald_pw, &
1253 cell, particle_set, xv0, rhs, rv0)
1254 resin = sqrt(sum(rv0*rv0))
1257 res = sqrt(sum(rv0*rv0))
1258 IF (res > 10._dp*resin)
EXIT
1259 IF (res < eps_diis)
EXIT
1261 now = mod(iv - 1, mdiis) + 1
1262 ndiis = min(iv, mdiis)
1263 xvec(1:ns, now) = xv0(1:ns)
1264 fvec(1:ns, now) = rv0(1:ns)
1266 vmat(now, i) = sum(fvec(:, now)*fvec(:, i))
1267 vmat(i, now) = vmat(now, i)
1269 IF (ndiis < sdiis)
THEN
1270 xv0(1:ns) = xv0(1:ns) + astep*rv0(1:ns)
1273 dvec(ndiis + 1) = 1.0_dp
1274 dmat(1:ndiis, 1:ndiis) = vmat(1:ndiis, 1:ndiis)
1275 dmat(ndiis + 1, 1:ndiis) = 1.0_dp
1276 dmat(1:ndiis, ndiis + 1) = 1.0_dp
1277 dmat(ndiis + 1, ndiis + 1) = 0.0_dp
1278 CALL invmat(dmat(1:ndiis + 1, 1:ndiis + 1), info)
1279 dvec(1:ndiis + 1) = matmul(dmat(1:ndiis + 1, 1:ndiis + 1), dvec(1:ndiis + 1))
1280 xv0(1:ns) = matmul(xvec(1:ns, 1:ndiis), dvec(1:ndiis))
1281 xv0(1:ns) = xv0(1:ns) + matmul(fvec(1:ns, 1:ndiis), dvec(1:ndiis))
1284 CALL get_energy_gradient(eeqn, eeq_mat, mmat, ewald_env, ewald_pw, &
1285 cell, particle_set, xv0, rhs, rv0)
1287 charges(1:natom) = xv0(1:natom)
1290 IF (res > eps_diis) ierror = 1
1292 DEALLOCATE (xvec, fvec, bvec)
1293 DEALLOCATE (vmat, dmat, dvec)
1294 DEALLOCATE (xv0, rv0)
1301 IF (ierror == 1)
THEN
1302 WRITE (iunit,
'(A)')
" EEQ| PBC solver failed to converge "
1304 WRITE (iunit,
'(A,T50,I4,T61,E20.5)')
" EEQ| PBC solver: iterations/accuracy ", iv, res
1306 WRITE (iunit,
'(A,T67,F14.3)')
" EEQ| Iterative PBC solver: time[s]", te - ti
1308 CALL timestop(handle)
1310 END SUBROUTINE pbc_solver
1329 SUBROUTINE fpbc_solver(charges, lambda, eeq_energy, eeq_mat, particle_set, &
1330 kind_of, cell, chia, gam, gab, qtot, ewald_env, ewald_pw, iounit)
1332 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: charges
1333 REAL(kind=
dp),
INTENT(INOUT) :: lambda, eeq_energy
1335 TYPE(
particle_type),
DIMENSION(:),
INTENT(IN) :: particle_set
1336 INTEGER,
DIMENSION(:),
INTENT(IN) :: kind_of
1338 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: chia, gam
1339 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: gab
1340 REAL(kind=
dp),
INTENT(IN) :: qtot
1343 INTEGER,
INTENT(IN),
OPTIONAL :: iounit
1345 CHARACTER(len=*),
PARAMETER :: routinen =
'fpbc_solver'
1347 INTEGER :: ewald_type, handle, ia, iac, iar, ic, &
1348 ikind, ir, iunit, ix, iy, iz, jkind, &
1349 natom, ncloc, ncvloc, nkind, nrloc, &
1351 INTEGER,
DIMENSION(3) :: cvec, ncell, periodic
1352 INTEGER,
DIMENSION(:),
POINTER :: cind, cvind, rind, rvind
1353 REAL(kind=
dp) :: ad, alpha, deth, dr, grc1, rcut, rmax, &
1355 REAL(kind=
dp),
DIMENSION(3) :: ri, rij, rijl, rj
1356 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
1357 REAL(kind=
dp),
DIMENSION(:),
POINTER :: pval, xval
1362 CALL timeset(routinen, handle)
1366 IF (
PRESENT(iounit)) iunit = iounit
1368 natom =
SIZE(particle_set)
1372 CALL get_cell(cell=cell, deth=deth)
1373 CALL ewald_env_get(ewald_env, alpha=alpha, rcut=rcut, ewald_type=ewald_type)
1376 CALL cp_abort(__location__,
"Only SPME Ewald method available with EEQ.")
1381 CALL get_cell(cell, h=hmat, periodic=periodic)
1385 IF (periodic(1) == 0) ncell(1) = 0
1386 IF (periodic(2) == 0) ncell(2) = 0
1387 IF (periodic(3) == 0) ncell(3) = 0
1389 CALL cp_fm_get_info(eeq_mat, matrix_struct=mat_struct, para_env=para_env)
1391 CALL cp_fm_get_info(eeq_mat, nrow_local=nrloc, ncol_local=ncloc, &
1392 row_indices=rind, col_indices=cind)
1394 nrow_global=ns, ncol_global=1)
1396 CALL cp_fm_get_info(rhs_vec, nrow_local=nrvloc, ncol_local=ncvloc, &
1397 row_indices=rvind, col_indices=cvind)
1402 IF (iar <= natom)
THEN
1403 ikind = kind_of(iar)
1404 ri(1:3) = particle_set(iar)%r(1:3)
1408 IF (iac > natom .AND. iar > natom)
THEN
1409 eeq_mat%local_data(ir, ic) = 0.0_dp
1411 ELSE IF ((iac > natom) .OR. (iar > natom))
THEN
1412 eeq_mat%local_data(ir, ic) = 1.0_dp
1415 jkind = kind_of(iac)
1416 rj(1:3) = particle_set(iac)%r(1:3)
1417 rij(1:3) = ri(1:3) - rj(1:3)
1418 rij =
pbc(rij, cell)
1419 DO ix = -ncell(1), ncell(1)
1420 DO iy = -ncell(2), ncell(2)
1421 DO iz = -ncell(3), ncell(3)
1423 rijl = rij + matmul(hmat, cvec)
1424 dr = sqrt(sum(rijl**2))
1425 IF (dr > rmax) cycle
1426 IF (iar == iac .AND. dr < 0.0001_dp)
THEN
1427 grc1 = gam(ikind) + 2.0_dp*gab(ikind, ikind)*
oorootpi - ad
1429 grc1 = erf(gab(ikind, jkind)*dr)/dr - erf(alpha*dr)/dr
1431 eeq_mat%local_data(ir, ic) = eeq_mat%local_data(ir, ic) + grc1
1438 ALLOCATE (xval(natom), pval(natom))
1442 CALL apply_potential(ewald_env, ewald_pw, cell, particle_set, xval, pval)
1446 IF (iar /= ia) cycle
1449 IF (iac > natom) cycle
1450 eeq_mat%local_data(ir, ic) = eeq_mat%local_data(ir, ic) + pval(iac)
1454 DEALLOCATE (xval, pval)
1462 IF (ia > natom)
THEN
1467 rhs_vec%local_data(ir, ic) = xr
1480 IF (ia <= natom)
THEN
1481 xr = rhs_vec%local_data(ir, ic)
1484 lambda = rhs_vec%local_data(ir, ic)
1488 CALL para_env%sum(lambda)
1489 CALL para_env%sum(charges)
1492 eeq_energy = 0.5*sum(charges(1:natom)*chia(1:natom)) - 0.5_dp*lambda*qtot
1499 WRITE (iunit,
'(A,T67,F14.3)')
" EEQ| Direct PBC solver: time[s]", te - ti
1501 CALL timestop(handle)
1503 END SUBROUTINE fpbc_solver
1514 SUBROUTINE apply_potential(ewald_env, ewald_pw, cell, particle_set, charges, potential)
1518 TYPE(
particle_type),
DIMENSION(:),
INTENT(IN) :: particle_set
1519 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN),
POINTER :: charges
1520 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: potential
1526 CALL spme_potential(ewald_env, ewald_pw, cell, particle_set, charges, &
1527 particle_set, potential)
1528 CALL para_env%sum(potential)
1530 END SUBROUTINE apply_potential
1545 SUBROUTINE get_energy_gradient(eeqn, fm_mat, mmat, ewald_env, ewald_pw, &
1546 cell, particle_set, charges, rhs, potential)
1547 REAL(kind=
dp),
INTENT(INOUT) :: eeqn
1552 TYPE(
particle_type),
DIMENSION(:),
INTENT(IN) :: particle_set
1553 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN),
POINTER :: charges
1554 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: rhs
1555 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: potential
1558 REAL(kind=
dp) :: lambda
1559 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: mvec
1566 CALL spme_potential(ewald_env, ewald_pw, cell, particle_set, charges(1:na), &
1567 particle_set, potential(1:na))
1568 CALL para_env%sum(potential(1:na))
1569 CALL cp_fm_matvec(fm_mat, charges, potential, alpha=1.0_dp, beta=1.0_dp)
1570 eeqn = 0.5_dp*sum(charges(1:na)*potential(1:na)) + sum(charges(1:na)*rhs(1:na))
1571 potential(1:ns) = potential(1:ns) + rhs(1:ns)
1573 CALL cp_fm_matvec(mmat, potential, mvec, alpha=-1.0_dp, beta=0.0_dp)
1574 lambda = -sum(mvec(1:na))/real(na, kind=
dp)
1575 potential(1:na) = mvec(1:na) + lambda
1578 END SUBROUTINE get_energy_gradient
1588 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: charges
1589 REAL(kind=
dp),
INTENT(OUT) :: ef_energy
1591 COMPLEX(KIND=dp) :: zdeta
1592 COMPLEX(KIND=dp),
DIMENSION(3) :: zi
1593 INTEGER :: ia, idir, natom
1595 REAL(kind=
dp) :: kr, omega, q
1596 REAL(kind=
dp),
DIMENSION(3) :: ci, dfilter, fieldpol, fpolvec, kvec, &
1598 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
1603 CALL get_qs_env(qs_env, natom=natom, dft_control=dft_control)
1604 CALL get_qs_env(qs_env, cell=cell, particle_set=particle_set)
1606 IF (dft_control%apply_period_efield)
THEN
1607 dfield = dft_control%period_efield%displacement_field
1609 IF (
ALLOCATED(dft_control%period_efield%strength_list))
THEN
1610 cpabort(
"use of strength_list not implemented for eeq_efield_energy")
1613 fieldpol = dft_control%period_efield%polarisation
1614 fieldpol = fieldpol/sqrt(dot_product(fieldpol, fieldpol))
1615 fieldpol = -fieldpol*dft_control%period_efield%strength
1616 hmat = cell%hmat(:, :)/
twopi
1618 fpolvec(idir) = fieldpol(1)*hmat(1, idir) + fieldpol(2)*hmat(2, idir) &
1619 + fieldpol(3)*hmat(3, idir)
1622 zi(:) = cmplx(1._dp, 0._dp,
dp)
1625 ria = particle_set(ia)%r
1626 ria =
pbc(ria, cell)
1628 kvec(:) =
twopi*cell%h_inv(idir, :)
1629 kr = sum(kvec(:)*ria(:))
1630 zdeta = cmplx(cos(kr), sin(kr), kind=
dp)**q
1631 zi(idir) = zi(idir)*zdeta
1636 dfilter(1:3) = dft_control%period_efield%d_filter(1:3)
1638 ci = matmul(hmat, qi)/omega
1641 ef_energy = ef_energy + dfilter(idir)*(fieldpol(idir) - 2._dp*
twopi*ci(idir))**2
1643 ef_energy = -0.25_dp*omega/
twopi*ef_energy
1645 ef_energy = sum(fpolvec(:)*qi(:))
1648 ELSE IF (dft_control%apply_efield)
THEN
1650 fieldpol = dft_control%efield_fields(1)%efield%polarisation* &
1651 dft_control%efield_fields(1)%efield%strength
1655 ria = particle_set(ia)%r
1656 ria =
pbc(ria, cell)
1658 ef_energy = ef_energy - q*sum(fieldpol*ria)
1662 cpabort(
"apply field")
1674 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: efr
1676 INTEGER :: ia, idir, natom
1679 REAL(kind=
dp),
DIMENSION(3) :: fieldpol, fpolvec, kvec, ria
1680 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
1685 CALL get_qs_env(qs_env, natom=natom, dft_control=dft_control)
1686 CALL get_qs_env(qs_env, cell=cell, particle_set=particle_set)
1688 IF (dft_control%apply_period_efield)
THEN
1689 dfield = dft_control%period_efield%displacement_field
1691 IF (
ALLOCATED(dft_control%period_efield%strength_list))
THEN
1692 cpabort(
"use of strength_list not implemented for eeq_efield_pot")
1695 fieldpol = dft_control%period_efield%polarisation
1696 fieldpol = fieldpol/sqrt(dot_product(fieldpol, fieldpol))
1697 fieldpol = -fieldpol*dft_control%period_efield%strength
1698 hmat = cell%hmat(:, :)/
twopi
1700 fpolvec(idir) = fieldpol(1)*hmat(1, idir) + fieldpol(2)*hmat(2, idir) &
1701 + fieldpol(3)*hmat(3, idir)
1710 ria = particle_set(ia)%r
1711 ria =
pbc(ria, cell)
1713 kvec(:) =
twopi*cell%h_inv(idir, :)
1714 kr = sum(kvec(:)*ria(:))
1715 efr(ia) = efr(ia) + kr*fpolvec(idir)
1720 ELSE IF (dft_control%apply_efield)
THEN
1722 fieldpol = dft_control%efield_fields(1)%efield%polarisation* &
1723 dft_control%efield_fields(1)%efield%strength
1726 ria = particle_set(ia)%r
1727 ria =
pbc(ria, cell)
1728 efr(ia) = -sum(fieldpol*ria)
1732 cpabort(
"apply field")
1745 SUBROUTINE eeq_dfield_pot(charges, dft_control, particle_set, cell, efr)
1746 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: charges
1748 TYPE(
particle_type),
DIMENSION(:),
INTENT(IN) :: particle_set
1750 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: efr
1752 COMPLEX(KIND=dp) :: zdeta
1753 COMPLEX(KIND=dp),
DIMENSION(3) :: zi
1754 INTEGER :: ia, idir, natom
1755 REAL(kind=
dp) :: kr, omega, q
1756 REAL(kind=
dp),
DIMENSION(3) :: ci, dfilter, fieldpol, kvec, qi, ria
1757 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
1759 natom =
SIZE(particle_set)
1761 IF (
ALLOCATED(dft_control%period_efield%strength_list))
THEN
1762 cpabort(
"use of strength_list not implemented for eeq_dfield_pot")
1765 dfilter(1:3) = dft_control%period_efield%d_filter(1:3)
1766 fieldpol = dft_control%period_efield%polarisation
1767 fieldpol = fieldpol/sqrt(dot_product(fieldpol, fieldpol))
1768 fieldpol = -fieldpol*dft_control%period_efield%strength
1769 hmat = cell%hmat(:, :)/
twopi
1772 zi(:) = cmplx(1._dp, 0._dp,
dp)
1775 ria = particle_set(ia)%r
1776 ria =
pbc(ria, cell)
1778 kvec(:) =
twopi*cell%h_inv(idir, :)
1779 kr = sum(kvec(:)*ria(:))
1780 zdeta = cmplx(cos(kr), sin(kr), kind=
dp)**q
1781 zi(idir) = zi(idir)*zdeta
1785 ci = matmul(hmat, qi)/omega
1786 ci = dfilter*(fieldpol - 2._dp*
twopi*ci)
1788 ria = particle_set(ia)%r
1789 ria =
pbc(ria, cell)
1790 efr(ia) = efr(ia) - sum(ci*ria)
1793 END SUBROUTINE eeq_dfield_pot
1803 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: charges, qlag
1805 INTEGER :: atom_a, ia, iatom, ikind, natom, nkind
1806 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: atom_of_kind
1808 REAL(kind=
dp),
DIMENSION(3) :: fieldpol
1818 dft_control=dft_control, &
1819 cell=cell, particle_set=particle_set, &
1820 nkind=nkind, natom=natom, &
1821 para_env=para_env, &
1822 local_particles=local_particles)
1824 fieldpol = dft_control%efield_fields(1)%efield%polarisation* &
1825 dft_control%efield_fields(1)%efield%strength
1827 CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set)
1832 force(ikind)%efield = 0.0_dp
1833 DO ia = 1, local_particles%n_el(ikind)
1834 iatom = local_particles%list(ikind)%array(ia)
1835 q = charges(iatom) - qlag(iatom)
1836 atom_a = atom_of_kind(iatom)
1837 force(ikind)%efield(1:3, atom_a) = -q*fieldpol(1:3)
1839 CALL para_env%sum(force(ikind)%efield)
1852 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: charges, qlag
1854 INTEGER :: atom_a, ia, iatom, ikind, natom, nkind
1855 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: atom_of_kind
1856 LOGICAL :: dfield, use_virial
1858 REAL(kind=
dp),
DIMENSION(3) :: fa, fieldpol, ria
1859 REAL(kind=
dp),
DIMENSION(3, 3) :: pve
1870 dft_control=dft_control, &
1871 cell=cell, particle_set=particle_set, &
1873 nkind=nkind, natom=natom, &
1874 para_env=para_env, &
1875 local_particles=local_particles)
1877 dfield = dft_control%period_efield%displacement_field
1878 cpassert(.NOT. dfield)
1880 IF (
ALLOCATED(dft_control%period_efield%strength_list))
THEN
1881 cpabort(
"use of strength_list not implemented for eeq_efield_force_periodic")
1884 fieldpol = dft_control%period_efield%polarisation
1885 fieldpol = fieldpol/sqrt(dot_product(fieldpol, fieldpol))
1886 fieldpol = -fieldpol*dft_control%period_efield%strength
1888 use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
1890 CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set)
1896 force(ikind)%efield = 0.0_dp
1897 DO ia = 1, local_particles%n_el(ikind)
1898 iatom = local_particles%list(ikind)%array(ia)
1899 q = charges(iatom) - qlag(iatom)
1900 fa(1:3) = q*fieldpol(1:3)
1901 atom_a = atom_of_kind(iatom)
1902 force(ikind)%efield(1:3, atom_a) = fa
1903 IF (use_virial)
THEN
1904 ria = particle_set(ia)%r
1905 ria =
pbc(ria, cell)
1910 CALL para_env%sum(force(ikind)%efield)
1912 virial%pv_virial = virial%pv_virial + pve
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind_set(atomic_kind_set, atom_of_kind, kind_of, natom_of_kind, maxatom, natom, nshell, fist_potential_present, shell_present, shell_adiabatic, shell_check_distance, damping_present)
Get attributes of an atomic kind set.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
Holds information on atomic properties.
Handles all functions related to the CELL.
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
real(kind=dp) function, public plane_distance(h, k, l, cell)
Calculate the distance between two lattice planes as defined by a triple of Miller indices (hkl).
methods related to the blacs parallel environment
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_solve(matrix_a, general_a)
computes the the solution to A*b=A_general using lu decomposition
subroutine, public cp_fm_invert(matrix_a, matrix_inverse, det_a, eps_svd, eigval)
Inverts a cp_fm_type matrix, optionally returning the determinant of the input matrix.
subroutine, public cp_fm_matvec(amat, xv, yv, alpha, beta)
Calculates yv = alpha*amat*xv + beta*yv where amat: fm matrix xv : vector replicated yv : vector repl...
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_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
various routines to log and control the output. The idea is that decisions about where to log should ...
recursive integer function, public cp_logger_get_default_unit_nr(logger, local, skip_not_ionode)
asks the default unit number of the given logger. try to use cp_logger_get_unit_nr
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 medium_print_level
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
stores a mapping of 2D info (e.g. matrix) on a 2D processor distribution (i.e. blacs grid) where cpus...
EEQ data from different sources.
subroutine, public get_eeq_data(za, model, chi, eta, kcn, rad)
...
Calculation of charge equilibration method.
subroutine, public eeq_efield_energy(qs_env, charges, ef_energy)
...
subroutine, public eeq_charges(qs_env, charges, eeq_sparam, eeq_model, enshift_type, exclude, cn_max)
...
subroutine, public eeq_efield_force_periodic(qs_env, charges, qlag)
...
subroutine, public eeq_efield_pot(qs_env, efr)
...
subroutine, public eeq_print(qs_env, iounit, print_level, ext)
...
subroutine, public eeq_forces(qs_env, charges, dcharges, gradient, stress, eeq_sparam, eeq_model, enshift_type, response_only, exclude, cn_max)
...
subroutine, public eeq_efield_force_loc(qs_env, charges, qlag)
...
subroutine, public eeq_solver(charges, lambda, eeq_energy, particle_set, kind_of, cell, chia, gam, gab, para_env, blacs_env, dft_control, eeq_sparam, totalcharge, ewald, ewald_env, ewald_pw, iounit)
...
subroutine, public ewald_env_set(ewald_env, ewald_type, alpha, epsilon, eps_pol, gmax, ns_max, precs, o_spline, para_env, poisson_section, interaction_cutoffs, cell_hmat)
Purpose: Set the EWALD environment.
subroutine, public ewald_env_create(ewald_env, para_env)
allocates and intitializes a ewald_env
subroutine, public ewald_env_release(ewald_env)
releases the given ewald_env (see doc/ReferenceCounting.html)
subroutine, public read_ewald_section_tb(ewald_env, ewald_section, hmat, silent, pset)
Purpose: read the EWALD section for TB methods.
subroutine, public ewald_env_get(ewald_env, ewald_type, alpha, eps_pol, epsilon, gmax, ns_max, o_spline, group, para_env, poisson_section, precs, rcut, do_multipoles, max_multipole, do_ipol, max_ipol_iter, interaction_cutoffs, cell_hmat)
Purpose: Get the EWALD environment.
subroutine, public ewald_pw_release(ewald_pw)
releases the memory used by the ewald_pw
subroutine, public ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section)
creates the structure ewald_pw_type
Defines the basic variable types.
integer, parameter, public dp
Machine interface based on Fortran 2003 and POSIX.
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 oorootpi
real(kind=dp), parameter, public twopi
Collection of simple mathematical functions and subroutines.
subroutine, public invmat(a, info)
returns inverse of matrix using the lapack routines DGETRF and DGETRI
Interface to the message passing library MPI.
Define the data structure for the molecule information.
Define the data structure for the particle information.
Definition of physical constants:
real(kind=dp), parameter, public bohr
functions related to the poisson solver on regular grids
integer, parameter, public do_ewald_spme
Coordination number routines for dispersion pairpotentials.
subroutine, public cnumber_release(cnumbers, dcnum, derivatives)
...
subroutine, public cnumber_init(qs_env, cnumbers, dcnum, ftype, derivatives, disp_env)
...
Definition of disperson types for DFT calculations.
subroutine, public qs_dispersion_release(dispersion_env)
...
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, mimic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, cneo_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, monovalent, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
Define the neighbor list data types and the corresponding functionality.
subroutine, public release_neighbor_list_sets(nlists)
releases an array of neighbor_list_sets
subroutine, public neighbor_list_iterator_create(iterator_set, nl, search, nthread)
Neighbor list iterator functions.
subroutine, public neighbor_list_iterator_release(iterator_set)
...
integer function, public neighbor_list_iterate(iterator_set, mepos)
...
subroutine, public get_iterator_info(iterator_set, mepos, ikind, jkind, nkind, ilist, nlist, inode, nnode, iatom, jatom, r, cell)
...
Generate the atomic neighbor lists.
subroutine, public atom2d_cleanup(atom2d)
free the internals of atom2d
subroutine, public pair_radius_setup(present_a, present_b, radius_a, radius_b, pair_radius, prmin)
...
subroutine, public build_neighbor_lists(ab_list, particle_set, atom, cell, pair_radius, subcells, mic, symmetric, molecular, subset_of_mol, current_subset, operator_type, nlname, atomb_to_keep)
Build simple pair neighbor lists.
subroutine, public atom2d_build(atom2d, distribution_1d, distribution_2d, atomic_kind_set, molecule_set, molecule_only, particle_set)
Build some distribution structure of atoms, refactored from build_qs_neighbor_lists.
Calculate the electrostatic energy by the Smooth Particle Ewald method.
subroutine, public spme_forces(ewald_env, ewald_pw, box, particle_set_a, charges_a, particle_set_b, charges_b, forces_b)
Calculate the forces on particles B for the electrostatic interaction betrween particles A and B.
subroutine, public spme_potential(ewald_env, ewald_pw, box, particle_set_a, charges_a, particle_set_b, potential)
Calculate the electrostatic potential from particles A (charge A) at positions of particles B.
subroutine, public spme_virial(ewald_env, ewald_pw, particle_set, box, mcharge, virial)
Internal Virial for 1/2 [rho||rho] (rho=mcharge)
pure subroutine, public virial_pair_force(pv_virial, f0, force, rab)
Computes the contribution to the stress tensor from two-body pair-wise forces.
Provides all information about an atomic kind.
type for the atomic properties
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...
structure to store local (to a processor) ordered lists of integers.
distributes pairs on a 2d grid of processors
to build arrays of pointers
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.