G+Smo  23.12.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsKnotVector< T > Class Template Reference

Detailed Description

template<typename T>
class gismo::gsKnotVector< T >

Class for representing a knot vector.

Consists of a vector of non-decreasing knots repeated according to their multiplicities and a vector of sums of multiplicities up to the corresponding unique knot. The repeated knots are stored in m_repKnots, whereas the multiplicity sum is stored in m_multSum.

Iterators

There are three means of iterating through the knots:

  • iterator, which is just an STL iterator over the repeated knots;
  • uiterator, which is similar but iterates over the unique knots (i.e., without multiplicities) and provides additional functionality;
  • smart_iterator, which iterates over the repeated knots and provides additional functionality.

Terminology

  • repeated knots is the non-decreasing sequence, where each knot is repeated according to its multiplicity.
  • unique knots is the increasing sequence of all the distinct values from repeated knots.
  • starting knot is the degree-th knot from the beginning, i.e., m_repKnots[m_deg].
  • ending knot is the degree-th knot from the end, i.e., m_repKnots[ m_repKnots.size() - m_deg - 1 ] or, equivalently, m_repKnots.end()[-m_deg-1].
  • domain is the interval between the starting knot and ending knot.
  • knot intervalis the interval between two consecutive unique knots. It is closed from left and open from right, except for the interval [*(domainUEnd()-1), *domainUEnd()], which is considered closed from both sides.

Compatibility notes

Earlier versions contained several functions that have been removed.

  • uiterator findElement( const T u ) const : use uFind( const T u ) instead;
  • unsigned findspan( T u ) const : use iFind(u) - begin() instead;
  • iterator findspanIter( T u ) const : use iFind(u) instead;
  • int findElementIndex(T u) const : use uFind(u).uIndex() instead;
  • unsigned Uniquefindspan (T u) const : use uFind(u).uIndex() instead;
  • gsMatrix<unsigned,1> * findspan (const gsMatrix<T,1> & u) const : call iFind(u(0,i)) - begin() for each i from 0 to u.cols().
  • int numKnotSpans() const : use uSize() - 1 instead;
  • unsigned spans() const : use uSize() - 1 instead;
+ Inheritance diagram for gsKnotVector< T >:
+ Collaboration diagram for gsKnotVector< T >:

Public Types

typedef std::vector< T >
::const_iterator 
iterator
 
typedef
internal::gsKnotIterator< T > 
smart_iterator
 
typedef
internal::gsUKnotIterator< T > 
uiterator
 

Public Member Functions

void addConstant (T amount)
 Adds amount to all the knots.
 
void addConstant (T start, T amount)
 Adds amount to all the knots, starting at knot start.
 
void affineTransformTo (T newBeg, T newEnd)
 
template<typename iterType >
void append (iterType begin, iterType end)
 Adds the knots between begin and end to the knot vector.
 
gsAsConstMatrix< T > asMatrix () const
 Returns the knots as a matrix of size 1 x size()
 
at (const size_t &i) const
 Returns the value of the i - th knot (counted with repetitions).
 
iterator begin () const
 Returns iterator pointing to the beginning of the repeated knots.
 
iterator beginAt (mult_t upos) const
 
virtual gsMatrix< T > boundingBox ()
 
virtual knotContainer breaks () const
 Returns unique knots of the domain (i.e., including the endpoints of the domain).
 
void centers_into (gsMatrix< T > &result) const
 
bool check () const
 Checks whether the knot vector is in a consistent state.
 
gsKnotVectorclone () const
 Returns a pointer to a copy.
 
std::vector< T > coarsen (index_t knotRemove=1, index_t knotSkip=1, mult_t mul=-1)
 
const T * data () const
 Returns a pointer to the beginning of the vector of knots.
 
int deduceDegree () const
 
int degree () const
 Returns the degree of the knot vector.
 
void degreeDecrease (int const &i=1)
 Inverse of degreeIncrease.
 
void degreeElevate (const short_t &i=1)
 
