TABLE OF CONTENTS


ABINIT/redgr [ Functions ]

[ Top ] [ Functions ]

NAME

 redgr

FUNCTION

 Compute reduced gradients of a real function on the usual unshifted
 fft grid. The gradient directions are the along the primitive
 reciprocal lattice vectors.
 The input function is intended to be a single spin component of
 the valence charge density, the valence + core charge densities
 or the first-order core charge density for use in frozen wf
 elastic tensor calculations within the GGA.

NOTES

 Closely linked to xcden, but limited to Q=0, real charge densities,
 and unshifted grids.

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DRH, 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

  mpi_enreg=informations about MPI parallelization
  nfft=(effective) number of FFT grid points (for this processor)
  ngfft(18)=contain all needed information about 3D FFT,
   see ~abinit/doc/variables/vargs.htm#ngfft
  frin(nfft)=real space input function

OUTPUT

  frredgr(nfft,3)= reduced gradient of input function (same units as frin)

PARENTS

      dfpt_eltfrxc

CHILDREN

      fourdp

SOURCE

 44 #if defined HAVE_CONFIG_H
 45 #include "config.h"
 46 #endif
 47 
 48 #include "abi_common.h"
 49 
 50 
 51 subroutine redgr (frin,frredgr,mpi_enreg,nfft,ngfft,paral_kgb)
 52 
 53  use defs_basis
 54  use defs_abitypes
 55  use m_profiling_abi
 56 
 57 !This section has been created automatically by the script Abilint (TD).
 58 !Do not modify the following lines by hand.
 59 #undef ABI_FUNC
 60 #define ABI_FUNC 'redgr'
 61  use interfaces_53_ffts
 62 !End of the abilint section
 63 
 64  implicit none
 65 
 66 !Arguments ------------------------------------
 67 !scalars
 68  integer,intent(in) :: nfft,paral_kgb
 69  type(MPI_type),intent(in) :: mpi_enreg
 70 !arrays
 71  integer,intent(in) :: ngfft(18)
 72  real(dp),intent(in) :: frin(nfft)
 73  real(dp),intent(out) :: frredgr(nfft,3)
 74 
 75 !Local variables-------------------------------
 76 !scalars
 77  integer :: cplex_tmp,i1,i2,i3,id,idir,ifft,ig,ii,ing,n1,n2,n3
 78 !arrays
 79  real(dp),allocatable :: gg(:,:),wkcmpx(:,:),work(:),workgr(:,:)
 80 
 81 ! *************************************************************************
 82 
 83 !Only real arrays are treated
 84  cplex_tmp=1
 85 
 86 !Keep local copy of fft dimensions
 87  n1=ngfft(1) ; n2=ngfft(2) ; n3=ngfft(3)
 88 
 89 !In order to speed the routine, precompute the components of g, including 2pi factor
 90  ABI_ALLOCATE(gg,(max(n1,n2,n3),3))
 91  do ii=1,3
 92    id=ngfft(ii)/2+2
 93    do ing=1,ngfft(ii)
 94      ig=ing-(ing/id)*ngfft(ii)-1
 95      gg(ing,ii)=two_pi*ig
 96    end do
 97 !  Note that the G <-> -G symmetry must be maintained
 98    if(mod(ngfft(ii),2)==0)gg(ngfft(ii)/2+1,ii)=zero
 99  end do
100 
101  ABI_ALLOCATE(wkcmpx,(2,nfft))
102  ABI_ALLOCATE(work,(nfft))
103  ABI_ALLOCATE(workgr,(2,nfft))
104 
105 !Obtain rho(G) in wkcmpx from input rho(r)
106  work(:)=frin(:)
107 
108  call fourdp(cplex_tmp,wkcmpx,work,-1,mpi_enreg,nfft,ngfft,paral_kgb,0)
109 
110 !Gradient calculation for three reduced components in turn.
111 !Code duplicated to remove logic from loops.
112  do idir=1,3
113    if(idir==1) then
114 !$OMP PARALLEL DO PRIVATE(ifft)
115      do i3=1,n3
116        ifft=(i3-1)*n1*n2
117        do i2=1,n2
118          do i1=1,n1
119            ifft=ifft+1
120 !          Multiply by i 2pi G(idir)
121            workgr(2,ifft)= gg(i1,idir)*wkcmpx(1,ifft)
122            workgr(1,ifft)=-gg(i1,idir)*wkcmpx(2,ifft)
123          end do
124        end do
125      end do
126    else if(idir==2) then
127 !$OMP PARALLEL DO PRIVATE(ifft)
128      do i3=1,n3
129        ifft=(i3-1)*n1*n2
130        do i2=1,n2
131          do i1=1,n1
132            ifft=ifft+1
133 !          Multiply by i 2pi G(idir)
134            workgr(2,ifft)= gg(i2,idir)*wkcmpx(1,ifft)
135            workgr(1,ifft)=-gg(i2,idir)*wkcmpx(2,ifft)
136          end do
137        end do
138      end do
139    else
140 !$OMP PARALLEL DO PRIVATE(ifft)
141      do i3=1,n3
142        ifft=(i3-1)*n1*n2
143        do i2=1,n2
144          do i1=1,n1
145            ifft=ifft+1
146 !          Multiply by i 2pi G(idir)
147            workgr(2,ifft)= gg(i3,idir)*wkcmpx(1,ifft)
148            workgr(1,ifft)=-gg(i3,idir)*wkcmpx(2,ifft)
149          end do
150        end do
151      end do
152    end if !idir
153 
154    call fourdp(cplex_tmp,workgr,work,1,mpi_enreg,nfft,ngfft,paral_kgb,0)
155 
156 !$OMP PARALLEL DO 
157    do ifft=1,nfft
158      frredgr(ifft,idir)=work(ifft)
159    end do
160 
161  end do !idir
162 
163  ABI_DEALLOCATE(gg)
164  ABI_DEALLOCATE(wkcmpx)
165  ABI_DEALLOCATE(work)
166  ABI_DEALLOCATE(workgr)
167 
168 end subroutine redgr