The Elk Code
plot2d.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: plot2d
8 ! !INTERFACE:
9 subroutine plot2d(tproj,fnum,nf,rfmt,rfir)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! tproj : .true. if nf=3 and the vector function should be projected onto the
14 ! 2D plotting plane axes (in,logical)
15 ! fnum : plot file number (in,integer)
16 ! nf : number of functions (in,integer)
17 ! rfmt : real muffin-tin function (in,real(npmtmax,natmtot,nf))
18 ! rfir : real intersitial function (in,real(ngtot,nf))
19 ! !DESCRIPTION:
20 ! Produces a 2D plot of the real functions contained in arrays {\tt rfmt} and
21 ! {\tt rfir} on the parallelogram defined by the corner vertices in the global
22 ! array {\tt vclp2d}. See routine {\tt rfpts}.
23 !
24 ! !REVISION HISTORY:
25 ! Created June 2003 (JKD)
26 !EOP
27 !BOC
28 implicit none
29 ! arguments
30 logical, intent(in) :: tproj
31 integer, intent(in) :: fnum,nf
32 real(8), intent(in) :: rfmt(npmtmax,natmtot,nf),rfir(ngtot,nf)
33 ! local variables
34 integer np,jf,ip
35 real(8) vpnl(3)
36 ! allocatable arrays
37 real(8), allocatable :: vpl(:,:),vppc(:,:),fp(:,:)
38 if ((nf < 1).or.(nf > 4)) then
39  write(*,*)
40  write(*,'("Error(plot2d): invalid number of functions : ",I8)') nf
41  write(*,*)
42  stop
43 end if
44 ! allocate local arrays
45 np=np2d(1)*np2d(2)
46 allocate(vpl(3,np),vppc(2,np),fp(np,nf))
47 ! generate the 2D plotting points
48 call plotpt2d(avec,ainv,vpnl,vpl,vppc)
49 ! evaluate the real functions at the grid points
50 do jf=1,nf
51  call rfpts(np,vpl,rfmt(:,:,jf),rfir(:,jf),fp(:,jf))
52 end do
53 ! project the vector function onto the 2D plotting plane axes if required
54 if (tproj.and.(nf == 3)) then
55  call proj2d(np,fp)
56 end if
57 ! write the functions to file
58 write(fnum,'(2I6," : grid size")') np2d(:)
59 do ip=1,np
60  write(fnum,'(6G18.10)') vppc(1,ip),vppc(2,ip),(fp(ip,jf),jf=1,nf)
61 end do
62 deallocate(vpl,vppc,fp)
63 end subroutine
64 !EOC
65 
subroutine plotpt2d(cvec, cinv, vpnl, vpl, vppc)
Definition: plotpt2d.f90:7
real(8), dimension(3, 3) ainv
Definition: modmain.f90:14
subroutine rfpts(np, vpl, rfmt, rfir, fp)
Definition: rfpts.f90:7
subroutine plot2d(tproj, fnum, nf, rfmt, rfir)
Definition: plot2d.f90:10
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
subroutine proj2d(np, fp)
Definition: proj2d.f90:7
integer, dimension(2) np2d
Definition: modmain.f90:1130