void degreeIncrease (int const &i=1)
 
void degreeReduce (const short_t &i)
 Converse to degreeElevate.
 
std::string detail () const
 Return a string with detailed information on the knot vector.
 
void difference (const gsKnotVector< T > &other, std::vector< T > &result) const
 Computes the difference between this knot-vector and other. More...
 
virtual short_t dim () const
 dimension of the domain
 
iterator domainBegin () const
 Returns an iterator pointing to the starting knot of the domain.
 
iterator domainEnd () const
 Returns an iterator pointing to the end-knot of the domain.
 
smart_iterator domainSBegin () const
 
smart_iterator domainSEnd () const
 
uiterator domainUBegin () const
 Returns a unique iterator pointing to the starting knot of the domain.
 
uiterator domainUEnd () const
 Returns a unique iterator pointing to the ending knot of the domain.
 
virtual gsMatrix< T > elements ()
 Returns a list of elements.
 
iterator end () const
 Returns iterator pointing past the end of the repeated knots.
 
iterator endAt (mult_t upos) const
 
void erase (const mult_t first, const mult_t last)
 Removes the knots in the range [first,last)
 
first () const
 Get the first knot.
 
unsigned firstKnotIndex (const size_t &i) const
 
const knotContainer & get () const
 Get a reference to the underlying std::vector of knots.
 
void getUniformRefinementKnots (mult_t knotsPerSpan, knotContainer &result, mult_t mult=1) const
 
gsMatrix< T > greville () const
 
greville (int i) const
 
void greville_into (gsMatrix< T > &result) const
 
 gsKnotVector ()
 Empty constructor sets the degree to -1 and leaves the knots empty.
 
 gsKnotVector (knotContainer knots, short_t degree=-1)
 
 gsKnotVector (short_t degree)
 Sets the degree and leaves the knots uninitialized.
 
 gsKnotVector (T first, T last, unsigned interior, mult_t mult_ends=1, mult_t mult_interior=1, short_t degree=-1)
 
template<typename iterType >
 gsKnotVector (short_t deg, const iterType begOfKnots, const iterType endOfKnots)
 
 gsKnotVector (const knotContainer &uKnots, int degree, int regularity)
 
bool has (T knot) const
 
iterator iFind (const T u) const
 Returns an iterator to the last occurrence of the knot at the beginning of the knot interval containing u. Note that if u == *domainEnd(), it returns the iterator domainEnd() - 1. Cf. knot interval.
 
bool includes (const gsKnotVector< T > &other) const
 Returns true if every knot of other (counted with repetitions) is found within this knot-vector. Also returns true if other is empty. More...
 
void increaseMultiplicity (const mult_t i=1, bool boundary=false)
 
bool inDomain (const T u) const
 Checks, whether the given value is inside the domain.
 
void initClamped (int degree, unsigned numKnots=2, unsigned mult_interior=1)
 
void initClamped (T u0, T u1, int degree, unsigned interior=0, unsigned mult_interior=1)
 Resets the knot vector so that it is uniform and has clamped endknots.
 
void initGraded (unsigned numKnots, int degree, T grading=0.5, unsigned mult_interior=1)
 
void initGraded (T u0, T u1, unsigned interior, int degree, T grading, unsigned mult_interior=1)
 Resets the knot vector so that its knots are graded.
 
void initUniform (T first, T last, unsigned interior, unsigned mult_ends, unsigned mult_interior, short_t degree=-1)
 
void initUniform (unsigned numKnots, unsigned mult_ends, unsigned mult_interior=1, short_t degree=-1)
 
void insert (T knot, mult_t mult=1)
 Inserts knot knot into the knot vector with multiplicity mult.
 
template<typename iterType >
void insert (iterType ibeg, iterType iend)
 
void insert (const knotContainer &knots, int mult=1)
 Insert knots into the knot vector. More...
 
bool isOpen () const
 
bool isUniform (T tol=1e-9) const
 Checks whether the knot vector is uniform.
 
