The Elk Code
plotpt1d.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: plotpt1d
8 ! !INTERFACE:
9 subroutine plotpt1d(cvec,nv,np,vvl,vpl,dv,dp)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! cvec : matrix of (reciprocal) lattice vectors stored column-wise
12 ! (in,real(3,3))
13 ! nv : number of vertices (in,integer)
14 ! np : number of connecting points (in,integer)
15 ! vvl : vertex vectors in lattice coordinates (in,real(3,nv))
16 ! vpl : connecting point vectors in lattice coordinates (out,real(3,np))
17 ! dv : cummulative distance to each vertex (out,real(nv))
18 ! dp : cummulative distance to each connecting point (out,real(np))
19 ! !DESCRIPTION:
20 ! Generates a set of points which interpolate between a given set of vertices.
21 ! Vertex points are supplied in lattice coordinates in the array {\tt vvl} and
22 ! converted to Cartesian coordinates with the matrix {\tt cvec}. Interpolating
23 ! points are stored in the array {\tt vpl}. The cummulative distances to the
24 ! vertices and points along the path are stored in arrays {\tt dv} and
25 ! {\tt dp}, respectively.
26 !
27 ! !REVISION HISTORY:
28 ! Created June 2003 (JKD)
29 ! Improved September 2007 (JKD)
30 ! Improved again, July 2010 (T. McQueen and JKD)
31 !EOP
32 !BOC
33 implicit none
34 ! arguments
35 real(8), intent(in) :: cvec(3,3)
36 integer, intent(in) :: nv,np
37 real(8), intent(in) :: vvl(3,nv)
38 real(8), intent(out) :: vpl(3,np),dv(nv),dp(np)
39 ! local variables
40 integer i,j,k,m,n
41 real(8) vl(3),vc(3)
42 real(8) dt,f,t1
43 ! alloctable arrays
44 real(8), allocatable :: seg(:)
45 if (nv < 1) then
46  write(*,*)
47  write(*,'("Error(plotpt1d): nv < 1 : ",I0)') nv
48  write(*,*)
49  stop
50 end if
51 if (np < nv) then
52  write(*,*)
53  write(*,'("Error(plotpt1d): np < nv :",2(X,I0))') np,nv
54  write(*,*)
55  stop
56 end if
57 ! special case of 1 vertex
58 if (nv == 1) then
59  dv(1)=0.d0
60  dp(:)=0.d0
61  do i=1,np
62  vpl(:,i)=vvl(:,1)
63  end do
64  return
65 end if
66 allocate(seg(nv))
67 ! find the length of each segment and total distance
68 dt=0.d0
69 do i=1,nv-1
70  dv(i)=dt
71  vl(:)=vvl(:,i+1)-vvl(:,i)
72  call r3mv(cvec,vl,vc)
73  seg(i)=sqrt(vc(1)**2+vc(2)**2+vc(3)**2)
74  dt=dt+seg(i)
75 end do
76 dv(nv)=dt
77 ! add small amount to total distance to avoid 0/0 condition
78 dt=dt+1.d-8
79 ! number of points to use between vertices
80 n=np-nv
81 ! construct the interpolating path
82 k=0
83 do i=1,nv-1
84  t1=dble(n)*seg(i)/dt
85  m=nint(t1)
86  if ((m > n).or.(i == (nv-1))) m=n
87  do j=1,m+1
88  k=k+1
89  f=dble(j-1)/dble(m+1)
90  dp(k)=dv(i)+f*seg(i)
91  vpl(:,k)=vvl(:,i)*(1.d0-f)+vvl(:,i+1)*f
92  end do
93  dt=dt-seg(i)
94  n=n-m
95 end do
96 dp(np)=dv(nv)
97 vpl(:,np)=vvl(:,nv)
98 deallocate(seg)
99 end subroutine
100 !EOC
101 
subroutine plotpt1d(cvec, nv, np, vvl, vpl, dv, dp)
Definition: plotpt1d.f90:10
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10