The Elk Code
eveqnfvr.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2011 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 !BOP
7 ! !ROUTINE: eveqnfvr
8 ! !INTERFACE:
9 subroutine eveqnfvr(nmatp,ngp,vpc,h_,o_,evalfv,evecfv)
10 ! !USES:
11 use modmain
12 use, intrinsic :: iso_c_binding
13 ! !INPUT/OUTPUT PARAMETERS:
14 ! nmatp : order of overlap and Hamiltonian matrices (in,integer)
15 ! ngp : number of G+p-vectors (in,integer)
16 ! vpc : p-vector in Cartesian coordinates (in,real(3))
17 ! h,o : Hamiltonian and overlap matrices in upper triangular form
18 ! (in,complex(*))
19 ! evalfv : first-variational eigenvalues (out,real(nstfv))
20 ! evecfv : first-variational eigenvectors (out,complex(nmatmax,nstfv))
21 ! !DESCRIPTION:
22 ! This routine solves the first-variational eigenvalue equation for the
23 ! special case when inversion symmetry is present. In this case the
24 ! Hamiltonian and overlap matrices can be made real by using appropriate
25 ! linear combinations of the local-orbitals for atoms related by inversion
26 ! symmetry. These are derived from the effect of parity and complex
27 ! conjugation on the spherical harmonics: $P Y_{lm}=(-1)^l Y_{lm}$ and
28 ! $(Y_{lm})^*=(-1)^m Y_{l-m}$.
29 !
30 ! !REVISION HISTORY:
31 ! Created May 2011 (JKD)
32 !EOP
33 !BOC
34 implicit none
35 ! arguments
36 integer, intent(in) :: nmatp,ngp
37 real(8), intent(in) :: vpc(3)
38 real(8), target :: h_(*),o_(*)
39 real(8), intent(out) :: evalfv(nstfv)
40 complex(8), intent(out) :: evecfv(nmatmax,nstfv)
41 ! local variables
42 integer is,ia,ja,jas,ilo
43 integer n2,i,j,k,l,m
44 integer i1,i2,j1,j2
45 integer k1,k2,k3,k4
46 integer l1,l2,m1,m2
47 real(8) v(3),s1,t1,t2,t3,t4
48 complex(8) z1
49 ! automatic arrays
50 logical tr(nlotot),tp(nlotot)
51 integer idx(nlotot),s(nlotot)
52 complex(8) zp(nlotot)
53 ! allocatable arrays and pointers
54 real(8), allocatable :: rh(:)
55 real(8), pointer, contiguous :: ro(:),rv(:,:)
56 complex(8), pointer, contiguous :: h(:),o(:)
57 ! reuse already allocated memory
58 n2=nmatp**2
59 call c_f_pointer(c_loc(h_),h,shape=[n2])
60 call c_f_pointer(c_loc(o_),o,shape=[n2])
61 call c_f_pointer(c_loc(h_),ro,shape=[n2])
62 call c_f_pointer(c_loc(h_(n2+1)),rv,shape=[nmatp,nstfv])
63 tp(1:nlotot)=.false.
64 i=0
65 do is=1,nspecies
66  do ia=1,natoms(is)
67 ! symmetry equivalent atom, mapped with inversion
68  ja=ieqatom(ia,is,2)
69  jas=idxas(ja,is)
70 ! residual phase factor
71  v(1:3)=atposc(1:3,ia,is)+atposc(1:3,ja,is)
72  t1=0.5d0*(vpc(1)*v(1)+vpc(2)*v(2)+vpc(3)*v(3))
73  z1=cmplx(cos(t1),sin(t1),8)
74  do ilo=1,nlorb(is)
75  l=lorbl(ilo,is)
76  do m=-l,l
77  i=i+1
78 ! index to conjugate local-orbital in symmetry equivalent atom
79  idx(i)=idxlo(l*(l+1)-m+1,ilo,jas)
80  if (ia /= ja) then
81 ! sign of parity and conjugation operators
82  if (mod(l+m,2) == 0) then
83  s(i)=1
84  else
85  s(i)=-1
86  end if
87  if (ia < ja) then
88 ! if ia < ja use the real part of the sum of matrix elements
89  tr(i)=.true.
90  else if (ia > ja) then
91 ! if ia > ja use the imaginary part of the difference of matrix elements
92  s(i)=-s(i)
93  tr(i)=.false.
94  end if
95  else
96 ! if ia = ja then use real function when l even and imaginary when l is odd
97  if (mod(m,2) == 0) then
98  s(i)=1
99  else
100  s(i)=-1
101  end if
102 ! new function should be real if symmetric or imaginary if antisymmetric
103  if (mod(l,2) == 0) then
104 ! l even
105  if (m >= 0) then
106  tr(i)=.true.
107  else
108  s(i)=-s(i)
109  tr(i)=.false.
110  end if
111  else
112 ! l odd
113  if (m >= 0) then
114  tr(i)=.false.
115  else
116  s(i)=-s(i)
117  tr(i)=.true.
118  end if
119  end if
120  end if
121 ! phase factors if required
122  if (abs(t1) > 1.d-8) then
123  zp(i)=z1
124  tp(i)=.true.
125  end if
126  end do
127  end do
128  end do
129 end do
130 !---------------------------------!
131 ! real Hamiltonian matrix !
132 !---------------------------------!
133 allocate(rh(n2))
134 ! ⟨APW|H|APW⟩
135 do j=1,ngp
136  k=(j-1)*nmatp+1
137  rh(k:k+j-1)=dble(h(k:k+j-1))
138 end do
139 ! ⟨APW|H|lo⟩
140 do m1=1,nlotot
141  j1=(ngp+m1-1)*nmatp
142  j2=(ngp+idx(m1)-1)*nmatp
143  z1=zp(m1); s1=s(m1)
144  if (tp(m1)) then
145  if (tr(m1)) then
146  rh(j1+1:j1+ngp)=dble(h(j1+1:j1+ngp)*z1)+s1*dble(h(j2+1:j2+ngp)*z1)
147  else
148  rh(j1+1:j1+ngp)=aimag(h(j1+1:j1+ngp)*z1)+s1*aimag(h(j2+1:j2+ngp)*z1)
149  end if
150  else
151  if (tr(m1)) then
152  rh(j1+1:j1+ngp)=dble(h(j1+1:j1+ngp))+s1*dble(h(j2+1:j2+ngp))
153  else
154  rh(j1+1:j1+ngp)=aimag(h(j1+1:j1+ngp))+s1*aimag(h(j2+1:j2+ngp))
155  end if
156  end if
157 end do
158 ! ⟨lo|H|lo⟩
159 do m1=1,nlotot
160  m2=idx(m1)
161  do l1=1,m1
162  l2=idx(l1)
163  k1=map(l1,m1); k2=map(l1,m2); k3=map(l2,m1); k4=map(l2,m2)
164  if ((tr(l1).and.tr(m1)).or.((.not.tr(l1)).and.(.not.tr(m1)))) then
165  rh(k1)=h(k1)%re+s(m1)*h(k2)%re+s(l1)*(h(k3)%re+s(m1)*h(k4)%re)
166  else
167  t2=h(k2)%im
168  if (l1 > m2) t2=-t2
169  t3=h(k3)%im
170  if (l2 > m1) t3=-t3
171  t4=h(k4)%im
172  if (l2 > m2) t4=-t4
173  rh(k1)=h(k1)%im+s(m1)*t2+s(l1)*(t3+s(m1)*t4)
174  if (.not.tr(l1)) rh(k1)=-rh(k1)
175  end if
176  end do
177 end do
178 !-----------------------------!
179 ! real overlap matrix !
180 !-----------------------------!
181 ! ⟨APW|O|APW⟩
182 do j=1,ngp
183  k=(j-1)*nmatp+1
184  ro(k:k+j-1)=dble(o(k:k+j-1))
185 end do
186 ! ⟨APW|O|lo⟩
187 do m1=1,nlotot
188  j1=(ngp+m1-1)*nmatp
189  j2=(ngp+idx(m1)-1)*nmatp
190  z1=zp(m1); s1=s(m1)
191  if (tp(m1)) then
192  if (tr(m1)) then
193  ro(j1+1:j1+ngp)=dble(o(j1+1:j1+ngp)*z1)+s1*dble(o(j2+1:j2+ngp)*z1)
194  else
195  ro(j1+1:j1+ngp)=aimag(o(j1+1:j1+ngp)*z1)+s1*aimag(o(j2+1:j2+ngp)*z1)
196  end if
197  else
198  if (tr(m1)) then
199  ro(j1+1:j1+ngp)=dble(o(j1+1:j1+ngp))+s1*dble(o(j2+1:j2+ngp))
200  else
201  ro(j1+1:j1+ngp)=aimag(o(j1+1:j1+ngp))+s1*aimag(o(j2+1:j2+ngp))
202  end if
203  end if
204 end do
205 ! ⟨lo|O|lo⟩
206 do m1=1,nlotot
207  m2=idx(m1)
208  do l1=1,m1
209  l2=idx(l1)
210  k1=map(l1,m1); k2=map(l1,m2); k3=map(l2,m1); k4=map(l2,m2)
211  if ((tr(l1).and.tr(m1)).or.((.not.tr(l1)).and.(.not.tr(m1)))) then
212  ro(k1)=o(k1)%re+s(m1)*o(k2)%re+s(l1)*(o(k3)%re+s(m1)*o(k4)%re)
213  else
214  t2=o(k2)%im
215  if (l1 > m2) t2=-t2
216  t3=o(k3)%im
217  if (l2 > m1) t3=-t3
218  t4=o(k4)%im
219  if (l2 > m2) t4=-t4
220  ro(k1)=o(k1)%im+s(m1)*t2+s(l1)*(t3+s(m1)*t4)
221  if (.not.tr(l1)) ro(k1)=-ro(k1)
222  end if
223  end do
224 end do
225 ! solve the generalised eigenvalue problem for real symmetric matrices
226 call dsygvxi(nmatp,nstfv,nmatp,rh,ro,evalfv,nmatp,rv,n2,o_)
227 ! reconstruct the complex eigenvectors
228 do j=1,nstfv
229  evecfv(1:ngp,j)=rv(1:ngp,j)
230  evecfv(ngp+1:nmatp,j)=0.d0
231  do l1=1,nlotot
232  i1=ngp+l1
233  i2=ngp+idx(l1)
234  t1=rv(i1,j)
235  if (tr(l1)) then
236  evecfv(i1,j)=evecfv(i1,j)+t1
237  evecfv(i2,j)=evecfv(i2,j)+s(l1)*t1
238  else
239  evecfv(i1,j)%im=evecfv(i1,j)%im-t1
240  evecfv(i2,j)%im=evecfv(i2,j)%im-s(l1)*t1
241  end if
242  end do
243  do l1=1,nlotot
244  if (tp(l1)) then
245  i1=ngp+l1
246  evecfv(i1,j)=evecfv(i1,j)*zp(l1)
247  end if
248  end do
249 end do
250 deallocate(rh)
251 
252 contains
253 
254 elemental integer function map(i,j)
255 implicit none
256 ! arguments
257 integer, intent(in) :: i,j
258 ! map from local-orbital indices to position in matrix
259 map=merge(ngp+i+(ngp+j-1)*nmatp,ngp+j+(ngp+i-1)*nmatp,i <= j)
260 end function
261 
262 end subroutine
263 !EOC
264 
integer, dimension(maxspecies) nlorb
Definition: modmain.f90:784
integer, dimension(:,:,:), allocatable idxlo
Definition: modmain.f90:851
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
subroutine dsygvxi(n, m, ld1, a, b, w, ld2, z, lwork, work)
Definition: dsygvxi.f90:7
integer, dimension(:,:,:), allocatable ieqatom
Definition: modmain.f90:368
subroutine eveqnfvr(nmatp, ngp, vpc, h_, o_, evalfv, evecfv)
Definition: eveqnfvr.f90:10
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer nspecies
Definition: modmain.f90:34
elemental integer function map(i, j)
Definition: eveqnfvr.f90:255
integer, dimension(maxlorb, maxspecies) lorbl
Definition: modmain.f90:794
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition: modmain.f90:54