The Elk Code
 
Loading...
Searching...
No Matches
genlofr.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 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!BOP
7! !ROUTINE: genlofr
8! !INTERFACE:
9subroutine genlofr
10! !USES:
11use modmain
12use modomp
13! !DESCRIPTION:
14! Generates the local-orbital radial functions. This is done by integrating
15! the scalar relativistic Schr\"{o}dinger equation (or its energy deriatives)
16! at the current linearisation energies using the spherical part of the
17! Kohn-Sham potential. For each local-orbital, a linear combination of
18! {\tt lorbord} radial functions is constructed such that its radial
19! derivatives up to order ${\tt lorbord}-1$ are zero at the muffin-tin radius.
20! This function is normalised and the radial Hamiltonian applied to it. The
21! results are stored in the global array {\tt lofr}.
22!
23! !REVISION HISTORY:
24! Created March 2003 (JKD)
25! Copied to equivalent atoms, February 2010 (A. Kozhevnikov and JKD)
26!EOP
27!BOC
28implicit none
29! local variables
30integer is,ia,ja,ias,jas
31integer nr,nri,iro,ir,i,j
32integer i0,i1,nn,l,info
33integer ilo,jlo,io,jo,nthd
34real(8) e,t1
35! automatic arrays
36integer ipiv(nplorb)
37real(8) vr(nrmtmax),p0(nrmtmax,lorbordmax),p1(nrmtmax)
38real(8) q0(nrmtmax),q1(nrmtmax),ep0(nrmtmax,lorbordmax)
39real(8) p0c(nrmtmax,nlomax),ep0c(nrmtmax,nlomax)
40real(8) a(nplorb,nplorb),b(nplorb)
41! external functions
42real(8), external :: polynm
43! loop over all species and atoms
44call holdthd(natmtot,nthd)
45!$OMP PARALLEL DO DEFAULT(SHARED) &
46!$OMP PRIVATE(ipiv,vr,p0,p1,q0,q1) &
47!$OMP PRIVATE(ep0,p0c,ep0c,a,b) &
48!$OMP PRIVATE(is,ia,nr,nri,iro,i0,i1) &
49!$OMP PRIVATE(i,j,ilo,jlo,l,io,jo) &
50!$OMP PRIVATE(e,nn,t1,ir,info,ja,jas) &
51!$OMP SCHEDULE(DYNAMIC) &
52!$OMP NUM_THREADS(nthd)
53do ias=1,natmtot
54 is=idxis(ias)
55 ia=idxia(ias)
56! perform calculation for only the first equivalent atom
57 if (.not.tfeqat(ia,is)) cycle
58 nr=nrmt(is)
59 nri=nrmti(is)
60 iro=nri+1
61! use spherical part of potential
62 i1=lmmaxi*(nri-1)+1
63 vr(1:nri)=vsmt(1:i1:lmmaxi,ias)*y00
64 i0=i1+lmmaxi
65 i1=lmmaxo*(nr-iro)+i0
66 vr(iro:nr)=vsmt(i0:i1:lmmaxo,ias)*y00
67! loop over local-orbitals in ascending energy order
68 do i=1,nlorb(is)
69 ilo=idxelo(i,is)
70 l=lorbl(ilo,is)
71 do jo=1,lorbord(ilo,is)
72! linearisation energy accounting for energy derivative
73 e=lorbe(jo,ilo,ias)+lorbdm(jo,ilo,is)*delorb
74! integrate the radial Schrodinger equation
75 call rschrodint(solsc,l,e,nr,rlmt(:,1,is),vr,nn,p0(:,jo),p1,q0,q1)
76! divide the radial wavefunction by r
77 p0(1:nr,jo)=p0(1:nr,jo)*rlmt(1:nr,-1,is)
78! multiply by the linearisation energy
79 ep0(1:nr,jo)=e*p0(1:nr,jo)
80! set up the matrix of radial derivatives
81 a(1,jo)=p0(nr,jo)
82 ir=nr-nplorb+1
83 do io=2,lorbord(ilo,is)
84 a(io,jo)=polynm(io-1,nplorb,rlmt(ir,1,is),p0(ir,jo),rmt(is))
85 end do
86 end do
87! set up the target vector
88 b(1:nplorb)=0.d0
89 b(lorbord(ilo,is))=1.d0
90! solve the linear system
91 call dgesv(lorbord(ilo,is),1,a,nplorb,ipiv,b,nplorb,info)
92 if (info /= 0) call errstp
93! generate linear combination of radial functions
94 p0c(1:nr,ilo)=0.d0
95 ep0c(1:nr,ilo)=0.d0
96 do io=1,lorbord(ilo,is)
97 t1=b(io)
98 p0c(1:nr,ilo)=p0c(1:nr,ilo)+t1*p0(1:nr,io)
99 ep0c(1:nr,ilo)=ep0c(1:nr,ilo)+t1*ep0(1:nr,io)
100 end do
101! normalise radial functions
102 t1=sum(wr2mt(1:nr,is)*p0c(1:nr,ilo)**2)
103 t1=1.d0/sqrt(abs(t1))
104 p0c(1:nr,ilo)=t1*p0c(1:nr,ilo)
105 ep0c(1:nr,ilo)=t1*ep0c(1:nr,ilo)
106! subtract linear combination of previous local-orbitals with same l
107 do j=1,i-1
108 jlo=idxelo(j,is)
109 if (lorbl(jlo,is) == l) then
110 t1=-sum(wr2mt(1:nr,is)*p0c(1:nr,ilo)*p0c(1:nr,jlo))
111 p0c(1:nr,ilo)=p0c(1:nr,ilo)+t1*p0c(1:nr,jlo)
112 ep0c(1:nr,ilo)=ep0c(1:nr,ilo)+t1*ep0c(1:nr,jlo)
113 end if
114 end do
115! normalise radial functions
116 if (i > 1) then
117 t1=abs(sum(wr2mt(1:nr,is)*p0c(1:nr,ilo)**2))
118 if (t1 < 1.d-25) call errstp
119 t1=1.d0/sqrt(t1)
120 p0c(1:nr,ilo)=t1*p0c(1:nr,ilo)
121 ep0c(1:nr,ilo)=t1*ep0c(1:nr,ilo)
122 end if
123! store in global array
124 lofr(1:nr,1,ilo,ias)=p0c(1:nr,ilo)
125 lofr(1:nr,2,ilo,ias)=ep0c(1:nr,ilo)
126 end do
127! copy to equivalent atoms
128 do ja=1,natoms(is)
129 if (eqatoms(ia,ja,is).and.(ia /= ja)) then
130 jas=idxas(ja,is)
131 do ilo=1,nlorb(is)
132 lofr(1:nr,1:2,ilo,jas)=lofr(1:nr,1:2,ilo,ias)
133 end do
134 end if
135 end do
136! end loop over atoms and species
137end do
138!$OMP END PARALLEL DO
139call freethd(nthd)
140! make single-precision copy if required
141if (tfr_sp) then
142 do ias=1,natmtot
143 is=idxis(ias)
144 do ilo=1,nlorb(is)
145 lofr_sp(1:nrcmt(is),ilo,ias)=lofr(1:nrmt(is):lradstp,1,ilo,ias)
146 end do
147 end do
148end if
149return
150
151contains
152
153subroutine errstp
154implicit none
155write(*,*)
156write(*,'("Error(genlofr): degenerate local-orbital radial functions")')
157write(*,'(" for species ",I4)') is
158write(*,'(" atom ",I4)') ia
159write(*,'(" and local-orbital ",I4)') ilo
160write(*,*)
161stop
162end subroutine
163
164end subroutine
165!EOC
166
subroutine errstp
Definition genlofr.f90:154
subroutine genlofr
Definition genlofr.f90:10
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), parameter y00
Definition modmain.f90:1233
logical, dimension(:,:), allocatable tfeqat
Definition modmain.f90:372
real(8) delorb
Definition modmain.f90:802
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
real(8), dimension(maxspecies) rmt
Definition modmain.f90:162
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
integer lmmaxi
Definition modmain.f90:207
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(maxatoms *maxspecies) idxia
Definition modmain.f90:45
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer lradstp
Definition modmain.f90:171
integer, dimension(maxlorbord, maxlorb, maxspecies) lorbdm
Definition modmain.f90:810
real(4), dimension(:,:,:), allocatable lofr_sp
Definition modmain.f90:816
integer natmtot
Definition modmain.f90:40
integer, dimension(maxlorb, maxspecies) lorbord
Definition modmain.f90:792
logical tfr_sp
Definition modmain.f90:818
integer, dimension(maxlorb, maxspecies) idxelo
Definition modmain.f90:806
real(8), dimension(:,:), pointer, contiguous vsmt
Definition modmain.f90:649
real(8) solsc
Definition modmain.f90:1252
logical, dimension(:,:,:), allocatable eqatoms
Definition modmain.f90:370
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
integer, dimension(maxspecies) nlorb
Definition modmain.f90:786
integer lmmaxo
Definition modmain.f90:203
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
real(8), dimension(:,:,:,:), allocatable lofr
Definition modmain.f90:814
integer, dimension(maxlorb, maxspecies) lorbl
Definition modmain.f90:796
real(8), dimension(:,:,:), allocatable lorbe
Definition modmain.f90:808
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
pure real(8) function polynm(m, np, xa, ya, x)
Definition polynm.f90:10
pure subroutine rschrodint(sol, l, e, nr, r, vr, nn, p0, p1, q0, q1)