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