TABLE OF CONTENTS


ABINIT/fourier_interpol [ Functions ]

[ Top ] [ Functions ]

NAME

 fourier_interpol

FUNCTION

  Perform a Fourier interpolation. Just a wrapper for transgrid, the table giving the correspondence 
  between the coarse and the mesh FFT grids are constructed inside the routine. This allows to specify an 
  arbitrary FFT mesh to be used for the interpolation. Besides the routine works also in 
  case of NC calculations since it does not require Pawfgr.

COPYRIGHT

  Copyright (C) 2007-2018 ABINIT group (MG)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

 cplex=1 if rhor[f] is real, 2 if rhor[f] is complex
 MPI_enreg<MPI_type>=Information about MPI parallelization
 nspden=number of spin-density components
 nfft_in =number of points in the input FFT box (WARNING no FFT parallelism)
 nfft_out=number of points in the output FFT box
 ngfft_in(18)=all needed information about 3D FFT, for the input grid
 ngfft_out(18) =all needed information about 3D FFT, for the output grid
 optin= 0: input density/potential is taken from rhor_in(:,nspden)
        1: input density/potential is taken from rhog_in(:)     (ispden=1)
                                             and rhor_in(:,2:4) (ispden=2,3,4)
 optout= 0: output density/potential is given in r space in rhor_out(:,nspden)
         1: output density/potential is given in r space in rhor_out(:,nspden)
                                          and in g space in rhog_out(:)
 ngfft_asked(18)=All info on the required FFT mesh.
 paral_kgb=Flag related to the kpoint-band-fft parallelism (TODO not implemented)

OUTPUT

  rhor_out(cplex*nfft_out,nspden)=output density/potential in r space on the required FFT mesh.
  if optout=1:
   rhog_out(2,nfftc)=Fourier transform of output density/potential on the coarse grid

PARENTS

      m_dvdb,m_ioarr,m_qparticles

CHILDREN

      indgrid,pawfgr_destroy,transgrid

SOURCE

 48 #if defined HAVE_CONFIG_H
 49 #include "config.h"
 50 #endif
 51 
 52 #include "abi_common.h"
 53 
 54 subroutine fourier_interpol(cplex,nspden,optin,optout,nfft_in,ngfft_in,nfft_out,ngfft_out,&
 55 & paral_kgb,MPI_enreg,rhor_in,rhor_out,rhog_in,rhog_out)
 56 
 57  use defs_basis
 58  use defs_abitypes
 59  use m_profiling_abi
 60  use m_errors
 61 
 62  use m_pawfgr, only : pawfgr_type, pawfgr_destroy, indgrid
 63 
 64 !This section has been created automatically by the script Abilint (TD).
 65 !Do not modify the following lines by hand.
 66 #undef ABI_FUNC
 67 #define ABI_FUNC 'fourier_interpol'
 68  use interfaces_65_paw, except_this_one => fourier_interpol
 69 !End of the abilint section
 70 
 71  implicit none
 72 
 73 !Arguments ------------------------------------
 74 !scalars
 75  integer,intent(in) :: cplex,nspden,optin,optout
 76  integer,intent(in) :: nfft_in,nfft_out,paral_kgb
 77  type(MPI_type),intent(in) :: MPI_enreg
 78 !arrays
 79  integer,intent(in) :: ngfft_in(18),ngfft_out(18)
 80  real(dp),intent(inout) :: rhor_in(cplex*nfft_in,nspden)
 81  real(dp),intent(inout) :: rhog_in(2,nfft_in)
 82  real(dp),intent(out) :: rhor_out(cplex*nfft_out,nspden)
 83  real(dp),intent(out) :: rhog_out(2,nfft_out)
 84 
 85 !Local variables ---------------------------------------
 86 !scalars
 87  integer :: nfftf,nfftc,nfftc_tot,nfftf_tot,optgrid
 88  logical :: ltest
 89  type(Pawfgr_type) :: Pawfgr
 90 !arrays
 91  integer :: ngfftc(18),ngfftf(18)
 92 
 93 ! *************************************************************************
 94 
 95 !=== FFT parallelism not implemented ===
 96  ABI_CHECK(paral_kgb==0,'paral_kgb/=0 not implemented')
 97 
 98  ltest= ALL( ngfft_in(7:) == ngfft_out(7:) )
 99  ABI_CHECK(ltest,'ngfftf_in(7:18)/=ngfftf_out(7:18)')
100 
101 !================================
102 !=== Which one is the coarse? ===
103 !================================
104  if (nfft_out>=nfft_in) then 
105    ! From coarse to fine grid. If meshes are equivalent, call transgrid anyway because of optout, optin.
106    nfftf    =nfft_out 
107    ngfftf(:)=ngfft_out(:)
108    nfftf_tot =PRODUCT(ngfft_out(1:3))
109 
110    nfftc    =nfft_in 
111    ngfftc(:)=ngfft_in(:)
112    nfftc_tot =PRODUCT(ngfft_in (1:3))
113 
114    Pawfgr%usefinegrid=1 
115    if (ALL(ngfft_in(1:3)==ngfft_out(1:3))) Pawfgr%usefinegrid=0 
116    optgrid=1
117 
118  else 
119    ! From fine towards coarse.
120    nfftf    =nfft_in 
121    ngfftf(:)=ngfft_in(:)
122    nfftf_tot =PRODUCT(ngfft_in(1:3))
123 
124    nfftc    =nfft_out 
125    ngfftc(:)=ngfft_out(:)
126    nfftc_tot =PRODUCT(ngfft_out (1:3))
127    Pawfgr%usefinegrid=1 
128    optgrid=-1
129  end if
130 
131  ABI_MALLOC(Pawfgr%coatofin,(nfftc_tot))
132  ABI_MALLOC(Pawfgr%fintocoa,(nfftf_tot))
133 
134  call indgrid(Pawfgr%coatofin,Pawfgr%fintocoa,nfftc_tot,nfftf_tot,ngfftc,ngfftf)
135 
136  Pawfgr%mgfft =MAXVAL (ngfftf(1:3))
137  Pawfgr%nfft  =PRODUCT(ngfftf(1:3)) ! no FFT parallelism!
138  Pawfgr%ngfft =ngfftf
139 
140  Pawfgr%mgfftc =MAXVAL (ngfftc(1:3))
141  Pawfgr%nfftc  =PRODUCT(ngfftc(1:3)) ! no FFT parallelism!
142  Pawfgr%ngfftc =ngfftc
143 
144  if (optgrid==1) then 
145    call transgrid(cplex,MPI_enreg,nspden,optgrid,optin,optout,paral_kgb,Pawfgr,rhog_in ,rhog_out,rhor_in ,rhor_out)
146  else
147    call transgrid(cplex,MPI_enreg,nspden,optgrid,optin,optout,paral_kgb,Pawfgr,rhog_out,rhog_in ,rhor_out,rhor_in)
148  end if
149 
150  call pawfgr_destroy(Pawfgr)
151  
152 end subroutine fourier_interpol