TABLE OF CONTENTS


ABINIT/transgrid [ Functions ]

[ Top ] [ Functions ]

NAME

 transgrid

FUNCTION

 Convert a given density (or potential) from the coarse to the fine rectangular grid and vice versa
 Used in PAW calculations

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (FJ, 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=1 if rhor[f] is real, 2 if rhor[f] is complex
  mpi_enreg=information about MPI parallelization
  nspden=number of spin-density components
  optgrid=+1 to go from the coarse grid towards the fine grid
          -1 to go from the fine grid towards the coarse grid
  optin= 0: input density/potential is taken from rhor(:,nspden)
         1: input density/potential is taken from rhog(:)     (ispden=1)
                                              and rhor(:,2:4) (ispden=2,3,4)
  optout= 0: output density/potential is given in r space in rhor(:,nspden)
          1: output density/potential is given in r space in rhor(:,nspden)
                                           and in g space in rhog(:)
  pawfgr <type(paw_fgr_type)>=fine rectangular grid parameters
    %nfftc=number of points in the coarse FFT box
    %nfft =number of points in the fine FFT box
    %ngfftc(18)=all needed information about 3D FFT, for the coarse grid
    %ngfft(18) =all needed information about 3D FFT, for the fine grid
    %coatofin(nfftc)=Index of the points of the coarse grid on the fine grid
    %fintocoa(nfft) =Index of the points of the fine grid on the coarse grid
    %usefinegrid= 1 if a fine FFT grid is used (0 otherwise)
  if optgrid=+1 and optin=1:
    rhog(2,nfftc)=Fourier transform of input density/potential on the coarse grid
  if optgrid=-1 and optin=1:
    rhogf(2,nfftf)=Fourier transform of input density/potential on the fine grid
  if optgrid=+1
    rhor(cplex*nfftc,nspden)=input density/potential in r space on the coarse grid
  if optgrid=-1:
    rhorf(cplex*nfftf,nspden)=input density/potential in r space on the fine grid

OUTPUT

  if optgrid=-1 and optout=1:
    rhog(2,nfftc)=Fourier transform of output density/potential on the coarse grid
  if optgrid=+1 and optout=1:
    rhogf(2,nfftf)=Fourier transform of output density/potential on the fine grid
  if optgrid=-1
    rhor(cplex*nfftc,nspden)=output density/potential in r space on the coarse grid
  if optgrid=+1:
    rhorf(cplex*nfftf,nspden)=output density/potential in r space on the fine grid

PARENTS

      dfpt_looppert,energy,fourier_interpol,getgh1c,gstate,ks_ddiago,m_io_kss
      pawmkrho,respfn,scfcv,vtorho,vtorhorec

CHILDREN

      fourdp,indirect_parallel_fourier,zerosym

SOURCE

 65 #if defined HAVE_CONFIG_H
 66 #include "config.h"
 67 #endif
 68 
 69 #include "abi_common.h"
 70 
 71 subroutine transgrid(cplex,mpi_enreg,nspden,optgrid,optin,optout,paral_kgb,pawfgr,rhog,rhogf,rhor,rhorf)
 72 
 73  use m_profiling_abi
 74  use defs_basis
 75  use defs_abitypes
 76  use m_errors
 77 
 78  use m_pawfgr, only : pawfgr_type
 79  use m_fft,    only : zerosym
 80 
 81 !This section has been created automatically by the script Abilint (TD).
 82 !Do not modify the following lines by hand.
 83 #undef ABI_FUNC
 84 #define ABI_FUNC 'transgrid'
 85  use interfaces_52_fft_mpi_noabirule
 86  use interfaces_53_ffts
 87 !End of the abilint section
 88 
 89  implicit none
 90 
 91 !Arguments ---------------------------------------------
 92 !scalars
 93  integer,intent(in) :: cplex,nspden,optgrid,optin,optout,paral_kgb
 94  type(MPI_type),intent(in) :: mpi_enreg
 95  type(pawfgr_type),intent(in) :: pawfgr
 96 !arrays
 97  real(dp),intent(inout) :: rhog(2,pawfgr%nfftc),rhogf(2,pawfgr%nfft)
 98  real(dp),intent(inout) :: rhor(cplex*pawfgr%nfftc,nspden),rhorf(cplex*pawfgr%nfft,nspden)
 99 
100 !Local variables ---------------------------------------
101 !scalars
102  integer :: i1,ispden,nfftc,nfftctot,nfftf,nfftftot
103  character(len=500) :: msg
104 !arrays
105  integer :: ngfftc(18),ngfftf(18)
106  real(dp),allocatable :: vectg(:,:),work(:,:),workfft(:)
107 
108 ! *************************************************************************
109 
110  DBG_ENTER("COLL")
111 
112 !Tests
113  if(pawfgr%nfft<pawfgr%nfftc) then
114    write(msg,'(a,2(i0,1x))')' nfft (fine grid) must be >= nfft (coarse grid) while: ',pawfgr%nfft, pawfgr%nfftc
115    MSG_ERROR(msg)
116  end if
117 
118 !Store FFT dimensions
119  nfftc=pawfgr%nfftc;ngfftc(:)=pawfgr%ngfftc(:);nfftctot=ngfftc(1)*ngfftc(2)*ngfftc(3)
120  nfftf=pawfgr%nfft ;ngfftf(:)=pawfgr%ngfft (:);nfftftot=ngfftf(1)*ngfftf(2)*ngfftf(3)
121 
122 !If no fine FFT grid is used, this is only a simple transfer
123  if (pawfgr%usefinegrid==0) then
124    if (optgrid==1) then
125      rhorf=rhor
126      if (optout==1.and.optin==1) rhogf=rhog
127      if (optout==1.and.optin/=1) then
128        ABI_ALLOCATE(workfft,(cplex*nfftc))
129        workfft(:)=rhor(:,1)
130        call fourdp(cplex,rhogf,workfft,-1,mpi_enreg,nfftc,ngfftc,paral_kgb,0)
131        ABI_DEALLOCATE(workfft)
132      end if
133    end if
134    if (optgrid==-1) then
135      rhor=rhorf
136      if (optout==1.and.optin==1) rhog=rhogf
137      if (optout==1.and.optin/=1) then
138        ABI_ALLOCATE(workfft,(cplex*nfftc))
139        workfft(:)=rhorf(:,1)
140        call fourdp(cplex,rhog,workfft,-1,mpi_enreg,nfftc,ngfftc,paral_kgb,0)
141        ABI_DEALLOCATE(workfft)
142      end if
143    end if
144    return
145  end if
146 
147 !====== FROM THE COARSE GRID TOWARDS THE FINE GRID =============
148 !===============================================================
149 !Calculate the FT of rhor to have it in the g space on the coarse grid
150 !Transfer the FT of rhor on the coarse grid towards the fine grid
151 !Then calculate the FT back to get rhorf on the fine grid
152  if (optgrid==1) then
153 
154    ABI_ALLOCATE(work,(2,nfftc))
155 
156 !  First spin component
157 !  --------------------------------------------------------------
158    if (optout==0) then
159 !    if optout=0, rhog on the fine grid is temporary (in vectg)
160      ABI_ALLOCATE(vectg,(2,nfftf))
161      vectg(:,:)=zero
162      if (optin==1) then
163        call zerosym(rhog,2,ngfftc(1),ngfftc(2),ngfftc(3),&
164 &       comm_fft=mpi_enreg%comm_fft,distribfft=mpi_enreg%distribfft)
165        if(mpi_enreg%nproc_fft > 1 .and. mpi_enreg%paral_kgb==1) then
166          call indirect_parallel_Fourier&
167 &         (pawfgr%coatofin,vectg,mpi_enreg,ngfftf,ngfftc,nfftf,nfftc,paral_kgb,rhog,nfftctot)
168        else
169          do i1=1,nfftc
170            vectg(:,pawfgr%coatofin(i1))=rhog(:,i1)
171          end do
172        end if
173      else
174        ABI_ALLOCATE(workfft,(cplex*nfftc))
175        workfft(:)=rhor(:,1)
176        call fourdp(cplex,work,workfft,-1,mpi_enreg,nfftc,ngfftc,paral_kgb,0)
177        ABI_DEALLOCATE(workfft)
178        call zerosym(work,2,ngfftc(1),ngfftc(2),ngfftc(3),&
179 &       comm_fft=mpi_enreg%comm_fft,distribfft=mpi_enreg%distribfft)
180        if(mpi_enreg%nproc_fft > 1 .and. mpi_enreg%paral_kgb==1) then
181          call indirect_parallel_Fourier&
182 &         (pawfgr%coatofin,vectg,mpi_enreg,ngfftf,ngfftc,nfftf,nfftc,paral_kgb,work,nfftctot)
183        else
184          do i1=1,nfftc
185            vectg(:,pawfgr%coatofin(i1))=work(:,i1)
186          end do
187        end if
188      end if
189 !    call zerosym(vectg,2,ngfftf(1),ngfftf(2),ngfftf(3),&
190 !    &        comm_fft=mpi_enreg%comm_fft,distribfft=mpi_enreg%distribfft)
191      ABI_ALLOCATE(workfft,(cplex*nfftf))
192      call fourdp(cplex,vectg,workfft,1,mpi_enreg,nfftf,ngfftf,paral_kgb,0)
193      rhorf(:,1)=workfft(:)
194      ABI_DEALLOCATE(workfft)
195      ABI_DEALLOCATE(vectg)
196    else
197 !    if optout=1, rhog on the fine grid is saved
198      call zerosym(rhog,2,ngfftc(1),ngfftc(2),ngfftc(3),&
199 &     comm_fft=mpi_enreg%comm_fft,distribfft=mpi_enreg%distribfft)
200      rhogf(:,:)=zero
201      if (optin==1) then
202        if(mpi_enreg%nproc_fft > 1 .and. mpi_enreg%paral_kgb==1) then
203          call indirect_parallel_Fourier&
204 &         (pawfgr%coatofin,rhogf,mpi_enreg,ngfftf,ngfftc,nfftf,nfftc,paral_kgb,rhog,nfftctot)
205        else
206          do i1=1,nfftc
207            rhogf(:,pawfgr%coatofin(i1))=rhog(:,i1)
208          end do
209        end if
210      else
211        ABI_ALLOCATE(workfft,(cplex*nfftc))
212        workfft(:)=rhor(:,1)
213        call fourdp(cplex,work,workfft,-1,mpi_enreg,nfftc,ngfftc,paral_kgb,0)
214        ABI_DEALLOCATE(workfft)
215        call zerosym(work,2,ngfftc(1),ngfftc(2),ngfftc(3),&
216 &       comm_fft=mpi_enreg%comm_fft,distribfft=mpi_enreg%distribfft)
217        if(mpi_enreg%nproc_fft > 1 .and. mpi_enreg%paral_kgb==1) then
218          call indirect_parallel_Fourier&
219 &         (pawfgr%coatofin,rhogf,mpi_enreg,ngfftf,ngfftc,nfftf,nfftc,paral_kgb,work,nfftctot)
220        else
221          do i1=1,nfftc
222            rhogf(:,pawfgr%coatofin(i1))=work(:,i1)
223          end do
224        end if
225      end if
226 !    call zerosym(rhogf,2,ngfftf(1),ngfftf(2),ngfftf(3),&
227 !    &        comm_fft=mpi_enreg%comm_fft,distribfft=mpi_enreg%distribfft)
228      ABI_ALLOCATE(workfft,(cplex*nfftf))
229      call fourdp(cplex,rhogf,workfft,1,mpi_enreg,nfftf,ngfftf,paral_kgb,0)
230      rhorf(:,1)=workfft(:)
231      ABI_DEALLOCATE(workfft)
232    end if
233 
234 !  Additional spin components
235 !  ----------------------------------------------------
236    if (nspden>=2) then
237      ABI_ALLOCATE(vectg,(2,nfftf))
238      do ispden=2,nspden
239        vectg(:,:)=zero
240        ABI_ALLOCATE(workfft,(cplex*nfftc))
241        workfft(:)=rhor(:,ispden)
242        call fourdp(cplex,work,workfft,-1,mpi_enreg,nfftc,ngfftc,paral_kgb,0)
243        ABI_DEALLOCATE(workfft)
244        call zerosym(work,2,ngfftc(1),ngfftc(2),ngfftc(3),&
245 &       comm_fft=mpi_enreg%comm_fft,distribfft=mpi_enreg%distribfft)
246        if(mpi_enreg%nproc_fft > 1 .and. mpi_enreg%paral_kgb==1) then
247          call indirect_parallel_Fourier&
248 &         (pawfgr%coatofin,vectg,mpi_enreg,ngfftf,ngfftc,nfftf,nfftc,paral_kgb,work,nfftctot)
249        else
250          do i1=1,nfftc
251            vectg(:,pawfgr%coatofin(i1))=work(:,i1)
252          end do
253        end if
254 !      call zerosym(vectg,2,ngfftf(1),ngfftf(2),ngfftf(3),&
255 !      &          comm_fft=mpi_enreg%comm_fft,distribfft=mpi_enreg%distribfft)
256        ABI_ALLOCATE(workfft,(cplex*nfftf))
257        call fourdp(cplex,vectg,workfft,1,mpi_enreg,nfftf,ngfftf,paral_kgb,0)
258        rhorf(:,ispden)=workfft(:)
259        ABI_DEALLOCATE(workfft)
260      end do
261      ABI_DEALLOCATE(vectg)
262    end if
263 
264    ABI_DEALLOCATE(work)
265 
266 
267 !  ====== FROM THE FINE GRID TOWARDS THE COARSE GRID =============
268 !  ==============================================================
269 !  Calculate the FT of rhorf to have it in the g space on the fine grid
270 !  Transfer the FT of rhorf on the fine grid towards the coarse grid
271 !  Then calculate the FT back to get rhor on the coarse grid
272  else if (optgrid==-1) then
273 
274    ABI_ALLOCATE(work,(2,nfftf))
275 
276 !  First spin component
277 !  --------------------------------------------------------------
278    if (optout==0) then
279 !    if optout=0, rhog on the fine grid is temporary (in vectg)
280      ABI_ALLOCATE(vectg,(2,nfftc))
281      vectg(:,:)=zero
282      if (optin==1) then
283        do i1=1,nfftf
284          if (pawfgr%fintocoa(i1)/=0) vectg(:,pawfgr%fintocoa(i1))=rhogf(:,i1)
285        end do
286      else
287        ABI_ALLOCATE(workfft,(cplex*nfftf))
288        workfft(:)=rhorf(:,1)
289        call fourdp(cplex,work,workfft,-1,mpi_enreg,nfftf,ngfftf,paral_kgb,0)
290        ABI_DEALLOCATE(workfft)
291        if(mpi_enreg%nproc_fft > 1 .and. mpi_enreg%paral_kgb==1) then
292          call indirect_parallel_Fourier&
293 &         (pawfgr%fintocoa,vectg,mpi_enreg,ngfftc,ngfftf,nfftc,nfftf,paral_kgb,work,nfftftot)
294        else
295          do i1=1,nfftf
296            if (pawfgr%fintocoa(i1)/=0) vectg(:,pawfgr%fintocoa(i1))=work(:,i1)
297          end do
298        end if
299      end if
300      call zerosym(vectg,2,ngfftc(1),ngfftc(2),ngfftc(3),&
301 &     comm_fft=mpi_enreg%comm_fft,distribfft=mpi_enreg%distribfft)
302      ABI_ALLOCATE(workfft,(cplex*nfftc))
303      call fourdp(cplex,vectg,workfft,1,mpi_enreg,nfftc,ngfftc,paral_kgb,0)
304      rhor(:,1)=workfft(:)
305      ABI_DEALLOCATE(workfft)
306      ABI_DEALLOCATE(vectg)
307    else
308 !    if optout=1, rhog on the fine grid is saved
309      rhog(:,:)=zero
310      if (optin==1) then
311        do i1=1,nfftf
312          if (pawfgr%fintocoa(i1)/=0) rhog(:,pawfgr%fintocoa(i1))=rhogf(:,i1)
313        end do
314      else
315        ABI_ALLOCATE(workfft,(cplex*nfftf))
316        workfft(:)=rhorf(:,1)
317        call fourdp(cplex,work,workfft,-1,mpi_enreg,nfftf,ngfftf,paral_kgb,0)
318        ABI_DEALLOCATE(workfft)
319        if(mpi_enreg%nproc_fft > 1 .and. mpi_enreg%paral_kgb==1) then
320          call indirect_parallel_Fourier&
321 &         (pawfgr%fintocoa,rhog,mpi_enreg,ngfftc,ngfftf,nfftc,nfftf,paral_kgb,work,nfftftot)
322        else
323          do i1=1,nfftf
324            if (pawfgr%fintocoa(i1)/=0) rhog(:,pawfgr%fintocoa(i1))=work(:,i1)
325          end do
326        end if
327      end if
328      call zerosym(rhog,2,ngfftc(1),ngfftc(2),ngfftc(3),&
329 &     comm_fft=mpi_enreg%comm_fft,distribfft=mpi_enreg%distribfft)
330      ABI_ALLOCATE(workfft,(cplex*nfftc))
331      call fourdp(cplex,rhog,workfft,1,mpi_enreg,nfftc,ngfftc,paral_kgb,0)
332      rhor(:,1)=workfft(:)
333      ABI_DEALLOCATE(workfft)
334    end if
335 
336 !  Additional spin components
337 !  ----------------------------------------------------
338    if (nspden>=2) then
339      ABI_ALLOCATE(vectg,(2,nfftc))
340      do ispden=2,nspden
341        vectg(:,:)=zero
342        ABI_ALLOCATE(workfft,(cplex*nfftf))
343        workfft(:)=rhorf(:,ispden)
344        call fourdp(cplex,work,workfft,-1,mpi_enreg,nfftf,ngfftf,paral_kgb,0)
345        ABI_DEALLOCATE(workfft)
346        if(mpi_enreg%nproc_fft > 1 .and. mpi_enreg%paral_kgb==1) then
347          call indirect_parallel_Fourier&
348 &         (pawfgr%fintocoa,vectg,mpi_enreg,ngfftc,ngfftf,nfftc,nfftf,paral_kgb,work,nfftftot)
349        else
350          do i1=1,nfftf
351            if (pawfgr%fintocoa(i1)/=0) vectg(:,pawfgr%fintocoa(i1))=work(:,i1)
352          end do
353        end if
354        call zerosym(vectg,2,ngfftc(1),ngfftc(2),ngfftc(3),&
355 &       comm_fft=mpi_enreg%comm_fft,distribfft=mpi_enreg%distribfft)
356        ABI_ALLOCATE(workfft,(cplex*nfftc))
357        call fourdp(cplex,vectg,workfft,1,mpi_enreg,nfftc,ngfftc,paral_kgb,0)
358        rhor(:,ispden)=workfft(:)
359        ABI_DEALLOCATE(workfft)
360      end do
361      ABI_DEALLOCATE(vectg)
362    end if
363 
364    ABI_DEALLOCATE(work)
365 
366  end if
367 
368  DBG_EXIT("COLL")
369 
370 end subroutine transgrid