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,n1,i1,i2,i3,j1,j2,j3
28 real(8) v1(3),v2(3),v3(3)
29 ! allocatable arrays
30 integer, allocatable :: idx(:),ivg_(:,:)
31 real(8), allocatable :: vgc_(:,:),gc_(:)
32 if (gmaxvr < 0.d0) gmaxvr=abs(gmaxvr)*gkmax
33 ! ensure |G| cut-off is at least twice |G+k| cut-off
34 gmaxvr=max(gmaxvr,2.d0*gkmax+epslat)
35 ! find the G-vector grid sizes
37 ! number of complex FFT elements for real-complex transforms
38 n1=ngridg(1)/2+1
39 nfgrz=n1*ngridg(2)*ngridg(3)
40 ! allocate local arrays
41 allocate(idx(ngtot),ivg_(3,ngtot))
42 allocate(vgc_(3,ngtot),gc_(ngtot))
43 ig=0
44 do i1=intgv(1,1),intgv(2,1)
45  v1(1:3)=dble(i1)*bvec(1:3,1)
46  do i2=intgv(1,2),intgv(2,2)
47  v2(1:3)=v1(1:3)+dble(i2)*bvec(1:3,2)
48  do i3=intgv(1,3),intgv(2,3)
49  v3(1:3)=v2(1:3)+dble(i3)*bvec(1:3,3)
50  ig=ig+1
51 ! map from G-vector to (i1,i2,i3) index
52  ivg_(1,ig)=i1; ivg_(2,ig)=i2; ivg_(3,ig)=i3
53 ! G-vector in Cartesian coordinates
54  vgc_(1:3,ig)=v3(1:3)
55 ! length of each G-vector
56  gc_(ig)=sqrt(v3(1)**2+v3(2)**2+v3(3)**2)
57  end do
58  end do
59 end do
60 ! sort by vector length
61 call sortidx(ngtot,gc_,idx)
62 ! allocate global G-vector arrays
63 if (allocated(ivg)) deallocate(ivg)
64 allocate(ivg(3,ngtot))
65 if (allocated(vgc)) deallocate(vgc)
66 allocate(vgc(3,ngtot))
67 if (allocated(gc)) deallocate(gc)
68 allocate(gc(ngtot))
69 ! reorder arrays and store globally
70 ivg(1:3,1:ngtot)=ivg_(1:3,idx(1:ngtot))
71 vgc(1:3,1:ngtot)=vgc_(1:3,idx(1:ngtot))
72 gc(1:ngtot)=gc_(idx(1:ngtot))
73 deallocate(idx,ivg_,vgc_,gc_)
74 ! find the number of vectors with |G| <= gmaxvr
75 ngvec=1
76 do ig=2,ngtot
77  if (gc(ig) > gmaxvr) then
78  ngvec=ig-1
79  exit
80  end if
81 end do
82 ! generate index arrays
83 if (allocated(ivgig)) deallocate(ivgig)
84 allocate(ivgig(intgv(1,1):intgv(2,1), &
85  intgv(1,2):intgv(2,2), &
86  intgv(1,3):intgv(2,3)))
87 if (allocated(igfft)) deallocate(igfft)
88 allocate(igfft(ngtot))
89 if (allocated(igrzf)) deallocate(igrzf)
90 allocate(igrzf(nfgrz))
91 do ig=1,ngtot
92  i1=ivg(1,ig); i2=ivg(2,ig); i3=ivg(3,ig)
93 ! map from (i1,i2,i3) to G-vector index
94  ivgig(i1,i2,i3)=ig
95 ! Fourier transform index
96  j1=merge(i1,ngridg(1)+i1,i1 >= 0)
97  j2=merge(i2,ngridg(2)+i2,i2 >= 0)
98  j3=merge(i3,ngridg(3)+i3,i3 >= 0)
99  igfft(ig)=j3*ngridg(2)*ngridg(1)+j2*ngridg(1)+j1+1
100 ! map from real-complex FFT index to G-vector index
101  if (i1 >= 0) igrzf(j3*ngridg(2)*n1+j2*n1+j1+1)=ig
102 end do
103 end subroutine
104 !EOC
105 
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