The Elk Code
 
Loading...
Searching...
No Matches
genppts.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2007 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: genppts
8! !INTERFACE:
9subroutine genppts(tfbz,nsym,sym,ngridp,npptnr,epslat,bvec,boxl,nppt,ipvip, &
10 ipvipnr,ivp,vpl,vpc,wppt,wpptnr)
11! !INPUT/OUTPUT PARAMETERS:
12! tfbz : .true. if vpl and vpc should be mapped to the first Brillouin zone
13! (in,logical)
14! nsym : number of point group symmetries used for reduction, set to 1 for
15! no reduction (in,integer)
16! sym : symmetry matrices in lattice coordinates (in,integer(3,3,*))
17! ngridp : p-point grid sizes (in,integer(3))
18! npptnr : number of non-reduced p-points: ngridp(1)*ngridp(2)*ngridp(3)
19! (in,integer)
20! epslat : tolerance for determining identical vectors (in,real)
21! bvec : reciprocal lattice vectors (in,real(3,3))
22! boxl : corners of box containing p-points in lattice coordinates, the
23! zeroth vector is the origin (in,real(3,0:3))
24! nppt : total number of p-points (out,integer)
25! ipvip : map from (i1,i2,i3) to reduced p-point index
26! (out,integer(0:ngridp(1)-1,0:ngridp(2)-1,0:ngridp(3)-1))
27! ipvipnr : map from (i1,i2,i3) to non-reduced p-point index
28! (out,integer(0:ngridp(1)-1,0:ngridp(2)-1,0:ngridp(3)-1))
29! ivp : integer coordinates of the p-points
30! (out,integer(3,ngridp(1)*ngridp(2)*ngridp(3)))
31! vpl : lattice coordinates of each p-point
32! (out,real(3,ngridp(1)*ngridp(2)*ngridp(3)))
33! vpc : Cartesian coordinates of each p-point
34! (out,real(3,ngridp(1)*ngridp(2)*ngridp(3)))
35! wppt : weights of each reduced p-point
36! (out,real(ngridp(1)*ngridp(2)*ngridp(3)))
37! wpptnr : weight of each non-reduced p-point (out,real)
38! !DESCRIPTION:
39! This routine is used for generating $k$-point or $q$-point sets. Since these
40! are stored in global arrays, the points passed to this and other routines
41! are referred to as $p$-points. In lattice coordinates, the ${\bf p}$ vectors
42! are given by
43! $$ {\bf p}=\left(\begin{matrix} & & \\
44! {\bf B}_2-{\bf B}_1 & {\bf B}_3-{\bf B}_1 & {\bf B}_4-{\bf B}_1 \\
45! & & \end{matrix}\right)
46! \left(\begin{matrix}i_1/n_1 \\ i_2/n_2 \\ i_3/n_3 \end{matrix}\right)
47! +{\bf B}_1 $$
48! where $i_j$ runs from 0 to $n_j-1$, and the ${\bf B}$ vectors define the
49! corners of a box with ${\bf B}_1$ as the origin. If {\tt tfbz} is
50! {\tt .true.} then each {\tt vpl} vector is mapped to the first Brillouin
51! zone. If {\tt tfbz} is {\tt .false.} and then the coordinates of each
52! {\tt vpl} are mapped to the $[0,1)$ interval. The $p$-point weights are
53! stored in {\tt wppt} and the array {\tt ipvip} contains the map from the
54! integer coordinates to the reduced index.
55!
56! !REVISION HISTORY:
57! Created August 2002 (JKD)
58! Updated April 2007 (JKD)
59! Added mapping to the first Brillouin zone, September 2008 (JKD)
60! Made independent of modmain, February 2010 (JKD)
61!EOP
62!BOC
63implicit none
64! arguments
65logical, intent(in) :: tfbz
66integer, intent(in) :: nsym,sym(3,3,*),ngridp(3),npptnr
67real(8), intent(in) :: epslat,bvec(3,3),boxl(3,0:3)
68integer, intent(out) :: nppt
69integer, intent(out) :: ipvip(0:ngridp(1)-1,0:ngridp(2)-1,0:ngridp(3)-1)
70integer, intent(out) :: ipvipnr(0:ngridp(1)-1,0:ngridp(2)-1,0:ngridp(3)-1)
71integer, intent(out) :: ivp(3,npptnr)
72real(8), intent(out) :: vpl(3,npptnr),vpc(3,npptnr)
73real(8), intent(out) :: wppt(npptnr),wpptnr
74! local variables
75integer i1,i2,i3,i
76integer isym,ip,jp
77real(8) v1(3),v2(3),v3(3)
78real(8) b(3,3),t1
79if ((ngridp(1) < 1).or.(ngridp(2) < 1).or.(ngridp(3) < 1)) then
80 write(*,*)
81 write(*,'("Error(genppts): invalid ngridp : ",3I8)') ngridp
82 write(*,*)
83 stop
84end if
85if (npptnr /= ngridp(1)*ngridp(2)*ngridp(3)) then
86 write(*,*)
87 write(*,'("Error(genppts): mismatched npptnr and ngridp : ",4I8)') npptnr, &
88 ngridp
89 write(*,*)
90 stop
91end if
92! box vector matrix
93b(1:3,1)=boxl(1:3,1)-boxl(1:3,0)
94b(1:3,2)=boxl(1:3,2)-boxl(1:3,0)
95b(1:3,3)=boxl(1:3,3)-boxl(1:3,0)
96! weight of each non-reduced p-point
97wpptnr=1.d0/dble(ngridp(1)*ngridp(2)*ngridp(3))
98ip=0
99jp=npptnr
100do i3=0,ngridp(3)-1
101 v1(3)=dble(i3)/dble(ngridp(3))
102 do i2=0,ngridp(2)-1
103 v1(2)=dble(i2)/dble(ngridp(2))
104 do i1=0,ngridp(1)-1
105 v1(1)=dble(i1)/dble(ngridp(1))
106 call r3mv(b,v1,v2)
107 v2(1:3)=v2(1:3)+boxl(1:3,0)
108! map vector components to [0,1)
109 call r3frac(epslat,v2)
110 if (nsym > 1) then
111! determine if this point is equivalent to one already in the set
112 do isym=1,nsym
113 call i3mtrv(sym(:,:,isym),v2,v3)
114 call r3frac(epslat,v3)
115 do i=1,ip
116 t1=abs(vpl(1,i)-v3(1))+abs(vpl(2,i)-v3(2))+abs(vpl(3,i)-v3(3))
117 if (t1 < epslat) then
118! equivalent p-point found so add to existing weight
119 ipvip(i1,i2,i3)=i
120 wppt(i)=wppt(i)+wpptnr
121! add new point to back of set
122 ipvipnr(i1,i2,i3)=jp
123 ivp(1,jp)=i1; ivp(2,jp)=i2; ivp(3,jp)=i3
124 vpl(1:3,jp)=v2(1:3)
125 wppt(jp)=0.d0
126 jp=jp-1
127 goto 10
128 end if
129 end do
130 end do
131 end if
132! add new point to front of set
133 ip=ip+1
134 ipvip(i1,i2,i3)=ip
135 ipvipnr(i1,i2,i3)=ip
136 ivp(1,ip)=i1; ivp(2,ip)=i2; ivp(3,ip)=i3
137 vpl(1:3,ip)=v2(1:3)
138 wppt(ip)=wpptnr
13910 continue
140 end do
141 end do
142end do
143nppt=ip
144do ip=1,npptnr
145! map vpl to the first Brillouin zone if required
146 if (tfbz) call vecfbz(epslat,bvec,vpl(:,ip))
147! determine the Cartesian coordinates of the p-points
148 call r3mv(bvec,vpl(:,ip),vpc(:,ip))
149end do
150return
151
152contains
153
154pure subroutine i3mtrv(a,x,y)
155implicit none
156! arguments
157integer, intent(in) :: a(3,3)
158real(8), intent(in) :: x(3)
159real(8), intent(out) :: y(3)
160y(1)=a(1,1)*x(1)+a(2,1)*x(2)+a(3,1)*x(3)
161y(2)=a(1,2)*x(1)+a(2,2)*x(2)+a(3,2)*x(3)
162y(3)=a(1,3)*x(1)+a(2,3)*x(2)+a(3,3)*x(3)
163end subroutine
164
165end subroutine
166!EOC
167
subroutine genppts(tfbz, nsym, sym, ngridp, npptnr, epslat, bvec, boxl, nppt, ipvip, ipvipnr, ivp, vpl, vpc, wppt, wpptnr)
Definition genppts.f90:11
pure subroutine i3mtrv(a, x, y)
Definition genppts.f90:155
pure subroutine r3frac(eps, v)
Definition r3frac.f90:10
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10
subroutine vecfbz(eps, bvec, vpl)
Definition vecfbz.f90:10