The Elk Code
 
Loading...
Searching...
No Matches
gridsize.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 Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: gridsize
8! !INTERFACE:
9subroutine gridsize(avec,gmaxvr,npfft,ngridg,ngtot,intgv)
10! !INPUT/OUTPUT PARAMETERS:
11! avec : lattice vectors (in,real(3,3))
12! gmaxvr : G-vector cut-off (in,real)
13! npfft : number of prime factors to use for the FFT (in,integer)
14! ngridg : G-vector grid sizes (out,integer(3))
15! ngtot : total number of G-vectors (out,integer)
16! intgv : integer grid intervals for each direction (out,integer(2,3))
17! !DESCRIPTION:
18! Finds the ${\bf G}$-vector grid which completely contains the vectors with
19! $G<G_{\rm max}$ and is compatible with the FFT routine. The optimal sizes
20! are given by
21! $$ n_i=\frac{G_{\rm max}|{\bf a}_i|}{\pi}+1, $$
22! where ${\bf a}_i$ is the $i$th lattice vector.
23!
24! !REVISION HISTORY:
25! Created July 2003 (JKD)
26! Removed modmain and added arguments, September 2012 (JKD)
27!EOP
28!BOC
29implicit none
30! arguments
31real(8), intent(in) :: avec(3,3),gmaxvr
32integer, intent(in) :: npfft
33integer, intent(out) :: ngridg(3),ngtot,intgv(2,3)
34! local variables
35real(8), parameter :: pi=3.1415926535897932385d0
36! find optimal grid size for potential and density
37ngridg(1:3)=int(gmaxvr*sqrt(avec(1,1:3)**2+avec(2,1:3)**2+avec(3,1:3)**2)/pi)+1
38! find next largest FFT-compatible grid size
39call nfftifc(npfft,3,ngridg)
40! total number of points in grid
41ngtot=ngridg(1)*ngridg(2)*ngridg(3)
42! determine integer ranges for grid
43intgv(1,1:3)=ngridg(1:3)/2-ngridg(1:3)+1
44intgv(2,1:3)=ngridg(1:3)/2
45end subroutine
46!EOC
47
subroutine gridsize(avec, gmaxvr, npfft, ngridg, ngtot, intgv)
Definition gridsize.f90:10
subroutine nfftifc(np, nd, n)
Definition nfftifc.f90:10