TABLE OF CONTENTS
ABINIT/ionion_surface [ Functions ]
NAME
ionion_surface
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
ionicenergyandforces,xred2xcart
SOURCE
38 #if defined HAVE_CONFIG_H 39 #include "config.h" 40 #endif 41 42 #include "abi_common.h" 43 44 subroutine ionion_surface(dtset, eew, grewtn, me, nproc, rprimd, wvl, wvl_den, xred) 45 46 use defs_basis 47 use defs_abitypes 48 use defs_wvltypes 49 use m_profiling_abi 50 #if defined HAVE_BIGDFT 51 use BigDFT_API, only: IonicEnergyandForces 52 #endif 53 54 use m_geometry, only : xred2xcart 55 56 !This section has been created automatically by the script Abilint (TD). 57 !Do not modify the following lines by hand. 58 #undef ABI_FUNC 59 #define ABI_FUNC 'ionion_surface' 60 !End of the abilint section 61 62 implicit none 63 64 !Arguments ------------------------------------ 65 !scalars 66 integer, intent(in) :: me, nproc 67 real(dp),intent(out) :: eew 68 type(dataset_type),intent(in) :: dtset 69 type(wvl_internal_type), intent(in) :: wvl 70 type(wvl_denspot_type), intent(inout) :: wvl_den 71 !arrays 72 real(dp),intent(in) :: rprimd(3,3) 73 real(dp),intent(in) :: xred(3,dtset%natom) 74 real(dp),intent(out) :: grewtn(3,dtset%natom) 75 76 !Local variables------------------------------- 77 !scalars 78 integer :: dispersion, iatom, igeo 79 real(dp) :: psoffset 80 !arrays 81 real(dp),allocatable :: xcart(:,:) 82 real(dp),pointer :: grew_cart(:,:),fdisp(:,:) 83 #if defined HAVE_BIGDFT 84 real(dp) :: edisp 85 real(dp) :: ewaldstr(6) 86 #endif 87 88 ! ************************************************************************* 89 90 !Store xcart for each atom 91 ABI_ALLOCATE(xcart,(3, dtset%natom)) 92 call xred2xcart(dtset%natom, rprimd, xcart, xred) 93 94 nullify(fdisp) 95 nullify(grew_cart) 96 dispersion = 0 97 psoffset = 0._dp 98 #if defined HAVE_BIGDFT 99 call IonicEnergyandForces(me, nproc, wvl_den%denspot%dpbox,& 100 & wvl%atoms, dtset%efield, xcart, & 101 & eew, grew_cart, dispersion, edisp, fdisp,& 102 & ewaldstr,wvl%Glr%d%n1,wvl%Glr%d%n2,wvl%Glr%d%n3,& 103 & wvl_den%denspot%V_ext, wvl_den%denspot%pkernel,psoffset) 104 105 if (associated(fdisp)) then 106 ABI_DEALLOCATE(fdisp) 107 end if 108 #endif 109 110 ABI_DEALLOCATE(xcart) 111 112 !Transform cartesian gradients to reduced gradients. 113 do iatom = 1, dtset%natom, 1 114 do igeo = 1, 3, 1 115 grewtn(igeo, iatom) = -rprimd(1, igeo) * grew_cart(1, iatom) - & 116 & rprimd(2, igeo) * grew_cart(2, iatom) - & 117 & rprimd(3, igeo) * grew_cart(3, iatom) 118 end do 119 end do 120 if (associated(grew_cart)) then 121 ABI_DEALLOCATE(grew_cart) 122 end if 123 124 #if !defined HAVE_BIGDFT 125 if (.false.) write(std_out,*) me,nproc,wvl%h(1),wvl_den%symObj 126 #endif 127 128 end subroutine ionion_surface