The Elk Code
gengvec.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2012 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: gengvec
8 ! !INTERFACE:
9 subroutine gengvec
10 ! !USES:
11 use modmain
12 ! !DESCRIPTION:
13 ! Generates a set of ${\bf G}$-vectors used for the Fourier transform of the
14 ! charge density and potential and sorts them according to length. Integers
15 ! corresponding to the vectors in lattice coordinates are stored, as well as
16 ! the map from these integer coordinates to the ${\bf G}$-vector index. A map
17 ! from the ${\bf G}$-vector set to the standard FFT array structure is also
18 ! generated. Finally, the number of ${\bf G}$-vectors with magnitude less than
19 ! {\tt gmaxvr} is determined.
20 !
21 ! !REVISION HISTORY:
22 ! Created October 2002 (JKD)
23 !EOP
24 !BOC
25 implicit none
26 ! local variables
27 integer ig,jg,ifg,n1
28 integer i1,i2,i3,j1,j2,j3
29 real(8) v1(3),v2(3),v3(3)
30 ! allocatable arrays
31 integer, allocatable :: idx(:),ivg_(:,:)
32 real(8), allocatable :: vgc_(:,:),gc_(:)
33 if (gmaxvr < 0.d0) gmaxvr=abs(gmaxvr)*gkmax
34 ! ensure |G| cut-off is at least twice |G+k| cut-off
35 gmaxvr=max(gmaxvr,2.d0*gkmax+epslat)
36 ! find the G-vector grid sizes
38 ! number of complex FFT elements for real-complex transforms
39 n1=ngridg(1)/2+1
40 nfgrz=n1*ngridg(2)*ngridg(3)
41 ! allocate local arrays
42 allocate(idx(ngtot),ivg_(3,ngtot))
43 allocate(vgc_(3,ngtot),gc_(ngtot))
44 ig=0
45 do i1=intgv(1,1),intgv(2,1)
46  v1(1:3)=dble(i1)*bvec(1:3,1)
47  do i2=intgv(1,2),intgv(2,2)
48  v2(1:3)=v1(1:3)+dble(i2)*bvec(1:3,2)
49  do i3=intgv(1,3),intgv(2,3)
50  v3(1:3)=v2(1:3)+dble(i3)*bvec(1:3,3)
51  ig=ig+1
52 ! map from G-vector to (i1,i2,i3) index
53  ivg_(1,ig)=i1; ivg_(2,ig)=i2; ivg_(3,ig)=i3
54 ! G-vector in Cartesian coordinates
55  vgc_(1:3,ig)=v3(1:3)
56 ! length of each G-vector
57  gc_(ig)=sqrt(v3(1)**2+v3(2)**2+v3(3)**2)
58  end do
59  end do
60 end do
61 ! sort by vector length
62 call sortidx(ngtot,gc_,idx)
63 ! find the number of vectors with |G| <= gmaxvr
64 ngvec=1
65 do ig=2,ngtot
66  jg=idx(ig)
67  if (gc_(jg) > gmaxvr) then
68  ngvec=ig-1
69  exit
70  end if
71 end do
72 ! allocate global G-vector arrays
73 if (allocated(ivg)) deallocate(ivg)
74 allocate(ivg(3,ngtot))
75 if (allocated(vgc)) deallocate(vgc)
76 allocate(vgc(3,ngtot))
77 if (allocated(gc)) deallocate(gc)
78 allocate(gc(ngtot))
79 ! reorder arrays and store globally
80 do ig=1,ngtot
81  jg=idx(ig)
82  ivg(1:3,ig)=ivg_(1:3,jg)
83  vgc(1:3,ig)=vgc_(1:3,jg)
84  gc(ig)=gc_(jg)
85 end do
86 deallocate(idx,ivg_,vgc_,gc_)
87 ! generate index arrays
88 if (allocated(ivgig)) deallocate(ivgig)
89 allocate(ivgig(intgv(1,1):intgv(2,1), &
90  intgv(1,2):intgv(2,2), &
91  intgv(1,3):intgv(2,3)))
92 ivgig(:,:,:)=0
93 if (allocated(igfft)) deallocate(igfft)
94 allocate(igfft(ngtot))
95 if (allocated(igrzf)) deallocate(igrzf)
96 allocate(igrzf(nfgrz))
97 do ig=1,ngtot
98  i1=ivg(1,ig); i2=ivg(2,ig); i3=ivg(3,ig)
99 ! map from (i1,i2,i3) to G-vector index
100  ivgig(i1,i2,i3)=ig
101 ! Fourier transform index
102  if (i1 >= 0) then
103  j1=i1
104  else
105  j1=ngridg(1)+i1
106  end if
107  if (i2 >= 0) then
108  j2=i2
109  else
110  j2=ngridg(2)+i2
111  end if
112  if (i3 >= 0) then
113  j3=i3
114  else
115  j3=ngridg(3)+i3
116  end if
117  igfft(ig)=j3*ngridg(2)*ngridg(1)+j2*ngridg(1)+j1+1
118 ! map from real-complex FFT index to G-vector index
119  if (i1 >= 0) then
120  ifg=j3*ngridg(2)*n1+j2*n1+j1+1
121  igrzf(ifg)=ig
122  end if
123 end do
124 end subroutine
125 !EOC
126 
integer, dimension(3) ngridg
Definition: modmain.f90:386
integer ngtot
Definition: modmain.f90:390
integer npfftg
Definition: modmain.f90:404
integer, dimension(:,:,:), allocatable ivgig
Definition: modmain.f90:402
integer nfgrz
Definition: modmain.f90:412
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
integer, dimension(:), allocatable igrzf
Definition: modmain.f90:416
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
subroutine gridsize(avec, gmaxvr, npfft, ngridg, ngtot, intgv)
Definition: gridsize.f90:10
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
integer ngvec
Definition: modmain.f90:396
subroutine gengvec
Definition: gengvec.f90:10
pure subroutine sortidx(n, x, idx)
Definition: sortidx.f90:10
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
integer, dimension(:,:), allocatable ivg
Definition: modmain.f90:400
integer, dimension(2, 3) intgv
Definition: modmain.f90:394
real(8) epslat
Definition: modmain.f90:24
real(8), dimension(:), allocatable gc
Definition: modmain.f90:422
real(8) gkmax
Definition: modmain.f90:495
real(8) gmaxvr
Definition: modmain.f90:384