TABLE OF CONTENTS


ABINIT/defs_fftdata [ Modules ]

[ Top ] [ Modules ]

NAME

 defs_fftdata

FUNCTION

 This module contains definitions for a number of parameters
 used to define fft grids within Abinit

COPYRIGHT

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

NOTES

 1) Most of those number are choosen arbitrarily but any modification
 should be done keeping in mind the overall coherency.
 2) ifftdata cannot be constructed directly in one big reshape instruction,
 because the F90 line is too long for some compiler

SOURCE

 23 #if defined HAVE_CONFIG_H
 24 #include "config.h"
 25 #endif
 26 
 27 module defs_fftdata
 28 
 29   use defs_basis 
 30 
 31   implicit none
 32 
 33   integer, parameter :: mdata=7
 34   integer, parameter :: ndata=231
 35   integer, parameter :: mg=65536   ! biggest value in the latter tables (ifftdata and ifftsizes)
 36 
 37   ! The factors 2, 7 and higher than 8 are forbidden
 38   ! 8 is to be the favored factor
 39   ! The factor 6 is only allowed in the first place!
 40   ! These data have to be coherent with the one in getng.F90
 41   integer, parameter :: ifftdata1_40(280)   = (/     &
 42        3,   3, 1, 1, 1, 1, 1,       4,   4, 1, 1, 1, 1, 1,    5,   5, 1, 1, 1, 1, 1,      6,   6, 1, 1, 1, 1, 1,   &
 43        8,   8, 1, 1, 1, 1, 1,       9,   3, 3, 1, 1, 1, 1,   10,   5, 2, 1, 1, 1, 1,     12,   4, 3, 1, 1, 1, 1,   &
 44        15,   5, 3, 1, 1, 1, 1,     16,   4, 4, 1, 1, 1, 1,   18,   6, 3, 1, 1, 1, 1,     20,   5, 4, 1, 1, 1, 1,   &
 45        24,   8, 3, 1, 1, 1, 1,     25,   5, 5, 1, 1, 1, 1,   27,   3, 3, 3, 1, 1, 1,     30,   6, 5, 1, 1, 1, 1,   &
 46        32,   8, 4, 1, 1, 1, 1,     36,   4, 3, 3, 1, 1, 1,   40,   8, 5, 1, 1, 1, 1,     45,   5, 3, 3, 1, 1, 1,   &
 47        48,   4, 4, 3, 1, 1, 1,     50,   5, 5, 2, 1, 1, 1,   54,   6, 3, 3, 1, 1, 1,     60,   5, 4, 3, 1, 1, 1,   &
 48        64,   8, 8, 1, 1, 1, 1,      72,   8, 3, 3, 1, 1, 1,  75,   5, 5, 3, 1, 1, 1,     80,   5, 4, 4, 1, 1, 1,   &
 49        81,   3, 3, 3, 3, 1, 1,      90,   6, 5, 3, 1, 1, 1,  96,   8, 4, 3, 1, 1, 1,     100,   5, 5, 4, 1, 1, 1,  &
 50       108,   4, 3, 3, 3, 1, 1,     120,   8, 5, 3, 1, 1, 1, 125,   5, 5, 5, 1, 1, 1,     128,   8, 4, 4, 1, 1, 1,  &
 51       135,   5, 3, 3, 3, 1, 1,     144,   6, 8, 3, 1, 1, 1, 150,   6, 5, 5, 1, 1, 1,     160,   8, 5, 4, 1, 1, 1   /)
 52   integer, parameter :: ifftdata41_80(280)   = (/     &
 53       162,   6, 3, 3, 3, 1, 1,     180,   5, 4, 3, 3, 1, 1,    192,   6, 8, 4, 1, 1, 1,     200,   8, 5, 5, 1, 1, 1,  &
 54       216,   8, 3, 3, 3, 1, 1,     225,   5, 5, 3, 3, 1, 1,    240,   6, 8, 5, 1, 1, 1,     243,   3, 3, 3, 3, 3, 1,  &
 55       250,   5, 5, 5, 2, 1, 1,     256,   8, 8, 4, 1, 1, 1,    270,   6, 5, 3, 3, 1, 1,     288,   8, 4, 3, 3, 1, 1,  &
 56       300,   5, 5, 4, 3, 1, 1,     320,   5, 4, 4, 4, 1, 1,    324,   4, 3, 3, 3, 3, 1,     360,   8, 5, 3, 3, 1, 1,  &
 57       375,   5, 5, 5, 3, 1, 1,     384,   8, 4, 4, 3, 1, 1,    400,   5, 5, 4, 4, 1, 1,     405,   5, 3, 3, 3, 3, 1,  &
 58       432,   4, 4, 3, 3, 3, 1,     450,   6, 5, 5, 3, 1, 1,    480,   8, 5, 4, 3, 1, 1,     486,   6, 3, 3, 3, 3, 1,  &
 59       500,   5, 5, 5, 4, 1, 1,     512,   8, 8, 8, 1, 1, 1,    540,   5, 4, 3, 3, 3, 1,     576,   4, 4, 4, 3, 3, 1,  &
 60       600,   8, 5, 5, 3, 1, 1,     625,   5, 5, 5, 5, 1, 1,    640,   8, 5, 4, 4, 1, 1,     648,   8, 3, 3, 3, 3, 1,  &
 61       675,   5, 5, 3, 3, 3, 1,     720,   5, 4, 4, 3, 3, 1,    729,   3, 3, 3, 3, 3, 3,     750,   6, 5, 5, 5, 1, 1,  &
 62       768,   4, 4, 4, 4, 3, 1,     800,   8, 5, 5, 4, 1, 1,     810,   6, 5, 3, 3, 3, 1,    864,   8, 4, 3, 3, 3, 1  /)
 63   integer, parameter :: ifftdata81_120(280)   = (/     &
 64       900,   5, 5, 4, 3, 3, 1,     960,   5, 4, 4, 4, 3, 1,     972,   4, 3, 3, 3, 3, 3,   1000,   8, 5, 5, 5, 1, 1,  &
 65      1024,   4, 4, 4, 4, 4, 1,    1080,   6, 5, 4, 3, 3, 1,    1125,   5, 5, 5, 3, 3, 1,   1152,   6, 4, 4, 4, 3, 1,  &
 66      1200,   6, 8, 5, 5, 1, 1,    1215,   5, 3, 3, 3, 3, 3 ,   1250,   5, 5, 5, 5, 2, 1,   1280 ,   8, 8, 5, 4, 1, 1,   &
 67      1296,   6, 8, 3, 3, 3, 1,    1350,   6, 5, 5, 3, 3, 1,   1440,   6, 5, 4, 4, 3, 1,    1458,   6, 3, 3, 3, 3, 3,   &
 68      1500,   5, 5, 5, 4, 3, 1,    1536,   6, 8, 8, 4, 1, 1,   1600,   8, 8, 5, 5, 1, 1,    1620,   5, 4, 3, 3, 3, 3,   &
 69      1728,   6, 8, 4, 3, 3, 1,    1800,   6, 5, 5, 4, 3, 1,   1875,   5, 5, 5, 5, 3, 1,    1920,   6, 5, 4, 4, 4, 1,   &
 70      1944,   6, 4, 3, 3, 3, 3,    2000,   5, 5, 5, 4, 4, 1,   2025,   5, 5, 3, 3, 3, 3,    2048,   8, 4, 4, 4, 4, 1,   &
 71      2160,   6, 8, 5, 3, 3, 1,    2250,   6, 5, 5, 5, 3, 1,   2304,   6, 8, 4, 4, 3, 1,    2400,   6, 5, 5, 4, 4, 1,   &
 72      2430,   6, 5, 3, 3, 3, 3,    2500,   5, 5, 5, 5, 4, 1,   2560,   8, 5, 4, 4, 4, 1,    2592,   6, 4, 4, 3, 3, 3,   &
 73      2700,   5, 5, 4, 3, 3, 3,    2880,   6, 8, 5, 4, 3, 1,   3000,   6, 5, 5, 5, 4, 1,    3072,   6, 8, 4, 4, 4, 1   /)
 74   integer, parameter :: ifftdata121_160(280)   = (/     &
 75      3125,   5, 5, 5, 5, 5, 1,    3200,   8, 5, 5, 4, 4, 1,   3240,   6, 5, 4, 3, 3, 3,    3375,   5, 5, 5, 3, 3, 3,   &
 76      3456,   6, 4, 4, 4, 3, 3,    3600,   6, 8, 5, 5, 3, 1,   3750,   6, 5, 5, 5, 5, 1,    3840,   6, 8, 5, 4, 4, 1,   &
 77      3888,   6, 8, 3, 3, 3, 3,    4000,   8, 5, 5, 5, 4, 1,   4050,   6, 5, 5, 3, 3, 3,    4096,   8, 8, 4, 4, 4, 1,   &
 78      4320,   6, 5, 4, 4, 3, 3,    4500,   5, 5, 5, 4, 3, 3,   4608,   6, 8, 8, 4, 3, 1,    4800,   6, 8, 5, 5, 4, 1,   &
 79      5000,   8, 5, 5, 5, 5, 1,    5120,   8, 8, 5, 4, 4, 1,   5184,   6, 8, 4, 3, 3, 3,    5400,   6, 5, 5, 4, 3, 3,   &
 80      5625,   5, 5, 5, 5, 3, 3,    5760,   6, 8, 8, 5, 3, 1,   6000,   6, 8, 5, 5, 5, 1,    6144,   6, 8, 8, 4, 4, 1,   &
 81      6400,   8, 8, 5, 5, 4, 1,    6480,   6, 8, 5, 3, 3, 3,   6750,   6, 5, 5, 5, 3, 3,    6912,   6, 8, 4, 4, 3, 3,   &
 82      7200,   6, 5, 5, 4, 4, 3,    7500,   5, 5, 5, 5, 4, 3,   7680,   6, 8, 8, 5, 4, 1,    7776,   8, 6, 6, 3, 3, 3,   &
 83      8000,   8, 8, 5, 5, 5, 1,    8100,   6, 6, 5, 5, 3, 3,   8192,   8, 8, 8, 4, 4, 1,   08640,   8, 8, 5, 3, 3, 3,   &
 84     09000, 8, 5, 5, 5, 3, 3,   09216 , 8, 8, 8, 6, 3, 1,     09375, 5, 5, 5, 5, 5, 3,     09600,   8, 8, 6, 5, 5, 1   /)
 85   integer, parameter :: ifftdata161_200(280)   = (/ &
 86     09720, 6, 6, 6, 5, 3, 3,   10000 , 5, 5, 5, 5, 4, 4,    10240, 8, 8, 8, 5, 4, 1,     10368, 8, 8, 6, 3, 3, 3,  &
 87     10800, 8, 6, 5, 5, 3, 3,   11250 , 6, 5, 5, 5, 5, 3,    11520, 8, 8, 6, 6, 5, 1,     11664, 6, 6, 6, 6, 3, 3,  &
 88     12000, 8, 5, 5, 5, 4, 3,   12288, 8, 8, 8, 8, 3, 1,     12500, 5, 5, 5, 5, 5, 4,     12800, 8, 8, 8, 5, 5, 1,  &
 89     12960, 8, 6, 6, 5, 3, 3,   13500, 6, 6, 5, 5, 5, 3,     13824, 8, 8, 8, 3, 3, 3,     14400, 8, 8, 5, 5, 3, 3,  &
 90     15000, 8, 5, 5, 5, 5, 3,   15360, 8, 8, 8, 6, 5, 1,     15552, 8, 6, 6, 6, 3, 3,     15625, 5, 5, 5, 5, 5, 5,  &
 91     16000, 8, 5, 5, 5, 4, 4,   16200, 6, 6, 6, 5, 5, 3,     16384, 8, 8, 8, 8, 4, 1,     17280, 8, 8, 6, 5, 3, 3,  &
 92     18000, 8, 6, 5, 5, 5, 3,   18432, 8, 8, 8, 6, 6, 1,     18750, 6, 5, 5, 5, 5, 5,     19200, 8, 8, 5, 5, 4, 3,  &
 93     19440, 6, 6, 6, 6, 5, 3,   20000, 8, 5, 5, 5, 5, 4,     20480, 8, 8, 8, 8, 5, 1,     20736, 8, 8, 6, 6, 3, 3,  &
 94     21600, 8, 6, 6, 5, 5, 3,   22500, 6, 6, 5, 5, 5, 5,     23040, 8, 8, 8, 5, 3, 3,     23328, 6, 6, 6, 6, 6, 3,  &
 95     24000, 8, 8, 5, 5, 5, 3,   24576, 8, 8, 8, 8, 6, 1,     25000, 8, 5, 5, 5, 5, 5,     25600, 8, 8, 5, 5, 4, 4   /)
 96   integer, parameter :: ifftdata201_231(217)   = (/ &
 97     25920, 8, 6, 6, 6, 5, 3,   27000, 6, 6, 6, 5, 5, 5,     27648, 8, 8, 8, 6, 3, 3,     28800, 8, 8, 6, 5, 5, 3,  &
 98     30000, 8, 6, 5, 5, 5, 5,   30720, 8, 8, 8, 5, 4, 3,     31104, 8, 6, 6, 6, 6, 3,     32000, 8, 8, 5, 5, 5, 4,  &
 99     32400, 6, 6, 6, 6, 5, 5,   32768, 8, 8, 8, 8, 8, 1,     34560, 8, 8, 6, 6, 5, 3,     36000, 8, 6, 6, 5, 5, 5,  &
100     36864, 8, 8, 8, 8, 3, 3,   38400, 8, 8, 8, 5, 5, 3,     38880, 6, 6, 6, 6, 6, 5,     40000, 8, 8, 5, 5, 5, 5,  &
101     40960, 8, 8, 8, 5, 4, 4,   41472, 8, 8, 6, 6, 6, 3,     43200, 8, 6, 6, 6, 5, 5,     46080, 8, 8, 8, 6, 5, 3,  &
102     46656, 6, 6, 6, 6, 6, 6,   48000, 8, 8, 6, 5, 5, 5,     49152, 8, 8, 8, 8, 4, 3,     51200, 8, 8, 8, 5, 5, 4,  &
103     51840, 8, 6, 6, 6, 6, 5,   55296, 8, 8, 8, 6, 6, 3,     57600, 8, 8, 6, 6, 5, 5,     61440, 8, 8, 8, 8, 5, 3,  &
104     62208, 8, 6, 6, 6, 6, 6,   64000, 8, 8, 8, 5, 5, 5,     65536, 8, 8, 8, 8, 4, 4   /)
105   integer, parameter :: ifftdata(mdata,ndata)= reshape(  (/ &
106   &  ifftdata1_40, &
107   &  ifftdata41_80, &
108   &  ifftdata81_120, &
109   &  ifftdata121_160, &
110   &  ifftdata161_200, &
111   &  ifftdata201_231  &
112   &     /), (/mdata,ndata/) )
113 
114   integer, parameter :: ifftsizes(ndata)= (/   &
115        3,    4,   5,     6,    8,    9,    10,   12,   15,   16,   18, &
116        20,   24,   25,   27,   30,   32,   36,   40,   45,   48, &
117        50,   54,   60,   64,   72,   75,   80,   81,   90,   96,  100, &
118        108,  120,  125,  128,  135,  144,  150,  160,  162,  180, &
119        192,  200,  216,  225,  240,  243,  250,  256,  270,  288,  300, &
120        320,  324,  360,  375,  384,  400,  405,  432,  450,  480, &
121        486,  500,  512,  540,  576,  600,  625,  640,  648,  675, &
122        720,  729,  750,  768,  800,  810,  864,  900,  960,  972, &
123        1000, 1024, 1080, 1125, 1152, 1200, 1215, 1250, 1280, 1296, 1350,&
124        1440, 1458, 1500, 1536, 1600, 1620, 1728, 1800, 1875, 1920,&
125        1944, 2000, 2025, 2048, 2160, 2250, 2304, 2400, 2430, 2500,&
126        2560, 2592, 2700, 2880, 3000, 3072, 3125, 3200, 3240, 3375,&
127        3456, 3600, 3750, 3840, 3888, 4000, 4050, 4096, 4320, 4500,&
128        4608, 4800, 5000, 5120, 5184, 5400, 5625, 5760, 6000, 6144,&
129        6400, 6480, 6750, 6912, 7200, 7500, 7680, 7776, 8000, 8100, 8192,&
130        08640, 09000, 09216, 09375, 09600, 09720, 10000, 10240, 10368, 10800, &
131        11250, 11520, 11664, 12000, 12288, 12500, 12800, 12960, 13500, 13824, &
132        14400, 15000, 15360, 15552, 15625, 16000, 16200, 16384, 17280, 18000, &
133        18432, 18750, 19200, 19440, 20000, 20480, 20736, 21600, 22500, 23040, &
134        23328, 24000, 24576, 25000, 25600, 25920, 27000, 27648, 28800, 30000, &
135        30720, 31104, 32000, 32400, 32768, 34560, 36000, 36864, 38400, 38880, &
136        40000, 40960, 41472, 43200, 46080, 46656, 48000, 49152, 51200, 51840, &
137        55296, 57600, 61440, 62208, 64000, 65536 &
138        /)
139 
140 CONTAINS  !===========================================================

