The Elk Code
Loading...
Searching...
No Matches
gengvc.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2017 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3
! This file is distributed under the terms of the GNU General Public License.
4
! See the file COPYING for license details.
5
6
subroutine
gengvc
7
use
modmain
8
implicit none
9
! local variables
10
integer
ig,ifg,n1
11
integer
i1,i2,i3,j1,j2,j3
12
real
(8) gk2
13
! find optimal grid size for |G| < 2 gkmax
14
gk2=2.d0*
gkmax
+
epslat
15
ngdgc
(1:3)=int(gk2*sqrt(
avec
(1,1:3)**2+
avec
(2,1:3)**2+
avec
(3,1:3)**2)/
pi
)+1
16
! find next largest FFT-compatible grid size
17
call
nfftifc
(
npfftgc
,3,
ngdgc
)
18
! total number of points in coarse grid
19
ngtc
=
ngdgc
(1)*
ngdgc
(2)*
ngdgc
(3)
20
! find the number of vectors with |G| < 2 gkmax
21
ngvc
=
ngvec
22
do
ig=2,
ngvec
23
if
(
gc
(ig) > gk2)
then
24
ngvc
=ig-1
25
exit
26
end if
27
end do
28
! number of complex FFT elements for real-complex transforms
29
n1=
ngdgc
(1)/2+1
30
nfgrzc
=n1*
ngdgc
(2)*
ngdgc
(3)
31
! Fourier transform index
32
if
(
allocated
(
igfc
))
deallocate
(
igfc
)
33
allocate
(
igfc
(
ngvc
))
34
if
(
allocated
(
igrzfc
))
deallocate
(
igrzfc
)
35
allocate
(
igrzfc
(
nfgrzc
))
36
igrzfc
(1:
nfgrzc
)=
ngvc
+1
37
do
ig=1,
ngvc
38
i1=
ivg
(1,ig)
39
i2=
ivg
(2,ig)
40
i3=
ivg
(3,ig)
41
if
(i1 >= 0)
then
42
j1=i1
43
else
44
j1=
ngdgc
(1)+i1
45
end if
46
if
(i2 >= 0)
then
47
j2=i2
48
else
49
j2=
ngdgc
(2)+i2
50
end if
51
if
(i3 >= 0)
then
52
j3=i3
53
else
54
j3=
ngdgc
(3)+i3
55
end if
56
igfc
(ig)=j3*
ngdgc
(2)*
ngdgc
(1)+j2*
ngdgc
(1)+j1+1
57
! map from real-complex FFT index to G-vector index
58
if
(i1 >= 0)
then
59
ifg=j3*
ngdgc
(2)*n1+j2*n1+j1+1
60
igrzfc
(ifg)=ig
61
end if
62
end do
63
end subroutine
64
gengvc
subroutine gengvc
Definition
gengvc.f90:7
modmain
Definition
modmain.f90:6
modmain::nfgrzc
integer nfgrzc
Definition
modmain.f90:414
modmain::pi
real(8), parameter pi
Definition
modmain.f90:1229
modmain::ngdgc
integer, dimension(3) ngdgc
Definition
modmain.f90:388
modmain::ngvec
integer ngvec
Definition
modmain.f90:396
modmain::gkmax
real(8) gkmax
Definition
modmain.f90:495
modmain::igfc
integer, dimension(:), allocatable igfc
Definition
modmain.f90:410
modmain::avec
real(8), dimension(3, 3) avec
Definition
modmain.f90:12
modmain::ngtc
integer ngtc
Definition
modmain.f90:392
modmain::epslat
real(8) epslat
Definition
modmain.f90:24
modmain::ivg
integer, dimension(:,:), allocatable ivg
Definition
modmain.f90:400
modmain::igrzfc
integer, dimension(:), allocatable igrzfc
Definition
modmain.f90:418
modmain::ngvc
integer ngvc
Definition
modmain.f90:398
modmain::gc
real(8), dimension(:), allocatable gc
Definition
modmain.f90:422
modmain::npfftgc
integer npfftgc
Definition
modmain.f90:408
nfftifc
subroutine nfftifc(np, nd, n)
Definition
nfftifc.f90:10
gengvc.f90
Generated by
1.9.8