Go to the documentation of this file.
40 #ifndef NURBS3DCurve_H
41 #define NURBS3DCurve_H
77 label nrmOrientation_;
83 label sgn(
const scalar val)
const;
85 scalar abs(
const scalar val)
const;
87 label mod(
const label
x,
const label interval)
const;
105 const label spacingCorrInterval,
106 const scalar tolerance
205 const label lenAcc = 25,
206 const label maxIter = 10,
207 const label spacingCorrInterval = -1,
208 const scalar tolerance = 1.
e-5
221 const vector& targetPoint,
222 const label maxIter = 1000,
223 const scalar tolerance = 1.
e-13
230 const vector& targetPoint,
231 const scalar initGuess,
232 const label maxIter = 1000,
233 const scalar tolerance = 1.
e-13
239 const vector nrm2D(
const scalar zVal,
const scalar u)
const;
245 const label nPts = 100,
246 const label lenAcc = 25,
247 const label maxIter = 10,
248 const label spacingCorrInterval=-1,
249 const scalar tolerance = 1.
e-5
358 return nrmOrientation_;
364 return givenInitNrm_;
void buildCurve()
Build curve.
scalar lengthDerivativeU(const scalar uStart, const scalar uEnd, const label nPts)
vector curveDerivativeUU(const scalar u) const
Curve second derivative wrt u at point ui.
const vector nrm2D(const scalar zVal, const scalar u) const
A class for handling words, derived from Foam::string.
A class for handling file names.
scalarList genEquidistant(const label nPts=100, const label lenAcc=25, const label maxIter=10, const label spacingCorrInterval=-1, const scalar tolerance=1.e-5)
NURBS3DCurve(const NURBSbasis &basis, const List< vector > &CPs, const List< scalar > &weights, const scalarField &u, const label nPts, const word name="NURBS3DCurve")
NURBSbasis function. Used to construct NURBS curves, surfaces and volumes.
void setCPs(const List< vector > &CPs)
Set CPs.
const vector nrm3D(const vector &refTan, const scalar u) const
Find the normal to the curve, with the option of forcing a z-plane.
vector curveDerivativeU(const scalar u) const
Curve derivative wrt u at point ui.
vector curvePoint(const scalar u) const
Curve point cartesian coordinates at ui.
void setName(const word &name)
Set name.
void makeEquidistant(const label lenAcc=25, const label maxIter=10, const label spacingCorrInterval=-1, const scalar tolerance=1.e-5)
Make curve points equidistant in cartesian space.
scalar curveDerivativeCP(const scalar u, const label CPI)
void setNrm3DOrientation(const vector &givenNrm, const vector &givenTan)
void write()
Write curve to file.
scalar findClosestCurvePoint(const vector &targetPoint, const label maxIter=1000, const scalar tolerance=1.e-13)
const NURBSbasis & getBasisFunction() const
Get basis function.
const vector & givenInitNrm() const
Return the initial normal given to compare to the Curve's normals.
bool checkRange(const scalar u, const label CPI, const label degree) const
Check if given parametric coordinate u and CP are linked.
void setNrm2DOrientation(const vector &givenNrm, const scalar zVal)
scalar length() const
Calculate length for the entire curve length.
const List< scalar > & getParametricCoordinates() const
Get parametric coordinates.
const dimensionedScalar e
Elementary charge.
~NURBS3DCurve()=default
Destructor.
const word & getName() const
Get name.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
void setWeights(const List< scalar > &weights)
Set weights.
void invert()
Invert CPs order and re-build curve.
label nrmOrientation() const
Return the nrm sgn relation to the u=0 nrm.
const List< vector > & getCPs() const
Get CPs.
void flipNrmOrientation()
Flip the orientation of the nrm.
vector curveDerivativeWeight(const scalar u, const label CPI)
Curve derivative wrt CPII at point u.
void insertKnot(const scalarField &oldKnots, const scalar uBar, const label kInsert)
Insert a knot by re-computing the control points.
const List< scalar > & getWeights() const
Get weights.