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