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