The Elk Code
 
Loading...
Searching...
No Matches
plotu1d.f90
Go to the documentation of this file.
1
2! Copyright (C) 2019 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine plotu1d(fnum1,fnum2,nf,zfmt,zfir)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: fnum1,fnum2,nf
11complex(8), intent(in) :: zfmt(npcmtmax,natmtot,nf,nfqrz)
12complex(8), intent(in) :: zfir(ngtot,nf,nfqrz)
13! local variables
14integer jf,ip,iv
15real(8) fmin,fmax,t1
16! allocatable arrays
17real(8), allocatable :: fp(:,:)
18if ((nf < 1).or.(nf > 4)) then
19 write(*,*)
20 write(*,'("Error(plotu1d): invalid number of functions : ",I8)') nf
21 write(*,*)
22 stop
23end if
24allocate(fp(npp1d,nf))
25! connect the 1D plotting vertices
27! evaluate function at each point
28call plotulr(npp1d,vplp1d,nf,zfmt,zfir,fp)
29do ip=ip01d,npp1d
30! write the point distances and function to file
31 write(fnum1,'(5G18.10)') dpp1d(ip),(fp(ip,jf),jf=1,nf)
32end do
33! write the vertex location lines
34fmin=minval(fp(:,:))
35fmax=maxval(fp(:,:))
36t1=0.5d0*(fmax-fmin)
37fmin=fmin-t1
38fmax=fmax+t1
39do iv=1,nvp1d
40 write(fnum2,'(2G18.10)') dvp1d(iv),fmin
41 write(fnum2,'(2G18.10)') dvp1d(iv),fmax
42 write(fnum2,*)
43end do
44deallocate(fp)
45end subroutine
46
integer npp1d
Definition modmain.f90:1113
integer nvp1d
Definition modmain.f90:1111
real(8), dimension(:,:), allocatable vvlp1d
Definition modmain.f90:1117
integer ip01d
Definition modmain.f90:1115
real(8), dimension(3, 3) avec
Definition modmain.f90:12
real(8), dimension(:,:), allocatable vplp1d
Definition modmain.f90:1121
real(8), dimension(:), allocatable dvp1d
Definition modmain.f90:1119
real(8), dimension(:), allocatable dpp1d
Definition modmain.f90:1123
subroutine plotpt1d(cvec, nv, np, vvl, vpl, dv, dp)
Definition plotpt1d.f90:10
subroutine plotu1d(fnum1, fnum2, nf, zfmt, zfir)
Definition plotu1d.f90:7
subroutine plotulr(np, vpl, nf, zfmt, zfir, fp)
Definition plotulr.f90:7