Go to the documentation of this file.
32 Foam::scalar Foam::CatmullRomSpline::derivative
64 const point derivativePoint
71 2 * (2*e0 - 5*
p0 + 4*p1 - e1)
72 +
mu * 3 * (-e0 + 3*
p0 - 3*p1 + e1)
76 return mag(derivativePoint);
101 else if (
mu > 1 - SMALL)
123 else if (
segment > nSegments())
176 (2*e0 - 5*
p0 + 4*p1 - e1)
177 +
mu*(-e0 + 3*
p0 - 3*p1 + e1)
186 const solveScalar xi[5]=
188 -0.9061798459386639927976,
189 -0.5384693101056830910363,
191 0.5384693101056830910363,
192 0.9061798459386639927976
194 const solveScalar wi[5]=
196 0.2369268850561890875143,
197 0.4786286704993664680413,
198 0.5688888888888888888889,
199 0.4786286704993664680413,
200 0.2369268850561890875143
205 for (
int i=0;i<5;i++)
207 sum+=wi[i]*derivative(
segment,(xi[i]+1.0)/2.0)/2.0;
A series of straight line segments, which can also be interpreted as a series of control points for s...
const dimensionedScalar mu
Atomic mass unit.
const pointField & points() const noexcept
Return const-access to the control-points.
point position(const scalar lambda) const
The point position corresponding to the curve parameter.
label nSegments() const noexcept
The number of line segments.
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
CatmullRomSpline(const pointField &knots, const bool notImplementedClosed=false)
Construct from components.
An ordered pair of two objects of type <T> with first() and second() elements.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const volScalarField & p0
vector point
Point is a vector.
scalar length() const
The length of the curve.