The Elk Code
 
Loading...
Searching...
No Matches
genapwfr.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: genapwfr
8! !INTERFACE:
9subroutine genapwfr
10! !USES:
11use modmain
12use modomp
13! !DESCRIPTION:
14! Generates the APW radial functions. This is done by integrating the scalar
15! relativistic Schr\"{o}dinger equation (or its energy deriatives) at the
16! current linearisation energies using the spherical part of the Kohn-Sham
17! potential. The number of radial functions at each $l$-value is given by the
18! variable {\tt apword} (at the muffin-tin boundary, the APW functions have
19! continuous derivatives up to order ${\tt apword}-1$). Within each $l$, these
20! functions are orthonormalised with the Gram-Schmidt method. The radial
21! Hamiltonian is applied to the orthonormalised functions and the results are
22! stored in the global array {\tt apwfr}.
23!
24! !REVISION HISTORY:
25! Created March 2003 (JKD)
26! Copied to equivalent atoms, February 2010 (A. Kozhevnikov and JKD)
27!EOP
28!BOC
29implicit none
30! local variables
31integer is,ia,ja,ias,jas
32integer nr,nri,iro,nthd
33integer i0,i1,nn,l,io,jo
34real(8) e,t1
35! automatic arrays
36real(8) vr(nrmtmax),p0(nrmtmax,apwordmax),p1(nrmtmax),p1s(apwordmax)
37real(8) q0(nrmtmax),q1(nrmtmax),ep0(nrmtmax,apwordmax)
38! loop over all species and atoms
39call holdthd(natmtot,nthd)
40!$OMP PARALLEL DO DEFAULT(SHARED) &
41!$OMP PRIVATE(vr,p0,p1,p1s,q0,q1,ep0) &
42!$OMP PRIVATE(is,ia,nr,nri,iro,i0,i1) &
43!$OMP PRIVATE(l,io,jo,e,nn,t1,ja,jas) &
44!$OMP SCHEDULE(DYNAMIC) &
45!$OMP NUM_THREADS(nthd)
46do ias=1,natmtot
47 is=idxis(ias)
48 ia=idxia(ias)
49! perform calculation for only the first equivalent atom
50 if (.not.tfeqat(ia,is)) cycle
51 nr=nrmt(is)
52 nri=nrmti(is)
53 iro=nri+1
54! use spherical part of potential
55 i1=lmmaxi*(nri-1)+1
56 vr(1:nri)=vsmt(1:i1:lmmaxi,ias)*y00
57 i0=i1+lmmaxi
58 i1=lmmaxo*(nr-iro)+i0
59 vr(iro:nr)=vsmt(i0:i1:lmmaxo,ias)*y00
60 do l=0,lmaxapw
61 do io=1,apword(l,is)
62! linearisation energy accounting for energy derivative
63 e=apwe(io,l,ias)+apwdm(io,l,is)*deapw
64! integrate the radial Schrodinger equation
65 call rschrodint(solsc,l,e,nr,rlmt(:,1,is),vr,nn,p0(:,io),p1,q0,q1)
66! divide the radial wavefunction by r
67 p0(1:nr,io)=p0(1:nr,io)*rlmt(1:nr,-1,is)
68! multiply by the linearisation energy
69 ep0(1:nr,io)=e*p0(1:nr,io)
70! normalise radial functions
71 t1=sum(wr2mt(1:nr,is)*p0(1:nr,io)**2)
72 t1=1.d0/sqrt(abs(t1))
73 p0(1:nr,io)=t1*p0(1:nr,io)
74 p1s(io)=t1*p1(nr)
75 ep0(1:nr,io)=t1*ep0(1:nr,io)
76! subtract linear combination of previous vectors
77 do jo=1,io-1
78 t1=-sum(wr2mt(1:nr,is)*p0(1:nr,io)*p0(1:nr,jo))
79 p0(1:nr,io)=p0(1:nr,io)+t1*p0(1:nr,jo)
80 p1s(io)=p1s(io)+t1*p1s(jo)
81 ep0(1:nr,io)=ep0(1:nr,io)+t1*ep0(1:nr,jo)
82 end do
83! normalise radial functions again
84 if (io > 1) then
85 t1=sum(wr2mt(1:nr,is)*p0(1:nr,io)**2)
86 if (t1 < 1.d-25) then
87 write(*,*)
88 write(*,'("Error(genapwfr): degenerate APW radial functions")')
89 write(*,'(" for species ",I4)') is
90 write(*,'(" atom ",I4)') ia
91 write(*,'(" angular momentum ",I4)') l
92 write(*,'(" and order ",I4)') io
93 write(*,*)
94 stop
95 end if
96 t1=1.d0/sqrt(t1)
97 p0(1:nr,io)=t1*p0(1:nr,io)
98 p1s(io)=t1*p1s(io)
99 ep0(1:nr,io)=t1*ep0(1:nr,io)
100 end if
101! store in global array
102 apwfr(1:nr,1,io,l,ias)=p0(1:nr,io)
103 apwfr(1:nr,2,io,l,ias)=ep0(1:nr,io)
104! derivative at the muffin-tin surface multiplied by R_MT²/2
105 apwdfr(io,l,ias)=(p1s(io)-p0(nr,io))*rmt(is)/2.d0
106 end do
107 end do
108! copy to equivalent atoms
109 do ja=1,natoms(is)
110 if (eqatoms(ia,ja,is).and.(ia /= ja)) then
111 jas=idxas(ja,is)
112 do l=0,lmaxapw
113 do io=1,apword(l,is)
114 apwfr(1:nr,1:2,io,l,jas)=apwfr(1:nr,1:2,io,l,ias)
115 apwdfr(io,l,jas)=apwdfr(io,l,ias)
116 end do
117 end do
118 end if
119 end do
120! end loop over atoms and species
121end do
122!$OMP END PARALLEL DO
123call freethd(nthd)
124! make single-precision copy if required
125if (tfr_sp) then
126 do ias=1,natmtot
127 is=idxis(ias)
128 do l=0,lmaxapw
129 do io=1,apword(l,is)
130 apwfr_sp(1:nrcmt(is),io,l,ias)=apwfr(1:nrmt(is):lradstp,1,io,l,ias)
131 end do
132 end do
133 end do
134end if
135end subroutine
136!EOC
137
subroutine genapwfr
Definition genapwfr.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), dimension(:,:,:), allocatable apwe
Definition modmain.f90:768
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
real(4), dimension(:,:,:,:), allocatable apwfr_sp
Definition modmain.f90:776
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(0:maxlapw, maxspecies) apword
Definition modmain.f90:758
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 natmtot
Definition modmain.f90:40
integer, dimension(maxapword, 0:maxlapw, maxspecies) apwdm
Definition modmain.f90:770
integer lmaxapw
Definition modmain.f90:197
real(8), dimension(:,:,:), allocatable apwdfr
Definition modmain.f90:778
real(8), dimension(:,:,:,:,:), allocatable apwfr
Definition modmain.f90:774
logical tfr_sp
Definition modmain.f90:818
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 lmmaxo
Definition modmain.f90:203
real(8) deapw
Definition modmain.f90:764
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
pure subroutine rschrodint(sol, l, e, nr, r, vr, nn, p0, p1, q0, q1)