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