Namespace Scine::Utils::BSplines

namespace BSplines
class BSpline
#include <BSpline.h>

Contains the standard implementation of a basis-spline curve. Evaluation and getter methods are defined so that knot vectors and control points of derivatives are calculated when they are needed.

class BSplineBasis
#include <BSplineBasis.h>

Evaluate the B-Spline basis functions that are defined by the knot vector

class BSplineTools
#include <BSplineTools.h>

Contains methods that can be used by several B-spline classes.

class Coefficients
#include <Coefficients.h>

Class for the B-Spline coefficients N for some given derivativeOrder and u value. Since many of the coefficients are equal to zero, just the non-zero coefficients are saved.

class ControlPolygonGenerator : public Scine::Utils::BSplines::Generator
#include <ControlPolygonGenerator.h>

Generates a B-spline curve from a set of data points that are used as control points.

class Exception : public runtime_error
#include <Exceptions.h>

Base class for exceptions related to the BSplines library.

Subclassed by Scine::Utils::BSplines::IncompatibleKnotVectorsForDimensionalMerge, Scine::Utils::BSplines::InvalidSplitDimension, Scine::Utils::BSplines::UnavailableDerivativeException

class FixedEndsPenalizedLeastSquares
#include <FixedEndsPenalizedLeastSquares.h>

Generate control points for a B-spline fitted from data points at given coordinates for a known knot vector.

See

FixedEndsPenalizedLeastSquaresGenerator, FixedEndsFixedKnotSplineGenerator

class FixedEndsPenalizedLeastSquaresGenerator : public Scine::Utils::BSplines::Generator
#include <FixedEndsPenalizedLeastSquaresGenerator.h>

Generates a B-Spline curve with a specified number of control points approximating the data points with fixed ends. „Fixed ends” means that the first and last control points are constrained to coincide with the first and last data point of the set. (see the NURBS book by Piegl 1997) Thus, all but the first and last control points (n-1 in total) are optimized. Still, the resulting B-Splines curve passes through the first and last control point since the first and last basis functions are equal to unity at the ends of the domain (N(u=0)=1, N(u=1)=1, therefore C(u=0)=P_0=R_0 and C(u=1)=P_last=R_last). Penalization can be turned on by specifying a lambda value > 0 and difference orders can be specified by kappa. (see the NURBS book by Piegl 1997)

class Generator
#include <Generator.h>

Basis class for BSpline generators. NB: We do not take a “const Eigen::MatrixXd&” since this would create a temporary when using it with Eigen::VectorXd.

Subclassed by Scine::Utils::BSplines::ControlPolygonGenerator, Scine::Utils::BSplines::FixedEndsPenalizedLeastSquaresGenerator, Scine::Utils::BSplines::InterpolationGenerator, Scine::Utils::BSplines::LooseEndsPenalizedLeastSquaresGenerator

class InterpolationGenerator : public Scine::Utils::BSplines::Generator
#include <InterpolationGenerator.h>

Generates a B-spline curve interpolating the data points. (see the NURBS book by Piegl 1997)

class KnotInserter
#include <KnotInserter.h>

Implementation of the knot insertion algorithm by Böhm. doi: 10.1016/0010-4485(80)90154-2

class LinearInterpolator
#include <LinearInterpolator.h>

This class generates a bspline interpolating between two given position collections.

class LooseEndsPenalizedLeastSquares
#include <LooseEndsPenalizedLeastSquares.h>

Generate control points for a B-spline fitted from data points at given coordinates for a known knot vector.

See

LooseEndsPenalizedLeastSquaresGenerator, LooseEndsFixedKnotSplineGenerator

class LooseEndsPenalizedLeastSquaresGenerator : public Scine::Utils::BSplines::Generator
#include <LooseEndsPenalizedLeastSquaresGenerator.h>

Generates a B-Spline curve with a specified number of control points approximating the dat a points with loose ends. „Loose ends” means that the first and last control points do not need to coincide with the first and last data point of the set. Thus, all n+1 control points are optimized. Still, the resulting B-Spline curve passes through the first and last control point since the first and last basis functions are equal to unity at the ends of the domain (N(u=0)=1, N(u=1)=1, therefore C(u=0)=P0 and C(u=1)=Plast). Penalization can be turned on by specifying a lambda value > 0 and difference orders can be specified by kappa.

