The Elk Code
Loading...
Searching...
No Matches
nfftifc.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3
! This file is distributed under the terms of the GNU General Public License.
4
! See the file COPYING for license details.
5
6
!BOP
7
! !ROUTINE: nfftifc
8
! !INTERFACE:
9
subroutine
nfftifc
(np,nd,n)
10
! !INPUT/OUTPUT PARAMETERS:
11
! np : number of allowed primes (in,integer)
12
! nd : number of dimensions (in,integer)
13
! n : required/available grid size (inout,integer(nd))
14
! !DESCRIPTION:
15
! Interface to the grid requirements of the fast Fourier transform routine.
16
! Most routines restrict $n$ to specific prime factorisations. This routine
17
! returns the next largest grid size allowed by the FFT routine.
18
!
19
! !REVISION HISTORY:
20
! Created October 2002 (JKD)
21
!EOP
22
!BOC
23
implicit none
24
! arguments
25
integer
,
intent(in)
:: np,nd
26
integer
,
intent(inout)
:: n(nd)
27
! local variables
28
integer
i,j,l
29
integer
,
parameter
:: p(10)=[2,3,5,7,11,13,17,19,23,29]
30
if
((np < 1).or.(np > 10))
then
31
write
(*,*)
32
write
(*,
'("Error(nfftifc): np out of range : ",I8)'
) np
33
write
(*,*)
34
stop
35
end if
36
if
(nd < 1)
then
37
write
(*,*)
38
write
(*,
'("Error(nfftifc): nd < 1 : ",I8)'
) nd
39
write
(*,*)
40
stop
41
end if
42
! loop over dimensions
43
do
l=1,nd
44
if
(n(l) < 1)
then
45
write
(*,*)
46
write
(*,
'("Error(nfftifc): n < 1 : ",I8)'
) n(l)
47
write
(*,
'(" for dimension : ",I4)'
) l
48
write
(*,*)
49
stop
50
end if
51
do
52
i=n(l)
53
do
j=1,np
54
do
while
(mod(i,p(j)) == 0)
55
i=i/p(j)
56
end do
57
end do
58
if
(i == 1)
exit
59
n(l)=n(l)+1
60
end do
61
end do
62
end subroutine
63
!EOC
64
nfftifc
subroutine nfftifc(np, nd, n)
Definition
nfftifc.f90:10
nfftifc.f90
Generated by
1.9.8