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
6
subroutine
emdplot3d
(emds)
7
use
modmain
8
use
modpw
9
use
modomp
10
implicit none
11
! arguments
12
real
(4),
intent(in)
:: emds(nhkmax,nkpt)
13
! local variables
14
integer
np,ip,nthd
15
real
(8) v1(3),t1
16
! allocatable arrays
17
real
(8),
allocatable
:: vpl(:,:)
18
! external functions
19
real
(8),
external
:: rfhkintp
20
! total number of plot points
21
np=
np3d
(1)*
np3d
(2)*
np3d
(3)
22
! generate the 3D plotting points
23
allocate
(vpl(3,np))
24
call
plotpt3d
(vpl)
25
open
(50,file=
'EMD3D.OUT'
,form=
'FORMATTED'
)
26
write
(50,
'(3I6," : grid size")'
)
np3d
(:)
27
call
holdthd
(np,nthd)
28
!$OMP PARALLEL DEFAULT(SHARED) &
29
!$OMP PRIVATE(t1,v1) &
30
!$OMP NUM_THREADS(nthd)
31
!$OMP DO ORDERED
32
do
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
38
end do
39
!$OMP END DO
40
!$OMP END PARALLEL
41
call
freethd
(nthd)
42
close
(50)
43
deallocate
(vpl)
44
end subroutine
45
emdplot3d
subroutine emdplot3d(emds)
Definition
emdplot3d.f90:7
modmain
Definition
modmain.f90:6
modmain::bvec
real(8), dimension(3, 3) bvec
Definition
modmain.f90:16
modmain::np3d
integer, dimension(3) np3d
Definition
modmain.f90:1131
modomp
Definition
modomp.f90:6
modomp::holdthd
subroutine holdthd(nloop, nthd)
Definition
modomp.f90:78
modomp::freethd
subroutine freethd(nthd)
Definition
modomp.f90:106
modpw
Definition
modpw.f90:6
plotpt3d
pure subroutine plotpt3d(vpl)
Definition
plotpt3d.f90:7
r3mv
pure subroutine r3mv(a, x, y)
Definition
r3mv.f90:10
rfhkintp
real(8) function rfhkintp(vhpl, rfhk)
Definition
rfhkintp.f90:7
emdplot3d.f90
Generated by
1.9.8