9subroutine match(ngp,vgpc,gpc,sfacgp,apwalm)
50integer,
intent(in) :: ngp
51real(8),
intent(in) :: vgpc(3,ngkmax),gpc(ngkmax)
52complex(8),
intent(in) :: sfacgp(ngkmax,natmtot)
53complex(8),
intent(out) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot)
55integer is,ia,ias,omax,ord
56integer l,lma,lmb,lm,io,jo
57integer nr,ir,igp,i,info
61integer ipiv(apwordmax)
62real(8) djl(0:lmaxapw,apwordmax,ngp)
63complex(8) a(apwordmax,apwordmax),ylmgp(lmmaxapw,ngp)
64complex(8) b(apwordmax,ngp*(2*lmaxapw+1))
66real(8),
external :: polynm
69 call genylmv(.true.,lmaxapw,vgpc(:,igp),ylmgp(:,igp))
76 omax=maxval(
apword(1:lmaxapw,is))
81 call sbessel(lmaxapw,t1,djl(:,1,igp))
86 t1=t0/
apwfr(nr,1,1,l,ias)
87 lma=l**2+1; lmb=lma+2*l
89 z1=t1*djl(l,1,igp)*sfacgp(igp,ias)
90 apwalm(igp,1,lma:lmb,ias)=z1*conjg(ylmgp(lma:lmb,igp))
101 call sbessel(lmaxapw,t1,djl(:,1,igp))
104 call sbesseldm(io-1,lmaxapw,t1,djl(:,io,igp))
106 djl(0:lmaxapw,io,igp)=t2*djl(0:lmaxapw,io,igp)
117 a(1,jo)=
apwfr(nr,1,jo,l,ias)
119 a(io,jo)=
polynm(io-1,
npapw,
rsp(ir,is),
apwfr(ir,1,jo,l,ias),
rmt(is))
122 lma=l**2+1; lmb=lma+2*l
126 z1=t0*sfacgp(igp,ias)
129 z2=z1*conjg(ylmgp(lm,igp))
130 b(1:ord,i)=djl(l,1:ord,igp)*z2
134 call zgesv(ord,i,a,apwordmax,ipiv,b,apwordmax,info)
137 write(*,
'("Error(match): could not find APW matching coefficients")')
138 write(*,
'(" for species ",I4," and atom ",I4)') is,ia
139 write(*,
'(" ZGESV returned INFO = ",I8)') info
147 apwalm(igp,1:ord,lm,ias)=b(1:ord,i)
integer, dimension(maxspecies) nrmt
integer, dimension(maxspecies) natoms
real(8), dimension(maxspecies) rmt
integer, dimension(maxatoms, maxspecies) idxas
integer, dimension(0:maxlapw, maxspecies) apword
real(8), dimension(:,:), allocatable rsp
real(8), dimension(:,:,:,:,:), allocatable apwfr