The Elk Code
 
Loading...
Searching...
No Matches
emdplot3d.f90
Go to the documentation of this file.
1
2! Copyright (C) 2015 D. Ernsting, S. Dugdale and J. K. Dewhurst.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine emdplot3d(emds)
7use modmain
8use modpw
9use modomp
10implicit none
11! arguments
12real(4), intent(in) :: emds(nhkmax,nkpt)
13! local variables
14integer np,ip,nthd
15real(8) v1(3),t1
16! allocatable arrays
17real(8), allocatable :: vpl(:,:)
18! external functions
19real(8), external :: rfhkintp
20! total number of plot points
21np=np3d(1)*np3d(2)*np3d(3)
22! generate the 3D plotting points
23allocate(vpl(3,np))
24call plotpt3d(vpl)
25open(50,file='EMD3D.OUT',form='FORMATTED')
26write(50,'(3I6," : grid size")') np3d(:)
27call holdthd(np,nthd)
28!$OMP PARALLEL DEFAULT(SHARED) &
29!$OMP PRIVATE(t1,v1) &
30!$OMP NUM_THREADS(nthd)
31!$OMP DO ORDERED
32do ip=1,np
33 t1=rfhkintp(vpl(:,ip),emds)
34 call r3mv(bvec,vpl(:,ip),v1)
35!$OMP ORDERED
36 write(50,'(4G18.10)') v1(:),t1
37!$OMP END ORDERED
38end do
39!$OMP END DO
40!$OMP END PARALLEL
41call freethd(nthd)
42close(50)
43deallocate(vpl)
44end subroutine
45
subroutine emdplot3d(emds)
Definition emdplot3d.f90:7
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
integer, dimension(3) np3d
Definition modmain.f90:1131
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
Definition modpw.f90:6
pure subroutine plotpt3d(vpl)
Definition plotpt3d.f90:7
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10
real(8) function rfhkintp(vhpl, rfhk)
Definition rfhkintp.f90:7