The Elk Code
gendsocfr.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2013 J. K. Dewhurst, S. Sharma and E. K. U. Gross
3 ! This file is distributed under the terms of the GNU Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 subroutine gendsocfr
7 use modmain
8 use modphonon
9 implicit none
10 integer is,ias,i
11 integer nr,nri,ir,irc
12 real(8) cso
13 complex(8) z1
14 ! allocatable arrays
15 real(8), allocatable :: vr1(:),vr2(:)
16 real(8), allocatable :: dvr1(:),dvr2(:)
17 if (.not.spinorb) return
18 ! coefficient of spin-orbit coupling
19 cso=y00*socscf/(4.d0*solsc**2)
20 allocate(vr1(nrmtmax),vr2(nrmtmax))
21 allocate(dvr1(nrmtmax),dvr2(nrmtmax))
22 do ias=1,natmtot
23  is=idxis(ias)
24  nr=nrmt(is)
25  nri=nrmti(is)
26  i=1
27  do ir=1,nri
28  vr1(ir)=dble(dvsmt(i,ias))
29  vr2(ir)=aimag(dvsmt(i,ias))
30  i=i+lmmaxi
31  end do
32  do ir=nri+1,nr
33  vr1(ir)=dble(dvsmt(i,ias))
34  vr2(ir)=aimag(dvsmt(i,ias))
35  i=i+lmmaxo
36  end do
37  call splined(nr,wcrmt(:,:,is),vr1,dvr1)
38  call splined(nr,wcrmt(:,:,is),vr2,dvr2)
39  do ir=1,nr,lradstp
40  irc=(ir-1)/lradstp+1
41  z1=cmplx(dvr1(ir),dvr2(ir),8)
42  dsocfr(irc,ias)=(cso/rsp(ir,is))*z1
43  end do
44 end do
45 deallocate(vr1,vr2,dvr1,dvr2)
46 end subroutine
47 
real(8) socscf
Definition: modmain.f90:232
integer lmmaxo
Definition: modmain.f90:203
complex(8), dimension(:,:), allocatable dsocfr
Definition: modphonon.f90:118
subroutine gendsocfr
Definition: gendsocfr.f90:7
complex(8), dimension(:,:), pointer, contiguous dvsmt
Definition: modphonon.f90:110
integer lradstp
Definition: modmain.f90:171
real(8) solsc
Definition: modmain.f90:1253
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
logical spinorb
Definition: modmain.f90:230
pure subroutine splined(n, wc, f, df)
Definition: splined.f90:7
real(8), dimension(:,:), allocatable rsp
Definition: modmain.f90:135
integer natmtot
Definition: modmain.f90:40
real(8), dimension(:,:,:), allocatable wcrmt
Definition: modmain.f90:187
integer nrmtmax
Definition: modmain.f90:152
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
real(8), parameter y00
Definition: modmain.f90:1236
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150