Octopus
ion_interaction.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2021 N. Tancogne-Dejean
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
23 use atom_oct_m
24 use comm_oct_m
25 use debug_oct_m
26 use global_oct_m
33 use mpi_oct_m
36 use parser_oct_m
38 use ps_oct_m
40 use space_oct_m
43
44 implicit none
45
46 private
47 public :: &
55
57 real(real64) :: alpha
58 type(distributed_t) :: dist
59 end type ion_interaction_t
60
61 integer, parameter :: &
62 ION_COMPONENT_REAL = 1, &
66
67contains
68
69 subroutine ion_interaction_init(this, namespace, space, natoms)
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
74
75 push_sub(ion_interaction_init)
76
77 !%Variable EwaldAlpha
78 !%Type float
79 !%Default 0.21
80 !%Section Hamiltonian
81 !%Description
82 !% The value 'Alpha' that controls the splitting of the Coulomb
83 !% interaction in the Ewald sum used to calculation the ion-ion
84 !% interaction for periodic systems. This value affects the speed
85 !% of the calculation, normally users do not need to modify it.
86 !%End
87 call parse_variable(namespace, 'EwaldAlpha', 0.21_real64, this%alpha)
88
89 call distributed_nullify(this%dist, natoms)
90
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.')
95 call messages_warning(namespace=namespace)
96 end if
97
99 end subroutine ion_interaction_init
100
101 subroutine ion_interaction_init_parallelization(this, natoms, mc)
102 type(ion_interaction_t), intent(inout) :: this
103 integer, intent(in) :: natoms
104 type(multicomm_t), intent(in) :: mc
105
107
108 !As the code below is not parallelized with any of k-point, states nor domain
109 !we can safely parallelize it over atoms
110 if (debug%info) then
111 call distributed_init(this%dist, natoms, mc%master_comm, "Ions")
112 else
113 call distributed_init(this%dist, natoms, mc%master_comm)
114 end if
118
119 subroutine ion_interaction_end(this)
120 type(ion_interaction_t), intent(inout) :: this
121
122 push_sub(ion_interaction_end)
123
124 this%alpha = -m_one
125
126 call distributed_end(this%dist)
127
128 pop_sub(ion_interaction_end)
129 end subroutine ion_interaction_end
130
135 subroutine ion_interaction_calculate(this, space, latt, atom, natoms, pos, lsize, energy, force, &
136 energy_components, force_components)
137 type(ion_interaction_t), intent(inout) :: this
138 class(space_t), intent(in) :: space
139 type(lattice_vectors_t), intent(in) :: latt
140 type(atom_t), intent(in) :: atom(:)
141 integer, intent(in) :: natoms
142 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
143 real(real64), intent(in) :: lsize(:)
144 real(real64), intent(out) :: energy
145 real(real64), intent(out) :: force(:, :)
146 real(real64), optional, intent(out) :: energy_components(:)
147 real(real64), optional, intent(out) :: force_components(:, :, :)
148
151 call profiling_in("ION_ION_INTERACTION")
152
153 if (present(energy_components)) then
154 assert(ubound(energy_components, dim = 1) == ion_num_components)
155 energy_components = m_zero
156 end if
157
158 if (present(force_components)) then
159 assert(all(ubound(force_components) == (/space%dim, natoms, ion_num_components/)))
160 force_components = m_zero
161 end if
163 if (space%is_periodic() .and. any_species_is_jellium_sphere(atom)) then
164 call messages_not_implemented('No periodic implementation of ion-ion energy for the jellium sphere')
165 end if
166
167 if (space%is_periodic()) then
168 if (all_species_are_jellium_slab(atom)) then
169 energy = jellium_slab_energy_periodic(space, atom, lsize)
170 force = 0._real64
171 else
172 call ion_interaction_periodic(this, space, latt, atom, natoms, pos, energy, force, energy_components, force_components)
173 end if
174 else
175 call ion_interaction_finite(this%dist, space, atom, pos, lsize, energy, force)
176 energy = energy + jellium_self_energy_finite(this%dist, latt, atom, lsize)
177 end if
178
179 call profiling_out("ION_ION_INTERACTION")
181
182 end subroutine ion_interaction_calculate
183
189 function jellium_slab_energy_periodic(space, atom, lsize) result(energy)
190 class(space_t), intent(in) :: space
191 type(atom_t), intent(in) :: atom(:)
192 real(real64), intent(in) :: lsize(:)
193 real(real64) :: energy
195 real(real64) :: area
196
197 ! Implementation assumes a single atom
198 assert(size(atom) == 1)
199 ! This is only allowed if periodic dim = 2. In that case the lattice volume is in fact an area.
200 assert(space%periodic_dim == 2)
201
202 select type(spec => atom(1)%species)
203 type is (jellium_slab_t)
204 area = lsize(1) * lsize(2) * m_four
205 energy = m_pi * spec%get_density(lsize) **2 * area * spec%thickness()**3 / m_three
206 class default
207 assert(.false.)
208 end select
209
211
225 function jellium_self_energy_finite(dist, latt, atom, lsize) result(energy)
226 type(distributed_t), intent(in) :: dist
227 type(lattice_vectors_t), intent(in) :: latt
228 type(atom_t), intent(in) :: atom(:)
229 real(real64), intent(in) :: lsize(:)
230 real(real64) :: energy
231
232 real(real64) :: zi
233 integer :: iatom
234 logical :: lattice_is_orthogonal
235 class(species_t), pointer :: spec
236
238
239 energy = 0._real64
240 lattice_is_orthogonal = .not. latt%nonorthogonal
241
242 do iatom = dist%start, dist%end
243 spec => atom(iatom)%species
244 zi = spec%get_zval()
245
246 select type(spec)
247 type is (jellium_sphere_t)
248 energy = energy + (m_three / m_five) * zi**2 / spec%radius()
249 ! The part depending on the simulation sphere is neglected
250
251 type is (jellium_slab_t)
252 ! Jellium slab energy only implemented for orthogonal cells.
253 ! One would need to replace (lsize(1) * lsize(2)) * spec%thickness()) with the triple product
254 assert(lattice_is_orthogonal)
255 energy = energy + m_pi * zi**2 / (m_four * lsize(1)*lsize(2)) * spec%thickness() / m_three
256 ! The part depending on the simulation box transverse dimension is neglected
257 end select
258 nullify(spec)
259 enddo
260
261 call comm_allreduce(dist%mpi_grp, energy)
262
264
265 end function jellium_self_energy_finite
266
268 subroutine ion_interaction_finite(dist, space, atom, pos, lsize, energy, force)
269 type(distributed_t), intent(in) :: dist
270 class(space_t), intent(in) :: space
271 type(atom_t), intent(in) :: atom(:)
272 real(real64), intent(in) :: pos(:,:)
273 real(real64), intent(in) :: lsize(:)
274 real(real64), intent(out) :: energy
275 real(real64), intent(out) :: force(:, :)
276
277 class(species_t), pointer :: species_i, species_j
278 real(real64) :: r(space%dim), f(space%dim)
279 real(real64) :: r_mag
280 real(real64) :: u_e
281 real(real64) :: zi, zj
282 integer :: iatom, jatom, natoms
283
284 push_sub(ion_interaction_finite)
285
286 natoms = size(atom)
287 energy = m_zero
288 force(1:space%dim, 1:natoms) = m_zero
289
290 do iatom = dist%start, dist%end
291 species_i => atom(iatom)%species
292 zi = species_i%get_zval()
293
294 do jatom = iatom + 1, natoms
295 species_j => atom(jatom)%species
296 zj = species_j%get_zval()
297
298 r = pos(:, iatom) - pos(:, jatom)
299 r_mag = norm2(r)
300 u_e = zi * zj / r_mag
301
302 energy = energy + u_e
303 f(1:space%dim) = (u_e / r_mag**2) * r(1:space%dim)
304 force(1:space%dim, iatom) = force(1:space%dim, iatom) + f(1:space%dim)
305 force(1:space%dim, jatom) = force(1:space%dim, jatom) - f(1:space%dim)
306 end do
307 end do
308
309 call comm_allreduce(dist%mpi_grp, energy)
310 call comm_allreduce(dist%mpi_grp, force)
311
312 nullify(species_i, species_j)
313
315
316 end subroutine ion_interaction_finite
317
319 subroutine ion_interaction_periodic(this, space, latt, atom, natoms, pos, energy, force, &
320 energy_components, force_components)
321 type(ion_interaction_t), intent(in) :: this
322 class(space_t), intent(in) :: space
323 type(lattice_vectors_t), intent(in) :: latt
324 type(atom_t), intent(in) :: atom(:)
325 integer, intent(in) :: natoms
326 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
327 real(real64), intent(out) :: energy
328 real(real64), intent(out) :: force(:, :)
329 real(real64), optional, intent(out) :: energy_components(:)
330 real(real64), optional, intent(out) :: force_components(:, :, :)
331
332 real(real64) :: ereal, efourier, epseudo, eself
333 real(real64) :: charge
334
336
337 energy = m_zero
338 force(1:space%dim, 1:natoms) = m_zero
339
340 call ewald_short(this%dist, space, latt, atom, pos, this%alpha, ereal, force)
341 if (present(force_components)) then
342 force_components(1:space%dim, 1:natoms, ion_component_real) = force(1:space%dim, 1:natoms)
343 end if
344
345 call ewald_self_interaction(this%dist, atom, this%alpha, eself, charge)
346
347 call profiling_in("EWALD_LONG")
348 select case (space%periodic_dim)
349 case (1)
350 ! Warning added in init routine, such that it is not displayed per SCF step
351 efourier = m_zero
352 ! Do not confuse the user and set to zero all the other components
353 ereal = m_zero
354 eself = m_zero
355 force = m_zero
356 epseudo = m_zero
357 case (2)
358 ! The energy contribution of the long range part of the pseudo is
359 ! not correctly accounted for in systems periodic in 1D or 2D, however
360 ! this term should not appear here anyway. See Issue #950.
361 epseudo = m_zero
362 call ewald_long_2d(this, space, latt, atom, natoms, pos, efourier, force)
363 case (3)
364 call ewald_long_3d(this, space, latt, atom, natoms, pos, efourier, force, charge)
365 !TODO(Alex/Nicolas) Issue #950. Refactor: Move G=0 correction from ion-ion energy to pseudopotential energy
366 call pseudopotential_correction_3d(this%dist, latt, atom, charge, epseudo)
367 end select
368 call profiling_out("EWALD_LONG")
369
370 if (present(energy_components)) then
371 energy_components(ion_component_real) = ereal
372 energy_components(ion_component_self) = eself
373 energy_components(ion_component_fourier) = efourier
374 end if
375
376 if (present(force_components)) then
377 ! This is dependent on the order in which the force terms are computed
378 force_components(1:space%dim, 1:natoms, ion_component_fourier) = &
379 force(1:space%dim, 1:natoms) - force_components(1:space%dim, 1:natoms, ion_component_real)
380 end if
381
382 energy = ereal + efourier + eself + epseudo
383
385 end subroutine ion_interaction_periodic
386
407 subroutine ewald_short(dist, space, latt, atom, pos, alpha, ereal, force)
408 type(distributed_t), intent(in) :: dist
409 class(space_t), intent(in) :: space
410 type(lattice_vectors_t), intent(in) :: latt
411 type(atom_t), intent(in) :: atom(:)
412 real(real64), intent(in) :: pos(:, :)
413
414 real(real64), intent(in) :: alpha
415 real(real64), intent(out) :: ereal
416 real(real64), intent(inout) :: force(:, :)
417 ! Intent(inout) allows force contributions to be summed
418 integer :: iatom, jatom, icopy, natoms
419 real(real64) :: rnorm, xi(space%dim)
420 real(real64) :: force_real(space%dim)
421 real(real64) :: zi, zj
422 real(real64) :: erfc
423 real(real64) :: rcut
424 type(lattice_iterator_t) :: latt_iter
425
426 push_sub_with_profile(ewald_short)
427
428 ereal = m_zero
429 ! Origin of this value is not documented
430 rcut = 6.0_real64 / alpha
431 latt_iter = lattice_iterator_t(latt, rcut)
432 natoms = size(atom)
433
434 !$omp parallel default(shared) private(iatom, jatom, zi, zj, icopy, xi, rnorm, erfc, force_real) reduction(+:ereal, force)
435 do iatom = dist%start, dist%end
436 if (.not. atom(iatom)%species%represents_real_atom()) cycle
437 zi = atom(iatom)%species%get_zval()
438
439 ! Diagonal terms iatom == jatom for all cells, except T=(0,0,0)
440 !$omp do
441 do icopy = 1, latt_iter%n_cells
442 rnorm = norm2(latt_iter%get(icopy))
443 if (rnorm < r_min_atom_dist) cycle
444 if (rnorm > rcut) cycle
445 erfc = m_one - loct_erf(alpha * rnorm)
446 ereal = ereal + m_half * zi * zi * erfc /rnorm
447 enddo
448 !$omp end do
449
450 ! Upper triangle, for all replica cells
451 do jatom = iatom + 1, natoms
452 zj = atom(jatom)%species%get_zval()
453
454 !$omp do
455 do icopy = 1, latt_iter%n_cells
456 xi = pos(:, iatom) + latt_iter%get(icopy)
457 rnorm = norm2(xi - pos(:, jatom))
458 if (rnorm > rcut) cycle
459
460 erfc = m_one - loct_erf(alpha * rnorm)
461 force_real(:) = zj * zi * (xi - pos(:, jatom)) * &
462 (erfc / rnorm + m_two * alpha / sqrt(m_pi) *exp(-(alpha*rnorm)**2)) / rnorm**2
463
464 ! Factor 1/2 omitted as one is only summing over upper triangle
465 ereal = ereal + zj*zi*erfc/rnorm
466
467 ! Upper trianglar contribution
468 force(1:space%dim, jatom) = force(1:space%dim, jatom) - force_real
469
470 ! Lower triangular contribution
471 force(1:space%dim, iatom) = force(1:space%dim, iatom) + force_real
472 end do
473 !$omp end do
474
475 end do
476 end do
477 !$omp end parallel
478
479 call comm_allreduce(dist%mpi_grp, ereal)
480 call comm_allreduce(dist%mpi_grp, force)
481
482 pop_sub_with_profile(ewald_short)
483 end subroutine ewald_short
484
489 subroutine ewald_self_interaction(dist, atom, alpha, eself, charge)
490 type(distributed_t), intent(in) :: dist
491 type(atom_t), intent(in) :: atom(:)
492 real(real64), intent(in) :: alpha
493 real(real64), intent(out) :: eself
494 real(real64), intent(out) :: charge
495
496 integer :: iatom
497 real(real64) :: zi
498
499 push_sub(ewald_self_interaction)
501 eself = m_zero
502 charge = m_zero
503
504 do iatom = dist%start, dist%end
505 zi = atom(iatom)%species%get_zval()
506 charge = charge + zi
507 eself = eself - alpha / sqrt(m_pi) * zi**2
508 end do
509
510 call comm_allreduce(dist%mpi_grp, eself)
511 call comm_allreduce(dist%mpi_grp, charge)
512
514 end subroutine ewald_self_interaction
515
517 subroutine ewald_long_3d(this, space, latt, atom, natoms, pos, efourier, force, charge)
518 type(ion_interaction_t), intent(in) :: this
519 class(space_t), intent(in) :: space
520 type(lattice_vectors_t), intent(in) :: latt
521 type(atom_t), intent(in) :: atom(:)
522 integer, intent(in) :: natoms
523 real(real64), intent(in) :: pos(:,:)
524 real(real64), intent(inout) :: efourier
525 real(real64), intent(inout) :: force(:, :)
526 real(real64), intent(in) :: charge
527
528 real(real64) :: rcut, gmax_squared
529 integer :: iatom
530 integer :: ix, iy, iz, isph
531 real(real64) :: gvec(3), gred(3), gg2, gx
532 real(real64) :: factor
533 complex(real64) :: sumatoms, tmp(3), aa
534
535 complex(real64), allocatable :: phase(:)
536
537 push_sub(ewald_long_3d)
538
539 assert(space%dim == 3)
540 assert(space%periodic_dim == 3)
541
542 ! And the long-range part, using an Ewald sum
543 safe_allocate(phase(1:natoms))
544
545 ! get a converged value for the cutoff in g
546 rcut = sqrt(minval(sum(latt%klattice**2, dim=1)))
547
548 ! 9.5 is a constant that controls the range separation
549 isph = ceiling(9.5_real64*this%alpha/rcut)
550
551 ! First the G = 0 term (charge was calculated previously)
552 efourier = -m_pi*charge**2/(m_two*this%alpha**2*latt%rcell_volume)
553
554 ! Getting the G-shell cutoff
555 gmax_squared = isph**2 * minval(sum(latt%klattice**2, dim=1))
556
557 do ix = -isph, isph
558 do iy = -isph, isph
559 do iz = -isph, isph
560
561 gred = [ix, iy, iz]
562 call kpoints_to_absolute(latt, gred, gvec)
563 gg2 = dot_product(gvec, gvec)
564
565 ! g=0 must be removed from the sum
566 if (gg2 < m_epsilon .or. gg2 > gmax_squared*1.001_real64) cycle
567
568 gx = -0.25_real64*gg2/this%alpha**2
569
570 if (gx < -36.0_real64) cycle
571
572 factor = m_two*m_pi/latt%rcell_volume*exp(gx)/gg2
573
574 if (factor < epsilon(factor)) cycle
575
576 sumatoms = m_z0
577 !$omp parallel do private(iatom, gx, aa) reduction(+:sumatoms)
578 do iatom = 1, natoms
579 gx = sum(gvec*pos(:,iatom))
580 aa = atom(iatom)%species%get_zval()*cmplx(cos(gx), sin(gx), real64)
581 phase(iatom) = aa
582 sumatoms = sumatoms + aa
583 end do
584
585 efourier = efourier + real(factor*sumatoms*conjg(sumatoms), real64)
586
587 do iatom = 1, natoms
588 tmp = m_zi*gvec*phase(iatom)
589 force(1:space%dim, iatom) = force(1:space%dim, iatom) - factor*real(conjg(tmp)*sumatoms + tmp*conjg(sumatoms), real64)
590 end do
591
592 end do
593 end do
594 end do
595
596 safe_deallocate_a(phase)
597
598 pop_sub(ewald_long_3d)
599
600 end subroutine ewald_long_3d
601
603 subroutine ewald_long_2d(this, space, latt, atom, natoms, pos, efourier, force)
604 type(ion_interaction_t), intent(in) :: this
605 class(space_t), intent(in) :: space
606 type(lattice_vectors_t), intent(in) :: latt
607 type(atom_t), intent(in) :: atom(:)
608 integer, intent(in) :: natoms
609 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
610 real(real64), intent(inout) :: efourier
611 real(real64), intent(inout) :: force(:, :)
612
613 real(real64) :: rcut, gmax_squared
614 integer :: iatom, jatom
615 integer :: ix, iy, ix_max, iy_max
616 real(real64) :: gvec(space%dim), gg2, gx, gg_abs
617 real(real64) :: factor,factor1,factor2, coeff
618 real(real64) :: dz_max, dz_ij, erfc1, erfc2, tmp_erf
619 real(real64), allocatable :: force_tmp(:,:)
620 real(real64), parameter :: tol = 1e-10_real64
621
622 push_sub(ewald_long_2d)
623
624 assert(space%periodic_dim == 2)
625 assert(space%dim == 2 .or. space%dim == 3)
626
627 ! And the long-range part, using an Ewald sum
628
629 ! Searching maximum distance
630 if (space%dim == 3) then
631 dz_max = m_zero
632 do iatom = 1, natoms
633 do jatom = iatom + 1, natoms
634 dz_max = max(dz_max, abs(pos(3, iatom) - pos(3, jatom)))
635 end do
636 end do
637 else
638 ! For a 2D system, all atoms are on the plane, so the distance is zero
639 dz_max = m_zero
640 end if
641
642 !get a converged value for the cutoff in g
643 rcut = m_two*this%alpha*4.6_real64 + m_two*this%alpha**2*dz_max
644 if (dz_max > tol) then
645 do
646 if (rcut * dz_max >= m_max_exp_arg) exit !Maximum double precision number
647 erfc1 = m_one - loct_erf(this%alpha*dz_max + m_half*rcut/this%alpha)
648 if (erfc1 * exp(rcut*dz_max) < 1.e-10_real64) exit
649 rcut = rcut * 1.414_real64
650 end do
651 end if
652
653 ix_max = ceiling(rcut/norm2(latt%klattice(:, 1)))
654 iy_max = ceiling(rcut/norm2(latt%klattice(:, 2)))
655
656 safe_allocate(force_tmp(1:space%dim, 1:natoms))
657 force_tmp = m_zero
658
659 ! First the G = 0 term
660 efourier = m_zero
661 factor = m_pi/latt%rcell_volume
662 !$omp parallel do private(jatom, dz_ij, tmp_erf, factor1, factor2) reduction(+:efourier,force_tmp) &
663 !$omp& collapse(2)
664 do iatom = this%dist%start, this%dist%end
665 do jatom = 1, natoms
666 ! efourier
667 if (space%dim == 3) then
668 dz_ij = pos(3, iatom) - pos(3, jatom)
669 else
670 dz_ij = m_zero
671 end if
672
673 tmp_erf = loct_erf(this%alpha*dz_ij)
674 factor1 = dz_ij*tmp_erf
675 factor2 = exp(-(this%alpha*dz_ij)**2)/(this%alpha*sqrt(m_pi))
676
677 efourier = efourier - factor &
678 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() * (factor1 + factor2)
679
680 ! force
681 if (iatom == jatom)cycle
682 if (abs(tmp_erf) < m_epsilon) cycle
683
684 if (space%dim == 3) then
685 force_tmp(3, iatom) = force_tmp(3, iatom) - (- m_two*factor) &
686 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() * tmp_erf
687 end if
688
689 end do
690 end do
691
692 ! Getting the G-shell cutoff
693 gmax_squared = sum(ix_max*latt%klattice(:, 1)**2)
694 gmax_squared = min(gmax_squared, sum(iy_max*latt%klattice(:, 2)**2))
695
696 !$omp parallel do private(iy, gvec, gg2, gg_abs, factor, iatom, jatom, gx, dz_ij, erfc1, factor1, erfc2, factor2, coeff) &
697 !$omp& collapse(2) reduction(+:efourier, force_tmp)
698 do ix = -ix_max, ix_max
699 do iy = -iy_max, iy_max
700
701 gvec = ix*latt%klattice(:, 1) + iy*latt%klattice(:, 2)
702 gg2 = sum(gvec**2)
703
704 ! g=0 must be removed from the sum
705 if (gg2 < m_epsilon .or. gg2 > gmax_squared*1.001_real64) cycle
706 gg_abs = sqrt(gg2)
707 factor = m_half*m_pi/(latt%rcell_volume*gg_abs)
708
709 do iatom = this%dist%start, this%dist%end
710 do jatom = iatom, natoms
711 ! efourier
712 gx = sum(gvec(1:2) * (pos(1:2, iatom) - pos(1:2, jatom)))
713 gx = gvec(1)*(pos(1, iatom) - pos(1, jatom)) + gvec(2)*(pos(2, iatom) - pos(2, jatom))
714 if (space%dim == 3) then
715 dz_ij = pos(3, iatom) - pos(3, jatom)
716 else
717 dz_ij = m_zero
718 end if
719
720 erfc1 = m_one - loct_erf(this%alpha*dz_ij + m_half*gg_abs/this%alpha)
721 if (abs(erfc1) > m_epsilon) then
722 factor1 = exp(gg_abs*dz_ij)*erfc1
723 else
724 factor1 = m_zero
725 end if
726 erfc2 = m_one - loct_erf(-this%alpha*dz_ij + m_half*gg_abs/this%alpha)
727 if (abs(erfc2) > m_epsilon) then
728 factor2 = exp(-gg_abs*dz_ij)*erfc2
729 else
730 factor2 = m_zero
731 end if
732
733 if (iatom == jatom) then
734 coeff = m_one
735 else
736 coeff = m_two
737 end if
738
739 efourier = efourier &
740 + factor * coeff &
741 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
742 * cos(gx)* ( factor1 + factor2)
743
744 ! force
745 if (iatom == jatom) cycle
746
747 force_tmp(1:2, iatom) = force_tmp(1:2, iatom) &
748 + m_two * factor * gvec(1:2) &
749 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
750 *sin(gx)*(factor1 + factor2)
751
752 force_tmp(1:2, jatom) = force_tmp(1:2, jatom) &
753 - m_two * factor * gvec(1:2) &
754 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
755 *sin(gx)*(factor1 + factor2)
756
757 factor1 = gg_abs*erfc1 &
758 - m_two*this%alpha/sqrt(m_pi)*exp(-(this%alpha*dz_ij + m_half*gg_abs/this%alpha)**2)
759 if (abs(factor1) > m_epsilon) then
760 factor1 = factor1*exp(gg_abs*dz_ij)
761 else
762 factor1 = m_zero
763 end if
764
765 factor2 = gg_abs*erfc2 &
766 - m_two*this%alpha/sqrt(m_pi)*exp(-(-this%alpha*dz_ij + m_half*gg_abs/this%alpha)**2)
767 if (abs(factor2) > m_epsilon) then
768 factor2 = factor2*exp(-gg_abs*dz_ij)
769 else
770 factor2 = m_zero
771 end if
772
773 if (space%dim == 3) then
774 force_tmp(3, iatom) = force_tmp(3, iatom) &
775 - m_two*factor &
776 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
777 * cos(gx)* ( factor1 - factor2)
778 force_tmp(3, jatom) = force_tmp(3, jatom) &
779 + m_two*factor &
780 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
781 * cos(gx)* ( factor1 - factor2)
782 end if
783
784 end do
785 end do
786
787
788 end do
789 end do
790
791 call comm_allreduce(this%dist%mpi_grp, efourier)
792 call comm_allreduce(this%dist%mpi_grp, force_tmp)
793
794 force = force + force_tmp
795
796 safe_deallocate_a(force_tmp)
797
798 pop_sub(ewald_long_2d)
799 end subroutine ewald_long_2d
800
801 !TODO(Alex/Nicolas) Issue #950. Refactor: Move G=0 correction from ion-ion energy to pseudopotential energy
808 subroutine pseudopotential_correction_3d(dist, latt, atom, charge, epseudo)
809 type(distributed_t), intent(in) :: dist
810 type(lattice_vectors_t), intent(in) :: latt
811 type(atom_t), intent(in) :: atom(:)
812 real(real64), intent(out) :: epseudo
813
814 real(real64) :: zi
815 real(real64) :: charge
816 integer :: iatom
817
819
820 epseudo = m_zero
821 do iatom = dist%start, dist%end
822 select type(spec => atom(iatom)%species)
823 class is(pseudopotential_t)
824 zi = spec%get_zval()
825 epseudo = epseudo + m_pi *zi * &
826 (spec%ps%sigma_erf * sqrt(m_two))**2 / latt%rcell_volume * charge
827 end select
828 end do
829 call comm_allreduce(dist%mpi_grp, epseudo)
830
832
833 end subroutine pseudopotential_correction_3d
834
836 subroutine ion_interaction_stress(this, space, latt, atom, natoms, pos, stress_ii)
837 type(ion_interaction_t), intent(inout) :: this
838 class(space_t), intent(in) :: space
839 type(lattice_vectors_t), intent(in) :: latt
840 type(atom_t), intent(in) :: atom(:)
841 integer, intent(in) :: natoms
842 real(real64), intent(in) :: pos(:,:)
843 real(real64), intent(out) :: stress_ii(space%dim, space%dim)
844
845 real(real64) :: stress_short(1:space%dim, 1:space%dim), stress_Ewald(1:space%dim, 1:space%dim)
846
847 push_sub(ion_interaction_stress)
848
849 stress_ii = m_zero
850
851 ! Only implemented in the periodic case
852 assert(space%is_periodic())
853
854 ! Short range part in real space
855 call ion_interaction_stress_short(this, space, latt, atom, natoms, pos, stress_short)
856
857 ! Long range part in Fourier space
858 select case(space%periodic_dim)
859 case(3)
860 call ewald_3d_stress(this, space, latt, atom, natoms, pos, stress_ewald)
861 case(2)
862 call ewald_2d_stress(this, space, latt, atom, natoms, pos, stress_ewald)
863 case default
864 assert(.false.)
865 end select
866
867 stress_ii = stress_short + stress_ewald
868
870 end subroutine ion_interaction_stress
871
872 ! ---------------------------------------------------------
888
889 subroutine ion_interaction_stress_short(this, space, latt, atom, natoms, pos, stress_short)
890 type(ion_interaction_t), intent(inout) :: this
891 class(space_t), intent(in) :: space
892 type(lattice_vectors_t), intent(in) :: latt
893 type(atom_t), intent(in) :: atom(:)
894 integer, intent(in) :: natoms
895 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
896 real(real64), intent(out) :: stress_short(1:space%dim, 1:space%dim)
897
898 real(real64) :: xi(space%dim)
899 real(real64) :: r_ij, zi, zj, erfc, Hp, factor
900 integer :: iatom, jatom, icopy, idir, jdir
901 real(real64) :: alpha, rcut
902 type(lattice_iterator_t) :: latt_iter
903
905 call profiling_in("ION_ION_STRESS_SHORT")
906
907 ! Only implemented in the periodic case
908 assert(space%is_periodic())
909
910 alpha = this%alpha
911
912 ! See the code for the energy above to understand this parameter
913 rcut = 6.0_real64/alpha
914
915 ! the short-range part is calculated directly
916 stress_short = m_zero
917 latt_iter = lattice_iterator_t(latt, rcut)
918
919 do iatom = this%dist%start, this%dist%end
920 select type(spec => atom(iatom)%species)
921 class is(jellium_t)
922 cycle
923 end select
924 zi = atom(iatom)%species%get_zval()
925
926 do icopy = 1, latt_iter%n_cells
927 xi = pos(:, iatom) + latt_iter%get(icopy)
928
929 do jatom = 1, natoms
930 zj = atom(jatom)%species%get_zval()
931 r_ij = norm2(xi - pos(:, jatom))
932
933 if (r_ij < r_min_atom_dist) cycle
934
935 erfc = m_one - loct_erf(alpha*r_ij)
936 hp = -m_two/sqrt(m_pi)*exp(-(alpha*r_ij)**2) - erfc/(alpha*r_ij)
937 factor = m_half*zj*zi*alpha*hp
938 do idir = 1, space%periodic_dim
939 do jdir = 1, space%periodic_dim
940 stress_short(idir, jdir) = stress_short(idir, jdir) &
941 - factor*(xi(idir) - pos(idir, jatom))*(xi(jdir) - pos(jdir, jatom))/(r_ij**2)
942 end do
943 end do
944
945 end do
946 end do
947 end do
948
949 if (this%dist%parallel) then
950 call comm_allreduce(this%dist%mpi_grp, stress_short)
951 end if
952
953 stress_short = stress_short/latt%rcell_volume
954
955 call profiling_out("ION_ION_STRESS_SHORT")
956
958 end subroutine ion_interaction_stress_short
959
960
961
962 ! ---------------------------------------------------------
974 subroutine ewald_3d_stress(this, space, latt, atom, natoms, pos, stress_Ewald)
975 type(ion_interaction_t), intent(inout) :: this
976 class(space_t), intent(in) :: space
977 type(lattice_vectors_t), intent(in) :: latt
978 type(atom_t), intent(in) :: atom(:)
979 integer, intent(in) :: natoms
980 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
981 real(real64), intent(out) :: stress_ewald(3, 3)
983 real(real64) :: zi, rcut, gmax_squared
984 integer :: iatom
985 integer :: ix, iy, iz, isph, idim, idir, jdir
986 real(real64) :: gred(3), gvec(3), gg2, gx
987 real(real64) :: factor, charge, charge_sq
988 complex(real64) :: sumatoms, aa
989
990 call profiling_in("STRESS_3D_EWALD")
991 push_sub(ewald_3d_stress)
992
993 ! Currently this is only implemented for 3D
994 assert(space%dim == 3)
995 assert(space%periodic_dim == 3) ! Not working for mixed periodicity
996 ! (klattice along the non-periodic directions is wrong)
997 ! Anyway gg/gg2 is not correct for mixed periodicity
998
999 stress_ewald = m_zero
1000
1001 ! And the long-range part, using an Ewald sum
1002 charge = m_zero
1003 charge_sq = m_zero
1004 do iatom = 1, natoms
1005 zi = atom(iatom)%species%get_zval()
1006 charge = charge + zi
1007 charge_sq = charge_sq + zi**2
1008 end do
1009
1010 ! get a converged value for the cutoff in g
1011 rcut = huge(rcut)
1012 do idim = 1, space%periodic_dim
1013 rcut = min(rcut, sum(latt%klattice(1:space%periodic_dim, idim)**2))
1014 end do
1015
1016 rcut = sqrt(rcut)
1017
1018 isph = ceiling(9.5_real64*this%alpha/rcut)
1019
1020 ! Getting the G-shell cutoff
1021 gmax_squared = isph**2 * minval(sum(latt%klattice**2, dim=1))
1022
1023 do ix = -isph, isph
1024 do iy = -isph, isph
1025 do iz = -isph, isph
1026
1027 gred = [ix, iy, iz]
1028 call kpoints_to_absolute(latt, gred, gvec)
1029 gg2 = sum(gvec**2)
1030
1031 ! g=0 must be removed from the sum
1032 if (gg2 < m_epsilon .or. gg2 > gmax_squared*1.001_real64) cycle
1033
1034 gx = -0.25_real64*gg2/this%alpha**2
1035
1036 if (gx < -36.0_real64) cycle
1037
1038 factor = m_two*m_pi*exp(gx)/(latt%rcell_volume*gg2)
1039
1040 if (factor < epsilon(factor)) cycle
1041
1042 sumatoms = m_z0
1043
1044 do iatom = 1, natoms
1045 gx = sum(gvec*pos(:, iatom))
1046 aa = atom(iatom)%species%get_zval()*cmplx(cos(gx), sin(gx), real64)
1047 sumatoms = sumatoms + aa
1048 end do
1049
1050 factor = factor*abs(sumatoms)**2
1051
1052 do idir = 1, 3
1053 do jdir = 1, 3
1054 stress_ewald(idir, jdir) = stress_ewald(idir, jdir) &
1055 - m_two*factor*gvec(idir)*gvec(jdir)/gg2*(0.25_real64*gg2/this%alpha**2+m_one)
1056 end do
1057 stress_ewald(idir, idir) = stress_ewald(idir, idir) + factor
1058 end do
1059
1060 end do
1061 end do
1062 end do
1063
1064
1065 ! The G = 0 term of the Ewald summation
1066 factor = m_half*m_pi*charge**2/(latt%rcell_volume*this%alpha**2)
1067 do idir = 1,3
1068 stress_ewald(idir,idir) = stress_ewald(idir,idir) - factor
1069 end do
1070
1071 stress_ewald = stress_ewald / latt%rcell_volume
1072
1073
1074 call profiling_out("STRESS_3D_EWALD")
1075 pop_sub(ewald_3d_stress)
1076
1077 end subroutine ewald_3d_stress
1078
1079 ! ---------------------------------------------------------
1092 subroutine ewald_2d_stress(this, space, latt, atom, natoms, pos, stress_Ewald)
1093 type(ion_interaction_t), intent(inout) :: this
1094 type(space_t), intent(in) :: space
1095 type(lattice_vectors_t), intent(in) :: latt
1096 type(atom_t), intent(in) :: atom(:)
1097 integer, intent(in) :: natoms
1098 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
1099 real(real64), intent(out) :: stress_Ewald(3, 3)
1100
1101 real(real64) :: rcut, efourier
1102 integer :: iatom, jatom, idir, jdir
1103 integer :: ix, iy, ix_max, iy_max
1104 real(real64) :: gvec(3), gred(3), gg2, cos_gx, gg_abs, gmax_squared
1105 real(real64) :: factor,factor1,factor2, coeff, e_ewald
1106 real(real64) :: dz_max, z_ij, erfc1, erfc2, diff(3)
1107 real(real64), parameter :: tol = 1e-10_real64
1108
1109 push_sub(ewald_2d_stress)
1110
1111 assert(space%periodic_dim == 2)
1112 assert(space%dim == 3)
1113
1114 stress_ewald = m_zero
1115
1116 ! Searching maximum distance
1117 dz_max = m_zero
1118 do iatom = 1, natoms
1119 do jatom = iatom + 1, natoms
1120 dz_max = max(dz_max, abs(pos(3, iatom) - pos(3, jatom)))
1121 end do
1122 end do
1123
1124 !get a converged value for the cutoff in g
1125 ! Note: to understand these numbers, one needs to look into the energy routine for Ewald 2D
1126 rcut = m_two*this%alpha*4.6_real64 + m_two*this%alpha**2*dz_max
1127 if (dz_max > tol) then ! Else the code here does not work properly
1128 do
1129 if (rcut * dz_max >= m_max_exp_arg) exit !Maximum double precision number
1130 erfc1 = m_one - loct_erf(this%alpha*dz_max + m_half*rcut/this%alpha)
1131 if (erfc1 * exp(rcut*dz_max) < tol) exit
1132 rcut = rcut * 1.414_real64
1133 end do
1134 end if
1135
1136 ! First the G = 0 term
1137 efourier = m_zero
1138 factor = m_pi/latt%rcell_volume
1139 !$omp parallel do private(jatom, z_ij, factor1, factor2) reduction(+:efourier) collapse(2)
1140 do iatom = 1, natoms
1141 do jatom = 1, natoms
1142 z_ij = pos(3, iatom) - pos(3, jatom)
1143
1144 factor1 = z_ij * loct_erf(this%alpha*z_ij)
1145 factor2 = exp(-(this%alpha*z_ij)**2)/(this%alpha*sqrt(m_pi))
1146
1147 efourier = efourier - factor &
1148 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() * (factor1 + factor2)
1149 end do
1150 end do
1151
1152 ! Adding the G=0 term
1153 do idir = 1, 2
1154 stress_ewald(idir, idir) = efourier
1155 end do
1156
1157 ! Getting the G-shell cutoff
1158 ix_max = ceiling(rcut/norm2(latt%klattice(:, 1)))
1159 iy_max = ceiling(rcut/norm2(latt%klattice(:, 2)))
1160 gmax_squared = sum(ix_max*latt%klattice(:, 1)**2)
1161 gmax_squared = min(gmax_squared, sum(iy_max*latt%klattice(:, 2)**2))
1162
1163 !$omp parallel do private(iy, gvec, gg2, gg_abs, factor, iatom, jatom, diff, cos_gx, z_ij, idir, jdir, erfc1, factor1) &
1164 !$omp& private(erfc2, factor2, coeff, e_ewald) &
1165 !$omp& collapse(2) reduction(+:stress_Ewald)
1166 do ix = -ix_max, ix_max
1167 do iy = -iy_max, iy_max
1168
1169 gred = [ix, iy, 0]
1170 call kpoints_to_absolute(latt, gred, gvec)
1171 gg2 = dot_product(gvec,gvec)
1172
1173 ! g=0 must be removed from the sum
1174 if (gg2 < m_epsilon .or. gg2 > gmax_squared*1.001_real64) cycle
1175
1176 gg_abs = sqrt(gg2)
1177 factor = m_fourth*m_pi/(latt%rcell_volume*this%alpha*gg2)
1178
1179 do iatom = 1, natoms
1180 do jatom = iatom, natoms
1181 diff = pos(:, iatom) - pos(:, jatom)
1182 cos_gx = cos(sum(gvec(1:2) * diff(1:2)))
1183 z_ij = diff(3)
1184
1185 factor1 = screening_function_2d(this%alpha, z_ij, gg_abs, erfc1)
1186 factor2 = screening_function_2d(this%alpha,-z_ij, gg_abs, erfc2)
1187
1188 if (iatom == jatom) then
1189 coeff = m_one
1190 else
1191 coeff = m_two
1192 end if
1193
1194 do idir = 1, 2
1195 do jdir = 1, 2
1196 stress_ewald(idir, jdir) = stress_ewald(idir, jdir) &
1197 - factor*gvec(idir)*gvec(jdir) * cos_gx * (factor1 + factor2) * coeff&
1198 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval()
1199 end do
1200 end do
1201
1202 if (abs(erfc1) > m_epsilon) then
1203 factor1 = exp(-gg_abs*z_ij)*erfc1
1204 else
1205 factor1 = m_zero
1206 end if
1207 if (abs(erfc2) > m_epsilon) then
1208 factor2 = exp(gg_abs*z_ij)*erfc2
1209 else
1210 factor2 = m_zero
1211 end if
1212
1213 e_ewald = m_half * m_pi/latt%rcell_volume * coeff &
1214 * atom(iatom)%species%get_zval() * atom(jatom)%species%get_zval() &
1215 * cos_gx / gg_abs * (factor1 + factor2)
1216
1217 do idir = 1, 2
1218 stress_ewald(idir, idir) = stress_ewald(idir, idir) + e_ewald
1219 end do
1220
1221 end do !jatom
1222 end do !iatom
1223 end do !iy
1224 end do !ix
1225
1226 !call comm_allreduce(this%dist%mpi_grp, stress_Ewald)
1227
1228 stress_ewald = stress_ewald / latt%rcell_volume
1229
1230 pop_sub(ewald_2d_stress)
1231 end subroutine ewald_2d_stress
1232
1233 ! ---------------------------------------------------------
1235 real(real64) function screening_function_2d(alpha, z_ij, gg_abs, erfc) result(factor)
1236 real(real64), intent(in) :: alpha
1237 real(real64), intent(in) :: z_ij
1238 real(real64), intent(in) :: gg_abs
1239 real(real64), intent(out) :: erfc
1240
1241 real(real64) :: arg
1242
1243 arg = -alpha*z_ij + m_half*gg_abs/alpha
1244 erfc = m_one - loct_erf(arg)
1245 factor = m_two*alpha*(m_one/gg_abs + z_ij)*erfc - m_two/sqrt(m_pi)*exp(-arg**2)
1246 factor = factor*exp(-gg_abs*z_ij)
1247
1248 end function screening_function_2d
1249
1250 ! ---------------------------------------------------------
1251
1252 subroutine ion_interaction_test(space, latt, atom, natoms, pos, lsize, &
1253 namespace, mc)
1254 class(space_t), intent(in) :: space
1255 type(lattice_vectors_t), intent(in) :: latt
1256 type(atom_t), intent(in) :: atom(:)
1257 integer, intent(in) :: natoms
1258 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
1259 real(real64), intent(in) :: lsize(:)
1260 type(namespace_t), intent(in) :: namespace
1261 type(multicomm_t), intent(in) :: mc
1262
1263 type(ion_interaction_t) :: ion_interaction
1264 real(real64) :: energy
1265 real(real64), allocatable :: force(:, :), force_components(:, :, :)
1266 real(real64) :: energy_components(1:ION_NUM_COMPONENTS)
1267 integer :: iatom, idir
1268
1269 push_sub(ion_interaction_test)
1270
1271 call ion_interaction_init(ion_interaction, namespace, space, natoms)
1272 call ion_interaction_init_parallelization(ion_interaction, natoms, mc)
1273
1274 safe_allocate(force(1:space%dim, 1:natoms))
1275 safe_allocate(force_components(1:space%dim, 1:natoms, 1:ion_num_components))
1276
1277 call ion_interaction_calculate(ion_interaction, space, latt, atom, natoms, pos, lsize, energy, force, &
1278 energy_components = energy_components, force_components = force_components)
1279
1280 call messages_write('Ionic energy =')
1281 call messages_write(energy, fmt = '(f20.10)')
1282 call messages_info(namespace=namespace)
1283
1284 call messages_write('Real space energy =')
1285 call messages_write(energy_components(ion_component_real), fmt = '(f20.10)')
1286 call messages_info(namespace=namespace)
1287
1288 call messages_write('Self energy =')
1289 call messages_write(energy_components(ion_component_self), fmt = '(f20.10)')
1290 call messages_info(namespace=namespace)
1291
1292 call messages_write('Fourier energy =')
1293 call messages_write(energy_components(ion_component_fourier), fmt = '(f20.10)')
1294 call messages_info(namespace=namespace)
1295
1296 call messages_info(namespace=namespace)
1297
1298 do iatom = 1, natoms
1299 call messages_write('Ionic force atom')
1300 call messages_write(iatom)
1301 call messages_write(' =')
1302 do idir = 1, space%dim
1303 call messages_write(force(idir, iatom), fmt = '(f20.10)')
1304 end do
1305 call messages_info(namespace=namespace)
1306
1307 call messages_write('Real space force atom')
1308 call messages_write(iatom)
1309 call messages_write(' =')
1310 do idir = 1, space%dim
1311 call messages_write(force_components(idir, iatom, ion_component_real), fmt = '(f20.10)')
1312 end do
1313 call messages_info(namespace=namespace)
1314
1315 call messages_write('Fourier space force atom')
1316 call messages_write(iatom)
1317 call messages_write(' =')
1318 do idir = 1, space%dim
1319 call messages_write(force_components(idir, iatom, ion_component_fourier), fmt = '(f20.10)')
1320 end do
1321 call messages_info(namespace=namespace)
1322
1323 call messages_info(namespace=namespace)
1324 end do
1325
1326 safe_deallocate_a(force)
1327 safe_deallocate_a(force_components)
1329 call ion_interaction_end(ion_interaction)
1330
1331 pop_sub(ion_interaction_test)
1332 end subroutine ion_interaction_test
1333
1334end module ion_interaction_oct_m
1335
1336!! Local Variables:
1337!! mode: f90
1338!! coding: utf-8
1339!! End:
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.
Definition: atom.F90:292
pure logical function, public any_species_is_jellium_sphere(atom)
Check if any species is a jellium sphere.
Definition: atom.F90:309
type(debug_t), save, public debug
Definition: debug.F90:154
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
Definition: global.F90:189
real(real64), parameter, public m_max_exp_arg
Definition: global.F90:207
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_four
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
real(real64), parameter, public m_fourth
Definition: global.F90:196
complex(real64), parameter, public m_z0
Definition: global.F90:197
complex(real64), parameter, public m_zi
Definition: global.F90:201
real(real64), parameter, public r_min_atom_dist
Minimal distance between two distinguishable atoms.
Definition: global.F90:182
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
real(real64), parameter, public m_three
Definition: global.F90:190
real(real64), parameter, public m_five
Definition: global.F90:192
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, 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)
Definition: kpoints.F90:1031
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
Definition: ps.F90:114
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...
Definition: species.F90:143
int true(void)