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