CatmullRomSpline.H
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-2019 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 Class
25  Foam::CatmullRomSpline
26 
27 Description
28  An implementation of Catmull-Rom splines
29  (sometimes known as Overhauser splines).
30 
31  In this implementation, the end tangents are created automatically
32  by reflection.
33 
34  In matrix form, the \e local interpolation on the interval t=[0..1] is
35  described as follows:
36  \verbatim
37  P(t) = 1/2 * [ t^3 t^2 t 1 ] * [ -1 3 -3 1 ] * [ P-1 ]
38  [ 2 -5 4 -1 ] [ P0 ]
39  [ -1 0 1 0 ] [ P1 ]
40  [ 0 2 0 0 ] [ P2 ]
41  \endverbatim
42 
43  Where P-1 and P2 represent the neighbouring points or the extrapolated
44  end points. Simple reflection is used to automatically create the end
45  points.
46 
47  The spline is discretized based on the chord length of the individual
48  segments. In rare cases (sections with very high curvatures), the
49  resulting distribution may be sub-optimal.
50 
51  A future implementation could also handle closed splines.
52 
53 See also
54  http://www.algorithmist.net/catmullrom.html provides a nice
55  introduction
56 
57 SourceFiles
58  CatmullRomSpline.C
59 
60 \*---------------------------------------------------------------------------*/
61 
62 #ifndef CatmullRomSpline_H
63 #define CatmullRomSpline_H
64 
65 #include "polyLine.H"
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 namespace Foam
70 {
71 
72 /*---------------------------------------------------------------------------*\
73  Class CatmullRomSpline Declaration
74 \*---------------------------------------------------------------------------*/
75 
76 class CatmullRomSpline
77 :
78  public polyLine
79 {
80 public:
81 
82  // Constructors
83 
84  //- Construct from components
86  (
87  const pointField& knots,
88  const bool notImplementedClosed = false
89  );
90 
91  //- Disallow default bitwise copy construction
92  CatmullRomSpline(const CatmullRomSpline&) = delete;
93 
94 
95  // Member Functions
96 
97  //- Return the point position corresponding to the curve parameter
98  // 0 <= lambda <= 1
99  point position(const scalar lambda) const;
100 
101  //- Return the point position corresponding to the local parameter
102  // 0 <= lambda <= 1 on the given segment
103  point position(const label segment, const scalar lambda) const;
104 
105  //- Return the length of the curve
106  scalar length() const;
107 
108 
109  // Member Operators
110 
111  //- Disallow default bitwise assignment
112  void operator=(const CatmullRomSpline&) = delete;
113 };
114 
115 
116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117 
118 } // End namespace Foam
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 #endif
123 
124 // ************************************************************************* //
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
void operator=(const CatmullRomSpline &)=delete
Disallow default bitwise assignment.
point position(const scalar lambda) const
Return the point position corresponding to the curve parameter.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
A series of straight line segments, which can also be interpreted as a series of control points for s...
Definition: polyLine.H:53
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
scalar length() const
Return the length of the curve.
An implementation of Catmull-Rom splines (sometimes known as Overhauser splines). ...
Namespace for OpenFOAM.