![]() |
MFEM v4.9.0
Finite element discretization library
|
A vector of knots in one dimension, with B-spline basis functions of a prescribed order. More...
#include <nurbs.hpp>
Public Member Functions | |
| KnotVector ()=default | |
| Create an empty KnotVector. | |
| KnotVector (std::istream &input) | |
| Create a KnotVector by reading data from stream input. Two integers are read, for order and number of control points. | |
| KnotVector (int order, int NCP) | |
| Create a KnotVector with undefined knots (initialized to -1) of order order and number of control points NCP. | |
| KnotVector (int order, const Vector &intervals, const Array< int > &continuity) | |
| Create a KnotVector by passing in a degree, a Vector of interval lengths of length n, and a list of continuity of length n + 1. | |
| KnotVector (const KnotVector &kv) | |
| Copy constructor. | |
| KnotVector & | operator= (const KnotVector &kv) |
| int | GetNE () const |
| Return the number of elements, defined by distinct knots. | |
| int | GetNCP () const |
| Return the number of control points. | |
| int | GetOrder () const |
| Return the order. | |
| int | Size () const |
| Return the number of knots, including multiplicities. | |
| void | GetElements () |
| Count the number of elements. | |
| bool | isElement (int i) const |
| Return whether the knot index Order plus i is the beginning of an element. | |
| int | GetNKS () const |
| Return the number of control points minus the order. This is not the number of knot spans, but it gives the number of knots to be checked with isElement for non-empty knot spans (elements). | |
| real_t | getKnotLocation (real_t xi, int ni) const |
| Return the parameter for element reference coordinate xi in [0,1], for the element beginning at knot ni. | |
| int | findKnotSpan (real_t u) const |
| Return the index of the knot span containing parameter u. | |
| void | CalcShape (Vector &shape, int i, real_t xi) const |
| Calculate the nonvanishing shape function values in shape for the element corresponding to knot index i and element reference coordinate xi. | |
| void | CalcDShape (Vector &grad, int i, real_t xi) const |
| Calculate derivatives of the nonvanishing shape function values in grad for the element corresponding to knot index i and element reference coordinate xi. | |
| void | CalcDnShape (Vector &gradn, int n, int i, real_t xi) const |
| Calculate n-th derivatives (order n) of the nonvanishing shape function values in grad for the element corresponding to knot index i and element reference coordinate xi. | |
| void | CalcD2Shape (Vector &grad2, int i, real_t xi) const |
| Calculate second-order shape function derivatives, using CalcDnShape. | |
| void | FindMaxima (Array< int > &ks, Vector &xi, Vector &u) const |
| Gives the locations of the maxima of the KnotVector in reference space. The function gives the knot span ks, the coordinate in the knot span xi, and the coordinate of the maximum in parameter space u. | |
| void | FindInterpolant (Array< Vector * > &x, bool reuse_inverse=false) |
| Global curve interpolation through the points x (overwritten). x is an array with the length of the spatial dimension containing vectors with spatial coordinates. The control points of the interpolated curve are returned in x in the same form. | |
| void | Difference (const KnotVector &kv, Vector &diff) const |
| void | UniformRefinement (Vector &new_knots, int rf) const |
| Uniformly refine by factor rf, by inserting knots in each span. | |
| void | Refinement (Vector &new_knots, int rf) const |
| Refine with refinement factor rf. | |
| int | GetCoarseningFactor () const |
| Vector | GetFineKnots (const int cf) const |
| KnotVector * | DegreeElevate (int t) const |
| Return a new KnotVector with elevated degree by repeating the endpoints of the KnotVector. | |
| void | Flip () |
| Reverse the knots. | |
| void | Print (std::ostream &os) const |
| Print the order, number of control points, and knots. | |
| void | PrintFunctions (std::ostream &os, int samples=11) const |
| Prints the non-zero shape functions and their first and second derivatives associated with the KnotVector per element. Use GetElements() to count the elements before using this function. samples is the number of samples of the shape functions per element. | |
| ~KnotVector () | |
| Destroys KnotVector. | |
| real_t & | operator[] (int i) |
| Access function to knot i. | |
| const real_t & | operator[] (int i) const |
| Const access function to knot i. | |
| KnotVector * | FullyCoarsen () |
| Coarsen to a single element. | |
Public Attributes | |
| std::shared_ptr< SpacingFunction > | spacing |
| Function to define the distribution of knots for any number of knot spans. | |
| bool | coarse |
| Flag to indicate whether the KnotVector has been coarsened, which means it is ready for non-nested refinement. | |
| DenseMatrix | fact_AB |
| Array< int > | fact_ipiv |
| Banded matrix factorization. | |
| DenseMatrix | A_coll_inv |
| Row pivot indices. | |
Protected Attributes | |
| Vector | knot |
| Stores the values of all knots. | |
| int | Order |
| Order of the B-spline basis functions. | |
| int | NumOfControlPoints |
| Number of control points. | |
| int | NumOfElements |
| Number of elements, defined by distinct knots. | |
Static Protected Attributes | |
| static const int | MaxOrder = 10 |
A vector of knots in one dimension, with B-spline basis functions of a prescribed order.
|
default |
Create an empty KnotVector.
| mfem::KnotVector::KnotVector | ( | std::istream & | input | ) |
Create a KnotVector by reading data from stream input. Two integers are read, for order and number of control points.
| mfem::KnotVector::KnotVector | ( | int | order, |
| int | NCP ) |
Create a KnotVector with undefined knots (initialized to -1) of order order and number of control points NCP.
| mfem::KnotVector::KnotVector | ( | int | order, |
| const Vector & | intervals, | ||
| const Array< int > & | continuity ) |
Create a KnotVector by passing in a degree, a Vector of interval lengths of length n, and a list of continuity of length n + 1.
The intervals refer to spans between unique knot values (not counting zero-size intervals at repeated knots), and the continuity values should be >= -1 (discontinuous) and <= order-1 (maximally-smooth for the given polynomial degree). Periodicity is not supported.
|
inline |
|
inline |
Destroys KnotVector.
| KnotVector * mfem::KnotVector::DegreeElevate | ( | int | t | ) | const |
Return a new KnotVector with elevated degree by repeating the endpoints of the KnotVector.
| void mfem::KnotVector::Difference | ( | const KnotVector & | kv, |
| Vector & | diff ) const |
Set diff, comprised of knots in kv not contained in this KnotVector. kv must be of the same order as this KnotVector. The current implementation is not well defined, and the function may have undefined behavior, as diff may have unset entries at the end.
Global curve interpolation through the points x (overwritten). x is an array with the length of the spatial dimension containing vectors with spatial coordinates. The control points of the interpolated curve are returned in x in the same form.
The inverse of the collocation matrix, used in the interpolation, is stored for repeated calls and used if reuse_inverse is true. Reuse is valid only if this KnotVector has not changed since the initial call with reuse_inverse false.
| int mfem::KnotVector::findKnotSpan | ( | real_t | u | ) | const |
Gives the locations of the maxima of the KnotVector in reference space. The function gives the knot span ks, the coordinate in the knot span xi, and the coordinate of the maximum in parameter space u.
| KnotVector * mfem::KnotVector::FullyCoarsen | ( | ) |
| int mfem::KnotVector::GetCoarseningFactor | ( | ) | const |
| void mfem::KnotVector::GetElements | ( | ) |
| Vector mfem::KnotVector::GetFineKnots | ( | const int | cf | ) | const |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
| KnotVector & mfem::KnotVector::operator= | ( | const KnotVector & | kv | ) |
|
inline |
|
inline |
| void mfem::KnotVector::Print | ( | std::ostream & | os | ) | const |
Print the order, number of control points, and knots.
The output is formatted for writing a mesh to file. This function is called by NURBSPatch::Print.
| void mfem::KnotVector::PrintFunctions | ( | std::ostream & | os, |
| int | samples = 11 ) const |
Prints the non-zero shape functions and their first and second derivatives associated with the KnotVector per element. Use GetElements() to count the elements before using this function. samples is the number of samples of the shape functions per element.
| void mfem::KnotVector::Refinement | ( | Vector & | new_knots, |
| int | rf ) const |
|
inline |
| void mfem::KnotVector::UniformRefinement | ( | Vector & | new_knots, |
| int | rf ) const |
| DenseMatrix mfem::KnotVector::A_coll_inv |
| bool mfem::KnotVector::coarse |
Flag to indicate whether the KnotVector has been coarsened, which means it is ready for non-nested refinement.
| DenseMatrix mfem::KnotVector::fact_AB |
| Array<int> mfem::KnotVector::fact_ipiv |
|
protected |
|
protected |
|
protected |
|
protected |
| std::shared_ptr<SpacingFunction> mfem::KnotVector::spacing |