The Elk Code
 
Loading...
Searching...
No Matches
plotu3d.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 plotu3d(fnum,nf,zfmt,zfir)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: fnum,nf
11complex(8), intent(in) :: zfmt(npcmtmax,natmtot,nf,nfqrz)
12complex(8), intent(in) :: zfir(ngtot,nf,nfqrz)
13! local variables
14integer np,jf,ip
15real(8) v1(3)
16! allocatable arrays
17real(8), allocatable :: vpl(:,:),fp(:,:)
18if ((nf < 1).or.(nf > 4)) then
19 write(*,*)
20 write(*,'("Error(plotu3d): invalid number of functions : ",I8)') nf
21 write(*,*)
22 stop
23end if
24! total number of plot points
25np=np3d(1)*np3d(2)*np3d(3)
26! allocate local arrays
27allocate(vpl(3,np),fp(np,nf))
28! generate the 3D plotting points
29call plotpt3d(vpl)
30! evaluate the functions at the grid points
31call plotulr(np,vpl,nf,zfmt,zfir,fp)
32! write functions to file
33write(fnum,'(3I6," : grid size")') np3d(:)
34do ip=1,np
35 call r3mv(avec,vpl(:,ip),v1)
36 write(fnum,'(7G18.10)') v1(:),(fp(ip,jf),jf=1,nf)
37end do
38deallocate(vpl,fp)
39end subroutine
40
real(8), dimension(3, 3) avec
Definition modmain.f90:12
integer, dimension(3) np3d
Definition modmain.f90:1131
pure subroutine plotpt3d(vpl)
Definition plotpt3d.f90:7
subroutine plotu3d(fnum, nf, zfmt, zfir)
Definition plotu3d.f90:7
subroutine plotulr(np, vpl, nf, zfmt, zfir, fp)
Definition plotulr.f90:7
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10