The Elk Code
plotpt2d.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2015 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
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 plotpt2d(cvec,cinv,vpnl,vpl,vppc)
7 use modmain
8 implicit none
9 ! arguments
10 real(8), intent(in) :: cvec(3,3),cinv(3,3)
11 real(8), intent(out) :: vpnl(3)
12 real(8), intent(out) :: vpl(3,np2d(1)*np2d(2))
13 real(8), intent(out) :: vppc(2,np2d(1)*np2d(2))
14 ! local variables
15 integer ip,i1,i2
16 real(8) vl1(3),vl2(3)
17 real(8) vc1(3),vc2(3),vc3(3)
18 real(8) d1,d2,d12,t1,t2
19 vl1(:)=vclp2d(:,1)-vclp2d(:,0)
20 vl2(:)=vclp2d(:,2)-vclp2d(:,0)
21 call r3mv(cvec,vl1,vc1)
22 call r3mv(cvec,vl2,vc2)
23 d1=sqrt(vc1(1)**2+vc1(2)**2+vc1(3)**2)
24 d2=sqrt(vc2(1)**2+vc2(2)**2+vc2(3)**2)
25 if ((d1 < epslat).or.(d2 < epslat)) then
26  write(*,*)
27  write(*,'("Error(plotpt2d): zero length plotting vectors")')
28  write(*,*)
29  stop
30 end if
31 d12=(vc1(1)*vc2(1)+vc1(2)*vc2(2)+vc1(3)*vc2(3))/(d1*d2)
32 ! vector normal to plane
33 call r3cross(vc1,vc2,vc3)
34 t1=sqrt(vc3(1)**2+vc3(2)**2+vc3(3)**2)
35 if (t1 < epslat) then
36  write(*,*)
37  write(*,'("Error(plotpt2d): 2D plotting plane vectors are collinear")')
38  write(*,*)
39  stop
40 end if
41 vc3(:)=vc3(:)/t1
42 call r3mv(cinv,vc3,vpnl)
43 ip=0
44 do i2=0,np2d(2)-1
45  do i1=0,np2d(1)-1
46  ip=ip+1
47  t1=dble(i1)/dble(np2d(1))
48  t2=dble(i2)/dble(np2d(2))
49 ! plot points in 3D space
50  vpl(:,ip)=t1*vl1(:)+t2*vl2(:)+vclp2d(:,0)
51 ! plot points on the plane
52  vppc(1,ip)=t1*d1+t2*d2*d12
53  vppc(2,ip)=t2*d2*sqrt(abs(1.d0-d12**2))
54  end do
55 end do
56 end subroutine
57 
subroutine plotpt2d(cvec, cinv, vpnl, vpl, vppc)
Definition: plotpt2d.f90:7
real(8), dimension(3, 0:2) vclp2d
Definition: modmain.f90:1128
real(8) epslat
Definition: modmain.f90:24
pure subroutine r3cross(x, y, z)
Definition: r3cross.f90:10
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10