9 subroutine eveqnfvr(nmatp,ngp,vpc,h_,o_,evalfv,evecfv)
13 use,
intrinsic :: iso_c_binding
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)
43 integer is,ia,ja,jas,ilo
50 real(8) s1,t1,t2,t3,t4
53 logical tr(nlotot),tp(nlotot)
54 integer idx(nlotot),s(nlotot)
55 integer iwork(5*nmatp),ifail(nmatp)
59 real(8),
allocatable :: rh(:)
60 real(8),
pointer,
contiguous :: ro(:),rv(:,:)
61 complex(8),
pointer,
contiguous :: h(:),o(:)
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])
77 t1=0.5d0*(vpc(1)*v(1)+vpc(2)*v(2)+vpc(3)*v(3))
78 z1=cmplx(cos(t1),sin(t1),8)
84 idx(i)=
idxlo(l*(l+1)-m+1,ilo,jas)
87 if (mod(l+m,2) == 0)
then 95 else if (ia > ja)
then 102 if (mod(m,2) == 0)
then 108 if (mod(l,2) == 0)
then 127 if (abs(t1) > 1.d-8)
then 142 rh(k:k+j-1)=dble(h(k:k+j-1))
147 j2=(ngp+idx(m1)-1)*nmatp
151 rh(j1+1:j1+ngp)=dble(h(j1+1:j1+ngp)*z1)+s1*dble(h(j2+1:j2+ngp)*z1)
153 rh(j1+1:j1+ngp)=aimag(h(j1+1:j1+ngp)*z1)+s1*aimag(h(j2+1:j2+ngp)*z1)
157 rh(j1+1:j1+ngp)=dble(h(j1+1:j1+ngp))+s1*dble(h(j2+1:j2+ngp))
159 rh(j1+1:j1+ngp)=aimag(h(j1+1:j1+ngp))+s1*aimag(h(j2+1:j2+ngp))
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)))
178 rh(k1)=aimag(h(k1))+s(m1)*t2+s(l1)*(t3+s(m1)*t4)
179 if (.not.tr(l1)) rh(k1)=-rh(k1)
189 ro(k:k+j-1)=dble(o(k:k+j-1))
194 j2=(ngp+idx(m1)-1)*nmatp
198 ro(j1+1:j1+ngp)=dble(o(j1+1:j1+ngp)*z1)+s1*dble(o(j2+1:j2+ngp)*z1)
200 ro(j1+1:j1+ngp)=aimag(o(j1+1:j1+ngp)*z1)+s1*aimag(o(j2+1:j2+ngp)*z1)
204 ro(j1+1:j1+ngp)=dble(o(j1+1:j1+ngp))+s1*dble(o(j2+1:j2+ngp))
206 ro(j1+1:j1+ngp)=aimag(o(j1+1:j1+ngp))+s1*aimag(o(j2+1:j2+ngp))
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)))
225 ro(k1)=aimag(o(k1))+s(m1)*t2+s(l1)*(t3+s(m1)*t4)
226 if (.not.tr(l1)) ro(k1)=-ro(k1)
233 nts=mkl_set_num_threads_local(nthd)
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)
241 write(*,
'("Error(eveqnfvr): diagonalisation failed")')
242 write(*,
'(" DSYGVX returned INFO = ",I0)') info
243 if (info > nmatp)
then 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
252 evalfv(1:nstfv)=w(1:nstfv)
255 evecfv(1:ngp,j)=rv(1:ngp,j)
256 evecfv(ngp+1:nmatp,j)=0.d0
262 evecfv(i1,j)=evecfv(i1,j)+t1
263 evecfv(i2,j)=evecfv(i2,j)+s(l1)*t1
265 evecfv(i1,j)%im=evecfv(i1,j)%im-t1
266 evecfv(i2,j)%im=evecfv(i2,j)%im-s(l1)*t1
272 evecfv(i1,j)=evecfv(i1,j)*zp(l1)
280 elemental integer function map(i,j)
283 integer,
intent(in) :: i,j
286 map=ngp+i+(ngp+j-1)*nmatp
288 map=ngp+j+(ngp+i-1)*nmatp
integer, dimension(maxspecies) nlorb
integer, dimension(:,:,:), allocatable idxlo
integer, dimension(maxatoms, maxspecies) idxas
integer, dimension(:,:,:), allocatable ieqatom
subroutine eveqnfvr(nmatp, ngp, vpc, h_, o_, evalfv, evecfv)
integer, dimension(maxspecies) natoms
subroutine holdthd(nloop, nthd)
elemental integer function map(i, j)
integer, dimension(maxlorb, maxspecies) lorbl
real(8), dimension(3, maxatoms, maxspecies) atposc