The Elk Code
plot1d.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: plot1d
8 ! !INTERFACE:
9 subroutine plot1d(fnum1,fnum2,nf,rfmt,rfir)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! fnum1 : plot file number (in,integer)
14 ! fnum2 : vertex location file number (in,integer)
15 ! nf : number of functions (in,integer)
16 ! rfmt : real muffin-tin function (in,real(npmtmax,natmtot,nf))
17 ! rfir : real intersitial function (in,real(ngtot,nf))
18 ! !DESCRIPTION:
19 ! Produces a 1D plot of the real functions contained in arrays {\tt rfmt} and
20 ! {\tt rfir} along the lines connecting the vertices in the global array
21 ! {\tt vvlp1d}. See routine {\tt rfpts}.
22 !
23 ! !REVISION HISTORY:
24 ! Created June 2003 (JKD)
25 !EOP
26 !BOC
27 implicit none
28 ! arguments
29 integer, intent(in) :: fnum1,fnum2,nf
30 real(8), intent(in) :: rfmt(npmtmax,natmtot,nf),rfir(ngtot,nf)
31 ! local variables
32 integer jf,ip,iv
33 real(8) fmin,fmax,t1
34 ! allocatable arrays
35 real(8), allocatable :: fp(:,:)
36 if ((nf < 1).or.(nf > 4)) then
37  write(*,*)
38  write(*,'("Error(plot1d): invalid number of functions : ",I8)') nf
39  write(*,*)
40  stop
41 end if
42 allocate(fp(npp1d,nf))
43 ! connect the 1D plotting vertices
45 do jf=1,nf
46 ! evaluate real function at each point
47  call rfpts(npp1d,vplp1d,rfmt(:,:,jf),rfir(:,jf),fp(:,jf))
48 end do
49 do ip=ip01d,npp1d
50 ! write the point distances and function to file
51  write(fnum1,'(5G18.10)') dpp1d(ip),(fp(ip,jf),jf=1,nf)
52 end do
53 ! write the vertex location lines
54 fmin=minval(fp(:,:))
55 fmax=maxval(fp(:,:))
56 t1=0.5d0*(fmax-fmin)
57 fmin=fmin-t1
58 fmax=fmax+t1
59 do iv=1,nvp1d
60  write(fnum2,'(2G18.10)') dvp1d(iv),fmin
61  write(fnum2,'(2G18.10)') dvp1d(iv),fmax
62  write(fnum2,*)
63 end do
64 deallocate(fp)
65 end subroutine
66 !EOC
67 
real(8), dimension(:), allocatable dpp1d
Definition: modmain.f90:1126
subroutine rfpts(np, vpl, rfmt, rfir, fp)
Definition: rfpts.f90:7
subroutine plot1d(fnum1, fnum2, nf, rfmt, rfir)
Definition: plot1d.f90:10
real(8), dimension(:), allocatable dvp1d
Definition: modmain.f90:1122
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
subroutine plotpt1d(cvec, nv, np, vvl, vpl, dv, dp)
Definition: plotpt1d.f90:10
integer ip01d
Definition: modmain.f90:1118
real(8), dimension(:,:), allocatable vvlp1d
Definition: modmain.f90:1120
integer nvp1d
Definition: modmain.f90:1114
integer npp1d
Definition: modmain.f90:1116
real(8), dimension(:,:), allocatable vplp1d
Definition: modmain.f90:1124