TABLE OF CONTENTS
ABINIT/pawpolev [ Functions ]
NAME
pawpolev
FUNCTION
Compute the PAW term for polarization, named expected value term
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (FJ, PH) 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
comm_atom=--optional-- MPI communicator over atoms my_natom=number of atoms treated by current processor natom=number of atoms in cell. ntypat = number of atom types pawrhoij(natom) <type(pawrhoij_type)>= paw rhoij occupancies and related data pawtab(ntypat) <type(pawtab_type)>=paw tabulated starting data
OUTPUT
pelev(3)= electronic polarisation. expectation value term (PAW only)
PARENTS
berryphase_new
CHILDREN
timab,xmpi_sum
SOURCE
35 #if defined HAVE_CONFIG_H 36 #include "config.h" 37 #endif 38 39 #include "abi_common.h" 40 41 subroutine pawpolev(my_natom,natom,ntypat,pawrhoij,pawtab,pelev,& 42 & comm_atom) ! optional argument (parallelism) 43 44 use defs_basis 45 use m_profiling_abi 46 use m_errors 47 48 use m_xmpi, only : xmpi_sum 49 use m_pawtab, only : pawtab_type 50 use m_pawrhoij, only : pawrhoij_type 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 'pawpolev' 56 use interfaces_18_timing 57 !End of the abilint section 58 59 implicit none 60 61 !Arguments --------------------------------------------- 62 !scalars 63 integer,intent(in) :: my_natom,natom,ntypat 64 integer,optional,intent(in) :: comm_atom 65 !arrays 66 real(dp),intent(out) :: pelev(3) 67 type(pawrhoij_type),intent(in) :: pawrhoij(my_natom) 68 type(pawtab_type),intent(in) :: pawtab(ntypat) 69 70 71 !Local variables --------------------------------------- 72 !scalars 73 integer :: iatom,idir,ierr,irhoij,ispden,itypat,jrhoij,klmn 74 logical :: paral_atom 75 real(dp) :: c1,ro_dlt 76 !arrays 77 integer,dimension(3) :: idirindx = (/4,2,3/) 78 real(dp) :: tsec(2) 79 80 81 ! ************************************************************************* 82 83 DBG_ENTER("COLL") 84 85 call timab(560,1,tsec) 86 87 !Check for parallelism over atoms 88 paral_atom=(present(comm_atom).and.(my_natom/=natom)) 89 90 !note that when vector r is expanded in real spherical harmonics, the factor 91 !sqrt(four_pi/three) appears, as in the following 92 !x = sqrt(four_pi/three)*r*S_{1,1} 93 !y = sqrt(four_pi/three)*r*S_{1,-1} 94 !z = sqrt(four_pi/three)*r*S_{1,0} 95 ! 96 !the moments pawtab()%qijl do not include such normalization factors 97 !see pawinit.F90 for their definition and computation 98 99 c1=sqrt(four_pi/three) 100 101 pelev=zero 102 do idir=1,3 103 do iatom=1,my_natom 104 itypat=pawrhoij(iatom)%itypat 105 do ispden=1,pawrhoij(iatom)%nspden 106 jrhoij=1 107 do irhoij=1,pawrhoij(iatom)%nrhoijsel 108 klmn=pawrhoij(iatom)%rhoijselect(irhoij) 109 ro_dlt=pawrhoij(iatom)%rhoijp(jrhoij,ispden)*pawtab(itypat)%dltij(klmn) 110 pelev(idir)=pelev(idir)+ro_dlt*c1*pawtab(itypat)%qijl(idirindx(idir),klmn) 111 jrhoij=jrhoij+pawrhoij(iatom)%cplex 112 end do 113 end do 114 end do 115 end do 116 117 if (paral_atom) then 118 call xmpi_sum(pelev,comm_atom,ierr) 119 end if 120 121 call timab(560,2,tsec) 122 123 DBG_EXIT("COLL") 124 125 end subroutine pawpolev