gsKnotVector knotIntersection (const gsKnotVector &b) const
 Computes the intersection of knot-vectors this and b.
 
unsigned knotsUntilSpan (const size_t &i) const
 
gsKnotVector knotUnion (const gsKnotVector &b) const
 Computes the union of knot-vectors this and b. More...
 
last () const
 Get the last knot.
 
unsigned lastKnotIndex (const size_t &i) const
 
mult_t maxInteriorMultiplicity () const
 Returns the maximum multiplicity in the interior.
 
maxIntervalLength () const
 Returns the maximum interval length of the knot sequence.
 
virtual gsMatrix< T > mesh ()
 Returns the mesh..
 
mult_t minInteriorMultiplicity () const
 Returns the minimum multiplicity in the interior.
 
minIntervalLength () const
 Returns the minimum interval length of the knot sequence.
 
mult_t multFirst () const
 Returns the multiplicity of the first knot.
 
multContainer multiplicities () const
 Returns vector of multiplicities of the knots.
 
mult_t multiplicity (T u) const
 
mult_t multiplicityIndex (mult_t i) const
 
mult_t multLast () const
 Returns the multiplicity of the last knot.
 
const mult_t * multSumData () const
 
size_t numElements () const
 Number of knot intervals inside domain.
 
index_t numLeftGhosts () const
 
index_t numRightGhosts () const
 
 operator const knotContainer & () const
 Cast to the full vector of knot values (with repetitions).
 
bool operator!= (const gsKnotVector< T > &other) const
 Compares with another knot vector.
 
const T & operator() (const mult_t i) const
 Provides the knot with unique index i.
 
bool operator== (const gsKnotVector< T > &other) const
 Compare with another knot vector.
 
const T & operator[] (const mult_t i) const
 Provides the i-th knot (numbered including repetitions).
 
std::ostream & print (std::ostream &os=gsInfo) const
 
reverse_iterator rbegin () const
 Returns reverse iterator pointing past the end of the repeated knots.
 
void reduceMultiplicity (const mult_t i=1, bool boundary=false)
 
template<typename iterType >
void refineSpans (iterType begin, iterType end, mult_t knotsPerSpan)
 Inserts knotsPerSpan knots between each two knots between begin and end.
 
void refineSpans (const std::vector< unsigned > &spanIndices, mult_t knotsPerSpan=1)
 Because of type compatibility, cf. the other version.
 
void refineSpans (const multContainer &spanIndices, mult_t knotsPerSpan=1)
 
void remove (uiterator uit, mult_t mult=1)
 
void remove (const T knot, mult_t mult=1)
 Decreases multiplicity of the knot knot by mult.
 
reverse_iterator rend () const
 Returns reverse iterator pointing to the beginning of the repeated knots.
 
void reverse ()
 Better directly use affineTransformTo.
 
reverse_smart_iterator rsbegin () const
 Returns the reverse smart iterator pointing past the end of the repeated knots.
 
reverse_smart_iterator rsend () const
 Returns the reverse smart iterator pointing to the beginning of the repeated knots.
 
smart_iterator sbegin () const
 Returns the smart iterator pointing to the beginning of the repeated knots.
 
smart_iterator send () const
 Returns the smart iterator pointing past the end of the repeated knots.
 
void set_degree (short_t p)
 Sets the degree to p.
 
size_t size () const
 Number of knots (including repetitions).
 
void supportIndex_into (const mult_t &i, gsMatrix< index_t > &result) const
 
void swap (gsKnotVector &other)
 Swaps with other knot vector.
 
void symDifference (const gsKnotVector< T > &other, std::vector< T > &result) const
 Computes the symmetric difference between this knot-vector and other. More...
 
void transform (T c, T d)
 See affineTransform().
 
void trimLeft (const mult_t numKnots)
 Removes the left-most numKnots from the knot-vector.
 
void trimRight (const mult_t numKnots)
 Removes the right-most numKnots from the knot-vector.
 
