|
Classes | |
class | array2 |
Functions | |
bool | isKnotVectorMonotonic (const std::vector< double > &knots) |
bool | close (double a, double b, double eps) |
int | findSpan (int degree, const std::vector< double > &knots, double u) |
double | bsplineOneBasis (int i, int deg, const std::vector< double > &U, double u) |
void | bsplineBasis (int deg, int span, const std::vector< double > &knots, double u, std::vector< double > &N) |
void | bsplineDerBasis (int deg, int span, const std::vector< double > &knots, double u, int nDers, std::vector< std::vector< double >> &ders) |
std::vector< Mdouble > | createUniformKnotVector (unsigned int numberOfControlPoints, unsigned int degree, bool clampedAtStart, bool clampedAtEnd) |
Creates a uniform (clamped) knot vector. More... | |
void | normalizeKnotVector (std::vector< Mdouble > &knots) |
Resets the knot vector to the interval [0, 1]. More... | |
void | extendKnotVector (std::vector< Mdouble > &knots, unsigned int degree, unsigned int numStart, unsigned int numEnd, bool forceBothEndsUniform=false) |
Extends the knot vector for when control points have been added at the start or end. More... | |
Vec3D | evaluate (Mdouble u, Mdouble v, std::vector< Mdouble > knotsU, std::vector< Mdouble > knotsV, std::vector< std::vector< Vec3D >> controlPoints, std::vector< std::vector< Mdouble >> weights) |
Evaluate point on a NURBS surface. More... | |
void NurbsUtils::bsplineBasis | ( | int | deg, |
int | span, | ||
const std::vector< double > & | knots, | ||
double | u, | ||
std::vector< double > & | N | ||
) |
Compute all non-zero B-spline basis functions
[in] | deg | Degree of the basis function. |
[in] | span | Index obtained from findSpan() corresponding the u and knots. |
[in] | knots | Knot vector corresponding to the basis functions. |
[in] | u | Parameter to evaluate the basis functions at. |
[in,out] | N | Values of (deg+1) non-zero basis functions. |
Referenced by NurbsSurface::evaluate(), and evaluate().
void NurbsUtils::bsplineDerBasis | ( | int | deg, |
int | span, | ||
const std::vector< double > & | knots, | ||
double | u, | ||
int | nDers, | ||
std::vector< std::vector< double >> & | ders | ||
) |
Compute all non-zero derivatives of B-spline basis functions
[in] | deg | Degree of the basis function. |
[in] | span | Index obtained from findSpan() corresponding the u and knots. |
[in] | knots | Knot vector corresponding to the basis functions. |
[in] | u | Parameter to evaluate the basis functions at. |
[in] | nDers | Number of derivatives to compute (nDers <= deg) |
[in,out] | ders | Values of non-zero derivatives of basis functions. |
References constants::i.
Referenced by NurbsSurface::evaluateDerivatives().
Compute a single B-spline basis function
[in] | i | The ith basis function to compute. |
[in] | deg | Degree of the basis function. |
[in] | U | Knot vector corresponding to the basis functions. |
[in] | u | Parameter to evaluate the basis functions at. |
References close(), and constants::i.
Referenced by bsplineOneBasis(), extendKnotVector(), findSpan(), and normalizeKnotVector().
std::vector< Mdouble > NurbsUtils::createUniformKnotVector | ( | unsigned int | numberOfControlPoints, |
unsigned int | degree, | ||
bool | clampedAtStart, | ||
bool | clampedAtEnd | ||
) |
Creates a uniform (clamped) knot vector.
numberOfControlPoints | The number of control points the knot vector should correspond to |
degree | The degree the knot vector should correspond to |
clampedAtStart | Whether or not the knot vector should be clamped at the start |
clampedAtEnd | Whether or not the knot vector should be clamped at the end |
References constants::degree, and constants::i.
Referenced by NurbsSurface::NurbsSurface().
Vec3D NurbsUtils::evaluate | ( | Mdouble | u, |
Mdouble | v, | ||
std::vector< Mdouble > | knotsU, | ||
std::vector< Mdouble > | knotsV, | ||
std::vector< std::vector< Vec3D >> | controlPoints, | ||
std::vector< std::vector< Mdouble >> | weights | ||
) |
Evaluate point on a NURBS surface.
u | Parameter to evaluate the surface at |
v | Parameter to evaluate the surface at |
knotsU | Knot vector in u-direction |
knotsV | Knot vector in v-direction |
controlPoints | 2D vector of control points |
weights | 2D vector of weights |
References bsplineBasis(), findSpan(), and Global_Physical_Variables::Nu.
Referenced by WearableNurbsWall::getVolumeUnderSurface(), WearableNurbsWall::getVolumeUnderSurfaceX(), main(), and NurbsSurface::set().
void NurbsUtils::extendKnotVector | ( | std::vector< Mdouble > & | knots, |
unsigned int | degree, | ||
unsigned int | numStart, | ||
unsigned int | numEnd, | ||
bool | forceBothEndsUniform = false |
||
) |
Extends the knot vector for when control points have been added at the start or end.
knots | The knot vector to be extended |
degree | The nurbs degree |
numStart | Number to be added at the start |
numEnd | Number to be added at the end |
forceBothEndsUniform | When only adding to one end, the other end can remain untouched or it can be forced to be made uniform |
References close(), constants::degree, constants::i, INFO, logger, and normalizeKnotVector().
Referenced by NurbsSurface::wrapAroundInU(), and NurbsSurface::wrapAroundInV().
Find the span of the given parameter in the knot vector.
[in] | degree | Degree of the curve. |
[in] | knots | Knot vector of the curve. |
[in] | u | Parameter value. |
References close(), constants::degree, and n.
Referenced by NurbsSurface::evaluate(), evaluate(), and NurbsSurface::evaluateDerivatives().
Referenced by NurbsSurface::set().
void NurbsUtils::normalizeKnotVector | ( | std::vector< Mdouble > & | knots | ) |
Resets the knot vector to the interval [0, 1].
knots | The knot vector to be normalized |
References close().
Referenced by extendKnotVector(), NurbsSurface::set(), and NurbsSurface::unclampKnots().