The Elk Code
 
Loading...
Searching...
No Matches
emdplot2d.f90
Go to the documentation of this file.
1
2! Copyright (C) 2014 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 emdplot2d(emds)
7use modmain
8use modpw
9use modomp
10implicit none
11! arguments
12real(4), intent(in) :: emds(nhkmax,nkpt)
13! local variables
14integer nh(3),np,ip,n,i,nthd
15real(8) vpnl(3),v1(3),t1
16! allocatable arrays
17real(8), allocatable :: vpl(:,:),vppc(:,:)
18real(8), allocatable :: x(:),wx(:),f(:)
19! external functions
20real(8), external :: rfhkintp
21! allocate local arrays
22np=np2d(1)*np2d(2)
23allocate(vpl(3,np),vppc(2,np))
24! generate the 2D plotting points
25call plotpt2d(bvec,binv,vpnl,vpl,vppc)
26! determine the number of integration points
27nh(:)=int(hkmax*sqrt(avec(1,:)**2+avec(2,:)**2+avec(3,:)**2)/pi)+1
28n=2*maxval(nh(:)*ngridk(:))
29allocate(x(n),wx(n))
30do i=1,n
31 t1=2.d0*dble(i-1)/dble(n-1)-1.d0
32 x(i)=t1*hkmax
33end do
34! determine the weights for spline integration
35call wsplint(n,x,wx)
36open(50,file='EMD2D.OUT',form='FORMATTED')
37write(50,'(2I6," : grid size")') np2d(:)
38! loop over plotting points in the 2D plane
39call holdthd(np,nthd)
40!$OMP PARALLEL DEFAULT(SHARED) &
41!$OMP PRIVATE(f,i,v1,t1) &
42!$OMP NUM_THREADS(nthd)
43allocate(f(n))
44!$OMP DO ORDERED
45do ip=1,np
46! integrate along normal to plane
47 do i=1,n
48 v1(:)=vpl(:,ip)+x(i)*vpnl(:)
49 f(i)=rfhkintp(v1,emds)
50 end do
51 t1=dot_product(wx(:),f(:))
52!$OMP ORDERED
53 write(50,'(3G18.10)') vppc(1,ip),vppc(2,ip),t1
54!$OMP END ORDERED
55end do
56!$OMP END DO
57deallocate(f)
58!$OMP END PARALLEL
59call freethd(nthd)
60close(50)
61deallocate(vpl,vppc,x,wx)
62end subroutine
63
subroutine emdplot2d(emds)
Definition emdplot2d.f90:7
integer, dimension(2) np2d
Definition modmain.f90:1127
real(8), parameter pi
Definition modmain.f90:1229
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
real(8), dimension(3, 3) binv
Definition modmain.f90:18
real(8), dimension(3, 3) avec
Definition modmain.f90:12
integer, dimension(3) ngridk
Definition modmain.f90:448
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
Definition modpw.f90:6
real(8) hkmax
Definition modpw.f90:38
subroutine plotpt2d(cvec, cinv, vpnl, vpl, vppc)
Definition plotpt2d.f90:7
real(8) function rfhkintp(vhpl, rfhk)
Definition rfhkintp.f90:7
subroutine wsplint(n, x, w)
Definition wsplint.f90:7