unsigned u_multiplicityIndex (size_t const &i) const
 
uiterator ubegin () const
 Returns unique iterator pointing to the beginning of the unique knots.
 
uiterator uend () const
 Returns unique iterator pointing past the end of the unique knots.
 
uiterator uFind (const T u) const
 Returns the uiterator pointing to the knot at the beginning of the knot interval containing u. Note that if u == *domainUEnd(), it returns the uiterator domainUEnd() - 1. Cf. knot interval.
 
uiterator uLowerBound (const T u) const
 Returns a uiterator pointing to the first knot which does not compare less than u. More...
 
void uniformRefine (mult_t numKnots=1, mult_t mult=1)
 
knotContainer unique () const
 Returns unique knots.
 
reverse_uiterator urbegin () const
 Returns reverse unique iterator pointing past the end of the unique knots.
 
reverse_uiterator urend () const
 Returns reverse unique iterator pointing to the beginning of the unique knots.
 
size_t uSize () const
 Number of unique knots (i.e., without repetitions).
 
uiterator uUpperBound (const T u) const
 Returns an iterator pointing to the first knot which compares greater than u. More...
 
uValue (const size_t &i) const
 Provides the knot with unique index i.
 

Static Public Member Functions

static bool isConsistent (const knotContainer &repKnots, const multContainer &multSums)
 Sanity check.
 

Member Typedef Documentation

typedef std::vector<T>::const_iterator iterator

This iterator is an iterator over the repeated knots. Can be obtained by calling [r]begin() and [r]end().

typedef internal::gsKnotIterator<T> smart_iterator

This iterator is an iterator over repeated knots with additional functions. Can be obtained by calling [r]sbegin() and [r]send().

typedef internal::gsUKnotIterator<T> uiterator

This iterator is an iterator over the unique knots. Has additional functions. Can be obtained by calling [r]ubegin() and [r]uend().

Constructor & Destructor Documentation

gsKnotVector ( knotContainer  knots,
short_t  degree = -1 
)
explicit

Constructs knot vector from the given knots (repeated according to multiplicities) and deduces the degree from the multiplicity of the endknots.

gsKnotVector ( first,
last,
unsigned  interior,
mult_t  mult_ends = 1,
mult_t  mult_interior = 1,
short_t  degree = -1 
)
Parameters
firststarting parameter
lastend parameter parameter
interiornumber of interior knots
mult_endsmultiplicity at the two end knots
mult_interiormultiplicity at the interior knots
degreeDegree of the knot. (default: -1)
gsKnotVector ( short_t  deg,
const iterType  begOfKnots,
const iterType  endOfKnots 
)
inline

Constructs knot vector from the given degree and iterators marking its endpoints.

Parameters
degDegree of the knot.
begOfKnotsiterator pointing to the beginning of the knots,
endOfKnotsiterator pointing past the end of the knots. The knots between begOfKnots and endOfKnots are assumed to be repeated according to their multiplicities and sorted.
gsKnotVector ( const knotContainer &  uKnots,
int  degree,
int  regularity 
)
Parameters
uKnotsunique knots (assumed to be sorted),
degree-> endknots have multiplicity degree + 1,
regularity-> internal knots have multiplicity degree - regularity

Member Function Documentation

void affineTransformTo ( newBeg,
newEnd 
)

Scales and translates the knot vector so that its new endpoints are newBeg and newEnd.

gsKnotVector< T >::iterator beginAt ( mult_t  upos) const

Returns an iterator pointing to the first appearance of the knot with unique index (ie. counted without repetitions, left ghosts mapped to negatives) equal to upos.

virtual gsMatrix<T> boundingBox ( )
inlinevirtualinherited

Returns a bounding box for the domain eg. This coincides to the domain in case of tensor-product domains

void centers_into ( gsMatrix< T > &  result) const

Writes the center points of each element (knot-span inside the domain) to result.

