The Elk Code
 
Loading...
Searching...
No Matches
plot3d.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2008 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: plot3d
8! !INTERFACE:
9subroutine plot3d(fnum,nf,rfmt,rfir)
10! !USES:
11use modmain
12! !INPUT/OUTPUT PARAMETERS:
13! fnum : plot file number (in,integer)
14! nf : number of functions (in,integer)
15! rfmt : real muffin-tin function (in,real(npmtmax,natmtot,nf))
16! rfir : real intersitial function (in,real(ngtot,nf))
17! !DESCRIPTION:
18! Produces a 3D plot of the real functions contained in arrays {\tt rfmt} and
19! {\tt rfir} in the parallelepiped defined by the corner vertices in the
20! global array {\tt vclp3d}. See routine {\tt rfarray}.
21!
22! !REVISION HISTORY:
23! Created June 2003 (JKD)
24! Modified, October 2008 (F. Bultmark, F. Cricchio, L. Nordstrom)
25!EOP
26!BOC
27implicit none
28! arguments
29integer, intent(in) :: fnum,nf
30real(8), intent(in) :: rfmt(npmtmax,natmtot,nf),rfir(ngtot,nf)
31! local variables
32integer np,jf,ip
33real(8) v1(3)
34! allocatable arrays
35real(8), allocatable :: vpl(:,:),fp(:,:)
36if ((nf < 1).or.(nf > 4)) then
37 write(*,*)
38 write(*,'("Error(plot3d): invalid number of functions : ",I8)') nf
39 write(*,*)
40 stop
41end if
42! total number of plot points
43np=np3d(1)*np3d(2)*np3d(3)
44! allocate local arrays
45allocate(vpl(3,np),fp(np,nf))
46! generate the 3D plotting points
47call plotpt3d(vpl)
48! evaluate the real functions at the grid points
49do jf=1,nf
50 call rfpts(np,vpl,rfmt(:,:,jf),rfir(:,jf),fp(:,jf))
51end do
52! write functions to file
53write(fnum,'(3I6," : grid size")') np3d(:)
54do ip=1,np
55 call r3mv(avec,vpl(:,ip),v1)
56 write(fnum,'(7G18.10)') v1(:),(fp(ip,jf),jf=1,nf)
57end do
58deallocate(vpl,fp)
59end subroutine
60!EOC
61
real(8), dimension(3, 3) avec
Definition modmain.f90:12
integer, dimension(3) np3d
Definition modmain.f90:1131
subroutine plot3d(fnum, nf, rfmt, rfir)
Definition plot3d.f90:10
pure subroutine plotpt3d(vpl)
Definition plotpt3d.f90:7
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10
subroutine rfpts(np, vpl, rfmt, rfir, fp)
Definition rfpts.f90:7