An implementation of Catmull-Rom splines (sometimes known as Overhauser splines). More...
Public Member Functions | |
CatmullRomSpline (const pointField &knots, const bool notImplementedClosed=false) | |
point | position (const scalar lambda) const |
point | position (const label segment, const scalar lambda) const |
scalar | length () const |
![]() | |
polyLine (const pointField &points, const bool notImplementedClosed=false) | |
polyLine (const point &start, const pointField &intermediate, const point &end, const bool notImplementedClosed=false) | |
const pointField & | points () const noexcept |
label | nSegments () const noexcept |
point | position (const scalar) const |
point | position (const label segment, const scalar) const |
scalar | length () const noexcept |
Additional Inherited Members | |
![]() | |
static tmp< pointField > | concat (const point &start, const pointField &intermediate, const point &end) |
![]() | |
void | calcParam () |
label | localParameter (scalar &lambda) const |
![]() | |
pointField | points_ |
scalar | lineLength_ |
scalarList | param_ |
An implementation of Catmull-Rom splines (sometimes known as Overhauser splines).
In this implementation, the end tangents are created automatically by reflection.
In matrix form, the local interpolation on the interval t=[0..1] is described as follows:
P(t) = 1/2 * [ t^3 t^2 t 1 ] * [ -1 3 -3 1 ] * [ P-1 ] [ 2 -5 4 -1 ] [ P0 ] [ -1 0 1 0 ] [ P1 ] [ 0 2 0 0 ] [ P2 ]
Where P-1 and P2 represent the neighbouring points or the extrapolated end points. Simple reflection is used to automatically create the end points.
The spline is discretized based on the chord length of the individual segments. In rare cases (sections with very high curvatures), the resulting distribution may be sub-optimal.
A future implementation could also handle closed splines.
Definition at line 74 of file CatmullRomSpline.H.
|
explicit |
Definition at line 76 of file CatmullRomSpline.C.
Foam::point position | ( | const scalar | lambda | ) | const |
Definition at line 87 of file CatmullRomSpline.C.
References lambda(), Foam::constant::physicoChemical::mu, and points.
Referenced by splineEdge::position().
Foam::point position | ( | const label | segment, |
const scalar | lambda | ||
) | const |
Definition at line 106 of file CatmullRomSpline.C.
References Foam::constant::physicoChemical::mu, p0, and points.
Foam::scalar length | ( | ) | const |
Definition at line 177 of file CatmullRomSpline.C.
References Foam::sum().
Referenced by splineEdge::length().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.