The Elk Code
 
Loading...
Searching...
No Matches
elfplot.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: elfplot
8! !INTERFACE:
9subroutine elfplot
10! !USES:
11use modmain
12! !DESCRIPTION:
13! Outputs the electron localisation function (ELF) for 1D, 2D or 3D plotting.
14! The spin-averaged ELF is given by
15! $$ f_{\rm ELF}({\bf r})=\frac{1}{1+[D({\bf r})/D^0({\bf r})]^2}, $$
16! where
17! $$ D({\bf r})=\frac{1}{2}\left(\tau({\bf r})-\frac{1}{4}
18! \frac{[\nabla n({\bf r})]^2}{n({\bf r})}\right) $$
19! and
20! $$ \tau({\bf r})=\sum_{i=1}^N \left|\nabla\Psi_i({\bf r})
21! \right|^2 $$
22! is the spin-averaged kinetic energy density from the spinor wavefunctions.
23! The function $D^0$ is the kinetic energy density for the homogeneous
24! electron gas evaluated for $n({\bf r})$:
25! $$ D^0({\bf r})=\frac{3}{5}(6\pi^2)^{2/3}\left(\frac{n({\bf r})}{2}
26! \right)^{5/3}. $$
27! The ELF is useful for the topological classification of bonding. See for
28! example T. Burnus, M. A. L. Marques and E. K. U. Gross [Phys. Rev. A 71,
29! 10501 (2005)].
30!
31! !REVISION HISTORY:
32! Created September 2003 (JKD)
33! Fixed bug found by F. Wagner (JKD)
34!EOP
35!BOC
36implicit none
37! local variables
38integer ik,is,ias
39integer nr,nri,ir
40integer ig,ifg,i
41real(8) r,t1,t2
42! allocatable arrays
43real(8), allocatable :: gwf2mt(:,:),gwf2ir(:)
44real(8), allocatable :: rfmt(:),grfmt(:,:)
45real(8), allocatable :: rfir(:),grfir(:)
46real(8), allocatable :: elfmt(:,:),elfir(:)
47complex(8), allocatable :: zfft1(:),zfft2(:)
48! initialise universal variables
49call init0
50call init1
51! allocate local arrays
52allocate(gwf2mt(npmtmax,natmtot),gwf2ir(ngtot))
53allocate(rfmt(npmtmax),grfmt(npmtmax,3))
54allocate(rfir(ngtot),grfir(ngtot))
55allocate(elfmt(npmtmax,natmtot),elfir(ngtot))
56allocate(zfft1(nfgrz),zfft2(nfgrz))
57! read density and potentials from file
58call readstate
59! generate the core wavefunctions and densities
60call gencore
61! find the new linearisation energies
62call linengy
63! generate the APW radial functions
64call genapwfr
65! generate the local-orbital radial functions
66call genlofr
67! get the occupation numbers from file
68call readoccsv
69! set the wavefunction gradient squared to zero
70gwf2mt(:,:)=0.d0
71gwf2ir(:)=0.d0
72do ik=1,nkpt
73! add to the valence wavefunction gradient squared
74 call gradwf2(ik,gwf2mt,gwf2ir)
75end do
76! convert muffin-tin gradient squared to spherical harmonics
77do ias=1,natmtot
78 is=idxis(ias)
79 call rfshtip(nrcmt(is),nrcmti(is),gwf2mt(:,ias))
80end do
81! symmetrise the wavefunction gradient squared
83 gwf2mt,gwf2ir)
84! convert back to spherical coordinates
85do ias=1,natmtot
86 is=idxis(ias)
87 call rbshtip(nrcmt(is),nrcmti(is),gwf2mt(:,ias))
88end do
89! convert from coarse to fine muffin-tin radial mesh
90call rfmtctof(gwf2mt)
91! convert from coarse to fine interstitial grid
92call rfirctof(gwf2ir,gwf2ir)
93! add core wavefunction gradient squared
94call gradwfcr2(gwf2mt)
95!------------------------!
96! muffin-tin ELF !
97!------------------------!
98do ias=1,natmtot
99 is=idxis(ias)
100 nr=nrmt(is)
101 nri=nrmti(is)
102! convert rho from spherical harmonics to spherical coordinates
103 call rbsht(nr,nri,rhomt(:,ias),rfmt)
104! compute the gradient of the density
105 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rhomt(:,ias),npmtmax,grfmt)
106! convert gradient to spherical coordinates
107 do i=1,3
108 call rbshtip(nr,nri,grfmt(:,i))
109 end do
110 do i=1,npmt(is)
111 r=abs(rfmt(i))
112! square of gradient of rho
113 t1=grfmt(i,1)**2+grfmt(i,2)**2+grfmt(i,3)**2
114! D for inhomogeneous density
115 t1=(1.d0/2.d0)*(gwf2mt(i,ias)-(1.d0/4.d0)*t1/r)
116! D0 for uniform electron gas
117 t2=(3.d0/5.d0)*((6.d0*pi**2)**(2.d0/3.d0))*(r/2.d0)**(5.d0/3.d0)
118! ELF function
119 elfmt(i,ias)=1.d0/(1.d0+(t1/t2)**2)
120 end do
121! convert ELF from spherical coordinates to spherical harmonics
122 call rfshtip(nr,nri,elfmt(:,ias))
123end do
124!--------------------------!
125! interstitial ELF !
126!--------------------------!
127! Fourier transform density to G-space
128call rzfftifc(3,ngridg,-1,rhoir,zfft1)
129! calculate the square of gradient of rho
130grfir(:)=0.d0
131do i=1,3
132 do ifg=1,nfgrz
133 ig=igrzf(ifg)
134 if (ig <= ngvc) then
135 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
136 else
137 zfft2(ifg)=0.d0
138 end if
139 end do
140! Fourier transform gradient to real-space
141 call rzfftifc(3,ngridg,1,rfir,zfft2)
142 do ir=1,ngtot
143 grfir(ir)=grfir(ir)+rfir(ir)**2
144 end do
145end do
146do ir=1,ngtot
147 r=abs(rhoir(ir))
148! D for inhomogeneous density
149 t1=(1.d0/2.d0)*(gwf2ir(ir)-(1.d0/4.d0)*grfir(ir)/r)
150! D0 for homogeneous electron gas
151 t2=(3.d0/5.d0)*((6.d0*pi**2)**(2.d0/3.d0))*(r/2.d0)**(5.d0/3.d0)
152! ELF function
153 elfir(ir)=1.d0/(1.d0+(t1/t2)**2)
154end do
155! plot the ELF to file
156select case(task)
157case(51)
158 open(50,file='ELF1D.OUT',form='FORMATTED')
159 open(51,file='ELFLINES.OUT',form='FORMATTED')
160 call plot1d(50,51,1,elfmt,elfir)
161 close(50)
162 close(51)
163 write(*,*)
164 write(*,'("Info(elfplot):")')
165 write(*,'(" 1D ELF plot written to ELF1D.OUT")')
166 write(*,'(" vertex location lines written to ELFLINES.OUT")')
167case(52)
168 open(50,file='ELF2D.OUT',form='FORMATTED')
169 call plot2d(.false.,50,1,elfmt,elfir)
170 close(50)
171 write(*,*)
172 write(*,'("Info(elfplot): 2D ELF plot written to ELF2D.OUT")')
173case(53)
174 open(50,file='ELF3D.OUT',form='FORMATTED')
175 call plot3d(50,1,elfmt,elfir)
176 close(50)
177 write(*,*)
178 write(*,'("Info(elfplot): 3D ELF plot written to ELF3D.OUT")')
179end select
180deallocate(gwf2mt,gwf2ir,rfmt,grfmt,rfir,grfir)
181deallocate(elfmt,elfir,zfft1,zfft2)
182end subroutine
183!EOC
184
subroutine elfplot
Definition elfplot.f90:10
subroutine genapwfr
Definition genapwfr.f90:10
subroutine gencore
Definition gencore.f90:10
subroutine genlofr
Definition genlofr.f90:10
subroutine gradrfmt(nr, nri, ri, wcr, rfmt, ld, grfmt)
Definition gradrfmt.f90:10
subroutine gradwf2(ik, gwf2mt, gwf2ir)
Definition gradwf2.f90:7
subroutine gradwfcr2(gwf2mt)
Definition gradwfcr2.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine linengy
Definition linengy.f90:10
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), dimension(:,:,:), allocatable wcrmt
Definition modmain.f90:187
integer nfgrzc
Definition modmain.f90:414
integer ngtot
Definition modmain.f90:390
integer, dimension(3) ngridg
Definition modmain.f90:386
real(8), parameter pi
Definition modmain.f90:1229
integer nfgrz
Definition modmain.f90:412
integer, dimension(3) ngdgc
Definition modmain.f90:388
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
real(8), dimension(:), pointer, contiguous rhoir
Definition modmain.f90:614
integer, dimension(:), allocatable igrzf
Definition modmain.f90:416
complex(8), parameter zi
Definition modmain.f90:1239
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer ngtc
Definition modmain.f90:392
integer nkpt
Definition modmain.f90:461
integer natmtot
Definition modmain.f90:40
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer task
Definition modmain.f90:1298
integer npmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
integer, dimension(:), allocatable igrzfc
Definition modmain.f90:418
integer ngvc
Definition modmain.f90:398
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
subroutine plot1d(fnum1, fnum2, nf, rfmt, rfir)
Definition plot1d.f90:10
subroutine plot2d(tproj, fnum, nf, rfmt, rfir)
Definition plot2d.f90:10
subroutine plot3d(fnum, nf, rfmt, rfir)
Definition plot3d.f90:10
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition rbsht.f90:7
subroutine rbshtip(nr, nri, rfmt)
Definition rbshtip.f90:7
subroutine readoccsv
Definition readoccsv.f90:7
subroutine readstate
Definition readstate.f90:10
subroutine rfirctof(rfirc, rfir)
Definition rfirctof.f90:10
subroutine rfmtctof(rfmt)
Definition rfmtctof.f90:10
subroutine rfshtip(nr, nri, rfmt)
Definition rfshtip.f90:7
subroutine symrf(nrmt_, nrmti_, npmt_, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld, rfmt, rfir)
Definition symrf.f90:11
subroutine rzfftifc(nd, n, sgn, r, z)