File: plot1d.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (70 lines) | stat: -rw-r--r-- 1,937 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: plot1d
! !INTERFACE:
subroutine plot1d(fnum1,fnum2,nf,rfmt,rfir)
! !USES:
use modmain
! !INPUT/OUTPUT PARAMETERS:
!   fnum1 : plot file number (in,integer)
!   fnum2 : vertex location file number (in,integer)
!   nf    : number of functions (in,integer)
!   rfmt  : real muffin-tin function (in,real(npmtmax,natmtot,nf))
!   rfir  : real intersitial function (in,real(ngtot,nf))
! !DESCRIPTION:
!   Produces a 1D plot of the real functions contained in arrays {\tt rfmt} and
!   {\tt rfir} along the lines connecting the vertices in the global array
!   {\tt vvlp1d}. See routine {\tt rfplot}.
!
! !REVISION HISTORY:
!   Created June 2003 (JKD)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: fnum1,fnum2,nf
real(8), intent(in) :: rfmt(npmtmax,natmtot,nf),rfir(ngtot,nf)
! local variables
integer i,ip,iv
real(8) fmin,fmax,t1
! allocatable arrays
real(8), allocatable :: fp(:,:)
if ((nf.lt.1).or.(nf.gt.4)) then
  write(*,*)
  write(*,'("Error(plot1d): invalid number of functions : ",I8)') nf
  write(*,*)
  stop
end if
allocate(fp(npp1d,nf))
! connect the 1D plotting vertices
call plotpt1d(avec,nvp1d,npp1d,vvlp1d,vplp1d,dvp1d,dpp1d)
do i=1,nf
! evaluate function at each point
  call rfplot(npp1d,vplp1d,rfmt(:,:,i),rfir(:,i),fp(:,i))
end do
fmin=fp(1,1)
fmax=fp(1,1)
do ip=1,npp1d
  do i=1,nf
    fmin=min(fmin,fp(ip,i))
    fmax=max(fmax,fp(ip,i))
  end do
! write the point distances and function to file
  write(fnum1,'(5G18.10)') dpp1d(ip),(fp(ip,i),i=1,nf)
end do
! write the vertex location lines
t1=0.5d0*(fmax-fmin)
do iv=1,nvp1d
  write(fnum2,'(2G18.10)') dvp1d(iv),fmax+t1
  write(fnum2,'(2G18.10)') dvp1d(iv),fmin-t1
  write(fnum2,'("     ")')
end do
deallocate(fp)
return
end subroutine
!EOC