defs_fftdata/size_goed_fft [ Functions ]

[ Top ] [ defs_fftdata ] [ Functions ]

NAME

 size_goed_fft

FUNCTION

 Calculate N, the size of the smallest allowed FFT that will encompass
 the 2M+1 points -M..+M. (Ideally, N=2M+1, but this may not be an allowable FFT).

INPUTS

  m_in=Input M value

OUTPUT

  n_out=Output N value
  ierr=Status error 

SOURCE

160 subroutine size_goed_fft(m_in,n_out,ierr)
161 
162  implicit none
163 
164 !Arguments ------------------------------------
165 !scalars
166  integer,intent(in) :: m_in
167  integer,intent(out) :: n_out,ierr
168 
169 !Local variables-------------------------------
170 !scalars
171  integer :: ii,nbest
172  character(len=500) :: msg
173 
174 ! *************************************************************************
175 
176    ierr = 0
177 
178    nbest=2*m_in+1
179 
180    if (nbest<2) then
181      write(msg,'(4a,i8)')ch10,&
182 &     ' size_goed_fft : BUG-',ch10,&
183 &     ' nbest = ',nbest
184      write(std_out,*)msg
185      ierr = 1 
186      RETURN
187    end if
188 
189    if (nbest>ifftsizes(ndata)) then
190      write(msg,'(4a,i8,2a)')ch10,&
191 &     ' size_goed_fft : ERROR-',ch10,&
192 &     ' nbest = ',nbest,ch10,&
193 &     ' is larger than any allowable FFT'
194      write(std_out,*)msg
195      ierr = 2
196      RETURN
197    end if
198 
199    do ii=ndata,1,-1
200      if (ifftsizes(ii)>=nbest) n_out=ifftsizes(ii)
201    end do
202 
203  end subroutine size_goed_fft
204 
205 end module defs_fftdata