The Elk Code
rhoplot.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: rhoplot
8 ! !INTERFACE:
9 subroutine rhoplot
10 ! !USES:
11 use modmain
12 ! !DESCRIPTION:
13 ! Outputs the charge density, read in from {\tt STATE.OUT}, for 1D, 2D or 3D
14 ! plotting.
15 !
16 ! !REVISION HISTORY:
17 ! Created June 2003 (JKD)
18 !EOP
19 !BOC
20 implicit none
21 ! initialise universal variables
22 call init0
23 ! read density from file
24 call readstate
25 ! write the density plot to file
26 select case(task)
27 case(31)
28  open(50,file='RHO1D.OUT',form='FORMATTED',action='WRITE')
29  open(51,file='RHOLINES.OUT',form='FORMATTED',action='WRITE')
30  call plot1d(50,51,1,rhomt,rhoir)
31  close(50)
32  close(51)
33  write(*,*)
34  write(*,'("Info(rhoplot):")')
35  write(*,'(" 1D density plot written to RHO1D.OUT")')
36  write(*,'(" vertex location lines written to RHOLINES.OUT")')
37 case(32)
38  open(50,file='RHO2D.OUT',form='FORMATTED',action='WRITE')
39  call plot2d(.false.,50,1,rhomt,rhoir)
40  close(50)
41  write(*,*)
42  write(*,'("Info(rhoplot): 2D density plot written to RHO2D.OUT")')
43 case(33)
44  open(50,file='RHO3D.OUT',form='FORMATTED',action='WRITE')
45  call plot3d(50,1,rhomt,rhoir)
46  close(50)
47  write(*,*)
48  write(*,'("Info(rhoplot): 3D density plot written to RHO3D.OUT")')
49 end select
50 end subroutine
51 !EOC
52 
integer task
Definition: modmain.f90:1299
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
subroutine plot3d(fnum, nf, rfmt, rfir)
Definition: plot3d.f90:10
subroutine plot1d(fnum1, fnum2, nf, rfmt, rfir)
Definition: plot1d.f90:10
subroutine plot2d(tproj, fnum, nf, rfmt, rfir)
Definition: plot2d.f90:10
subroutine rhoplot
Definition: rhoplot.f90:10
subroutine readstate
Definition: readstate.f90:10
subroutine init0
Definition: init0.f90:10