9subroutine eveqnfvr(nmatp,ngp,vpc,h_,o_,evalfv,evecfv)
13use,
intrinsic :: iso_c_binding
37integer,
intent(in) :: nmatp,ngp
38real(8),
intent(in) :: vpc(3)
39real(8),
target :: h_(*),o_(*)
40real(8),
intent(out) :: evalfv(nstfv)
41complex(8),
intent(out) :: evecfv(nmatmax,nstfv)
43integer is,ia,ja,jas,ilo
54logical tr(nlotot),tp(nlotot)
55integer idx(nlotot),s(nlotot)
56integer iwork(5*nmatp),ifail(nmatp)
60real(8),
allocatable :: rh(:)
61real(8),
pointer,
contiguous :: ro(:),rv(:,:)
62complex(8),
pointer,
contiguous :: h(:),o(:)
66call c_f_pointer(c_loc(h_),h,shape=[n2])
67call c_f_pointer(c_loc(o_),o,shape=[n2])
68call c_f_pointer(c_loc(h_),ro,shape=[n2])
69call c_f_pointer(c_loc(h_(n2+1)),rv,shape=[nmatp,nstfv])
79 t1=0.5d0*(vpc(1)*v(1)+vpc(2)*v(2)+vpc(3)*v(3))
80 z1=cmplx(cos(t1),sin(t1),8)
86 idx(i)=
idxlo(l*(l+1)-m+1,ilo,jas)
89 if (mod(l+m,2) == 0)
then
97 else if (ia > ja)
then
104 if (mod(m,2) == 0)
then
110 if (mod(l,2) == 0)
then
129 if (abs(t1) > 1.d-8)
then
144 rh(k:k+j-1)=dble(h(k:k+j-1))
149 j2=(ngp+idx(m1)-1)*nmatp
153 rh(j1+1:j1+ngp)=dble(h(j1+1:j1+ngp)*z1)+s1*dble(h(j2+1:j2+ngp)*z1)
155 rh(j1+1:j1+ngp)=aimag(h(j1+1:j1+ngp)*z1)+s1*aimag(h(j2+1:j2+ngp)*z1)
159 rh(j1+1:j1+ngp)=dble(h(j1+1:j1+ngp))+s1*dble(h(j2+1:j2+ngp))
161 rh(j1+1:j1+ngp)=aimag(h(j1+1:j1+ngp))+s1*aimag(h(j2+1:j2+ngp))
170 k1=
map(l1,m1); k2=
map(l1,m2); k3=
map(l2,m1); k4=
map(l2,m2)
171 if ((tr(l1).and.tr(m1)).or.((.not.tr(l1)).and.(.not.tr(m1))))
then
172 rh(k1)=dble(h(k1))+s(m1)*dble(h(k2))+s(l1)*(dble(h(k3))+s(m1)*dble(h(k4)))
180 rh(k1)=aimag(h(k1))+s(m1)*t2+s(l1)*(t3+s(m1)*t4)
181 if (.not.tr(l1)) rh(k1)=-rh(k1)
191 ro(k:k+j-1)=dble(o(k:k+j-1))
196 j2=(ngp+idx(m1)-1)*nmatp
200 ro(j1+1:j1+ngp)=dble(o(j1+1:j1+ngp)*z1)+s1*dble(o(j2+1:j2+ngp)*z1)
202 ro(j1+1:j1+ngp)=aimag(o(j1+1:j1+ngp)*z1)+s1*aimag(o(j2+1:j2+ngp)*z1)
206 ro(j1+1:j1+ngp)=dble(o(j1+1:j1+ngp))+s1*dble(o(j2+1:j2+ngp))
208 ro(j1+1:j1+ngp)=aimag(o(j1+1:j1+ngp))+s1*aimag(o(j2+1:j2+ngp))
217 k1=
map(l1,m1); k2=
map(l1,m2); k3=
map(l2,m1); k4=
map(l2,m2)
218 if ((tr(l1).and.tr(m1)).or.((.not.tr(l1)).and.(.not.tr(m1))))
then
219 ro(k1)=dble(o(k1))+s(m1)*dble(o(k2))+s(l1)*(dble(o(k3))+s(m1)*dble(o(k4)))
227 ro(k1)=aimag(o(k1))+s(m1)*t2+s(l1)*(t3+s(m1)*t4)
228 if (.not.tr(l1)) ro(k1)=-ro(k1)
235nts=mkl_set_num_threads_local(nthd)
237call dsygvx(1,
'V',
'I',
'U',nmatp,rh,nmatp,ro,nmatp,vl,vu,1,nstfv,
evaltol,m,w, &
238 rv,nmatp,o_,n2,iwork,ifail,info)
239nts=mkl_set_num_threads_local(0)
243 write(*,
'("Error(eveqnfvr): diagonalisation failed")')
244 write(*,
'(" DSYGVX returned INFO = ",I8)') info
245 if (info > nmatp)
then
247 write(*,
'(" The leading minor of the overlap matrix of order ",I8)') i
248 write(*,
'(" is not positive definite")')
249 write(*,
'(" Order of overlap matrix : ",I8)') nmatp
254evalfv(1:nstfv)=w(1:nstfv)
257 evecfv(1:ngp,j)=rv(1:ngp,j)
258 evecfv(ngp+1:nmatp,j)=0.d0
264 evecfv(i1,j)=evecfv(i1,j)+t1
265 evecfv(i2,j)=evecfv(i2,j)+s(l1)*t1
267 evecfv(i1,j)%im=evecfv(i1,j)%im-t1
268 evecfv(i2,j)%im=evecfv(i2,j)%im-s(l1)*t1
274 evecfv(i1,j)=evecfv(i1,j)*zp(l1)
285elemental integer function map(i,j)
288integer,
intent(in) :: i,j
291 map=ngp+i+(ngp+j-1)*nmatp
293 map=ngp+j+(ngp+i-1)*nmatp
real(8), dimension(3, maxatoms, maxspecies) atposc
integer, dimension(maxspecies) natoms
integer, dimension(:,:,:), allocatable ieqatom
integer, dimension(:,:,:), allocatable idxlo
integer, dimension(maxatoms, maxspecies) idxas
integer, dimension(maxspecies) nlorb
integer, dimension(maxlorb, maxspecies) lorbl