The Elk Code
 
Loading...
Searching...
No Matches
gencore.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: gencore
8! !INTERFACE:
9subroutine gencore
10! !USES:
11use modmain
12use modomp
13! !DESCRIPTION:
14! Computes the core radial wavefunctions, eigenvalues and densities. The
15! radial Dirac equation is solved in the spherical part of the Kohn-Sham
16! potential to which the atomic potential has been appended for
17! $r>R_{\rm MT}$. In the case of spin-polarised calculations, and when
18! {\tt spincore} is set to {\tt .true.}, the Dirac equation is solved in the
19! spin-up and -down potentials created from the Kohn-Sham scalar potential and
20! magnetic field magnitude, with the occupancy divided equally between up and
21! down. The up and down densities determined in this way are added to both the
22! scalar density and the magnetisation in the routine {\tt rhocore}. Note
23! that this procedure is a simple, but inexact, approach to solving the radial
24! Dirac equation in a magnetic field.
25!
26! !REVISION HISTORY:
27! Created April 2003 (JKD)
28! Added polarised cores, November 2009 (JKD)
29! Fixed race condition, February 2025 (JKD)
30!EOP
31!BOC
32implicit none
33! local variables
34integer ist,ispn,idm
35integer is,ia,ja,ias,jas
36integer nr,nri,nrs,nthd
37real(8) v(3),t1
38! automatic arrays
39real(8) vr(nrspmax),br(nrmtmax),fr(nrmtmax),eval(nstspmax)
40! external functions
41real(8), external :: rfmtint
42! loop over species and atoms
43call holdthd(natmtot,nthd)
44!$OMP PARALLEL DO DEFAULT(SHARED) &
45!$OMP PRIVATE(vr,br,fr,eval) &
46!$OMP PRIVATE(is,ia,nr,nri,nrs,idm) &
47!$OMP PRIVATE(v,t1,ispn,ist,ja,jas) &
48!$OMP SCHEDULE(DYNAMIC) &
49!$OMP NUM_THREADS(nthd)
50do ias=1,natmtot
51 is=idxis(ias)
52 ia=idxia(ias)
53! perform calculation for only the first equivalent atom
54 if (.not.tfeqat(ia,is)) cycle
55 nr=nrmt(is)
56 nri=nrmti(is)
57 nrs=nrsp(is)
58! Kohn-Sham magnetic field for spin-polarised core
59 if (spincore) then
60! compute the averaged direction of the magnetisation
61 do idm=1,ndmag
62 v(idm)=rfmtint(nr,nri,wr2mt(:,is),magmt(:,ias,idm))
63 end do
64! normalise
65 if (ncmag) then
66 t1=sqrt(v(1)**2+v(2)**2+v(3)**2)
67 else
68 t1=abs(v(1))
69 end if
70 if (t1 > 1.d-10) v(1:ndmag)=v(1:ndmag)/t1
71! determine the component of the field along the averaged direction
72 br(1:nr)=0.d0
73 do idm=1,ndmag
74! extract the spherical (l = m = 0) component of B_xc
75 call rfmtlm(1,nr,nri,bxcmt(:,ias,idm),fr)
76 t1=v(idm)*y00
77 br(1:nr)=br(1:nr)+t1*fr(1:nr)
78 end do
79 end if
80! loop over spin channels
81 do ispn=1,nspncr
82! use the spherical part of the crystal Kohn-Sham potential
83 call rfmtlm(1,nr,nri,vsmt(:,ias),vr)
84 vr(1:nr)=vr(1:nr)*y00
85! spin-up and -down potentials for polarised core
86 if (spincore) then
87 if (ispn == 1) then
88 vr(1:nr)=vr(1:nr)+br(1:nr)
89 else
90 vr(1:nr)=vr(1:nr)-br(1:nr)
91 end if
92 end if
93! append the Kohn-Sham potential from the atomic calculation for r > R_MT
94 t1=vr(nr)-vrsp(nr,is)
95 vr(nr+1:nrs)=vrsp(nr+1:nrs,is)+t1
96 rhocr(1:nr,ias,ispn)=0.d0
97 do ist=1,nstsp(is)
98 if (spcore(ist,is)) then
99! solve the Dirac equation
100 eval(ist)=evalcr(ist,ias)
101 call rdirac(solsc,nsp(ist,is),lsp(ist,is),ksp(ist,is),nrs,rsp(:,is),vr,&
102 eval(ist),rwfcr(:,1,ist,ias),rwfcr(:,2,ist,ias))
103 if (spincore) then
104! use the spin-averaged eigenvalue for the polarised core
105 if (ispn == 1) then
106 evalcr(ist,ias)=eval(ist)
107 else
108 evalcr(ist,ias)=0.5d0*(evalcr(ist,ias)+eval(ist))
109 end if
110 t1=0.5d0*occcr(ist,ias)
111 else
112 evalcr(ist,ias)=eval(ist)
113 t1=occcr(ist,ias)
114 end if
115! add to the core density
116 rhocr(1:nr,ias,ispn)=rhocr(1:nr,ias,ispn) &
117 +t1*(rwfcr(1:nr,1,ist,ias)**2+rwfcr(1:nr,2,ist,ias)**2)
118 end if
119 end do
120 rhocr(1:nr,ias,ispn)=rhocr(1:nr,ias,ispn)*rlmt(1:nr,-2,is)*y00
121! end loop over spin channels
122 end do
123! copy to equivalent atoms
124 do ja=1,natoms(is)
125 if (eqatoms(ia,ja,is).and.(ia /= ja)) then
126 jas=idxas(ja,is)
127 do ist=1,nstsp(is)
128 if (spcore(ist,is)) then
129 evalcr(ist,jas)=evalcr(ist,ias)
130 rwfcr(1:nrs,1:2,ist,jas)=rwfcr(1:nrs,1:2,ist,ias)
131 end if
132 end do
133 rhocr(1:nr,jas,1:nspncr)=rhocr(1:nr,ias,1:nspncr)
134 end if
135 end do
136! end loop over species and atoms
137end do
138!$OMP END PARALLEL DO
139call freethd(nthd)
140end subroutine
141!EOC
142
subroutine gencore
Definition gencore.f90:10
integer, dimension(maxstsp, maxspecies) lsp
Definition modmain.f90:123
real(8), dimension(:,:), allocatable evalcr
Definition modmain.f90:934
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), parameter y00
Definition modmain.f90:1233
logical, dimension(:,:), allocatable tfeqat
Definition modmain.f90:372
integer nspncr
Definition modmain.f90:942
real(8), dimension(:,:,:), allocatable bxcmt
Definition modmain.f90:636
logical ncmag
Definition modmain.f90:240
integer, dimension(maxstsp, maxspecies) nsp
Definition modmain.f90:121
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition modmain.f90:616
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
logical, dimension(maxstsp, maxspecies) spcore
Definition modmain.f90:127
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
logical spincore
Definition modmain.f90:940
real(8), dimension(:,:), allocatable rsp
Definition modmain.f90:135
integer, dimension(maxatoms *maxspecies) idxia
Definition modmain.f90:45
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer, dimension(maxspecies) nstsp
Definition modmain.f90:113
real(8), dimension(:,:,:,:), allocatable rwfcr
Definition modmain.f90:936
integer natmtot
Definition modmain.f90:40
integer, dimension(maxstsp, maxspecies) ksp
Definition modmain.f90:125
real(8), dimension(:,:), pointer, contiguous vsmt
Definition modmain.f90:649
real(8) solsc
Definition modmain.f90:1252
real(8), dimension(:,:), allocatable vrsp
Definition modmain.f90:139
logical, dimension(:,:,:), allocatable eqatoms
Definition modmain.f90:370
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
integer, dimension(maxspecies) nrsp
Definition modmain.f90:107
real(8), dimension(:,:,:), allocatable rhocr
Definition modmain.f90:938
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
integer ndmag
Definition modmain.f90:238
real(8), dimension(:,:), allocatable occcr
Definition modmain.f90:932
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine rdirac(sol, n, l, k, nr, r, vr, eval, g0, f0)
Definition rdirac.f90:10
pure real(8) function rfmtint(nr, nri, wr, rfmt)
Definition rfmtint.f90:7
pure subroutine rfmtlm(lm, nr, nri, rfmt, fr)
Definition rfmtlm.f90:7