The Elk Code
 
Loading...
Searching...
No Matches
emdplot1d.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 emdplot1d(emds)
7use modmain
8use modpw
9use modomp
10use modtest
11implicit none
12! arguments
13real(4), intent(in) :: emds(nhkmax,nkpt)
14! local variables
15integer nh(3),ip,n,i,j,nthd
16real(8) vl1(3),vl2(3),vl3(3)
17real(8) vc1(3),vc2(3),vc3(3),t1
18! allocatable arrays
19real(8), allocatable :: x(:),wx(:),f1(:),f2(:)
20! external functions
21real(8), external :: rfhkintp
22! generate the 1D plotting points: use only the first segment
24! compute two vectors orthogonal to each other and the plotting vector; these
25! are the directions to be used for integration
26vl1(:)=vvlp1d(:,2)-vvlp1d(:,1)
27call r3mv(bvec,vl1,vc1)
28t1=sqrt(vc1(1)**2+vc1(2)**2+vc1(3)**2)
29if (t1 < epslat) then
30 write(*,*)
31 write(*,'("Error(emdplot1d): zero length plotting vector")')
32 write(*,*)
33 stop
34end if
35vc1(:)=vc1(:)/t1
36i=1
37do j=2,3
38 if (abs(vc1(j)) < abs(vc1(i))) i=j
39end do
40vc2(:)=0.d0
41vc2(i)=1.d0
42t1=dot_product(vc1,vc2)
43vc2(:)=vc2(:)-t1*vc1(:)
44t1=sqrt(vc2(1)**2+vc2(2)**2+vc2(3)**2)
45vc2(:)=vc2(:)/t1
46call r3cross(vc1,vc2,vc3)
47! integration directions in lattice coordinates
48call r3mv(binv,vc2,vl2)
49call r3mv(binv,vc3,vl3)
50! determine the number of integration points
51nh(:)=int(hkmax*sqrt(avec(1,:)**2+avec(2,:)**2+avec(3,:)**2)/pi)+1
52n=2*maxval(nh(:)*ngridk(:))
53allocate(x(n),wx(n))
54do i=1,n
55 t1=2.d0*dble(i-1)/dble(n-1)-1.d0
56 x(i)=t1*hkmax
57end do
58! determine the weights for spline integration
59call wsplint(n,x,wx)
60open(50,file='EMD1D.OUT',form='FORMATTED')
61write(*,*)
62! loop over plotting points along 1D line
63call holdthd(npp1d,nthd)
64!$OMP PARALLEL DEFAULT(SHARED) &
65!$OMP PRIVATE(f1,f2,i,j,vl1,t1) &
66!$OMP NUM_THREADS(nthd)
67allocate(f1(n),f2(n))
68!$OMP DO ORDERED
69do ip=ip01d,npp1d
70 do i=1,n
71 do j=1,n
72 vl1(:)=vplp1d(:,ip)+x(i)*vl2(:)+x(j)*vl3(:)
73 f1(j)=rfhkintp(vl1,emds)
74 end do
75 f2(i)=dot_product(wx(:),f1(:))
76 end do
77 t1=dot_product(wx(:),f2(:))
78!$OMP ORDERED
79 write(*,'("Info(emdplot1d): done ",I6," of ",I6," points")') ip,npp1d
80 write(50,'(2G18.10)') dpp1d(ip),t1
81 flush(50)
82!$OMP END ORDERED
83! write to test file if required
84 if (test.and.(ip == 1)) then
85 call writetest(171,'integrated EMD',nv=n,tol=1.d-4,rva=f2)
86 end if
87end do
88!$OMP END DO
89deallocate(f1,f2)
90!$OMP END PARALLEL
91call freethd(nthd)
92close(50)
93deallocate(x,wx)
94end subroutine
95
subroutine emdplot1d(emds)
Definition emdplot1d.f90:7
integer npp1d
Definition modmain.f90:1113
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(:,:), allocatable vvlp1d
Definition modmain.f90:1117
integer ip01d
Definition modmain.f90:1115
real(8), dimension(3, 3) avec
Definition modmain.f90:12
real(8) epslat
Definition modmain.f90:24
integer, dimension(3) ngridk
Definition modmain.f90:448
real(8), dimension(:,:), allocatable vplp1d
Definition modmain.f90:1121
real(8), dimension(:), allocatable dvp1d
Definition modmain.f90:1119
real(8), dimension(:), allocatable dpp1d
Definition modmain.f90:1123
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
logical test
Definition modtest.f90:11
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
subroutine plotpt1d(cvec, nv, np, vvl, vpl, dv, dp)
Definition plotpt1d.f90:10
pure subroutine r3cross(x, y, z)
Definition r3cross.f90:10
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10
real(8) function rfhkintp(vhpl, rfhk)
Definition rfhkintp.f90:7
subroutine wsplint(n, x, w)
Definition wsplint.f90:7