The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine gendsocfr
7use modmain
8use modphonon
9implicit none
10integer is,ias,i
11integer nr,nri,ir,irc
12real(8) cso
13complex(8) z1
14! allocatable arrays
15real(8), allocatable :: vr1(:),vr2(:)
16real(8), allocatable :: dvr1(:),dvr2(:)
17if (.not.spinorb) return
18! coefficient of spin-orbit coupling
19cso=y00*socscf/(4.d0*solsc**2)
20allocate(vr1(nrmtmax),vr2(nrmtmax))
21allocate(dvr1(nrmtmax),dvr2(nrmtmax))
22do 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
44end do
45deallocate(vr1,vr2,dvr1,dvr2)
46end subroutine
47
subroutine gendsocfr
Definition gendsocfr.f90:7
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), parameter y00
Definition modmain.f90:1233
real(8), dimension(:,:,:), allocatable wcrmt
Definition modmain.f90:187
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer lmmaxi
Definition modmain.f90:207
logical spinorb
Definition modmain.f90:230
real(8), dimension(:,:), allocatable rsp
Definition modmain.f90:135
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer lradstp
Definition modmain.f90:171
integer natmtot
Definition modmain.f90:40
real(8) socscf
Definition modmain.f90:232
real(8) solsc
Definition modmain.f90:1252
integer nrmtmax
Definition modmain.f90:152
integer lmmaxo
Definition modmain.f90:203
complex(8), dimension(:,:), pointer, contiguous dvsmt
complex(8), dimension(:,:), allocatable dsocfr
pure subroutine splined(n, wc, f, df)
Definition splined.f90:7