The Elk Code
 
Loading...
Searching...
No Matches
rhocore.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: rhocore
8! !INTERFACE:
9subroutine rhocore
10! !USES:
11use modmain
12! !DESCRIPTION:
13! Adds the core density and magnetisation to the muffin-tin functions. Also
14! computes the amount of leakage of core charge from the muffin-tin spheres
15! into the interstitial.
16!
17! !REVISION HISTORY:
18! Created April 2003 (JKD)
19! Fixed core moment direction, October 2012 (M. Meinert)
20!EOP
21!BOC
22implicit none
23! local variables
24integer ispn,idm,is,ias
25integer nr,nri,iro,ir,i,i0,i1
26real(8) v(ndmag),sm,t1
27! external functions
28real(8), external :: rfmtint
29do ias=1,natmtot
30 is=idxis(ias)
31 nr=nrmt(is)
32 nri=nrmti(is)
33 iro=nri+1
34 sm=0.d0
35! loop over spin channels
36 do ispn=1,nspncr
37! add the core density to the muffin-tin density
38 i1=lmmaxi*(nri-1)+1
39 rhomt(1:i1:lmmaxi,ias)=rhomt(1:i1:lmmaxi,ias)+rhocr(1:nri,ias,ispn)
40 i0=i1+lmmaxi
41 i1=lmmaxo*(nr-iro)+i0
42 rhomt(i0:i1:lmmaxo,ias)=rhomt(i0:i1:lmmaxo,ias)+rhocr(iro:nr,ias,ispn)
43! compute the core charge inside the muffin-tins
44 sm=sm+fourpi*y00*sum(wr2mt(1:nr,is)*rhocr(1:nr,ias,ispn))
45 end do
46! core leakage charge
47 chgcrlk(ias)=chgcr(is)-sm
48! add to the magnetisation in the case of a spin-polarised core
49 if (spincore) then
50! compute the moment in the muffin-tin
51 do idm=1,ndmag
52 v(idm)=rfmtint(nr,nri,wr2mt(:,is),magmt(:,ias,idm))
53 end do
54! normalise
55 if (ncmag) then
56 t1=sqrt(v(1)**2+v(2)**2+v(3)**2)
57 else
58 t1=abs(v(1))
59 end if
60! add the core magnetisation to the total
61 if (t1 > 1.d-10) then
62 v(1:ndmag)=v(1:ndmag)/t1
63 i=1
64 do ir=1,nri
65 t1=rhocr(ir,ias,1)-rhocr(ir,ias,2)
66 magmt(i,ias,1:ndmag)=magmt(i,ias,1:ndmag)+t1*v(1:ndmag)
67 i=i+lmmaxi
68 end do
69 do ir=iro,nr
70 t1=rhocr(ir,ias,1)-rhocr(ir,ias,2)
71 magmt(i,ias,1:ndmag)=magmt(i,ias,1:ndmag)+t1*v(1:ndmag)
72 i=i+lmmaxo
73 end do
74 end if
75 end if
76end do
77end subroutine
78!EOC
79
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), parameter y00
Definition modmain.f90:1233
integer nspncr
Definition modmain.f90:942
logical ncmag
Definition modmain.f90:240
real(8), dimension(:), allocatable chgcrlk
Definition modmain.f90:720
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition modmain.f90:616
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer lmmaxi
Definition modmain.f90:207
logical spincore
Definition modmain.f90:940
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), dimension(maxspecies) chgcr
Definition modmain.f90:716
real(8), parameter fourpi
Definition modmain.f90:1231
integer natmtot
Definition modmain.f90:40
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
integer lmmaxo
Definition modmain.f90:203
real(8), dimension(:,:,:), allocatable rhocr
Definition modmain.f90:938
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
pure real(8) function rfmtint(nr, nri, wr, rfmt)
Definition rfmtint.f90:7
subroutine rhocore
Definition rhocore.f90:10