The Elk Code
 
Loading...
Searching...
No Matches
match.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2006 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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: match
8! !INTERFACE:
9subroutine match(ngp,vgpc,gpc,sfacgp,apwalm)
10! !USES:
11use modmain
12! !INPUT/OUTPUT PARAMETERS:
13! ngp : number of G+p-vectors (in,integer)
14! vgpc : G+p-vectors in Cartesian coordinates (in,real(3,ngkmax))
15! gpc : length of G+p-vectors (in,real(ngkmax))
16! sfacgp : structure factors of G+p-vectors (in,complex(ngkmax,natmtot))
17! apwalm : APW matching coefficients
18! (out,complex(ngkmax,apwordmax,lmmaxapw,natmtot))
19! !DESCRIPTION:
20! Computes the $({\bf G+p})$-dependent matching coefficients for the APW basis
21! functions. Inside muffin-tin $\alpha$, the APW functions are given by
22! $$ \phi^{\alpha}_{\bf G+p}({\bf r})=\sum_{l=0}^{l_{\rm max}}
23! \sum_{m=-l}^{l}\sum_{j=1}^{M^{\alpha}_l}A^{\alpha}_{jlm}({\bf G+p})
24! u^{\alpha}_{jl}(r)Y_{lm}(\hat{{\bf r}}), $$
25! where $A^{\alpha}_{jlm}({\bf G+p})$ is the matching coefficient,
26! $M^{\alpha}_l$ is the order of the APW and $u^{\alpha}_{jl}$ is the radial
27! function. In the interstitial region, an APW function is a plane wave,
28! $\exp(i({\bf G+p})\cdot{\bf r})/\sqrt{\Omega}$, where $\Omega$ is the unit
29! cell volume. Ensuring continuity up to the $(M^{\alpha}_l-1)$th derivative
30! across the muffin-tin boundary therefore requires that the matching
31! coefficients satisfy
32! $$ \sum_{j=1}^{M^{\alpha}_l}D_{ij}A^{\alpha}_{jlm}({\bf G+p})=b_i\;, $$
33! where
34! $$ D_{ij}=\left.\frac{d^{i-1}u^{\alpha}_{jl}(r)}{dr^{i-1}}
35! \right|_{r=R_{\alpha}} $$
36! and
37! $$ b_i=\frac{4\pi i^l}{\sqrt{\Omega}}|{\bf G+p}|^{i-1}j^{(i-1)}_l
38! (|{\bf G+p}|R_{\alpha})\exp(i({\bf G+p})\cdot{\bf r}_{\alpha})Y^*_{lm}
39! (\widehat{{\bf G+p}}), $$
40! with ${\bf r}_{\alpha}$ the atomic position and $R_{\alpha}$ the muffin-tin
41! radius. See routine {\tt wfmtfv}.
42!
43! !REVISION HISTORY:
44! Created April 2003 (JKD)
45! Fixed documentation, June 2006 (JKD)
46!EOP
47!BOC
48implicit none
49! arguments
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)
54! local variables
55integer is,ia,ias,omax,ord
56integer l,lma,lmb,lm,io,jo
57integer nr,ir,igp,i,info
58real(8) t0,t1,t2
59complex(8) z1,z2
60! automatic arrays
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))
65! external functions
66real(8), external :: polynm
67! compute the spherical harmonics of the G+p-vectors
68do igp=1,ngp
69 call genylmv(.true.,lmaxapw,vgpc(:,igp),ylmgp(:,igp))
70end do
71t0=1.d0/sqrt(omega)
72! loop over species
73do is=1,nspecies
74 nr=nrmt(is)
75! maximum APW order for this species
76 omax=maxval(apword(1:lmaxapw,is))
77! special case of omax=1
78 if (omax == 1) then
79 do igp=1,ngp
80 t1=gpc(igp)*rmt(is)
81 call sbessel(lmaxapw,t1,djl(:,1,igp))
82 end do
83 do ia=1,natoms(is)
84 ias=idxas(ia,is)
85 do l=0,lmaxapw
86 t1=t0/apwfr(nr,1,1,l,ias)
87 lma=l**2+1; lmb=lma+2*l
88 do igp=1,ngp
89 z1=t1*djl(l,1,igp)*sfacgp(igp,ias)
90 apwalm(igp,1,lma:lmb,ias)=z1*conjg(ylmgp(lma:lmb,igp))
91 end do
92 end do
93 end do
94 cycle
95 end if
96! starting point on radial mesh for fitting polynomial of order npapw
97 ir=nr-npapw+1
98! evaluate the spherical Bessel function derivatives for all G+p-vectors
99 do igp=1,ngp
100 t1=gpc(igp)*rmt(is)
101 call sbessel(lmaxapw,t1,djl(:,1,igp))
102 t2=1.d0
103 do io=2,omax
104 call sbesseldm(io-1,lmaxapw,t1,djl(:,io,igp))
105 t2=t2*gpc(igp)
106 djl(0:lmaxapw,io,igp)=t2*djl(0:lmaxapw,io,igp)
107 end do
108 end do
109! loop over atoms
110 do ia=1,natoms(is)
111 ias=idxas(ia,is)
112! begin loop over l
113 do l=0,lmaxapw
114 ord=apword(l,is)
115! set up matrix of derivatives
116 do jo=1,ord
117 a(1,jo)=apwfr(nr,1,jo,l,ias)
118 do io=2,ord
119 a(io,jo)=polynm(io-1,npapw,rsp(ir,is),apwfr(ir,1,jo,l,ias),rmt(is))
120 end do
121 end do
122 lma=l**2+1; lmb=lma+2*l
123! set up target vectors
124 i=0
125 do igp=1,ngp
126 z1=t0*sfacgp(igp,ias)
127 do lm=lma,lmb
128 i=i+1
129 z2=z1*conjg(ylmgp(lm,igp))
130 b(1:ord,i)=djl(l,1:ord,igp)*z2
131 end do
132 end do
133! solve the general complex linear systems
134 call zgesv(ord,i,a,apwordmax,ipiv,b,apwordmax,info)
135 if (info /= 0) then
136 write(*,*)
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
140 write(*,*)
141 stop
142 end if
143 i=0
144 do igp=1,ngp
145 do lm=lma,lmb
146 i=i+1
147 apwalm(igp,1:ord,lm,ias)=b(1:ord,i)
148 end do
149 end do
150! end loop over l
151 end do
152! end loops over atoms and species
153 end do
154end do
155end subroutine
156!EOC
157
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition genylmv.f90:10
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition match.f90:10
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
real(8), dimension(maxspecies) rmt
Definition modmain.f90:162
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
real(8) omega
Definition modmain.f90:20
integer, dimension(0:maxlapw, maxspecies) apword
Definition modmain.f90:758
real(8), dimension(:,:), allocatable rsp
Definition modmain.f90:135
real(8), dimension(:,:,:,:,:), allocatable apwfr
Definition modmain.f90:774
integer, parameter npapw
Definition modmain.f90:756
integer nspecies
Definition modmain.f90:34
pure real(8) function polynm(m, np, xa, ya, x)
Definition polynm.f90:10
subroutine sbessel(lmax, x, jl)
Definition sbessel.f90:10
subroutine sbesseldm(m, lmax, x, djl)
Definition sbesseldm.f90:10