58 type(distributed_t) :: dist
61 integer,
parameter :: &
62 ION_COMPONENT_REAL = 1, &
70 type(ion_interaction_t),
intent(out) :: this
71 type(namespace_t),
intent(in) :: namespace
72 class(space_t),
intent(in) :: space
73 integer,
intent(in) :: natoms
87 call parse_variable(namespace,
'EwaldAlpha', 0.21_real64, this%alpha)
91 if (space%periodic_dim == 1)
then
92 call messages_write(
'For systems that are periodic in 1D, the interaction between', new_line = .
true.)
93 call messages_write(
'ions is not implemented. This affects the calculation', new_line = .
true.)
94 call messages_write(
'of total energy and forces, so both are zeroed.')
102 type(ion_interaction_t),
intent(inout) :: this
103 integer,
intent(in) :: natoms
104 type(multicomm_t),
intent(in) :: mc
116 type(ion_interaction_t),
intent(inout) :: this
132 energy_components, force_components)
133 type(ion_interaction_t),
intent(inout) :: this
134 class(space_t),
intent(in) :: space
135 type(lattice_vectors_t),
intent(in) :: latt
136 type(atom_t),
intent(in) :: atom(:)
137 integer,
intent(in) :: natoms
138 real(real64),
intent(in) :: pos(1:space%dim,1:natoms)
139 real(real64),
intent(in) :: lsize(:)
140 real(real64),
intent(out) :: energy
141 real(real64),
intent(out) :: force(:, :)
142 real(real64),
optional,
intent(out) :: energy_components(:)
143 real(real64),
optional,
intent(out) :: force_components(:, :, :)
149 if (
present(energy_components))
then
154 if (
present(force_components))
then
163 if (space%is_periodic())
then
168 call ion_interaction_periodic(this, space, latt, atom, natoms, pos, energy, force, energy_components, force_components)
185 class(
space_t),
intent(in) :: space
187 type(
atom_t),
intent(in) :: atom(:)
188 real(real64),
intent(in) :: lsize(:)
189 real(real64) :: energy
192 assert(
size(atom) == 1)
194 assert(space%periodic_dim == 2)
196 select type(spec => atom(1)%species)
198 energy =
m_pi * spec%get_zval()**2 /latt%rcell_volume * (lsize(3) - spec%thickness() /
m_three)
221 type(
atom_t),
intent(in) :: atom(:)
222 real(real64),
intent(in) :: lsize(:)
223 real(real64) :: energy
227 logical :: lattice_is_orthogonal
233 lattice_is_orthogonal = .not. latt%nonorthogonal
235 do iatom = dist%start, dist%end
236 spec => atom(iatom)%species
247 assert(lattice_is_orthogonal)
248 energy = energy -
m_pi * zi**2 / (
m_four * lsize(1)*lsize(2)) * spec%thickness() /
m_three
263 class(
space_t),
intent(in) :: space
264 type(
atom_t),
intent(in) :: atom(:)
265 real(real64),
intent(in) :: pos(:,:)
266 real(real64),
intent(in) :: lsize(:)
267 real(real64),
intent(out) :: energy
268 real(real64),
intent(out) :: force(:, :)
270 class(
species_t),
pointer :: species_i, species_j
271 real(real64) :: r(space%dim), f(space%dim)
272 real(real64) :: r_mag
274 real(real64) :: zi, zj
275 integer :: iatom, jatom, natoms
281 force(1:space%dim, 1:natoms) =
m_zero
283 do iatom = dist%start, dist%end
284 species_i => atom(iatom)%species
285 zi = species_i%get_zval()
287 do jatom = iatom + 1, natoms
288 species_j => atom(jatom)%species
289 zj = species_j%get_zval()
291 r = pos(:, iatom) - pos(:, jatom)
293 u_e = zi * zj / r_mag
295 energy = energy + u_e
296 f(1:space%dim) = (u_e / r_mag**2) * r(1:space%dim)
297 force(1:space%dim, iatom) = force(1:space%dim, iatom) + f(1:space%dim)
298 force(1:space%dim, jatom) = force(1:space%dim, jatom) - f(1:space%dim)
305 nullify(species_i, species_j)
313 energy_components, force_components)
315 class(
space_t),
intent(in) :: space
317 type(
atom_t),
intent(in) :: atom(:)
318 integer,
intent(in) :: natoms
319 real(real64),
intent(in) :: pos(1:space%dim,1:natoms)
320 real(real64),
intent(out) :: energy
321 real(real64),
intent(out) :: force(:, :)
322 real(real64),
optional,
intent(out) :: energy_components(:)
323 real(real64),
optional,
intent(out) :: force_components(:, :, :)
325 real(real64) :: ereal, efourier, epseudo, eself
326 real(real64) :: charge
331 force(1:space%dim, 1:natoms) =
m_zero
333 call ewald_short(this%dist, space, latt, atom, pos, this%alpha, ereal, force)
334 if (
present(force_components))
then
335 force_components(1:space%dim, 1:natoms, ion_component_real) = force(1:space%dim, 1:natoms)
341 select case (space%periodic_dim)
355 call ewald_long_2d(this, space, latt, atom, natoms, pos, efourier, force)
357 call ewald_long_3d(this, space, latt, atom, natoms, pos, efourier, force, charge)
363 if (
present(energy_components))
then
364 energy_components(ion_component_real) = ereal
369 if (
present(force_components))
then
372 force(1:space%dim, 1:natoms) - force_components(1:space%dim, 1:natoms, ion_component_real)
375 energy = ereal + efourier + eself + epseudo
400 subroutine ewald_short(dist, space, latt, atom, pos, alpha, ereal, force)
402 class(
space_t),
intent(in) :: space
404 type(
atom_t),
intent(in) :: atom(:)
405 real(real64),
intent(in) :: pos(:, :)
407 real(real64),
intent(in) :: alpha
408 real(real64),
intent(out) :: ereal
409 real(real64),
intent(inout) :: force(:, :)
411 integer :: iatom, jatom, icopy, natoms
412 real(real64) :: rnorm, xi(space%dim)
413 real(real64) :: force_real(space%dim)
414 real(real64) :: zi, zj
423 rcut = 6.0_real64 / alpha
428 do iatom = dist%start, dist%end
429 if (.not. atom(iatom)%species%represents_real_atom()) cycle
430 zi = atom(iatom)%species%get_zval()
434 do icopy = 1, latt_iter%n_cells
435 rnorm = norm2(latt_iter%get(icopy))
437 if (rnorm > rcut) cycle
439 ereal = ereal +
m_half * zi * zi * erfc /rnorm
444 do jatom = iatom + 1, natoms
445 zj = atom(jatom)%species%get_zval()
448 do icopy = 1, latt_iter%n_cells
449 xi = pos(:, iatom) + latt_iter%get(icopy)
450 rnorm = norm2(xi - pos(:, jatom))
451 if (rnorm > rcut) cycle
454 force_real(:) = zj * zi * (xi - pos(:, jatom)) * &
455 (erfc / rnorm +
m_two * alpha /
sqrt(
m_pi) *
exp(-(alpha*rnorm)**2)) / rnorm**2
458 ereal = ereal + zj*zi*erfc/rnorm
461 force(1:space%dim, jatom) = force(1:space%dim, jatom) - force_real
464 force(1:space%dim, iatom) = force(1:space%dim, iatom) + force_real
484 type(
atom_t),
intent(in) :: atom(:)
485 real(real64),
intent(in) :: alpha
486 real(real64),
intent(out) :: eself
487 real(real64),
intent(out) :: charge
497 do iatom = dist%start, dist%end
498 zi = atom(iatom)%species%get_zval()
500 eself = eself - alpha /
sqrt(
m_pi) * zi**2
510 subroutine ewald_long_3d(this, space, latt, atom, natoms, pos, efourier, force, charge)
512 class(
space_t),
intent(in) :: space
514 type(
atom_t),
intent(in) :: atom(:)
515 integer,
intent(in) :: natoms
516 real(real64),
intent(in) :: pos(:,:)
517 real(real64),
intent(inout) :: efourier
518 real(real64),
intent(inout) :: force(:, :)
519 real(real64),
intent(in) :: charge
521 real(real64) :: rcut, gmax_squared
523 integer :: ix, iy, iz, isph
524 real(real64) :: gvec(3), gred(3), gg2, gx
525 real(real64) :: factor
526 complex(real64) :: sumatoms, tmp(3), aa
528 complex(real64),
allocatable :: phase(:)
532 assert(space%dim == 3)
533 assert(space%periodic_dim == 3)
536 safe_allocate(phase(1:natoms))
539 rcut =
sqrt(minval(sum(latt%klattice**2, dim=1)))
542 isph = ceiling(9.5_real64*this%alpha/rcut)
545 efourier = -
m_pi*charge**2/(
m_two*this%alpha**2*latt%rcell_volume)
548 gmax_squared = isph**2 * minval(sum(latt%klattice**2, dim=1))
556 gg2 = dot_product(gvec, gvec)
559 if (gg2 <
m_epsilon .or. gg2 > gmax_squared*1.001_real64) cycle
561 gx = -0.25_real64*gg2/this%alpha**2
563 if (gx < -36.0_real64) cycle
567 if (factor < epsilon(factor)) cycle
572 gx = sum(gvec*pos(:,iatom))
573 aa = atom(iatom)%species%get_zval()*cmplx(
cos(gx),
sin(gx), real64)
575 sumatoms = sumatoms + aa
578 efourier = efourier + real(factor*sumatoms*conjg(sumatoms), real64)
581 tmp =
m_zi*gvec*phase(iatom)
582 force(1:space%dim, iatom) = force(1:space%dim, iatom) - factor*real(conjg(tmp)*sumatoms + tmp*conjg(sumatoms), real64)
589 safe_deallocate_a(phase)
596 subroutine ewald_long_2d(this, space, latt, atom, natoms, pos, efourier, force)
598 class(
space_t),
intent(in) :: space
600 type(
atom_t),
intent(in) :: atom(:)
601 integer,
intent(in) :: natoms
602 real(real64),
intent(in) :: pos(1:space%dim,1:natoms)
603 real(real64),
intent(inout) :: efourier
604 real(real64),
intent(inout) :: force(:, :)
606 real(real64) :: rcut, gmax_squared
607 integer :: iatom, jatom
608 integer :: ix, iy, ix_max, iy_max
609 real(real64) :: gvec(space%dim), gg2, gx, gg_abs
610 real(real64) :: factor,factor1,factor2, coeff
611 real(real64) :: dz_max, dz_ij, erfc1, erfc2, tmp_erf
612 real(real64),
allocatable :: force_tmp(:,:)
613 real(real64),
parameter :: tol = 1e-10_real64
617 assert(space%periodic_dim == 2)
618 assert(space%dim == 2 .or. space%dim == 3)
623 if (space%dim == 3)
then
626 do jatom = iatom + 1, natoms
627 dz_max = max(dz_max, abs(pos(3, iatom) - pos(3, jatom)))
636 rcut =
m_two*this%alpha*4.6_real64 +
m_two*this%alpha**2*dz_max
637 if (dz_max > tol)
then
641 if (erfc1 *
exp(rcut*dz_max) < 1.e-10_real64)
exit
642 rcut = rcut * 1.414_real64
646 ix_max = ceiling(rcut/norm2(latt%klattice(:, 1)))
647 iy_max = ceiling(rcut/norm2(latt%klattice(:, 2)))
649 safe_allocate(force_tmp(1:space%dim, 1:natoms))
654 factor =
m_pi/latt%rcell_volume
657 do iatom = this%dist%start, this%dist%end
660 if (space%dim == 3)
then
661 dz_ij = pos(3, iatom) - pos(3, jatom)
666 tmp_erf =
loct_erf(this%alpha*dz_ij)
667 factor1 = dz_ij*tmp_erf
668 factor2 =
exp(-(this%alpha*dz_ij)**2)/(this%alpha*
sqrt(
m_pi))
670 efourier = efourier - factor &
671 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() * (factor1 + factor2)
674 if (iatom == jatom)cycle
677 if (space%dim == 3)
then
678 force_tmp(3, iatom) = force_tmp(3, iatom) - (-
m_two*factor) &
679 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() * tmp_erf
686 gmax_squared = sum(ix_max*latt%klattice(:, 1)**2)
687 gmax_squared = min(gmax_squared, sum(iy_max*latt%klattice(:, 2)**2))
691 do ix = -ix_max, ix_max
692 do iy = -iy_max, iy_max
694 gvec = ix*latt%klattice(:, 1) + iy*latt%klattice(:, 2)
698 if (gg2 <
m_epsilon .or. gg2 > gmax_squared*1.001_real64) cycle
700 factor =
m_half*
m_pi/(latt%rcell_volume*gg_abs)
702 do iatom = this%dist%start, this%dist%end
703 do jatom = iatom, natoms
705 gx = sum(gvec(1:2) * (pos(1:2, iatom) - pos(1:2, jatom)))
706 gx = gvec(1)*(pos(1, iatom) - pos(1, jatom)) + gvec(2)*(pos(2, iatom) - pos(2, jatom))
707 if (space%dim == 3)
then
708 dz_ij = pos(3, iatom) - pos(3, jatom)
715 factor1 =
exp(gg_abs*dz_ij)*erfc1
721 factor2 =
exp(-gg_abs*dz_ij)*erfc2
726 if (iatom == jatom)
then
732 efourier = efourier &
734 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
735 *
cos(gx)* ( factor1 + factor2)
738 if (iatom == jatom) cycle
740 force_tmp(1:2, iatom) = force_tmp(1:2, iatom) &
741 +
m_two * factor * gvec(1:2) &
742 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
743 *
sin(gx)*(factor1 + factor2)
745 force_tmp(1:2, jatom) = force_tmp(1:2, jatom) &
746 -
m_two * factor * gvec(1:2) &
747 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
748 *
sin(gx)*(factor1 + factor2)
750 factor1 = gg_abs*erfc1 &
753 factor1 = factor1*
exp(gg_abs*dz_ij)
758 factor2 = gg_abs*erfc2 &
761 factor2 = factor2*
exp(-gg_abs*dz_ij)
766 if (space%dim == 3)
then
767 force_tmp(3, iatom) = force_tmp(3, iatom) &
769 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
770 *
cos(gx)* ( factor1 - factor2)
771 force_tmp(3, jatom) = force_tmp(3, jatom) &
773 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
774 *
cos(gx)* ( factor1 - factor2)
787 force = force + force_tmp
789 safe_deallocate_a(force_tmp)
804 type(
atom_t),
intent(in) :: atom(:)
805 real(real64),
intent(out) :: epseudo
808 real(real64) :: charge
814 do iatom = dist%start, dist%end
815 select type(spec => atom(iatom)%species)
818 epseudo = epseudo +
m_pi *zi * &
819 (spec%ps%sigma_erf *
sqrt(
m_two))**2 / latt%rcell_volume * charge
831 class(
space_t),
intent(in) :: space
833 type(
atom_t),
intent(in) :: atom(:)
834 integer,
intent(in) :: natoms
835 real(real64),
intent(in) :: pos(:,:)
836 real(real64),
intent(out) :: stress_ii(space%dim, space%dim)
838 real(real64) :: stress_short(1:space%dim, 1:space%dim), stress_Ewald(1:space%dim, 1:space%dim)
845 assert(space%is_periodic())
851 select case(space%periodic_dim)
853 call ewald_3d_stress(this, space, latt, atom, natoms, pos, stress_ewald)
855 call ewald_2d_stress(this, space, latt, atom, natoms, pos, stress_ewald)
860 stress_ii = stress_short + stress_ewald
884 class(
space_t),
intent(in) :: space
886 type(
atom_t),
intent(in) :: atom(:)
887 integer,
intent(in) :: natoms
888 real(real64),
intent(in) :: pos(1:space%dim,1:natoms)
889 real(real64),
intent(out) :: stress_short(1:space%dim, 1:space%dim)
891 real(real64) :: xi(space%dim)
892 real(real64) :: r_ij, zi, zj, erfc, Hp, factor
893 integer :: iatom, jatom, icopy, idir, jdir
894 real(real64) :: alpha, rcut
901 assert(space%is_periodic())
906 rcut = 6.0_real64/alpha
912 do iatom = this%dist%start, this%dist%end
913 select type(spec => atom(iatom)%species)
917 zi = atom(iatom)%species%get_zval()
919 do icopy = 1, latt_iter%n_cells
920 xi = pos(:, iatom) + latt_iter%get(icopy)
923 zj = atom(jatom)%species%get_zval()
924 r_ij = norm2(xi - pos(:, jatom))
930 factor =
m_half*zj*zi*alpha*hp
931 do idir = 1, space%periodic_dim
932 do jdir = 1, space%periodic_dim
933 stress_short(idir, jdir) = stress_short(idir, jdir) &
934 - factor*(xi(idir) - pos(idir, jatom))*(xi(jdir) - pos(jdir, jatom))/(r_ij**2)
942 if (this%dist%parallel)
then
946 stress_short = stress_short/latt%rcell_volume
967 subroutine ewald_3d_stress(this, space, latt, atom, natoms, pos, stress_Ewald)
969 class(
space_t),
intent(in) :: space
971 type(
atom_t),
intent(in) :: atom(:)
972 integer,
intent(in) :: natoms
973 real(real64),
intent(in) :: pos(1:space%dim,1:natoms)
974 real(real64),
intent(out) :: stress_ewald(3, 3)
976 real(real64) :: zi, rcut, gmax_squared
978 integer :: ix, iy, iz, isph, idim, idir, jdir
979 real(real64) :: gred(3), gvec(3), gg2, gx
980 real(real64) :: factor, charge, charge_sq
981 complex(real64) :: sumatoms, aa
987 assert(space%dim == 3)
988 assert(space%periodic_dim == 3)
998 zi = atom(iatom)%species%get_zval()
1000 charge_sq = charge_sq + zi**2
1005 do idim = 1, space%periodic_dim
1006 rcut = min(rcut, sum(latt%klattice(1:space%periodic_dim, idim)**2))
1011 isph = ceiling(9.5_real64*this%alpha/rcut)
1014 gmax_squared = isph**2 * minval(sum(latt%klattice**2, dim=1))
1025 if (gg2 <
m_epsilon .or. gg2 > gmax_squared*1.001_real64) cycle
1027 gx = -0.25_real64*gg2/this%alpha**2
1029 if (gx < -36.0_real64) cycle
1033 if (factor < epsilon(factor)) cycle
1037 do iatom = 1, natoms
1038 gx = sum(gvec*pos(:, iatom))
1039 aa = atom(iatom)%species%get_zval()*cmplx(
cos(gx),
sin(gx), real64)
1040 sumatoms = sumatoms + aa
1043 factor = factor*abs(sumatoms)**2
1047 stress_ewald(idir, jdir) = stress_ewald(idir, jdir) &
1048 -
m_two*factor*gvec(idir)*gvec(jdir)/gg2*(0.25_real64*gg2/this%alpha**2+
m_one)
1050 stress_ewald(idir, idir) = stress_ewald(idir, idir) + factor
1059 factor =
m_half*
m_pi*charge**2/(latt%rcell_volume*this%alpha**2)
1061 stress_ewald(idir,idir) = stress_ewald(idir,idir) - factor
1064 stress_ewald = stress_ewald / latt%rcell_volume
1085 subroutine ewald_2d_stress(this, space, latt, atom, natoms, pos, stress_Ewald)
1087 type(
space_t),
intent(in) :: space
1089 type(
atom_t),
intent(in) :: atom(:)
1090 integer,
intent(in) :: natoms
1091 real(real64),
intent(in) :: pos(1:space%dim,1:natoms)
1092 real(real64),
intent(out) :: stress_Ewald(3, 3)
1094 real(real64) :: rcut, efourier
1095 integer :: iatom, jatom, idir, jdir
1096 integer :: ix, iy, ix_max, iy_max
1097 real(real64) :: gvec(3), gred(3), gg2, cos_gx, gg_abs, gmax_squared
1098 real(real64) :: factor,factor1,factor2, coeff, e_ewald
1099 real(real64) :: dz_max, z_ij, erfc1, erfc2, diff(3)
1100 real(real64),
parameter :: tol = 1e-10_real64
1104 assert(space%periodic_dim == 2)
1105 assert(space%dim == 3)
1111 do iatom = 1, natoms
1112 do jatom = iatom + 1, natoms
1113 dz_max = max(dz_max, abs(pos(3, iatom) - pos(3, jatom)))
1119 rcut =
m_two*this%alpha*4.6_real64 +
m_two*this%alpha**2*dz_max
1120 if (dz_max > tol)
then
1124 if (erfc1 *
exp(rcut*dz_max) < tol)
exit
1125 rcut = rcut * 1.414_real64
1131 factor =
m_pi/latt%rcell_volume
1133 do iatom = 1, natoms
1134 do jatom = 1, natoms
1135 z_ij = pos(3, iatom) - pos(3, jatom)
1137 factor1 = z_ij *
loct_erf(this%alpha*z_ij)
1138 factor2 =
exp(-(this%alpha*z_ij)**2)/(this%alpha*
sqrt(
m_pi))
1140 efourier = efourier - factor &
1141 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() * (factor1 + factor2)
1147 stress_ewald(idir, idir) = efourier
1151 ix_max = ceiling(rcut/norm2(latt%klattice(:, 1)))
1152 iy_max = ceiling(rcut/norm2(latt%klattice(:, 2)))
1153 gmax_squared = sum(ix_max*latt%klattice(:, 1)**2)
1154 gmax_squared = min(gmax_squared, sum(iy_max*latt%klattice(:, 2)**2))
1159 do ix = -ix_max, ix_max
1160 do iy = -iy_max, iy_max
1164 gg2 = dot_product(gvec,gvec)
1167 if (gg2 <
m_epsilon .or. gg2 > gmax_squared*1.001_real64) cycle
1170 factor =
m_fourth*
m_pi/(latt%rcell_volume*this%alpha*gg2)
1172 do iatom = 1, natoms
1173 do jatom = iatom, natoms
1174 diff = pos(:, iatom) - pos(:, jatom)
1175 cos_gx =
cos(sum(gvec(1:2) * diff(1:2)))
1181 if (iatom == jatom)
then
1189 stress_ewald(idir, jdir) = stress_ewald(idir, jdir) &
1190 - factor*gvec(idir)*gvec(jdir) * cos_gx * (factor1 + factor2) * coeff&
1191 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval()
1196 factor1 =
exp(-gg_abs*z_ij)*erfc1
1201 factor2 =
exp(gg_abs*z_ij)*erfc2
1206 e_ewald =
m_half *
m_pi/latt%rcell_volume * coeff &
1207 * atom(iatom)%species%get_zval() * atom(jatom)%species%get_zval() &
1208 * cos_gx / gg_abs * (factor1 + factor2)
1211 stress_ewald(idir, idir) = stress_ewald(idir, idir) + e_ewald
1221 stress_ewald = stress_ewald / latt%rcell_volume
1228 real(real64) function screening_function_2d(alpha, z_ij, gg_abs, erfc) result(factor)
1229 real(real64),
intent(in) :: alpha
1230 real(real64),
intent(in) :: z_ij
1231 real(real64),
intent(in) :: gg_abs
1232 real(real64),
intent(out) :: erfc
1236 arg = -alpha*z_ij +
m_half*gg_abs/alpha
1239 factor = factor*
exp(-gg_abs*z_ij)
1247 class(space_t),
intent(in) :: space
1248 type(lattice_vectors_t),
intent(in) :: latt
1249 type(atom_t),
intent(in) :: atom(:)
1250 integer,
intent(in) :: natoms
1251 real(real64),
intent(in) :: pos(1:space%dim,1:natoms)
1252 real(real64),
intent(in) :: lsize(:)
1253 type(namespace_t),
intent(in) :: namespace
1254 type(multicomm_t),
intent(in) :: mc
1257 real(real64) :: energy
1258 real(real64),
allocatable :: force(:, :), force_components(:, :, :)
1259 real(real64) :: energy_components(1:ION_NUM_COMPONENTS)
1260 integer :: iatom, idir
1267 safe_allocate(force(1:space%dim, 1:natoms))
1268 safe_allocate(force_components(1:space%dim, 1:natoms, 1:ion_num_components))
1271 energy_components = energy_components, force_components = force_components)
1273 call messages_write(
'Ionic energy =')
1274 call messages_write(energy, fmt =
'(f20.10)')
1275 call messages_info(namespace=namespace)
1277 call messages_write(
'Real space energy =')
1279 call messages_info(namespace=namespace)
1281 call messages_write(
'Self energy =')
1283 call messages_info(namespace=namespace)
1285 call messages_write(
'Fourier energy =')
1287 call messages_info(namespace=namespace)
1289 call messages_info(namespace=namespace)
1291 do iatom = 1, natoms
1292 call messages_write(
'Ionic force atom')
1293 call messages_write(iatom)
1294 call messages_write(
' =')
1295 do idir = 1, space%dim
1296 call messages_write(force(idir, iatom), fmt =
'(f20.10)')
1298 call messages_info(namespace=namespace)
1300 call messages_write(
'Real space force atom')
1301 call messages_write(iatom)
1302 call messages_write(
' =')
1303 do idir = 1, space%dim
1304 call messages_write(force_components(idir, iatom,
ion_component_real), fmt =
'(f20.10)')
1306 call messages_info(namespace=namespace)
1308 call messages_write(
'Fourier space force atom')
1309 call messages_write(iatom)
1310 call messages_write(
' =')
1311 do idir = 1, space%dim
1314 call messages_info(namespace=namespace)
1316 call messages_info(namespace=namespace)
1319 safe_deallocate_a(force)
1320 safe_deallocate_a(force_components)
double exp(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
pure logical function, public all_species_are_jellium_slab(atom)
Check if all species are jellium slab.
pure logical function, public any_species_is_jellium_sphere(atom)
Check if any species is a jellium sphere.
subroutine, public distributed_end(this)
subroutine, public distributed_nullify(this, total)
subroutine, public distributed_init(this, total, comm, tag, scalapack_compat)
Distribute N instances across M processes of communicator comm
real(real64), parameter, public m_two
real(real64), parameter, public m_max_exp_arg
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_fourth
complex(real64), parameter, public m_z0
complex(real64), parameter, public m_zi
real(real64), parameter, public r_min_atom_dist
Minimal distance between two distinguishable atoms.
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_three
real(real64), parameter, public m_five
real(real64) function screening_function_2d(alpha, z_ij, gg_abs, erfc)
Auxiliary function for the Ewald 2D stress.
subroutine, public ion_interaction_stress(this, space, latt, atom, natoms, pos, stress_ii)
Computes the contribution to the stress tensor the ion-ion energy.
subroutine, public ion_interaction_init_parallelization(this, natoms, mc)
integer, parameter ion_component_self
real(real64) function jellium_slab_energy_periodic(space, latt, atom, lsize)
Electrostatic energy of a periodic jellium slab.
subroutine, public ion_interaction_test(space, latt, atom, natoms, pos, lsize, namespace, mc)
subroutine ewald_long_2d(this, space, latt, atom, natoms, pos, efourier, force)
In-Chul Yeh and Max L. Berkowitz, J. Chem. Phys. 111, 3155 (1999).
subroutine ion_interaction_stress_short(this, space, latt, atom, natoms, pos, stress_short)
Computes the short-range contribution to the stress tensor the ion-ion energy.
subroutine ion_interaction_periodic(this, space, latt, atom, natoms, pos, energy, force, energy_components, force_components)
Total Ewald electrostatic energy and forces, for 1D, 2D and 3D systems.
real(real64) function jellium_self_energy_finite(dist, latt, atom, lsize)
Electrostatic self-interaction for jellium instances, with orthogonal cells.
subroutine, public ion_interaction_init(this, namespace, space, natoms)
subroutine ewald_short(dist, space, latt, atom, pos, alpha, ereal, force)
Short range component of the Ewald electrostatic energy and force.
subroutine pseudopotential_correction_3d(dist, latt, atom, charge, epseudo)
G=0 component of Ewald energy arising from the pseudopotentials, for 3D systems.
subroutine ewald_long_3d(this, space, latt, atom, natoms, pos, efourier, force, charge)
Computes the long-range part of the 3D Ewald summation.
integer, parameter ion_component_real
integer, parameter ion_num_components
subroutine ewald_3d_stress(this, space, latt, atom, natoms, pos, stress_Ewald)
Computes the contribution to the stress tensor from the 3D Ewald sum.
integer, parameter ion_component_fourier
subroutine ion_interaction_finite(dist, space, atom, pos, lsize, energy, force)
Electrostatic Ewald energy and forces for finite systems.
subroutine, public ion_interaction_end(this)
subroutine, public ion_interaction_calculate(this, space, latt, atom, natoms, pos, lsize, energy, force, energy_components, force_components)
Top level routine for computing electrostatic energies and forces between ions.
subroutine ewald_2d_stress(this, space, latt, atom, natoms, pos, stress_Ewald)
Computes the contribution to the stress tensor from the 2D Ewald sum.
subroutine ewald_self_interaction(dist, atom, alpha, eself, charge)
@ brief Ewald self-interaction energy
subroutine, public kpoints_to_absolute(latt, kin, kout)
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
This module handles the communicators for the various parallelization strategies.
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....
The following class implements a lattice iterator. It allows one to loop over all cells that are with...
An abstract class for species. Derived classes include jellium, all electron, and pseudopotential spe...