TABLE OF CONTENTS


ABINIT/wvl_nhatgrid [ Functions ]

[ Top ] [ Functions ]

NAME

 wvl_nhatgrid

FUNCTION

 Determine parts of the rectangular (fine) grid that are contained
 inside spheres around atoms (used to compute n_hat density).
 If corresponding option is selected, compute also g_l(r)*Y_lm(r)
 (and derivatives) on this grid (g_l=radial shape function).

COPYRIGHT

 Copyright (C) 2011-2018 ABINIT group (T Rangel)
 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

OUTPUT

  pawfgrtab(natom)%ifftsph(nfgd)=FFT index (fine grid) of a points in paw spheres around each atom
  pawfgrtab(natom)%nfgd= number of (fine grid) FFT points in paw spheres around atoms
  if (optgr0==1)
    pawfgrtab(natom)%gylm(nfgd,l_size**2)= g_l(r)*Y_lm(r) around each atom
  if (optgr1==1)
    pawfgrtab(natom)%gylmgr(3,nfgd,l_size**2)= derivatives of g_l(r)*Y_lm(r) wrt cart. coordinates
  if (optgr2==1)
    pawfgrtab(natom)%gylmgr2(6,nfgd,l_size**2)= second derivatives of g_l(r)*Y_lm(r) wrt cart. coordinates
  if (optrad==1)
    pawfgrtab(natom)%rfgd(3,nfgd)= coordinates of r-r_atom around each atom

PARENTS

      afterscfloop,scfcv

CHILDREN

      pawgylm,pawrfgd_wvl,timab,xred2xcart

SOURCE

 41 #if defined HAVE_CONFIG_H
 42 #include "config.h"
 43 #endif
 44 
 45 #include "abi_common.h"
 46 
 47 !PENDING: ADD PARALELLISM OVER ATOMS:
 48 !COPY NHATGRID
 49 
 50 subroutine wvl_nhatgrid(atindx1,geocode,h,i3s,natom,natom_tot,&
 51 & nattyp,ntypat,n1,n1i,n2,n2i,n3,n3pi,optcut,optgr0,optgr1,optgr2,optrad,&
 52 & pawfgrtab,pawtab,psppar,rprimd,shift,xred)
 53 
 54  use defs_basis
 55  use m_profiling_abi
 56  use m_errors
 57 
 58  use m_geometry,     only : xred2xcart
 59  use m_pawtab,       only : pawtab_type
 60  use m_pawfgrtab,    only : pawfgrtab_type
 61  use m_paw_finegrid, only : pawgylm, pawrfgd_wvl
 62 
 63 !This section has been created automatically by the script Abilint (TD).
 64 !Do not modify the following lines by hand.
 65 #undef ABI_FUNC
 66 #define ABI_FUNC 'wvl_nhatgrid'
 67  use interfaces_18_timing
 68 !End of the abilint section
 69 
 70  implicit none
 71 
 72 !Arguments ---------------------------------------------
 73 !scalars
 74  integer,intent(in) :: i3s,natom,natom_tot,ntypat,optcut,optgr0,optgr1,optgr2,optrad
 75  integer,intent(in) :: n1,n2,n3,n1i,n2i,n3pi,shift
 76  real(dp),intent(in) :: h(3)
 77  character(1),intent(in) :: geocode
 78 !integer,intent(in),optional :: mpi_comm_wvl
 79 !arrays
 80  integer,intent(in) :: atindx1(natom),nattyp(ntypat)
 81  real(dp),intent(in) :: psppar(0:4,0:6,ntypat),rprimd(3,3)
 82  real(dp),intent(inout) :: xred(3,natom)
 83  type(pawfgrtab_type),intent(inout) :: pawfgrtab(natom)
 84  type(pawtab_type),intent(in) :: pawtab(ntypat)
 85 
 86 !Local variables ------------------------------
 87 !scalars
 88 !buffer to be added at the end of the last dimension of an array to control bounds_check
 89  integer :: iat,iatm,iatom,iatom_tot,itypat,lm_size,nfgd
 90  real(dp) :: rloc,rshp,xcart(3,natom)
 91 !arrays
 92  integer,allocatable :: ifftsph_tmp(:)
 93  real(dp) :: hh(3) !fine grid spacing for wavelets
 94  real(dp) :: tsec(2)
 95  real(dp),allocatable :: rfgd_tmp(:,:)
 96 
 97 ! *************************************************************************
 98 
 99  DBG_ENTER("COLL")
