CatmullRomSpline.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "CatmullRomSpline.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
33 (
34  const pointField& knots,
35  const bool closed
36 )
37 :
38  polyLine(knots, closed)
39 {}
40 
41 
42 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
43 
45 {
46  // endpoints
47  if (mu < SMALL)
48  {
49  return points().first();
50  }
51  else if (mu > 1 - SMALL)
52  {
53  return points().last();
54  }
55 
56  scalar lambda = mu;
57  label segment = localParameter(lambda);
58  return position(segment, lambda);
59 }
60 
61 
63 (
64  const label segment,
65  const scalar mu
66 ) const
67 {
68  // out-of-bounds
69  if (segment < 0)
70  {
71  return points().first();
72  }
73  else if (segment > nSegments())
74  {
75  return points().last();
76  }
77 
78  const point& p0 = points()[segment];
79  const point& p1 = points()[segment+1];
80 
81  // special cases - no calculation needed
82  if (mu <= 0.0)
83  {
84  return p0;
85  }
86  else if (mu >= 1.0)
87  {
88  return p1;
89  }
90 
91 
92  // determine the end points
93  point e0;
94  point e1;
95 
96  if (segment == 0)
97  {
98  // end: simple reflection
99  e0 = 2*p0 - p1;
100  }
101  else
102  {
103  e0 = points()[segment-1];
104  }
105 
106  if (segment+1 == nSegments())
107  {
108  // end: simple reflection
109  e1 = 2*p1 - p0;
110  }
111  else
112  {
113  e1 = points()[segment+2];
114  }
115 
116 
117  return 0.5 *
118  (
119  ( 2*p0 )
120  + mu *
121  (
122  ( -e0 + p1 )
123  + mu *
124  (
125  ( 2*e0 - 5*p0 + 4*p1 - e1 )
126  + mu *
127  ( -e0 + 3*p0 - 3*p1 + e1 )
128  )
129  )
130  );
131 }
132 
133 
134 Foam::scalar Foam::CatmullRomSpline::length() const
135 {
137  return 1;
138 }
139 
140 
141 // ************************************************************************* //
Foam::polyLine
A series of straight line segments, which can also be interpreted as a series of control points for s...
Definition: polyLine.H:55
Foam::polyLine::points
const pointField & points() const
Return const-access to the control-points.
Definition: polyLine.C:112
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
CatmullRomSpline.H
Foam::CatmullRomSpline::position
point position(const scalar lambda) const
The point position corresponding to the curve parameter.
Definition: CatmullRomSpline.C:44
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:445
Foam::Field< vector >
Foam::segment
Pair< point > segment
Definition: distributedTriSurfaceMesh.H:75
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
Foam::CatmullRomSpline::CatmullRomSpline
CatmullRomSpline(const pointField &knots, const bool notImplementedClosed=false)
Construct from components.
Definition: CatmullRomSpline.C:33
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::Vector< scalar >
points
const pointField & points
Definition: gmvOutputHeader.H:1
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::polyLine::localParameter
label localParameter(scalar &lambda) const
Definition: polyLine.C:124
Foam::CatmullRomSpline::length
scalar length() const
The length of the curve.
Definition: CatmullRomSpline.C:134