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 Copyright (C) 2020-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::CatmullRomSpline
29
30Description
31 An implementation of Catmull-Rom splines
32 (sometimes known as Overhauser splines).
33
34 In this implementation, the end tangents are created automatically
35 by reflection.
36
37 In matrix form, the \e local interpolation on the interval t=[0..1] is
38 described as follows:
39 \verbatim
40 P(t) = 1/2 * [ t^3 t^2 t 1 ] * [ -1 3 -3 1 ] * [ P-1 ]
41 [ 2 -5 4 -1 ] [ P0 ]
42 [ -1 0 1 0 ] [ P1 ]
43 [ 0 2 0 0 ] [ P2 ]
44 \endverbatim
45
46 Where P-1 and P2 represent the neighbouring points or the extrapolated
47 end points. Simple reflection is used to automatically create the end
48 points.
49
50 The spline is discretized based on the chord length of the individual
51 segments. In rare cases (sections with very high curvatures), the
52 resulting distribution may be sub-optimal.
53
54 A future implementation could also handle closed splines.
55
56See also
57 http://www.algorithmist.net/catmullrom.html provides a nice
58 introduction
59
60SourceFiles
61 CatmullRomSpline.C
62
63\*---------------------------------------------------------------------------*/
64
65#ifndef CatmullRomSpline_H
66#define CatmullRomSpline_H
67
68#include "polyLine.H"
69
70// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71
72namespace Foam
73{
74
75/*---------------------------------------------------------------------------*\
76 Class CatmullRomSpline Declaration
77\*---------------------------------------------------------------------------*/
80:
81 public polyLine
82{
83
84 //- Derivative of spline
85 scalar derivative
86 (
87 const label segment,
88 const scalar mu
89 ) const;
90
91
92public:
93
94 // Constructors
95
96 //- Construct from components
97 explicit CatmullRomSpline
98 (
99 const pointField& knots,
100 const bool notImplementedClosed = false
101 );
102
103
104 // Member Functions
105
106 //- The point position corresponding to the curve parameter
107 // 0 <= lambda <= 1
108 point position(const scalar lambda) const;
109
110 //- The point position corresponding to the local parameter
111 // 0 <= lambda <= 1 on the given segment
112 point position(const label segment, const scalar lambda) const;
113
114 //- The length of the curve
115 scalar length() const;
116};
117
118
119// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120
121} // End namespace Foam
122
123// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
124
125#endif
126
127// ************************************************************************* //
An implementation of Catmull-Rom splines (sometimes known as Overhauser splines).
scalar length() const
The length of the curve.
point position(const scalar lambda) const
The point position corresponding to the curve parameter.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
A series of straight line segments, which can also be interpreted as a series of control points for s...
Definition: polyLine.H:56
const volScalarField & mu
Namespace for OpenFOAM.
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)