TABLE OF CONTENTS
ABINIT/pawmkrho [ Functions ]
NAME
pawmkrho
FUNCTION
PAW only: Build total pseudo (compensated) density (\tild_rho + \hat_rho) Build compensation charge density (\hat_rho) Build occupation matrix (packed storage)
COPYRIGHT
Copyright (C) 2010-2018 ABINIT group (MT) 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
cplex: if 1, real space 1-order functions on FFT grid are REAL, if 2, COMPLEX 1 for GS calculations gprimd(3,3)=dimensional primitive translations for reciprocal space(bohr^-1). indsym(4,nsym,natom)=indirect indexing array for atom labels ipert=index of perturbation if pawrhoij is a pertubed rhoij no meaning for ground-state calculations (should be 0) idir=direction of atomic displacement (in case of atomic displ. perturb.) mpi_enreg=information about MPI parallelization my_natom=number of atoms treated by current processor natom=number of atoms in cell nspden=number of spin-density components nsym=number of symmetry elements in space group ntypat=number of types of atoms in unit cell. paral_kgb=option for (kpt,g vectors,bands) parallelism pawang <type(pawang_type)>=angular mesh discretization and related data pawang_sym <type(pawang_type)>=angular data used for symmetrization optional parameter only needed for RF calculations pawfgr <type(paw_fgr_type)>=fine rectangular grid parameters pawfgrtab(natom) <type(pawfgrtab_type)>=atomic data given on fine rectangular grid pawprtvol=control print volume and debugging output for PAW pawrhoij0(natom) <type(pawrhoij_type)>= GS paw rhoij occupancies and related data (used only if ipert>0) optional parameter only needed for RF calculations pawtab(ntypat) <type(pawtab_type)>=paw tabulated starting data qphon(3)=wavevector of the phonon (RF only) rhopsg(2,pawfgr%nfftc)= pseudo density given on the coarse grid in reciprocal space rhopsr(pawfgr%nfftc,nspden)= pseudo density given on the coarse grid in real space rprimd(3,3)=real space primitive translations. symafm(nsym)=(anti)ferromagnetic part of symmetry operations symrec(3,3,nsym)=symmetries of group in terms of operations on reciprocal space primitive translations typat(natom)=type for each atom ucvol=volume of the unit cell xred(3,natom)= reduced atomic coordinates
OUTPUT
compch_fft=compensation charge inside spheres integrated over fine fft grid pawnhat(pawfgr%nfft,nspden)=compensation charge density on fine rectangular grid (optional argument) rhog(2,pawfgr%nfft)= compensated pseudo density given on the fine grid in reciprocal space This output is optional rhor(pawfgr%nfft,nspden)= compensated pseudo density given on the fine grid in real space
SIDE EFFECTS
pawrhoij(my_natom)= PAW occupancies At input : values at previous step in packed storage (pawrhoij()%rhoijp) At output: values (symmetrized) in packed storage (pawrhoij()%rhoijp) pawrhoij_unsym(:)= unsymmetrized PAW occupancies At input : values (unsymmetrized) in unpacked storage (pawrhoij()%rhoij_) At output: values in unpacked storage (pawrhoij()%rhoij_) are destroyed
NOTES
pawrhoij and pawrhoij_unsym can be identical (refer to the same pawrhoij datastructure). They should be different only if pawrhoij is distributed over atomic sites (in that case pawrhoij_unsym should not be distributed over atomic sites).
PARENTS
afterscfloop,dfpt_nstpaw,dfpt_rhofermi,dfpt_vtorho,scfcv,vtorho
CHILDREN
fourdp,pawmknhat,pawrhoij_copy,pawrhoij_free,pawrhoij_free_unpacked pawrhoij_nullify,symrhoij,timab,transgrid
SOURCE
83 #if defined HAVE_CONFIG_H 84 #include "config.h" 85 #endif 86 87 #include "abi_common.h" 88 89 subroutine pawmkrho(compch_fft,cplex,gprimd,idir,indsym,ipert,mpi_enreg,& 90 & my_natom,natom,nspden,nsym,ntypat,paral_kgb,pawang,pawfgr,pawfgrtab,pawprtvol,& 91 & pawrhoij,pawrhoij_unsym,& 92 & pawtab,qphon,rhopsg,rhopsr,rhor,rprimd,symafm,symrec,typat,ucvol,usewvl,xred,& 93 & pawang_sym,pawnhat,pawrhoij0,rhog) ! optional arguments 94 95 use defs_basis 96 use defs_abitypes 97 use m_profiling_abi 98 use m_errors 99 100 use m_pawang, only : pawang_type 101 use m_pawtab, only : pawtab_type 102 use m_pawfgrtab,only : pawfgrtab_type 103 use m_pawrhoij, only : pawrhoij_type, pawrhoij_copy, pawrhoij_free_unpacked, & 104 & pawrhoij_nullify, pawrhoij_free, symrhoij 105 use m_pawfgr, only : pawfgr_type 106 107 !This section has been created automatically by the script Abilint (TD). 108 !Do not modify the following lines by hand. 109 #undef ABI_FUNC 110 #define ABI_FUNC 'pawmkrho' 111 use interfaces_18_timing 112 use interfaces_53_ffts 113 use interfaces_65_paw, except_this_one => pawmkrho 114 !End of the abilint section 115 116 implicit none 117 118 !Arguments ------------------------------------ 119 !scalars 120 integer,intent(in) :: cplex,idir,ipert,my_natom,natom,nspden,nsym,ntypat,paral_kgb,pawprtvol 121 integer,intent(in) :: usewvl 122 real(dp),intent(in) :: ucvol 123 real(dp),intent(out) :: compch_fft 124 type(MPI_type),intent(in) :: mpi_enreg 125 type(pawang_type),intent(in) :: pawang 126 type(pawang_type),intent(in),optional :: pawang_sym 127 type(pawfgr_type),intent(in) :: pawfgr 128 !arrays 129 integer,intent(in) :: indsym(4,nsym,natom) 130 integer,intent(in) :: symafm(nsym),symrec(3,3,nsym),typat(natom) 131 real(dp),intent(in) :: gprimd(3,3),qphon(3),rprimd(3,3),xred(3,natom) 132 real(dp),intent(inout),target,optional :: pawnhat(cplex*pawfgr%nfft,nspden) !vz_i 133 real(dp),intent(inout) :: rhor(cplex*pawfgr%nfft,nspden) 134 real(dp),intent(out),optional :: rhog(2,pawfgr%nfft) 135 real(dp),intent(inout) :: rhopsg(2,pawfgr%nfftc),rhopsr(cplex*pawfgr%nfftc,nspden) 136 type(pawfgrtab_type),intent(inout) :: pawfgrtab(my_natom) 137 type(pawrhoij_type),intent(inout),target :: pawrhoij(:) 138 type(pawrhoij_type),intent(inout) :: pawrhoij_unsym(:) 139 type(pawrhoij_type),intent(in),target,optional :: pawrhoij0(my_natom) 140 type(pawtab_type),intent(in) :: pawtab(ntypat) 141 142 !Local variables------------------------------- 143 !scalars 144 integer :: choice,ider,izero,option 145 character(len=500) :: msg 146 !arrays 147 real(dp) :: tsec(2) 148 real(dp) :: rhodum(0,0,0) 149 real(dp),pointer :: pawnhat_ptr(:,:) 150 type(pawrhoij_type),pointer :: pawrhoij_ptr(:),pawrhoij0_ptr(:) 151 152 ! *********************************************************************** 153 154 DBG_ENTER("COLL") 155 156 call timab(556,1,tsec) 157 158 !Compatibility tests 159 if (size(pawrhoij_unsym)>0) then 160 if (pawrhoij_unsym(1)%use_rhoij_==0) then 161 msg=' rhoij_ field must be allocated in pawrhoij_unsym !' 162 MSG_BUG(msg) 163 end if 164 end if 165 if (ipert>0.and.(.not.present(pawrhoij0))) then 166 msg=' pawrhoij0 must be present when ipert>0 !' 167 MSG_BUG(msg) 168 end if 169 170 !Symetrize PAW occupation matrix and store it in packed storage 171 call timab(557,1,tsec) 172 option=1;choice=1 173 if (present(pawang_sym)) then 174 call symrhoij(pawrhoij,pawrhoij_unsym,choice,gprimd,indsym,ipert,natom,nsym,ntypat,& 175 & option,pawang_sym,pawprtvol,pawtab,rprimd,symafm,symrec,typat,& 176 & comm_atom=mpi_enreg%comm_atom,mpi_atmtab=mpi_enreg%my_atmtab,& 177 & qphon=qphon) 178 else 179 call symrhoij(pawrhoij,pawrhoij_unsym,choice,gprimd,indsym,ipert,natom,nsym,ntypat,& 180 & option,pawang,pawprtvol,pawtab,rprimd,symafm,symrec,typat,& 181 & comm_atom=mpi_enreg%comm_atom,mpi_atmtab=mpi_enreg%my_atmtab,& 182 & qphon=qphon) 183 end if 184 call pawrhoij_free_unpacked(pawrhoij_unsym) 185 call timab(557,2,tsec) 186 187 !In somes cases (parallelism), has to distribute the PAW occupation matrix 188 if (size(pawrhoij)==natom.and.(my_natom/=natom)) then 189 ABI_DATATYPE_ALLOCATE(pawrhoij_ptr,(my_natom)) 190 call pawrhoij_nullify(pawrhoij_ptr) 191 call pawrhoij_copy(pawrhoij,pawrhoij_ptr,& 192 & mpi_atmtab=mpi_enreg%my_atmtab,comm_atom=mpi_enreg%comm_atom, & 193 & keep_cplex=.false.,keep_itypat=.false.,keep_nspden=.false.) 194 else 195 pawrhoij_ptr=>pawrhoij 196 end if 197 198 !Compute compensation charge density 199 ider=0;izero=0 200 if (present(pawnhat)) then 201 pawnhat_ptr => pawnhat 202 else 203 ABI_ALLOCATE(pawnhat_ptr,(pawfgr%nfft,nspden)) 204 end if 205 if (present(pawrhoij0)) then 206 pawrhoij0_ptr => pawrhoij0 207 else 208 pawrhoij0_ptr => pawrhoij_ptr 209 end if 210 211 call pawmknhat(compch_fft,cplex,ider,idir,ipert,izero,gprimd,my_natom,natom,& 212 & pawfgr%nfft,pawfgr%ngfft,ider,nspden,ntypat,pawang,pawfgrtab,& 213 & rhodum,pawnhat_ptr,pawrhoij_ptr,pawrhoij0_ptr,pawtab,qphon,rprimd,ucvol,usewvl,xred,& 214 & comm_atom=mpi_enreg%comm_atom,mpi_atmtab=mpi_enreg%my_atmtab,& 215 & comm_fft=mpi_enreg%comm_fft,paral_kgb=paral_kgb,me_g0=mpi_enreg%me_g0,& 216 & distribfft=mpi_enreg%distribfft,mpi_comm_wvl=mpi_enreg%comm_wvl) 217 218 !Transfer pseudo density from coarse grid to fine grid 219 if(usewvl==0) then 220 call transgrid(cplex,mpi_enreg,nspden,+1,1,0,paral_kgb,pawfgr,rhopsg,rhodum,rhopsr,rhor) 221 end if 222 223 !Add pseudo density and compensation charge density (on fine grid) 224 rhor(:,:)=rhor(:,:)+pawnhat_ptr(:,:) 225 226 !Free temporary memory spaces 227 if (.not.present(pawnhat)) then 228 ABI_DEALLOCATE(pawnhat_ptr) 229 end if 230 if (size(pawrhoij)==natom.and.(my_natom/=natom)) then 231 call pawrhoij_free(pawrhoij_ptr) 232 ABI_DATATYPE_DEALLOCATE(pawrhoij_ptr) 233 end if 234 nullify(pawnhat_ptr) 235 nullify(pawrhoij_ptr) 236 237 !Compute compensated pseudo density in reciprocal space 238 if (present(rhog)) then 239 call fourdp(cplex,rhog,rhor(:,1),-1,mpi_enreg,pawfgr%nfft,pawfgr%ngfft,paral_kgb,0) 240 end if 241 242 call timab(556,2,tsec) 243 244 DBG_EXIT("COLL") 245 246 end subroutine pawmkrho