File: plot3d.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 (62 lines) | stat: -rw-r--r-- 1,753 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

! Copyright (C) 2002-2008 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: plot3d
! !INTERFACE:
subroutine plot3d(fnum,nf,rfmt,rfir)
! !USES:
use modmain
! !INPUT/OUTPUT PARAMETERS:
!   fnum : plot 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 3D plot of the real functions contained in arrays {\tt rfmt} and
!   {\tt rfir} in the parallelepiped defined by the corner vertices in the
!   global array {\tt vclp3d}. See routine {\tt rfarray}.
!
! !REVISION HISTORY:
!   Created June 2003 (JKD)
!   Modified, October 2008 (F. Bultmark, F. Cricchio, L. Nordstrom)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: fnum,nf
real(8), intent(in) :: rfmt(npmtmax,natmtot,nf),rfir(ngtot,nf)
! local variables
integer np,ip,i
real(8) v1(3)
! allocatable arrays
real(8), allocatable :: vpl(:,:),fp(:,:)
if ((nf.lt.1).or.(nf.gt.4)) then
  write(*,*)
  write(*,'("Error(plot3d): invalid number of functions : ",I8)') nf
  write(*,*)
  stop
end if
! total number of plot points
np=np3d(1)*np3d(2)*np3d(3)
! allocate local arrays
allocate(vpl(3,np),fp(np,nf))
! generate the 3D plotting points
call plotpt3d(vpl)
! evaluate the functions at the grid points
do i=1,nf
  call rfplot(np,vpl,rfmt(:,:,i),rfir(:,i),fp(:,i))
end do
! write functions to file
write(fnum,'(3I6," : grid size")') np3d(:)
do ip=1,np
  call r3mv(avec,vpl(:,ip),v1)
  write(fnum,'(7G18.10)') v1(:),(fp(ip,i),i=1,nf)
end do
deallocate(vpl,fp)
return
end subroutine
!EOC