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 subroutine gensocfr
7 use modmain
8 use modomp
9 implicit none
10 ! local variables
11 integer is,ias,nthd
12 integer nr,nri,ir,irc
13 real(8) cso,rm
14 ! automatic arrays
15 real(8) vr(nrmtmax),dvr(nrmtmax)
16 if (.not.spinorb) return
17 ! coefficient of spin-orbit coupling
18 cso=y00*socscf/(4.d0*solsc**2)
19 call holdthd(natmtot,nthd)
20 !$OMP PARALLEL DO DEFAULT(SHARED) &
21 !$OMP PRIVATE(vr,dvr,is,nr,nri) &
22 !$OMP PRIVATE(ir,irc,rm) &
23 !$OMP SCHEDULE(DYNAMIC) &
24 !$OMP NUM_THREADS(nthd)
25 do ias=1,natmtot
26  is=idxis(ias)
27  nr=nrmt(is)
28  nri=nrmti(is)
29 ! radial derivative of the spherical part of the Kohn-Sham potential
30  call rfmtlm(1,nr,nri,vsmt(:,ias),vr)
31  call splined(nr,wcrmt(:,:,is),vr,dvr)
32  do ir=1,nr,lradstp
33  irc=(ir-1)/lradstp+1
34  rm=1.d0-2.d0*cso*vr(ir)
35  socfr(irc,ias)=cso*dvr(ir)/(rsp(ir,is)*rm**2)
36  end do
37 end do
38 !$OMP END PARALLEL DO
39 call freethd(nthd)
40 
41 contains
42 
43 pure subroutine splined(n,wc,f,df)
44 implicit none
45 ! arguments
46 integer, intent(in) :: n
47 real(8), intent(in) :: wc(12,n),f(n)
48 real(8), intent(out) :: df(n)
49 ! local variables
50 integer i
51 real(8) f1,f2,f3,f4
52 f1=f(1); f2=f(2); f3=f(3); f4=f(4)
53 df(1)=wc(1,1)*f1+wc(2,1)*f2+wc(3,1)*f3+wc(4,1)*f4
54 df(2)=wc(1,2)*f1+wc(2,2)*f2+wc(3,2)*f3+wc(4,2)*f4
55 !$OMP SIMD LASTPRIVATE(f1,f2,f3,f4) SIMDLEN(8)
56 do i=3,n-2
57  f1=f(i-1); f2=f(i); f3=f(i+1); f4=f(i+2)
58  df(i)=wc(1,i)*f1+wc(2,i)*f2+wc(3,i)*f3+wc(4,i)*f4
59 end do
60 i=n-1
61 df(i)=wc(1,i)*f1+wc(2,i)*f2+wc(3,i)*f3+wc(4,i)*f4
62 df(n)=wc(1,n)*f1+wc(2,n)*f2+wc(3,n)*f3+wc(4,n)*f4
63 end subroutine
64 
65 end subroutine
66 
real(8) socscf
Definition: modmain.f90:232
Definition: modomp.f90:6
subroutine gensocfr
Definition: gensocfr.f90:7
pure subroutine rfmtlm(lm, nr, nri, rfmt, fr)
Definition: rfmtlm.f90:10
integer lradstp
Definition: modmain.f90:171
real(8) solsc
Definition: modmain.f90:1247
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:44
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:649
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
real(8), parameter y00
Definition: modmain.f90:1230
real(8), dimension(:,:), allocatable socfr
Definition: modmain.f90:670
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150