The Elk Code
rhosplot.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine rhosplot
7 use modmain
8 use modtddft
9 implicit none
10 ! local variables
11 integer is,ias,ispn
12 integer nr,nri,iro,ir,i
13 ! determine the static density and charge
14 call rhostatic
15 ! remove the core charge density
16 do ias=1,natmtot
17  is=idxis(ias)
18  nr=nrmt(is)
19  nri=nrmti(is)
20  iro=nri+1
21  do ispn=1,nspncr
22  i=1
23  do ir=1,nri
24  rhosmt(i,ias,:)=rhosmt(i,ias,:)-rhocr(ir,ias,ispn)
25  i=i+lmmaxi
26  end do
27  do ir=iro,nr
28  rhosmt(i,ias,:)=rhosmt(i,ias,:)-rhocr(ir,ias,ispn)
29  i=i+lmmaxo
30  end do
31  end do
32 end do
33 ! produce 1D plot of the static density
34 open(50,file='RHOS1D.OUT',form='FORMATTED',action='WRITE')
35 open(51,file='RHOSLINES.OUT',form='FORMATTED',action='WRITE')
36 call plot1d(50,51,3,rhosmt,rhosir)
37 close(50)
38 close(51)
39 write(*,*)
40 write(*,'("Info(rhosplot):")')
41 write(*,'(" 1D static density plot written to RHOS1D.OUT")')
42 write(*,'(" vertex location lines written to RHOSLINES.OUT")')
43 write(*,*)
44 write(*,'(" The core density is not included in the plot")')
45 end subroutine
46 
integer lmmaxo
Definition: modmain.f90:203
real(8), dimension(:,:,:), allocatable rhocr
Definition: modmain.f90:941
subroutine plot1d(fnum1, fnum2, nf, rfmt, rfir)
Definition: plot1d.f90:10
subroutine rhostatic
Definition: rhostatic.f90:7
integer nspncr
Definition: modmain.f90:945
real(8), dimension(:,:), allocatable rhosir
Definition: modtddft.f90:82
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
subroutine rhosplot
Definition: rhosplot.f90:7
integer natmtot
Definition: modmain.f90:40
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
real(8), dimension(:,:,:), allocatable rhosmt
Definition: modtddft.f90:82
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150