The Elk Code
 
Loading...
Searching...
No Matches
gengkvec.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 General Public License.
4! See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: gengkvec
8! !INTERFACE:
9pure subroutine gengkvec(ngv,ivg,vgc,vkl,vkc,gkmax,ngkmax,ngk,igkig,vgkl,vgkc, &
10 gkc)
11! !INPUT/OUTPUT PARAMETERS:
12! ngv : number of G-vectors (in,integer)
13! ivg : G-vector integer coordinates (in,integer(3,ngv))
14! vgc : G-vectors in Cartesian coordinates (in,real(3,ngv))
15! vkl : k-point vector in lattice coordinates (in,real(3))
16! vkc : k-point vector in Cartesian coordinates (in,real(3))
17! gkmax : G+k-vector cut-off (in,real)
18! ngkmax : maximum number of G+k-vectors (in,integer)
19! ngk : number of G+k-vectors returned (out,integer)
20! igkig : index from G+k-vectors to G-vectors (out,integer(ngkmax))
21! vgkl : G+k-vectors in lattice coordinates (out,real(3,ngkmax))
22! vgkc : G+k-vectors in Cartesian coordinates (out,real(3,ngkmax))
23! gkc : length of G+k-vectors (out,real(ngkmax))
24! !DESCRIPTION:
25! Generates a set of ${\bf G+k}$-vectors for the input $k$-point with length
26! less than {\tt gkmax}.
27!
28! !REVISION HISTORY:
29! Created April 2003 (JKD)
30! Removed spherical coordinate generation, May 2010 (JKD)
31! Removed modmain and added arguments, September 2012 (JKD)
32!EOP
33!BOC
34implicit none
35! arguments
36integer, intent(in) :: ngv,ivg(3,ngv)
37real(8), intent(in) :: vgc(3,ngv),vkl(3),vkc(3),gkmax
38integer, intent(in) :: ngkmax
39integer, intent(out) :: ngk,igkig(ngkmax)
40real(8), intent(out) :: vgkl(3,ngkmax),vgkc(3,ngkmax),gkc(ngkmax)
41! local variables
42integer ig
43real(8) v1,v2,v3,t0,t1
44t0=gkmax**2
45ngk=0
46do ig=1,ngv
47 v1=vgc(1,ig)+vkc(1)
48 v2=vgc(2,ig)+vkc(2)
49 v3=vgc(3,ig)+vkc(3)
50 t1=v1**2+v2**2+v3**2
51 if (t1 < t0) then
52 ngk=ngk+1
53! index to G-vector
54 igkig(ngk)=ig
55! G+k-vector in lattice coordinates
56 vgkl(1:3,ngk)=dble(ivg(1:3,ig))+vkl(1:3)
57! G+k-vector in Cartesian coordinates
58 vgkc(1,ngk)=v1; vgkc(2,ngk)=v2; vgkc(3,ngk)=v3
59! length of G+k-vector
60 gkc(ngk)=sqrt(t1)
61 if (ngk == ngkmax) return
62 end if
63end do
64end subroutine
65!EOC
66
pure subroutine gengkvec(ngv, ivg, vgc, vkl, vkc, gkmax, ngkmax, ngk, igkig, vgkl, vgkc, gkc)
Definition gengkvec.f90:11