100 
101 #if !defined HAVE_BIGDFT
102  BIGDFT_NOTENABLED_ERROR()
103 #endif
104 
105  call timab(559,1,tsec)
106 
107 !Set up parallelism for wvl
108 !for debug: use me_wvl=xmpi_comm_rank(MPI_COMM_WORLD)
109 !if (present(mpi_comm_wvl)) then
110 !me_wvl=xmpi_comm_rank(mpi_comm_wvl)
111 !nproc_wvl=xmpi_comm_size(mpi_comm_wvl)
112 !else
113 !me_wvl=0;nproc_wvl=1
114 !end if
115 !Pending: parallelism over atoms: see nhatgrid
116 
117  if (natom_tot<natom) then   ! This test has to be remove when natom_tot is used
118    MSG_BUG(' natom_tot<natom !')
119  end if
120 
121 !Fine grid
122  hh(:)=0.5d0*h(:)
123 
124 !Compute xcart from xred
125  call xred2xcart(natom,rprimd,xcart,xred)
126 
127 !Loop over types of atom
128  iatm=0
129  do itypat=1,ntypat
130 
131    rloc=psppar(0,0,itypat)
132    if (optcut==1) then
133      rshp=pawtab(itypat)%rpaw
134    else
135      rshp=pawtab(itypat)%rshp
136    end if
137 
138 !  Loop over atoms
139    do iat=1,nattyp(itypat)
140      iatm=iatm+1;iatom=atindx1(iatm)
141      iatom_tot=iatom; !if (paral_atom) iatom_tot=my_atmtab(iatom)
142      lm_size=pawfgrtab(iatom)%l_size**2
143 
144 !    Determine FFT points and r-R vectors around the atom
145      call pawrfgd_wvl(geocode,hh,ifftsph_tmp,i3s,n1,n1i,n2,n2i,n3,n3pi,nfgd,rshp,rloc,&
146 &     rfgd_tmp,shift,xcart(:,iatom_tot))
147 
148 !    Allocate arrays defining sphere (and related data) around current atom
149      if (allocated(pawfgrtab(iatom)%ifftsph)) then
150        ABI_DEALLOCATE(pawfgrtab(iatom)%ifftsph)
151      end if
152      ABI_ALLOCATE(pawfgrtab(iatom)%ifftsph,(nfgd))
153      pawfgrtab(iatom)%nfgd=nfgd
154      pawfgrtab(iatom)%ifftsph(1:nfgd)=ifftsph_tmp(1:nfgd)
155 
156      if (optrad==1) then
157        if (allocated(pawfgrtab(iatom)%rfgd)) then
158          ABI_DEALLOCATE(pawfgrtab(iatom)%rfgd)
159        end if
160        ABI_ALLOCATE(pawfgrtab(iatom)%rfgd,(3,nfgd))
161        pawfgrtab(iatom)%rfgd_allocated=1
162        pawfgrtab(iatom)%rfgd(1:3,1:nfgd)=rfgd_tmp(1:3,1:nfgd)
163      end if
164 
165      if (optgr0==1) then
166        if (allocated(pawfgrtab(iatom)%gylm)) then
167          ABI_DEALLOCATE(pawfgrtab(iatom)%gylm)
168        end if
169        ABI_ALLOCATE(pawfgrtab(iatom)%gylm,(nfgd,lm_size))
170        pawfgrtab(iatom)%gylm_allocated=1
171      end if
172 
173      if (optgr1==1) then
174        if (allocated(pawfgrtab(iatom)%gylmgr)) then
175          ABI_DEALLOCATE(pawfgrtab(iatom)%gylmgr)
176        end if
177        ABI_ALLOCATE(pawfgrtab(iatom)%gylmgr,(3,nfgd,lm_size))
178        pawfgrtab(iatom)%gylmgr_allocated=1
179      end if
180 
181      if (optgr2==1) then
182        if (allocated(pawfgrtab(iatom)%gylmgr2)) then
183          ABI_DEALLOCATE(pawfgrtab(iatom)%gylmgr2)
184        end if
185        ABI_ALLOCATE(pawfgrtab(iatom)%gylmgr2,(6,nfgd,lm_size))
186        pawfgrtab(iatom)%gylmgr2_allocated=1
187      end if
188 
189 !    Calculate g_l(r-R)*Y_lm(r-R) for each r around the atom R
190      if (optgr0+optgr1+optgr2>0) then
191        call pawgylm(pawfgrtab(iatom)%gylm,pawfgrtab(iatom)%gylmgr,pawfgrtab(iatom)%gylmgr2,&
192 &       lm_size,nfgd,optgr0,optgr1,optgr2,pawtab(itypat),rfgd_tmp(:,1:nfgd))
193      end if
194 
195 !    End loops over types/atoms
196      ABI_DEALLOCATE(ifftsph_tmp)
197      ABI_DEALLOCATE(rfgd_tmp)
198    end do
199  end do
200 
201  call timab(559,2,tsec)
202 
203  DBG_EXIT("COLL")
204 
205 end subroutine wvl_nhatgrid