std::vector< T > coarsen ( index_t  knotRemove = 1,
index_t  knotSkip = 1,
mult_t  mul = -1 
)

Coarsen the knot vector: for each group of knotRemove consecutive interior unique knots, reduce their multiplicity by mul (or remove completely, if mul = -1), then skip knotSkip consecutive unique knots and go on

Parameters
knotRemoveNumber of consecutive knots to remove
knotSkipNumber of consecutive knots which stay untouched in between
mulMultiplicity drop for each processed knot (or -1, full remove)
int deduceDegree ( ) const

Attempts to deduce the degree from the multiplicities of the endpoints (assuming clamped knots).

void degreeElevate ( const short_t i = 1)

Elevate the degree. I.e., increase the multiplicity of all the knots by one and increment the degree.

void degreeIncrease ( int const &  i = 1)
inline

Increase the degree keeping interior knots intact (add clamped knots only). Cf. degreeElevate.

void difference ( const gsKnotVector< T > &  other,
std::vector< T > &  result 
) const

Computes the difference between this knot-vector and other.

All knots in this which are not found in other (counted with repetitions), are stored in result

smart_iterator domainSBegin ( ) const
inline

Returns a smart iterator pointing to the starting knot of the domain

smart_iterator domainSEnd ( ) const
inline

Returns a smart iterator pointing to the ending knot of the domain

gsKnotVector< T >::iterator endAt ( mult_t  upos) const

Returns an iterator pointing one past the last appearance of the knot with cardinal index (ie. counted without repetitions, left ghosts mapped to negatives) equal to upos.

unsigned firstKnotIndex ( const size_t &  i) const
inline

Returns the first knot-index of cardinal index i (i.e. counted without repetitions)

void getUniformRefinementKnots ( mult_t  knotsPerSpan,
knotContainer &  result,
mult_t  mult = 1 
) const

Compute the new knots needed for uniform refinement with the given number of knots per span and return them in result. TODO: Think who is in charge of uniform refinement: basis or knot vector?

gsMatrix<T> greville ( ) const
inline

Returns the greville points of the B-splines defined on this knot vector.

T greville ( int  i) const

Returns Greville abscissa of the i - the B-spline defined on the knot vector.

void greville_into ( gsMatrix< T > &  result) const

Writes Greville abscissae of the B-splines defined on this knot vector to result.

bool has ( knot) const
inline

True iff the knot exists in the vector

Parameters
knotparameter value.
bool includes ( const gsKnotVector< T > &  other) const

Returns true if every knot of other (counted with repetitions) is found within this knot-vector. Also returns true if other is empty.

Compresses the knot-vector by making the knot-sequence strictly increasing

void increaseMultiplicity ( const mult_t  i = 1,
bool  boundary = false 
)

Increase the multiplicity of all the knots by i. If boundary is set to false, the boundary knots are not adjusted

void initClamped ( int  degree,
unsigned  numKnots = 2,
unsigned  mult_interior = 1 
)

Resets the knot vector so that it is uniform from 0 to 1 and the multiplicities of the endpoints are chosen according to the degree.

void initGraded ( unsigned  numKnots,
int  degree,
grading = 0.5,
unsigned  mult_interior = 1 
)
inline

Resets the knot vector so that its knots are graded towards 0 according to the grading parameter and the endpoints are 0 and 1.

void initUniform ( first,
last,
unsigned  interior,
unsigned  mult_ends,
unsigned  mult_interior,
short_t  degree = -1 
)

Resets the knot vector so that it has interior knots between first and last, each of them repeated mult_interior times (and the endpoints repeated mult_ends times) and the degree is equal to degree

void initUniform ( unsigned  numKnots,
unsigned  mult_ends,
unsigned  mult_interior = 1,
short_t  degree = - 1 
)
inline

Resets the knot vector so that it has numKnots knots between 0 and 1, each repeated mult_interior times (whereas 0 and 1 are repeated mult_ends times each) and the degree is equal to degree.

void insert ( iterType  ibeg,
iterType  iend 
)
inline