class MolecularSpline
#include <MolecularSpline.h>

Class for a molecular trajectory saved as a B-spline. In principle, just a struct of a BSpline and a ElementTypeCollection.

class Splitter
#include <Splitter.h>

Splits a B-spline curve at a give parameter.

namespace ControlPointDerivatives

Calculation of the derivatives of BSpline curves C (or of dC/du, d2C/du2, …) with respect to the control points.

Not only can the derivatives of the curve itself with respect to the control points be calculated, but also of dC/du, d2C/du2, etc. The degree of the curve to derive is given by “curveDerivative” in some functions below.

The functions exist in two variants: one to calculate the derivatives of the curve with respect to one control point, or with respect to all the control points.

Note that the derivative of some control point coordinate i is non-zero only for the coordinate i of the curve.

Functions

Eigen::VectorXd oneDerivative(const BSpline &spline, double u, int controlPointIndex, int curveDerivative)
Eigen::MatrixXd allDerivatives(const BSpline &spline, double u, int curveDerivative)
Eigen::VectorXd curveDerivatives(const BSpline &spline, double u, int controlPointIndex)

Get the derivatives of the curve at a given u with respect to some control point coordinates. returns a vector of dCk / dPk for k in [0, nDimensions). NB: all other dCi / dPj are zero.

Eigen::MatrixXd curveDerivatives(const BSpline &spline, double u)

Get the derivatives of the curve at a given u with respect to all control point coordinates. returns a matrix of dCk / dP^A_k for k in [0, nDimensions), and A being a given control point. NB: all other dCi / dPj are zero.

Eigen::VectorXd firstOrderCurveDerivatives(const BSpline &spline, double u, int controlPointIndex)

Get the derivatives of the curve tangent at a given u with respect to some control point coordinates. returns a vector of dC’k / dPk for k in [0, nDimensions). NB: all other dC’i / dPj are zero.

Eigen::MatrixXd firstOrderCurveDerivatives(const BSpline &spline, double u)

Get the derivatives of the curve tangent at a given u with respect to all control point coordinates. returns a matrix of dC’k / dP^A_k for k in [0, nDimensions), and A being a given control point. NB: all other dC’i / dPj are zero.

Eigen::VectorXd secondOrderCurveDerivatives(const BSpline &spline, double u, int controlPointIndex)

Get the derivatives of the curve tangent at a given u with respect to some control point coordinates. returns a vector of dC’k / dPk for k in [0, nDimensions). NB: all other dC’i / dPj are zero.

Eigen::MatrixXd secondOrderCurveDerivatives(const BSpline &spline, double u)

Get the derivatives of the curve tangent at a given u with respect to all control point coordinates. returns a matrix of dC’k / dP^A_k for k in [0, nDimensions), and A being a given control point. NB: all other dC’i / dPj are zero.

namespace GeneratorUtils

Utility functions for BSpline generation. Includes methods and attributes that can be used by all BSplineGenerators Methods for creating parameters and knotVectors from the NURBS book by Piegl 1997: generateParametersByEquallySpacedMethod() uses eq. (9.3) generateParametersByChordLengthMethod() uses eqs. (9.4)-(9.5) generateKnotVectorByDeBoorMethod() uses the eqs. (9.68-9.69) generateKnotVectorByEquallySpacedMethod() uses eq. (9.7) generateKnotVectorByKnotAveraging() uses eq. (9.8)

See

Generator

Functions

Eigen::VectorXd generateParametersByEquallySpacedMethod(const Eigen::MatrixXd &dataPoints)
Eigen::VectorXd generateParametersByCentripetalMethod(const Eigen::MatrixXd &dataPoints)
Eigen::VectorXd generateParametersByChordLengthMethod(const Eigen::MatrixXd &dataPoints)
Eigen::VectorXd generateKnotVectorByDeBoorsMethod(int splineDegree, int numberPolynomialSegments, const Eigen::VectorXd &uBar)
Eigen::VectorXd generateKnotVectorByUniformMethod(int splineDegree, int numberPolynomialSegments)
Eigen::VectorXd generateKnotVectorByKnotAveraging(int splineDegree, int numberPolynomialSegments, const Eigen::VectorXd &uBar)