File: emdplot3d.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 (47 lines) | stat: -rw-r--r-- 1,059 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

! Copyright (C) 2015 D. Ernsting, S. Dugdale and J. K. Dewhurst.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine emdplot3d(emds)
use modmain
use modpw
use modomp
implicit none
! arguments
real(4), intent(in) :: emds(nhkmax,nkpt)
! local variables
integer np,ip,nthd
real(8) v1(3),t1
! allocatable arrays
real(8), allocatable :: vpl(:,:)
! external functions
real(8) rfhkintp
external rfhkintp
! total number of plot points
np=np3d(1)*np3d(2)*np3d(3)
! generate the 3D plotting points
allocate(vpl(3,np))
call plotpt3d(vpl)
open(50,file='EMD3D.OUT',form='FORMATTED')
write(50,'(3I6," : grid size")') np3d(:)
call omp_hold(np,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(t1,v1) &
!$OMP NUM_THREADS(nthd)
!$OMP DO ORDERED
do ip=1,np
  t1=rfhkintp(vpl(:,ip),emds)
  call r3mv(bvec,vpl(:,ip),v1)
!$OMP ORDERED
  write(50,'(4G18.10)') v1(:),t1
!$OMP END ORDERED
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
close(50)
deallocate(vpl)
return
end subroutine