TABLE OF CONTENTS


ABINIT/wvl_setngfft [ Functions ]

[ Top ] [ Functions ]

NAME

 wvl_setngfft

FUNCTION

 When wavelets are used, the FFT grid is used to store potentials and
 density. The size of the grid takes into account the two resolution in wavelet
 description and also the distribution over processor in the parallel case.

 The FFT grid is not in strict terms an FFT grid but rather a real space grid.
 Its dimensions are not directly compatible with FFTs. This is not relevant
 when using the wavelet part of the code and in the Poisson solver the arrays
 are extended to match FFT dimensions internally. But for other parts of the
 code, this must be taken into account.

 see doc/variables/vargs.html#ngfft for details about ngfft

COPYRIGHT

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

SIDE EFFECTS

  mpi_enreg=informations about MPI parallelization (description of the
            density and potentials scatterring is allocated and updated).
  dtset <type(dataset_type)>=the FFT grid is changed.

PARENTS

      gstate,wvl_wfsinp_reformat

CHILDREN

      wrtout

SOURCE

 39 #if defined HAVE_CONFIG_H
 40 #include "config.h"
 41 #endif
 42 
 43 #include "abi_common.h"
 44 
 45 subroutine wvl_setngfft(me_wvl, mgfft, nfft, ngfft, nproc_wvl, n1i, n2i, n3i,n3d)
 46 
 47  use defs_basis
 48  use m_errors
 49 
 50 !This section has been created automatically by the script Abilint (TD).
 51 !Do not modify the following lines by hand.
 52 #undef ABI_FUNC
 53 #define ABI_FUNC 'wvl_setngfft'
 54  use interfaces_14_hidewrite
 55 !End of the abilint section
 56 
 57  implicit none
 58 
 59 !Arguments ------------------------------------
 60 !scalars
 61  integer, intent(out) :: mgfft, nfft
 62  integer, intent(in)  :: n1i, n2i, n3i,n3d, nproc_wvl, me_wvl
 63 !arrays
 64  integer, intent(out) :: ngfft(18)
 65 
 66 !Local variables-------------------------------
 67 !scalars
 68 #if defined HAVE_BIGDFT
 69  character(len=500) :: message
 70 #endif
 71 
 72 ! *************************************************************************
 73 
 74 #if defined HAVE_BIGDFT
 75  write(message, '(a,a,a,a)' ) ch10,&
 76 & ' wvl_setngfft : Changing the FFT grid definition.'
 77  call wrtout(std_out,message,'COLL')
 78 
 79 !Change nfft and ngfft
 80 !Now ngfft will use the density definition (since the potential size
 81 !is always smaller than the density one). ????
 82  ngfft(1) = n1i
 83  ngfft(2) = n2i
 84  ngfft(3) = n3i
 85 
 86  nfft = n1i*n2i*n3d
 87 !Set up fft array dimensions ngfft(4,5,6) to avoid cache conflicts
 88 !Code paste from getng()
 89  ngfft(4) = 2 * (ngfft(1) / 2) + 1
 90  ngfft(5) = 2 * (ngfft(2) / 2) + 1
 91  ngfft(6) = ngfft(3)
 92  if (nproc_wvl == 0) then
 93    ngfft(9)  = 0    ! paral_fft
 94    ngfft(10) = 1    ! nproc_fft
 95    ngfft(11) = 0    ! me_fft
 96    ngfft(12) = 0    ! n2proc
 97    ngfft(13) = 0    ! n3proc
 98  else
 99    ngfft(9)  = 1    ! paral_fft
100    ngfft(10) = nproc_wvl
101    ngfft(11) = me_wvl
102    ngfft(12) = ngfft(2)
103    ngfft(13) = n3d
104  end if
105 
106  write(message, '(a,3I12)' ) &
107 & '  | ngfft(1:3) is now:    ', ngfft(1:3)
108  call wrtout(std_out,message,'COLL')
109  write(message, '(a,3I12)' ) &
110 & '  | ngfft(4:6) is now:    ', ngfft(4:6)
111  call wrtout(std_out,message,'COLL')
112 
113 !Set mgfft
114  mgfft= max(ngfft(1), ngfft(2), ngfft(3))
115 
116 #else
117  BIGDFT_NOTENABLED_ERROR()
118  if (.false.) write(std_out,*) mgfft,nfft,n1i,n2i,n3i,n3d,nproc_wvl,me_wvl,ngfft(1)
119 #endif
120 
121 end subroutine wvl_setngfft