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 real(real64) :: charge, coeff
426
427 push_sub_with_profile(ewald_short)
428
429 ereal = m_zero
430 ! Origin of this value is not documented
431 rcut = 6.0_real64 / alpha
432 latt_iter = lattice_iterator_t(latt, rcut)
433 natoms = size(atom)
434
435 charge = m_zero
436 do iatom = dist%start, dist%end
437 if (.not. atom(iatom)%species%represents_real_atom()) cycle
438 zi = atom(iatom)%species%get_zval()
439 charge = charge + zi**2
440 end do
441
442 ! Diagonal terms iatom == jatom for all cells, except T=(0,0,0)
443 ! Note: Only half of the copies are needed, by symmetries
444 do icopy = 1, latt_iter%n_cells
445 rnorm = norm2(latt_iter%get(icopy))
446 if (rnorm < r_min_atom_dist) cycle
447 if (rnorm > rcut) cycle
448 erfc = loct_erfc(alpha * rnorm)
449 ereal = ereal + m_half * charge * erfc /rnorm
450 end do
451
452 coeff = m_two * alpha / sqrt(m_pi)
453
454 !$omp parallel default(shared) private(iatom, jatom, zi, zj, icopy, xi, rnorm, erfc, force_real, charge) reduction(+:ereal, force)
455 do iatom = dist%start, dist%end
456 if (.not. atom(iatom)%species%represents_real_atom()) cycle
457 zi = atom(iatom)%species%get_zval()
458
459 ! Upper triangle, for all replica cells
460 do jatom = iatom + 1, natoms
461 zj = atom(jatom)%species%get_zval()
462
463 charge = zi*zj
464
465 !$omp do
466 do icopy = 1, latt_iter%n_cells
467 xi = pos(:, iatom) + latt_iter%get(icopy)
468 rnorm = norm2(xi - pos(:, jatom))
469 if (rnorm > rcut) cycle
470
471 erfc = loct_erfc(alpha * rnorm) / rnorm
472
473 ! Factor 1/2 omitted as one is only summing over upper triangle
474 ereal = ereal + charge * erfc
475
476 force_real(:) = charge * (xi - pos(:, jatom)) * &
477 (erfc + coeff *exp(-(alpha*rnorm)**2)) / rnorm**2
478
479 ! Upper trianglar contribution
480 force(1:space%dim, jatom) = force(1:space%dim, jatom) - force_real
481
482 ! Lower triangular contribution
483 force(1:space%dim, iatom) = force(1:space%dim, iatom) + force_real
484 end do
485 !$omp end do
486
487 end do
488 end do
489 !$omp end parallel
490
491 call comm_allreduce(dist%mpi_grp, ereal)
492 call comm_allreduce(dist%mpi_grp, force)
493
494 pop_sub_with_profile(ewald_short)
495 end subroutine ewald_short
496
501 subroutine ewald_self_interaction(dist, atom, alpha, eself, charge)
502 type(distributed_t), intent(in) :: dist
503 type(atom_t), intent(in) :: atom(:)
504 real(real64), intent(in) :: alpha
505 real(real64), intent(out) :: eself
506 real(real64), intent(out) :: charge
507
508 integer :: iatom
509 real(real64) :: zi
510
511 push_sub(ewald_self_interaction)
512
513 eself = m_zero
514 charge = m_zero
515
516 do iatom = dist%start, dist%end
517 zi = atom(iatom)%species%get_zval()
518 charge = charge + zi
519 eself = eself - alpha / sqrt(m_pi) * zi**2
520 end do
521
522 call comm_allreduce(dist%mpi_grp, eself)
523 call comm_allreduce(dist%mpi_grp, charge)
524
526 end subroutine ewald_self_interaction
527
529 subroutine ewald_long_3d(this, space, latt, atom, natoms, pos, efourier, force, charge)
530 type(ion_interaction_t), intent(in) :: this
531 class(space_t), intent(in) :: space
532 type(lattice_vectors_t), intent(in) :: latt
533 type(atom_t), intent(in) :: atom(:)
534 integer, intent(in) :: natoms
535 real(real64), intent(in) :: pos(:,:)
536 real(real64), intent(inout) :: efourier
537 real(real64), intent(inout) :: force(:, :)
538 real(real64), intent(in) :: charge
539
540 real(real64) :: rcut, gmax_squared
541 integer :: iatom
542 integer :: ix, iy, iz, isph
543 real(real64) :: gvec(3), gred(3), gg2, gx
544 real(real64) :: factor
545 complex(real64) :: sumatoms, tmp(3), aa
546
547 complex(real64), allocatable :: phase(:)
548
549 push_sub(ewald_long_3d)
550
551 assert(space%dim == 3)
552 assert(space%periodic_dim == 3)
553
554 ! And the long-range part, using an Ewald sum
555 safe_allocate(phase(1:natoms))
556
557 ! get a converged value for the cutoff in g
558 rcut = sqrt(minval(sum(latt%klattice**2, dim=1)))
559
560 ! 9.5 is a constant that controls the range separation
561 isph = ceiling(9.5_real64*this%alpha/rcut)
562
563 ! First the G = 0 term (charge was calculated previously)
564 efourier = -m_pi*charge**2/(m_two*this%alpha**2*latt%rcell_volume)
565
566 ! Getting the G-shell cutoff
567 gmax_squared = isph**2 * minval(sum(latt%klattice**2, dim=1))
568
569 do ix = -isph, isph
570 do iy = -isph, isph
571 do iz = -isph, isph
572
573 ! Exploit k <-> -k symmetry
574 ! Only process one half of reciprocal space.
575 ! g=0 must also be removed from the sum
576 if (ix < 0) cycle
577 if (ix == 0 .and. iy < 0) cycle
578 if (ix == 0 .and. iy == 0 .and. iz <= 0) cycle
579
580 gred = [ix, iy, iz]
581 call kpoints_to_absolute(latt, gred, gvec)
582 gg2 = dot_product(gvec, gvec)
583
584 if (gg2 > gmax_squared*1.001_real64) cycle
585
586 gx = -0.25_real64*gg2/this%alpha**2
587
588 if (gx < -36.0_real64) cycle
589
590 ! We have used the k-> -k symmetry, hence the factor 4
591 factor = m_four*m_pi/latt%rcell_volume*exp(gx)/gg2
592
593 if (factor < epsilon(factor)) cycle
595 sumatoms = m_z0
596 !$omp parallel do private(iatom, gx, aa) reduction(+:sumatoms)
597 do iatom = 1, natoms
598 gx = sum(gvec*pos(:,iatom))
599 aa = atom(iatom)%species%get_zval()*cmplx(cos(gx), sin(gx), real64)
600 phase(iatom) = aa
601 sumatoms = sumatoms + aa
602 end do
603
604 efourier = efourier + factor * real(sumatoms*conjg(sumatoms), real64)
605
606 do iatom = 1, natoms
607 tmp = m_zi*gvec*phase(iatom)
608 force(1:space%dim, iatom) = force(1:space%dim, iatom) - factor*real(conjg(tmp)*sumatoms + tmp*conjg(sumatoms), real64)
609
610 end do
611
612 end do
613 end do
614 end do
615
616 safe_deallocate_a(phase)
617
618 pop_sub(ewald_long_3d)
619
620 end subroutine ewald_long_3d
621
625 subroutine ewald_long_2d(this, space, latt, atom, natoms, pos, efourier, force)
626 type(ion_interaction_t), intent(in) :: this
627 class(space_t), intent(in) :: space
628 type(lattice_vectors_t), intent(in) :: latt
629 type(atom_t), intent(in) :: atom(:)
630 integer, intent(in) :: natoms
631 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
632 real(real64), intent(inout) :: efourier
633 real(real64), intent(inout) :: force(:, :)
634
635 real(real64) :: rcut, gmax_squared
636 integer :: iatom, jatom
637 integer :: ix, iy, ix_max, iy_max
638 real(real64) :: gvec(space%dim), gg2, gx, gg_abs
639 real(real64) :: factor,factor1,factor2, coeff
640 real(real64) :: dz_max, dz_ij, erfc1, erfc2, tmp_erf
641 real(real64), allocatable :: force_tmp(:,:)
642 real(real64), parameter :: tol = 1e-10_real64
643
644 push_sub(ewald_long_2d)
645
646 assert(space%periodic_dim == 2)
647 assert(space%dim == 2 .or. space%dim == 3)
648
649 ! And the long-range part, using an Ewald sum
650
651 ! Searching maximum distance
652 if (space%dim == 3) then
653 dz_max = m_zero
654 do iatom = 1, natoms
655 do jatom = iatom + 1, natoms
656 dz_max = max(dz_max, abs(pos(3, iatom) - pos(3, jatom)))
657 end do
658 end do
659 else
660 ! For a 2D system, all atoms are on the plane, so the distance is zero
661 dz_max = m_zero
662 end if
663
664 !get a converged value for the cutoff in g
665 rcut = m_two*this%alpha*4.6_real64 + m_two*this%alpha**2*dz_max
666 if (dz_max > tol) then
667 do
668 if (rcut * dz_max >= m_max_exp_arg) exit !Maximum double precision number
669 erfc1 = m_one - loct_erf(this%alpha*dz_max + m_half*rcut/this%alpha)
670 if (erfc1 * exp(rcut*dz_max) < 1.e-10_real64) exit
671 rcut = rcut * 1.414_real64
672 end do
673 end if
674
675 ix_max = ceiling(rcut/norm2(latt%klattice(:, 1)))
676 iy_max = ceiling(rcut/norm2(latt%klattice(:, 2)))
677
678 safe_allocate(force_tmp(1:space%dim, 1:natoms))
679 force_tmp = m_zero
680
681 ! First the G = 0 term
682 efourier = m_zero
683 factor = m_pi/latt%rcell_volume
684 !$omp parallel do private(jatom, dz_ij, tmp_erf, factor1, factor2) reduction(+:efourier,force_tmp) &
685 !$omp& collapse(2)
686 do iatom = this%dist%start, this%dist%end
687 do jatom = 1, natoms
688 ! efourier
689 if (space%dim == 3) then
690 dz_ij = pos(3, iatom) - pos(3, jatom)
691 else
692 dz_ij = m_zero
693 end if
694
695 tmp_erf = loct_erf(this%alpha*dz_ij)
696 factor1 = dz_ij*tmp_erf
697 factor2 = exp(-(this%alpha*dz_ij)**2)/(this%alpha*sqrt(m_pi))
698
699 efourier = efourier - factor &
700 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() * (factor1 + factor2)
701
702 ! force
703 if (iatom == jatom)cycle
704 if (abs(tmp_erf) < m_epsilon) cycle
705
706 if (space%dim == 3) then
707 force_tmp(3, iatom) = force_tmp(3, iatom) - (- m_two*factor) &
708 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() * tmp_erf
709 end if
710
711 end do
712 end do
713
714 ! Getting the G-shell cutoff
715 gmax_squared = sum(ix_max*latt%klattice(:, 1)**2)
716 gmax_squared = min(gmax_squared, sum(iy_max*latt%klattice(:, 2)**2))
717
718 !$omp parallel do private(iy, gvec, gg2, gg_abs, factor, iatom, jatom, gx, dz_ij, erfc1, factor1, erfc2, factor2, coeff) &
719 !$omp& collapse(2) reduction(+:efourier, force_tmp)
720 do ix = -ix_max, ix_max
721 do iy = -iy_max, iy_max
722
723 gvec = ix*latt%klattice(:, 1) + iy*latt%klattice(:, 2)
724 gg2 = sum(gvec**2)
725
726 ! g=0 must be removed from the sum
727 if (gg2 < m_epsilon .or. gg2 > gmax_squared*1.001_real64) cycle
728 gg_abs = sqrt(gg2)
729 factor = m_half*m_pi/(latt%rcell_volume*gg_abs)
730
731 do iatom = this%dist%start, this%dist%end
732 do jatom = iatom, natoms
733 ! efourier
734 gx = sum(gvec(1:2) * (pos(1:2, iatom) - pos(1:2, jatom)))
735 gx = gvec(1)*(pos(1, iatom) - pos(1, jatom)) + gvec(2)*(pos(2, iatom) - pos(2, jatom))
736 if (space%dim == 3) then
737 dz_ij = pos(3, iatom) - pos(3, jatom)
738 else
739 dz_ij = m_zero
740 end if
741
742 erfc1 = m_one - loct_erf(this%alpha*dz_ij + m_half*gg_abs/this%alpha)
743 if (abs(erfc1) > m_epsilon) then
744 factor1 = exp(gg_abs*dz_ij)*erfc1
745 else
746 factor1 = m_zero
747 end if
748 erfc2 = m_one - loct_erf(-this%alpha*dz_ij + m_half*gg_abs/this%alpha)
749 if (abs(erfc2) > m_epsilon) then
750 factor2 = exp(-gg_abs*dz_ij)*erfc2
751 else
752 factor2 = m_zero
753 end if
754
755 if (iatom == jatom) then
756 coeff = m_one
757 else
758 coeff = m_two
759 end if
760
761 efourier = efourier &
762 + factor * coeff &
763 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
764 * cos(gx)* ( factor1 + factor2)
765
766 ! force
767 if (iatom == jatom) cycle
768
769 force_tmp(1:2, iatom) = force_tmp(1:2, iatom) &
770 + m_two * factor * gvec(1:2) &
771 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
772 *sin(gx)*(factor1 + factor2)
773
774 force_tmp(1:2, jatom) = force_tmp(1:2, jatom) &
775 - m_two * factor * gvec(1:2) &
776 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
777 *sin(gx)*(factor1 + factor2)
778
779 factor1 = gg_abs*erfc1 &
780 - m_two*this%alpha/sqrt(m_pi)*exp(-(this%alpha*dz_ij + m_half*gg_abs/this%alpha)**2)
781 if (abs(factor1) > m_epsilon) then
782 factor1 = factor1*exp(gg_abs*dz_ij)
783 else
784 factor1 = m_zero
785 end if
786
787 factor2 = gg_abs*erfc2 &
788 - m_two*this%alpha/sqrt(m_pi)*exp(-(-this%alpha*dz_ij + m_half*gg_abs/this%alpha)**2)
789 if (abs(factor2) > m_epsilon) then
790 factor2 = factor2*exp(-gg_abs*dz_ij)
791 else
792 factor2 = m_zero
793 end if
794
795 if (space%dim == 3) then
796 force_tmp(3, iatom) = force_tmp(3, iatom) &
797 - m_two*factor &
798 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
799 * cos(gx)* ( factor1 - factor2)
800 force_tmp(3, jatom) = force_tmp(3, jatom) &
801 + m_two*factor &
802 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() &
803 * cos(gx)* ( factor1 - factor2)
804 end if
805
806 end do
807 end do
808
809
810 end do
811 end do
812
813 call comm_allreduce(this%dist%mpi_grp, efourier)
814 call comm_allreduce(this%dist%mpi_grp, force_tmp)
815
816 force = force + force_tmp
817
818 safe_deallocate_a(force_tmp)
819
820 pop_sub(ewald_long_2d)
821 end subroutine ewald_long_2d
822
823 !TODO(Alex/Nicolas) Issue #950. Refactor: Move G=0 correction from ion-ion energy to pseudopotential energy
830 subroutine pseudopotential_correction_3d(dist, latt, atom, charge, epseudo)
831 type(distributed_t), intent(in) :: dist
832 type(lattice_vectors_t), intent(in) :: latt
833 type(atom_t), intent(in) :: atom(:)
834 real(real64), intent(out) :: epseudo
835
836 real(real64) :: zi
837 real(real64) :: charge
838 integer :: iatom
839
841
842 epseudo = m_zero
843 do iatom = dist%start, dist%end
844 select type(spec => atom(iatom)%species)
845 class is(pseudopotential_t)
846 zi = spec%get_zval()
847 epseudo = epseudo + m_pi *zi * &
848 (spec%ps%sigma_erf * sqrt(m_two))**2 / latt%rcell_volume * charge
849 end select
850 end do
851 call comm_allreduce(dist%mpi_grp, epseudo)
852
854
855 end subroutine pseudopotential_correction_3d
856
858 subroutine ion_interaction_stress(this, space, latt, atom, natoms, pos, stress_ii)
859 type(ion_interaction_t), intent(inout) :: this
860 class(space_t), intent(in) :: space
861 type(lattice_vectors_t), intent(in) :: latt
862 type(atom_t), intent(in) :: atom(:)
863 integer, intent(in) :: natoms
864 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
865 real(real64), intent(out) :: stress_ii(space%dim, space%dim)
866
867 real(real64) :: stress_short(1:space%dim, 1:space%dim), stress_Ewald(1:space%dim, 1:space%dim)
868
869 push_sub(ion_interaction_stress)
870
871 stress_ii = m_zero
872
873 ! Only implemented in the periodic case
874 assert(space%is_periodic())
875
876 ! Short range part in real space
877 call ion_interaction_stress_short(this, space, latt, atom, natoms, pos, stress_short)
878
879 ! Long range part in Fourier space
880 select case(space%periodic_dim)
881 case(3)
882 call ewald_3d_stress(this, space, latt, atom, natoms, pos, stress_ewald)
883 case(2)
884 call ewald_2d_stress(this, space, latt, atom, natoms, pos, stress_ewald)
885 case default
886 assert(.false.)
887 end select
888
889 stress_ii = stress_short + stress_ewald
890
892 end subroutine ion_interaction_stress
893
894 ! ---------------------------------------------------------
910
911 subroutine ion_interaction_stress_short(this, space, latt, atom, natoms, pos, stress_short)
912 type(ion_interaction_t), intent(inout) :: this
913 class(space_t), intent(in) :: space
914 type(lattice_vectors_t), intent(in) :: latt
915 type(atom_t), intent(in) :: atom(:)
916 integer, intent(in) :: natoms
917 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
918 real(real64), intent(out) :: stress_short(1:space%dim, 1:space%dim)
919
920 real(real64) :: xi(space%dim)
921 real(real64) :: r_ij, zi, zj, erfc, Hp, factor
922 integer :: iatom, jatom, icopy, idir, jdir
923 real(real64) :: alpha, rcut
924 type(lattice_iterator_t) :: latt_iter
925
927 call profiling_in("ION_ION_STRESS_SHORT")
928
929 ! Only implemented in the periodic case
930 assert(space%is_periodic())
931
932 alpha = this%alpha
933
934 ! See the code for the energy above to understand this parameter
935 rcut = 6.0_real64/alpha
936
937 ! the short-range part is calculated directly
938 stress_short = m_zero
939 latt_iter = lattice_iterator_t(latt, rcut)
940
941 do iatom = this%dist%start, this%dist%end
942 select type(spec => atom(iatom)%species)
943 class is(jellium_t)
944 cycle
945 end select
946 zi = atom(iatom)%species%get_zval()
947
948 do icopy = 1, latt_iter%n_cells
949 xi = pos(:, iatom) + latt_iter%get(icopy)
950
951 do jatom = 1, natoms
952 zj = atom(jatom)%species%get_zval()
953 r_ij = norm2(xi - pos(:, jatom))
954
955 if (r_ij < r_min_atom_dist) cycle
956
957 erfc = loct_erfc(alpha*r_ij)
958 hp = -m_two/sqrt(m_pi)*exp(-(alpha*r_ij)**2) - erfc/(alpha*r_ij)
959 factor = m_half*zj*zi*alpha*hp
960 do idir = 1, space%periodic_dim
961 do jdir = 1, space%periodic_dim
962 stress_short(idir, jdir) = stress_short(idir, jdir) &
963 - factor*(xi(idir) - pos(idir, jatom))*(xi(jdir) - pos(jdir, jatom))/(r_ij**2)
964 end do
965 end do
966
967 end do
968 end do
969 end do
970
971 if (this%dist%parallel) then
972 call comm_allreduce(this%dist%mpi_grp, stress_short)
973 end if
974
975 stress_short = stress_short/latt%rcell_volume
976
977 call profiling_out("ION_ION_STRESS_SHORT")
978
980 end subroutine ion_interaction_stress_short
981
982
983
984 ! ---------------------------------------------------------
996 subroutine ewald_3d_stress(this, space, latt, atom, natoms, pos, stress_Ewald)
997 type(ion_interaction_t), intent(inout) :: this
998 class(space_t), intent(in) :: space
999 type(lattice_vectors_t), intent(in) :: latt
1000 type(atom_t), intent(in) :: atom(:)
1001 integer, intent(in) :: natoms
1002 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
1003 real(real64), intent(out) :: stress_ewald(3, 3)
1005 real(real64) :: zi, rcut, gmax_squared
1006 integer :: iatom
1007 integer :: ix, iy, iz, isph, idim, idir, jdir
1008 real(real64) :: gred(3), gvec(3), gg2, gx
1009 real(real64) :: factor, charge, charge_sq, off_diagonal_weight
1010 complex(real64) :: sumatoms, aa
1011
1012 call profiling_in("STRESS_3D_EWALD")
1013 push_sub(ewald_3d_stress)
1014
1015 ! Currently this is only implemented for 3D
1016 assert(space%dim == 3)
1017 assert(space%periodic_dim == 3) ! Not working for mixed periodicity
1018 ! (klattice along the non-periodic directions is wrong)
1019 ! Anyway gg/gg2 is not correct for mixed periodicity
1020
1021 stress_ewald = m_zero
1022
1023 ! And the long-range part, using an Ewald sum
1024 charge = m_zero
1025 charge_sq = m_zero
1026 do iatom = 1, natoms
1027 zi = atom(iatom)%species%get_zval()
1028 charge = charge + zi
1029 charge_sq = charge_sq + zi**2
1030 end do
1031
1032 ! get a converged value for the cutoff in g
1033 rcut = huge(rcut)
1034 do idim = 1, space%periodic_dim
1035 rcut = min(rcut, sum(latt%klattice(1:space%periodic_dim, idim)**2))
1036 end do
1037
1038 rcut = sqrt(rcut)
1039
1040 isph = ceiling(9.5_real64*this%alpha/rcut)
1041
1042 ! Getting the G-shell cutoff
1043 gmax_squared = isph**2 * minval(sum(latt%klattice**2, dim=1))
1044
1045 do ix = -isph, isph
1046 do iy = -isph, isph
1047 do iz = -isph, isph
1048
1049 ! Exploit k <-> -k symmetry
1050 ! Only process one half of reciprocal space.
1051 ! g=0 must also be removed from the sum
1052 if (ix < 0) cycle
1053 if (ix == 0 .and. iy < 0) cycle
1054 if (ix == 0 .and. iy == 0 .and. iz <= 0) cycle
1055
1056 gred = [ix, iy, iz]
1057 call kpoints_to_absolute(latt, gred, gvec)
1058 gg2 = sum(gvec**2)
1059
1060 ! g=0 must be removed from the sum
1061 if (gg2 > gmax_squared*1.001_real64) cycle
1062
1063 gx = -0.25_real64*gg2/this%alpha**2
1064
1065 if (gx < -36.0_real64) cycle
1066
1067 ! We have used the k-> -k symmetry, hence the factor 4
1068 factor = m_four*m_pi*exp(gx)/(latt%rcell_volume*gg2)
1069
1070 if (factor < epsilon(factor)) cycle
1071
1072 sumatoms = m_z0
1073
1074 do iatom = 1, natoms
1075 gx = sum(gvec*pos(:, iatom))
1076 aa = atom(iatom)%species%get_zval()*cmplx(cos(gx), sin(gx), real64)
1077 sumatoms = sumatoms + aa
1078 end do
1079
1080 factor = factor*abs(sumatoms)**2
1081 off_diagonal_weight = - m_two*factor/gg2*(0.25_real64*gg2/this%alpha**2+m_one)
1082
1083 do idir = 1, 3
1084 do jdir = 1, 3
1085 stress_ewald(idir, jdir) = stress_ewald(idir, jdir) &
1086 + gvec(idir) * gvec(jdir) * off_diagonal_weight
1087 end do
1088 stress_ewald(idir, idir) = stress_ewald(idir, idir) + factor
1089 end do
1090
1091 end do
1092 end do
1093 end do
1094
1095
1096 ! The G = 0 term of the Ewald summation
1097 factor = m_half*m_pi*charge**2/(latt%rcell_volume*this%alpha**2)
1098 do idir = 1,3
1099 stress_ewald(idir,idir) = stress_ewald(idir,idir) - factor
1100 end do
1101
1102 stress_ewald = stress_ewald / latt%rcell_volume
1103
1104
1105 call profiling_out("STRESS_3D_EWALD")
1106 pop_sub(ewald_3d_stress)
1107
1108 end subroutine ewald_3d_stress
1109
1110 ! ---------------------------------------------------------
1123 subroutine ewald_2d_stress(this, space, latt, atom, natoms, pos, stress_Ewald)
1124 type(ion_interaction_t), intent(inout) :: this
1125 type(space_t), intent(in) :: space
1126 type(lattice_vectors_t), intent(in) :: latt
1127 type(atom_t), intent(in) :: atom(:)
1128 integer, intent(in) :: natoms
1129 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
1130 real(real64), intent(out) :: stress_Ewald(3, 3)
1131
1132 real(real64) :: rcut, efourier
1133 integer :: iatom, jatom, idir, jdir
1134 integer :: ix, iy, ix_max, iy_max
1135 real(real64) :: gvec(3), gred(3), gg2, cos_gx, gg_abs, gmax_squared
1136 real(real64) :: factor,factor1,factor2, coeff, e_ewald
1137 real(real64) :: dz_max, z_ij, erfc1, erfc2, diff(3)
1138 real(real64), parameter :: tol = 1e-10_real64
1139
1140 push_sub(ewald_2d_stress)
1141
1142 assert(space%periodic_dim == 2)
1143 assert(space%dim == 3)
1144
1145 stress_ewald = m_zero
1146
1147 ! Searching maximum distance
1148 dz_max = m_zero
1149 do iatom = 1, natoms
1150 do jatom = iatom + 1, natoms
1151 dz_max = max(dz_max, abs(pos(3, iatom) - pos(3, jatom)))
1152 end do
1153 end do
1154
1155 !get a converged value for the cutoff in g
1156 ! Note: to understand these numbers, one needs to look into the energy routine for Ewald 2D
1157 rcut = m_two*this%alpha*4.6_real64 + m_two*this%alpha**2*dz_max
1158 if (dz_max > tol) then ! Else the code here does not work properly
1159 do
1160 if (rcut * dz_max >= m_max_exp_arg) exit !Maximum double precision number
1161 erfc1 = m_one - loct_erf(this%alpha*dz_max + m_half*rcut/this%alpha)
1162 if (erfc1 * exp(rcut*dz_max) < tol) exit
1163 rcut = rcut * 1.414_real64
1164 end do
1165 end if
1166
1167 ! First the G = 0 term
1168 efourier = m_zero
1169 factor = m_pi/latt%rcell_volume
1170 !$omp parallel do private(jatom, z_ij, factor1, factor2) reduction(+:efourier) collapse(2)
1171 do iatom = 1, natoms
1172 do jatom = 1, natoms
1173 z_ij = pos(3, iatom) - pos(3, jatom)
1174
1175 factor1 = z_ij * loct_erf(this%alpha*z_ij)
1176 factor2 = exp(-(this%alpha*z_ij)**2)/(this%alpha*sqrt(m_pi))
1177
1178 efourier = efourier - factor &
1179 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval() * (factor1 + factor2)
1180 end do
1181 end do
1182
1183 ! Adding the G=0 term
1184 do idir = 1, 2
1185 stress_ewald(idir, idir) = efourier
1186 end do
1187
1188 ! Getting the G-shell cutoff
1189 ix_max = ceiling(rcut/norm2(latt%klattice(:, 1)))
1190 iy_max = ceiling(rcut/norm2(latt%klattice(:, 2)))
1191 gmax_squared = sum(ix_max*latt%klattice(:, 1)**2)
1192 gmax_squared = min(gmax_squared, sum(iy_max*latt%klattice(:, 2)**2))
1193
1194 !$omp parallel do private(iy, gvec, gg2, gg_abs, factor, iatom, jatom, diff, cos_gx, z_ij, idir, jdir, erfc1, factor1) &
1195 !$omp& private(erfc2, factor2, coeff, e_ewald) &
1196 !$omp& collapse(2) reduction(+:stress_Ewald)
1197 do ix = -ix_max, ix_max
1198 do iy = -iy_max, iy_max
1199
1200 gred = [ix, iy, 0]
1201 call kpoints_to_absolute(latt, gred, gvec)
1202 gg2 = dot_product(gvec,gvec)
1203
1204 ! g=0 must be removed from the sum
1205 if (gg2 < m_epsilon .or. gg2 > gmax_squared*1.001_real64) cycle
1206
1207 gg_abs = sqrt(gg2)
1208 factor = m_fourth*m_pi/(latt%rcell_volume*this%alpha*gg2)
1209
1210 do iatom = 1, natoms
1211 do jatom = iatom, natoms
1212 diff = pos(:, iatom) - pos(:, jatom)
1213 cos_gx = cos(sum(gvec(1:2) * diff(1:2)))
1214 z_ij = diff(3)
1215
1216 factor1 = screening_function_2d(this%alpha, z_ij, gg_abs, erfc1)
1217 factor2 = screening_function_2d(this%alpha,-z_ij, gg_abs, erfc2)
1218
1219 if (iatom == jatom) then
1220 coeff = m_one
1221 else
1222 coeff = m_two
1223 end if
1224
1225 do idir = 1, 2
1226 do jdir = 1, 2
1227 stress_ewald(idir, jdir) = stress_ewald(idir, jdir) &
1228 - factor*gvec(idir)*gvec(jdir) * cos_gx * (factor1 + factor2) * coeff&
1229 * atom(iatom)%species%get_zval()*atom(jatom)%species%get_zval()
1230 end do
1231 end do
1232
1233 if (abs(erfc1) > m_epsilon) then
1234 factor1 = exp(-gg_abs*z_ij)*erfc1
1235 else
1236 factor1 = m_zero
1237 end if
1238 if (abs(erfc2) > m_epsilon) then
1239 factor2 = exp(gg_abs*z_ij)*erfc2
1240 else
1241 factor2 = m_zero
1242 end if
1243
1244 e_ewald = m_half * m_pi/latt%rcell_volume * coeff &
1245 * atom(iatom)%species%get_zval() * atom(jatom)%species%get_zval() &
1246 * cos_gx / gg_abs * (factor1 + factor2)
1247
1248 do idir = 1, 2
1249 stress_ewald(idir, idir) = stress_ewald(idir, idir) + e_ewald
1250 end do
1251
1252 end do !jatom
1253 end do !iatom
1254 end do !iy
1255 end do !ix
1256
1257 !call comm_allreduce(this%dist%mpi_grp, stress_Ewald)
1258
1259 stress_ewald = stress_ewald / latt%rcell_volume
1260
1261 pop_sub(ewald_2d_stress)
1262 end subroutine ewald_2d_stress
1263
1264 ! ---------------------------------------------------------
1266 real(real64) function screening_function_2d(alpha, z_ij, gg_abs, erfc) result(factor)
1267 real(real64), intent(in) :: alpha
1268 real(real64), intent(in) :: z_ij
1269 real(real64), intent(in) :: gg_abs
1270 real(real64), intent(out) :: erfc
1271
1272 real(real64) :: arg
1273
1274 arg = -alpha*z_ij + m_half*gg_abs/alpha
1275 erfc = m_one - loct_erf(arg)
1276 factor = m_two*alpha*(m_one/gg_abs + z_ij)*erfc - m_two/sqrt(m_pi)*exp(-arg**2)
1277 factor = factor*exp(-gg_abs*z_ij)
1278
1279 end function screening_function_2d
1280
1281 ! ---------------------------------------------------------
1282
1283 subroutine ion_interaction_test(space, latt, atom, natoms, pos, lsize, &
1284 namespace, mc)
1285 class(space_t), intent(in) :: space
1286 type(lattice_vectors_t), intent(in) :: latt
1287 type(atom_t), intent(in) :: atom(:)
1288 integer, intent(in) :: natoms
1289 real(real64), intent(in) :: pos(1:space%dim,1:natoms)
1290 real(real64), intent(in) :: lsize(:)
1291 type(namespace_t), intent(in) :: namespace
1292 type(multicomm_t), intent(in) :: mc
1293
1294 type(ion_interaction_t) :: ion_interaction
1295 real(real64) :: energy
1296 real(real64), allocatable :: force(:, :), force_components(:, :, :)
1297 real(real64) :: energy_components(1:ION_NUM_COMPONENTS)
1298 integer :: iatom, idir
1299
1300 push_sub(ion_interaction_test)
1301
1302 call ion_interaction_init(ion_interaction, namespace, space, natoms)
1303 call ion_interaction_init_parallelization(ion_interaction, natoms, mc)
1304
1305 safe_allocate(force(1:space%dim, 1:natoms))
1306 safe_allocate(force_components(1:space%dim, 1:natoms, 1:ion_num_components))
1307
1308 call ion_interaction_calculate(ion_interaction, space, latt, atom, natoms, pos, lsize, energy, force, &
1309 energy_components = energy_components, force_components = force_components)
1310
1311 call messages_write('Ionic energy =')
1312 call messages_write(energy, fmt = '(f20.10)')
1313 call messages_info(namespace=namespace)
1314
1315 call messages_write('Real space energy =')
1316 call messages_write(energy_components(ion_component_real), fmt = '(f20.10)')
1317 call messages_info(namespace=namespace)
1318
1319 call messages_write('Self energy =')
1320 call messages_write(energy_components(ion_component_self), fmt = '(f20.10)')
1321 call messages_info(namespace=namespace)
1322
1323 call messages_write('Fourier energy =')
1324 call messages_write(energy_components(ion_component_fourier), fmt = '(f20.10)')
1325 call messages_info(namespace=namespace)
1326
1327 call messages_info(namespace=namespace)
1328
1329 do iatom = 1, natoms
1330 call messages_write('Ionic force atom')
1331 call messages_write(iatom)
1332 call messages_write(' =')
1333 do idir = 1, space%dim
1334 call messages_write(force(idir, iatom), fmt = '(f20.10)')
1335 end do
1336 call messages_info(namespace=namespace)
1337
1338 call messages_write('Real space force atom')
1339 call messages_write(iatom)
1340 call messages_write(' =')
1341 do idir = 1, space%dim
1342 call messages_write(force_components(idir, iatom, ion_component_real), fmt = '(f20.10)')
1343 end do
1344 call messages_info(namespace=namespace)
1345
1346 call messages_write('Fourier space force atom')
1347 call messages_write(iatom)
1348 call messages_write(' =')
1349 do idir = 1, space%dim
1350 call messages_write(force_components(idir, iatom, ion_component_fourier), fmt = '(f20.10)')
1351 end do
1352 call messages_info(namespace=namespace)
1353
1354 call messages_info(namespace=namespace)
1355 end do
1356
1357 safe_deallocate_a(force)
1358 safe_deallocate_a(force_components)
1360 call ion_interaction_end(ion_interaction)
1361
1362 pop_sub(ion_interaction_test)
1363 end subroutine ion_interaction_test
1364
1365end module ion_interaction_oct_m
1366
1367!! Local Variables:
1368!! mode: f90
1369!! coding: utf-8
1370!! 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:156
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:190
real(real64), parameter, public m_max_exp_arg
Definition: global.F90:208
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
real(real64), parameter, public m_fourth
Definition: global.F90:197
complex(real64), parameter, public m_z0
Definition: global.F90:198
complex(real64), parameter, public m_zi
Definition: global.F90:202
real(real64), parameter, public r_min_atom_dist
Minimal distance between two distinguishable atoms.
Definition: global.F90:183
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
real(real64), parameter, public m_three
Definition: global.F90:191
real(real64), parameter, public m_five
Definition: global.F90:193
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)
Computes the long-range part of the 2D Ewald summation.
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:1113
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:537
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
static double f(double w, void *p)
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)