The Elk Code
wavefcr.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2011 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 pure subroutine wavefcr(tsh,lrstp,is,ia,ist,m,ld,wfcr)
7 use modmain
8 implicit none
9 ! arguments
10 logical, intent(in) :: tsh
11 integer, intent(in) :: lrstp,is,ia,ist
12 ! pass in m-1/2
13 integer, intent(in) :: m
14 integer, intent(in) :: ld
15 complex(4), intent(out) :: wfcr(ld,2)
16 ! local variables
17 integer ias,nr,nri,ir
18 integer k,l,lm1,lm2
19 integer i,i1,i2
20 real(8) c1,c2,t0,t1,t2
21 l=lsp(ist,is)
22 k=ksp(ist,is)
23 if (l > lmaxo) then
24  wfcr(:,:)=0.e0
25  return
26 end if
27 ias=idxas(ia,is)
28 ! calculate the Clebsch-Gordon coefficients
29 t1=sqrt(dble(l+m+1)/dble(2*l+1))
30 t2=sqrt(dble(l-m)/dble(2*l+1))
31 if (k == l+1) then
32  c1=t1
33  c2=t2
34 else
35  c1=t2
36  c2=-t1
37 end if
38 if (abs(m) <= l) then
39  lm1=l*(l+1)+m+1
40 else
41  lm1=0
42 end if
43 if (abs(m+1) <= l) then
44  lm2=l*(l+1)+m+2
45 else
46  lm2=0
47 end if
48 nr=nrmt(is)
49 nri=nrmti(is)
50 ! zero the wavefunction
51 if (lrstp == 1) then
52  wfcr(1:npmt(is),1:2)=0.e0
53 else
54  wfcr(1:npcmt(is),1:2)=0.e0
55 end if
56 !----------------------------------!
57 ! inner part of muffin-tin !
58 !----------------------------------!
59 if (l > lmaxi) goto 10
60 if (tsh) then
61  i1=lm1
62  i2=lm2
63 else
64  i=0
65 end if
66 do ir=1,nri,lrstp
67 ! major component of radial wavefunction
68  t0=rwfcr(ir,1,ist,ias)*rlmt(ir,-1,is)
69  if (tsh) then
70  if (lm1 > 0) wfcr(i1,1)=t0*c1
71  if (lm2 > 0) wfcr(i2,2)=t0*c2
72  i1=i1+lmmaxi
73  i2=i2+lmmaxi
74  else
75  t1=t0*c1
76  t2=t0*c2
77  if (lm1 > 0) wfcr(i+1:i+lmmaxi,1)=t1*zbshti(1:lmmaxi,lm1)
78  if (lm2 > 0) wfcr(i+1:i+lmmaxi,2)=t2*zbshti(1:lmmaxi,lm2)
79  i=i+lmmaxi
80  end if
81 end do
82 !----------------------------------!
83 ! outer part of muffin-tin !
84 !----------------------------------!
85 10 continue
86 if (lrstp == 1) then
87  i=lmmaxi*nrmti(is)
88 else
89  i=lmmaxi*nrcmti(is)
90 end if
91 if (tsh) then
92  i1=i+lm1
93  i2=i+lm2
94 end if
95 do ir=nri+lrstp,nr,lrstp
96  t0=rwfcr(ir,1,ist,ias)*rlmt(ir,-1,is)
97  if (tsh) then
98  if (lm1 > 0) wfcr(i1,1)=t0*c1
99  if (lm2 > 0) wfcr(i2,2)=t0*c2
100  i1=i1+lmmaxo
101  i2=i2+lmmaxo
102  else
103  t1=t0*c1
104  t2=t0*c2
105  if (lm1 > 0) wfcr(i+1:i+lmmaxo,1)=t1*zbshto(1:lmmaxo,lm1)
106  if (lm2 > 0) wfcr(i+1:i+lmmaxo,2)=t2*zbshto(1:lmmaxo,lm2)
107  i=i+lmmaxo
108  end if
109 end do
110 end subroutine
111 
integer, dimension(maxstsp, maxspecies) ksp
Definition: modmain.f90:125
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
integer, dimension(maxstsp, maxspecies) lsp
Definition: modmain.f90:123
integer lmmaxo
Definition: modmain.f90:203
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
integer lmaxo
Definition: modmain.f90:201
pure subroutine wavefcr(tsh, lrstp, is, ia, ist, m, ld, wfcr)
Definition: wavefcr.f90:7
integer lmmaxi
Definition: modmain.f90:207
complex(8), dimension(:,:), allocatable zbshto
Definition: modmain.f90:577
real(8), dimension(:,:,:,:), allocatable rwfcr
Definition: modmain.f90:939
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
complex(8), dimension(:,:), allocatable zbshti
Definition: modmain.f90:573
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150
integer lmaxi
Definition: modmain.f90:205