TABLE OF CONTENTS
ABINIT/ionion_realspace [ Functions ]
NAME
ionion_realspace
FUNCTION
Compute the ion/ion interaction energies and forces in real space case. Use ewald() instead if computations are done in reciprocal space since it also includes the correction for the shift done in potentials calculations and includes replica interactions.
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
dtset <type(dataset_type)>=all input variables in this dataset rmet(3,3)=metric tensor in real space (bohr^2) xred(3,natom)=relative coords of atoms in unit cell (dimensionless) zion(ntypat)=charge on each type of atom (real number)
OUTPUT
eew=final ion/ion energy in hartrees grewtn(3,natom)=grads of ion/ion wrt xred(3,natom), hartrees.
PARENTS
setvtr
CHILDREN
xred2xcart
SOURCE
38 #if defined HAVE_CONFIG_H 39 #include "config.h" 40 #endif 41 42 #include "abi_common.h" 43 44 subroutine ionion_realSpace(dtset, eew, grewtn, rprimd, xred, zion) 45 46 use defs_basis 47 use defs_abitypes 48 use m_profiling_abi 49 50 use m_geometry, only : xred2xcart 51 52 !This section has been created automatically by the script Abilint (TD). 53 !Do not modify the following lines by hand. 54 #undef ABI_FUNC 55 #define ABI_FUNC 'ionion_realSpace' 56 !End of the abilint section 57 58 implicit none 59 60 !Arguments ------------------------------------ 61 !scalars 62 real(dp),intent(out) :: eew 63 type(dataset_type),intent(in) :: dtset 64 !arrays 65 real(dp),intent(in) :: rprimd(3,3),zion(dtset%ntypat) 66 real(dp),intent(in) :: xred(3,dtset%natom) 67 real(dp),intent(out) :: grewtn(3,dtset%natom) 68 69 !Local variables------------------------------- 70 !scalars 71 integer :: ia1,ia2,iatom,igeo 72 real(dp) :: r 73 !arrays 74 real(dp),allocatable :: grew_cart(:,:),xcart(:,:) 75 76 ! ************************************************************************* 77 78 !Store xcart for each atom 79 ABI_ALLOCATE(xcart,(3, dtset%natom)) 80 call xred2xcart(dtset%natom, rprimd, xcart, xred) 81 82 !Summing the interaction between ions. 83 eew = 0._dp 84 do ia1 = 1, dtset%natom, 1 85 do ia2 = ia1 + 1, dtset%natom, 1 86 r = sqrt((xcart(1, ia1) - xcart(1, ia2)) ** 2 + & 87 & (xcart(2, ia1) - xcart(2, ia2)) ** 2 + & 88 & (xcart(3, ia1) - xcart(3, ia2)) ** 2) 89 eew = eew + zion(dtset%typat(ia1)) * zion(dtset%typat(ia2)) / r 90 end do 91 end do 92 93 !Allocate temporary array to store cartesian gradients. 94 ABI_ALLOCATE(grew_cart,(3, dtset%natom)) 95 96 !Summing the forces for each atom 97 do ia1 = 1, dtset%natom, 1 98 grew_cart(:, ia1) = 0._dp 99 do ia2 = 1, dtset%natom, 1 100 if (ia1 /= ia2) then 101 r = (xcart(1, ia1) - xcart(1, ia2)) ** 2 + & 102 & (xcart(2, ia1) - xcart(2, ia2)) ** 2 + & 103 & (xcart(3, ia1) - xcart(3, ia2)) ** 2 104 do igeo = 1, 3, 1 105 grew_cart(igeo, ia1) = grew_cart(igeo, ia1) - (xcart(igeo, ia1) - xcart(igeo, ia2)) * & 106 & zion(dtset%typat(ia1)) * zion(dtset%typat(ia2)) / (r ** 1.5_dp) 107 end do 108 end if 109 end do 110 end do 111 112 ABI_DEALLOCATE(xcart) 113 114 !Transform cartesian gradients to reduced gradients. 115 do iatom = 1, dtset%natom, 1 116 do igeo = 1, 3, 1 117 grewtn(igeo, iatom) = rprimd(1, igeo) * grew_cart(1, iatom) + & 118 & rprimd(2, igeo) * grew_cart(2, iatom) + & 119 & rprimd(3, igeo) * grew_cart(3, iatom) 120 end do 121 end do 122 ABI_DEALLOCATE(grew_cart) 123 124 end subroutine ionion_realSpace