Inserts each knot between ibeg and iend with multiplicity 1. Knots can be repeated in the sequence but they must be sorted. The complexity is size() + (iend-ibeg) comparisons

void insert ( const knotContainer &  knots,
int  mult = 1 
)
inline

Insert knots into the knot vector.

Parameters
knotsparameter values of the new knots, stored in a std::vector
multmultiplicity of the new knots
bool isOpen ( ) const
inline

Returns true iff the knot is open (i.e., both endpoints have multiplicities equal to degree+1).

unsigned knotsUntilSpan ( const size_t &  i) const
inline

Returns the multiplicity sum at unique index i (i.e. counted without repetitions)

gsKnotVector< T > knotUnion ( const gsKnotVector< T > &  b) const

Computes the union of knot-vectors this and b.

Example: gsKnotVector<> testKv1(0, 1, 2, 2, 1); // 0 0 1/3 2/3 1 1 gsKnotVector<> testKv2(0, 1, 1, 2, 2); // 0 0 1/2 1/2 1 1 gsKnotVector<T> kv = kv1.knotUnion(kv2);

Results in: 0 0 1/3 1/2 1/2 2/3 1 1

unsigned lastKnotIndex ( const size_t &  i) const
inline

Returns the last knot-index of cardinal index i (i.e. counted without repetitions)

gsKnotVector< T >::mult_t multiplicity ( u) const

Returns the multiplicity of the knot-value u (or zero if the value is not a knot)

gsKnotVector< T >::mult_t multiplicityIndex ( mult_t  i) const

Returns the multiplicity of the knot number i (counter with repetitions).

const mult_t* multSumData ( ) const
inline

Returns a pointer to the beginning of the vector of running multiplicities.

index_t numLeftGhosts ( ) const
inline

Computes the number of left ghosts, i.e., of the knots to the left of the domain beginning.

index_t numRightGhosts ( ) const
inline

Computes the number of right ghosts, i.e., of the knots to the right of the domain end.

std::ostream & print ( std::ostream &  os = gsInfo) const
virtual

Print the knot vector to the given stream. TODO: Improve.

Implements gsDomain< T >.

void reduceMultiplicity ( const mult_t  i = 1,
bool  boundary = false 
)

Reduce the multiplicity of all knots by i. If boundary is set to false, the boundary knots are not adjusted

void refineSpans ( const multContainer &  spanIndices,
mult_t  knotsPerSpan = 1 
)

Inserts knotsPerSpan knots into the spans corresponding to indices listed in span Indices.

void remove ( uiterator  uit,
mult_t  mult = 1 
)

Decreases multiplicity of the knot pointed by uit by mult; removes it altogether, should its multiplicity be decreased to zero.

void supportIndex_into ( const mult_t &  i,
gsMatrix< index_t > &  result 
) const

Writes unique indices of the knots of the endpoints of i - th B-spline defined on this knot vector to result.

void symDifference ( const gsKnotVector< T > &  other,
std::vector< T > &  result 
) const

Computes the symmetric difference between this knot-vector and other.

All knots which do not exist in both knot-vectors, this and other (counted with repetitions), are stored in result

unsigned u_multiplicityIndex ( size_t const &  i) const
inline

Get the multiplicity of the unique knot indexed i

Parameters
iindex of the knot (without repetitions)
gsKnotVector< T >::uiterator uLowerBound ( const T  u) const

Returns a uiterator pointing to the first knot which does not compare less than u.

Unlike uUpperBound, the value pointed by the iterator returned by this function may also be equivalent to u, and not only greater.

void uniformRefine ( mult_t  numKnots = 1,
mult_t  mult = 1 
)

Inserts numKnots (and each of them mult - times) between each two knots.

gsKnotVector< T >::uiterator uUpperBound ( const T  u) const

Returns an iterator pointing to the first knot which compares greater than u.

Unlike uLowerBound, the value pointed by the iterator returned by this function cannot be equal to u, only greater.