The Elk Code
gensocfr.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2009 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 !BOP
7 ! !ROUTINE: gensocfr
8 ! !INTERFACE:
9 subroutine gensocfr
10 ! !USES:
11 use modmain
12 use modomp
13 ! !DESCRIPTION:
14 ! Calculates the radial part of the spin-orbit coupling term which is added to
15 ! the second-variational Hamiltonian. In a particular muffin-tin, this is
16 ! given by
17 ! $$ f_{\rm soc}(r)=\frac{1}{(2Mc)^2}\frac{1}{r}
18 ! \frac{\partial V_s}{\partial r}, $$
19 ! where
20 ! $$ M(r)=1+\frac{1}{2c^2}(E-V_s) $$
21 ! (with $E$ set to zero) and $V_s$ is the spherical part of the Kohn-Sham
22 ! potential. The term added to the Hamiltonian is then
23 ! $$ \hat{H}_{\rm soc}(r)=f_{\rm soc}\hat{\bf L}\cdot\boldsymbol{\sigma}. $$
24 ! See Koelling and Harmon, {\it J. Phys. C: Solid State Phys.} {\bf 10}, 3107
25 ! (1977) for details.
26 !
27 ! !REVISION HISTORY:
28 ! Created April 2009 (JKD)
29 !EOP
30 !BOC
31 implicit none
32 ! local variables
33 integer is,ias,nthd
34 integer nr,nri,ir,irc
35 real(8) cso,rm
36 ! automatic arrays
37 real(8) vr(nrmtmax),dvr(nrmtmax)
38 if (.not.spinorb) return
39 ! coefficient of spin-orbit coupling
40 cso=y00*socscf/(4.d0*solsc**2)
41 call holdthd(natmtot,nthd)
42 !$OMP PARALLEL DO DEFAULT(SHARED) &
43 !$OMP PRIVATE(vr,dvr,is,nr,nri) &
44 !$OMP PRIVATE(ir,irc,rm) &
45 !$OMP SCHEDULE(DYNAMIC) &
46 !$OMP NUM_THREADS(nthd)
47 do ias=1,natmtot
48  is=idxis(ias)
49  nr=nrmt(is)
50  nri=nrmti(is)
51 ! radial derivative of the spherical part of the Kohn-Sham potential
52  call rfmtlm(1,nr,nri,vsmt(:,ias),vr)
53  call splined(nr,wcrmt(:,:,is),vr,dvr)
54  do ir=1,nr,lradstp
55  irc=(ir-1)/lradstp+1
56  rm=1.d0-2.d0*cso*vr(ir)
57  socfr(irc,ias)=cso*dvr(ir)/(rsp(ir,is)*rm**2)
58  end do
59 end do
60 !$OMP END PARALLEL DO
61 call freethd(nthd)
62 
63 contains
64 
65 pure subroutine splined(n,wc,f,df)
66 implicit none
67 ! arguments
68 integer, intent(in) :: n
69 real(8), intent(in) :: wc(12,n),f(n)
70 real(8), intent(out) :: df(n)
71 ! local variables
72 integer i
73 df(1)=wc(1,1)*f(1)+wc(2,1)*f(2)+wc(3,1)*f(3)+wc(4,1)*f(4)
74 df(2)=wc(1,2)*f(1)+wc(2,2)*f(2)+wc(3,2)*f(3)+wc(4,2)*f(4)
75 do i=3,n-2
76  df(i)=wc(1,i)*f(i-1)+wc(2,i)*f(i)+wc(3,i)*f(i+1)+wc(4,i)*f(i+2)
77 end do
78 i=n-1
79 df(i)=wc(1,i)*f(n-3)+wc(2,i)*f(n-2)+wc(3,i)*f(n-1)+wc(4,i)*f(n)
80 df(n)=wc(1,n)*f(n-3)+wc(2,n)*f(n-2)+wc(3,n)*f(n-1)+wc(4,n)*f(n)
81 end subroutine
82 
83 end subroutine
84 !EOC
85 
real(8) socscf
Definition: modmain.f90:232
Definition: modomp.f90:6
subroutine gensocfr
Definition: gensocfr.f90:10
pure subroutine rfmtlm(lm, nr, nri, rfmt, fr)
Definition: rfmtlm.f90:10
integer lradstp
Definition: modmain.f90:171
real(8) solsc
Definition: modmain.f90:1240
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
logical spinorb
Definition: modmain.f90:230
pure subroutine splined(n, wc, f, df)
Definition: gensocfr.f90:66
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
real(8), dimension(:,:), allocatable rsp
Definition: modmain.f90:135
integer natmtot
Definition: modmain.f90:40
real(8), dimension(:,:,:), allocatable wcrmt
Definition: modmain.f90:187
real(8), dimension(:,:), pointer, contiguous vsmt
Definition: modmain.f90:647
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
real(8), parameter y00
Definition: modmain.f90:1223
real(8), dimension(:,:), allocatable socfr
Definition: modmain.f90:668
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150