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

Detailed Description

template<short_t d, class T>
class gismo::gsHTensorBasis< d, T >

Class representing a (scalar) hierarchical tensor basis of functions \( \mathbb R^d \to \mathbb R \).

The principal idea for constructing the hierarchical basis is as follows (in simplified version):

  1. Take a sequence of simple tensor-product bases \(B^0,\ B^1,\ldots, B^L\). Each of these bases \( B^\ell \) defines a level \(\ell\) of the hierarchy. Note that we assume that \( B^{k+1} \) is always a "finer" basis than \( B^k \).
  2. From each of these basis \( B^\ell \), select a set of basis functions in a very smart way. This gives you a set of basis functions \(S^\ell \subseteq B^\ell \) of level \(\ell\).
  3. Take the union of these sets \( H = \bigcup_{\ell = 0,\ldots,L} S^\ell \). This is your hierarchical basis \( H \) (assuming that you selected the sets of functions \( S^\ell \) in a smart and appropriate way).

Remark on the numbering of the basis functions of \( H \):

The functions in \( H \) have global indices \(0, \ldots, N\). The numbering is sorted by levels in the following sense. Let \(n^\ell\) be the number of basis functions selected from level \(\ell\) (i.e., \( n^\ell = | S^\ell |\)), then the global indices \(0,\ldots,n^0-1\) correspond to functions which are taken from \( B^0 \), indices \( n^0,\ldots,n^0+n^1\) to functions from \(B^1\) and so forth.

Template parameters

Parameters
dis the domain dimension
Tis the coefficient type
+ Inheritance diagram for gsHTensorBasis< d, T >:
+ Collaboration diagram for gsHTensorBasis< d, T >:

Public Types

typedef memory::shared_ptr
< gsHTensorBasis
Ptr
 Shared pointer for gsHTensorBasis.
 
typedef memory::unique_ptr
< gsHTensorBasis
uPtr
 Unique pointer for gsHTensorBasis.
 

Public Member Functions

gsMatrix< index_tactive (const gsMatrix< T > &u) const
 Returns the indices of active (nonzero) functions at points u, as a list of indices. More...
 
virtual void active_into (const gsMatrix< T > &u, gsMatrix< index_t > &result) const
 Returns the indices of active basis functions at points u, as a list of indices, in result. A function is said to be active in a point if this point lies in the closure of the function's support. More...
 
void activeBoundaryFunctionsOfLevel (const unsigned level, const boxSide &s, std::vector< bool > &actives) const
 
void addLevel (const gsTensorBSplineBasis< d, T > &next_basis)
 Adds a level, only if manual levels are activated.
 
gsMatrix< index_tallBoundary () const
 Returns the indices of the basis functions that are nonzero at the domain boundary.
 
virtual void anchor_into (index_t i, gsMatrix< T > &result) const
 Returns the anchor point for member i of the basis.
 
gsMatrix< T > anchors () const
 Returns the anchor points that represent the members of the basis. There is exactly one anchor point per basis function. More...
 
virtual void anchors_into (gsMatrix< T > &result) const
 
const gsBasis< T > & basis (const index_t k) const
 Helper which casts and returns the k-th piece of this function set as a gsBasis.
 
gsMatrix< index_tboundary (boxSide const &s) const
 Returns the indices of the basis functions that are nonzero at the domain boundary as single-column-matrix.
 
gsBasis< T >::uPtr boundaryBasis (boxSide const &s)
 Returns the boundary basis for side s.
 
virtual gsMatrix< index_tboundaryOffset (boxSide const &s, index_t offset) const
 
uPtr clone ()
 Clone methode. Produceds a deep copy inside a uPtr.
 
gsSparseMatrix< T > collocationMatrix (gsMatrix< T > const &u) const
 Computes the collocation matrix w.r.t. points u. More...
 
virtual gsBSplineBasis< T > & component (short_t i)
 The 1-d basis for the i-th parameter component at the highest level.
 
virtual const gsBSplineBasis< T > & component (short_t i) const
 The 1-d basis for the i-th parameter component at the highest level.
 
virtual uPtr componentBasis (boxComponent b) const
 Returns the basis that corresponds to the component.
 
virtual uPtr componentBasis_withIndices (boxComponent b, gsMatrix< index_t > &indices, bool noBoundary=true) const
 Returns the basis that corresponds to the component. More...
 
virtual void compute (const gsMatrix< T > &in, gsFuncData< T > &out) const
 Computes function data. More...
 
virtual void connectivity (const gsMatrix< T > &nodes, gsMesh< T > &mesh) const
 
virtual void connectivityAtAnchors (gsMesh< T > &mesh) const
 
virtual gsBasis::uPtr create () const
 Create an empty basis of the derived type and return a pointer to it.
 
virtual short_t degree (short_t i) const
 If the basis is a tensor product of (piecewise) polynomial bases, then this function returns the polynomial degree of the i-th component.
 
virtual void degreeDecrease (short_t const &i=1, short_t const dir=-1)
 Lower the degree of the basis by the given amount, preserving knots multiplicity.
 
virtual void degreeElevate (short_t const &i=1, short_t const dir=-1)
 Elevate the degree of the basis by the given amount, preserve smoothness.
 
virtual void degreeIncrease (short_t const &i=1, short_t const dir=-1)
 Elevate the degree of the basis by the given amount, preserve knots multiplicity.
 
virtual void degreeReduce (short_t const &i=1, short_t const dir=-1)
 Reduce the degree of the basis by the given amount, preserve smoothness.
 
gsMatrix< T > deriv (const gsMatrix< T > &u) const
 Evaluate the derivatives,. More...
 
gsMatrix< T > deriv2 (const gsMatrix< T > &u) const
 Evaluates the second derivatives of active (i.e., non-zero) basis at points u. More...
 
virtual std::string detail () const
 Prints the object as a string with extended details.
 
virtual short_t dim () const
 Returns the dimension of the parameter space.
 
virtual gsDomain< T > * domain () const
 
std::vector< std::vector
< std::vector< index_t > > > 
domainBoundariesIndices (std::vector< std::vector< std::vector< std::vector< index_t > > > > &result) const
 Gives polylines on the boundaries between different levels of the mesh. More...
 
std::vector< std::vector
< std::vector< index_t > > > 
domainBoundariesParams (std::vector< std::vector< std::vector< std::vector< T > > > > &result) const
 Gives polylines on the boundaries between different levels of the mesh. More...
 
virtual short_t domainDim () const =0
 Dimension of the (source) domain. More...
 
virtual size_t elementIndex (const gsVector< T > &u) const
 Returns an index for the element which contains point u.
 
gsMatrix< T > elementInSupportOf (index_t j) const
 Returns (the coordinates of) an element in the support of basis function j.
 
virtual void elevateContinuity (int const &i=1)
 Elevates the continuity of the basis along element boundaries.
 
gsMatrix< T > eval (const gsMatrix< T > &u) const
 Evaluate the function,. More...
 
std::vector< gsMatrix< T > > evalAllDers (const gsMatrix< T > &u, int n) const
 Evaluate all derivatives upto order n,. More...
 
void flatTensorIndexesToHierachicalIndexes (gsSortedVector< int > &indexes, const int level) const
 transformes a sortedVector indexes of flat tensor index of the bspline basis of level to hierachical indexes in place. If a flat tensor index is not found, it will transform to -1. More...
 
index_t flatTensorIndexOf (const index_t i) const
 Returns the tensor index of the function indexed i (in continued indices). More...
 
index_t flatTensorIndexOf (const index_t i, const index_t level) const
 Returns the tensor index of the function indexed i (in continued indices). More...
 
int flatTensorIndexToHierachicalIndex (index_t index, const int level) const
 takes a flat tensor index of the bspline basis of level and gives back the hierachical index. If a flat tensor index is not found, it will return -1. More...
 
gsBasisFun< T > function (index_t i) const
 Returns the i-th basis function as a gsFunction. More...
 
const std::vector< tensorBasis * > & getBases () const
 Returns the tensor B-spline space of all levels. More...
 
index_t getLevelAtIndex (const point &Pt) const
 Returns the level(s) at indexes in the parameter domain. More...
 
index_t getLevelAtPoint (const gsMatrix< T > &Pt) const
 Returns the level(s) at point(s) in the parameter domain. More...
 
void getLevelUniqueSpanAtPoints (const gsMatrix< T > &Pt, gsVector< index_t > &lvl, gsMatrix< index_t > &loIdx) const
 Returns the level(s) and knot span(s) at point(s) in the parameter domain. More...
 
virtual T getMaxCellLength () const
 Get the maximum mesh size, as expected for approximation error estimates.
 
virtual T getMinCellLength () const
 Get the minimum mesh size, as expected for inverse inequalities.
 
 gsHTensorBasis ()
 Default empty constructor.
 
 gsHTensorBasis (gsTensorBSplineBasis< d, T > const &tbasis, gsMatrix< T > const &boxes)
 gsHTensorBasis More...
 
 gsHTensorBasis (gsTensorBSplineBasis< d, T > const &tbasis, gsMatrix< T > const &boxes, const std::vector< index_t > &levels)
 gsHTensorBasis More...
 
 gsHTensorBasis (const gsHTensorBasis &o)
 Copy constructor.
 
virtual void increaseMultiplicity (index_t lvl, int dir, T knotValue, int mult=1)
 Increases the multiplicity of a knot with the value knotValue in level lvl in direction dir by mult. If knotValue is not currently in the given knot vector its not added. More...
 
virtual void increaseMultiplicity (index_t lvl, int dir, const std::vector< T > &knotValue, int mult=1)
 Increases the multiplicity of several knots with the value knotValue in level lvl in direction dir by mult. If knotValue is not currently in the given knot vector its not added. More...
 
virtual memory::unique_ptr
< gsGeometry< T > > 
interpolateAtAnchors (gsMatrix< T > const &vals) const
 Applies interpolation of values pts using the anchors as parameter points. May be reimplemented in derived classes with more efficient algorithms. (by default uses interpolateData(pts,vals)
 
memory::unique_ptr< gsGeometry
< T > > 
interpolateData (gsMatrix< T > const &vals, gsMatrix< T > const &pts) const
 Applies interpolation given the parameter values pts and values vals.
 
virtual bool isRational () const
 Returns false, since all bases that inherit from gsBasis are not rational.
 
knot (int lvl, int k, int i) const
 Returns the i-th knot in direction k at level lvl.
 
index_t levelOf (index_t i) const
 Returns the level of the function indexed i (in continued indices)
 
void makeCompressed ()
 Cleans the basis, removing any inactive levels.
 
gsBasis< T >::domainIter makeDomainIterator () const
 Create a domain iterator for the computational mesh of this basis, that points to the first element of the domain.
 
gsBasis< T >::domainIter makeDomainIterator (const boxSide &s) const
 Create a boundary domain iterator for the computational mesh this basis, that points to the first element on the boundary of the domain.
 
virtual memory::unique_ptr
< gsGeometry< T > > 
makeGeometry (gsMatrix< T > coefs) const =0
 Create a gsGeometry of proper type for this basis with the given coefficient matrix.
 
virtual memory::unique_ptr
< gsBasis< T > > 
makeNonRational () const
 
bool manualLevels () const
 Returns true if levels are assigned manually.
 
void matchWith (const boundaryInterface &bi, const gsBasis< T > &other, gsMatrix< index_t > &bndThis, gsMatrix< index_t > &bndOther) const
 Computes the indices of DoFs that match on the interface bi. The interface is assumed to be a common face between this patch and other. The output is two lists of indices bndThis and bndOther, with indices that match one-to-one on the boundary bi.
 
void matchWith (const boundaryInterface &bi, const gsBasis< T > &other, gsMatrix< index_t > &bndThis, gsMatrix< index_t > &bndOther, index_t offset) const
 Computes the indices of DoFs that match on the interface bi. The interface is assumed to be a common face between this patch and other, with an offset offset. The output is two lists of indices bndThis and bndOther, with indices that match one-to-one on the boundary bi. More...
 
short_t maxDegree () const
 If the basis is of polynomial or piecewise polynomial type, then this function returns the maximum polynomial degree.
 
unsigned maxLevel () const
 Returns the level in which the indices are stored internally.
 
short_t minDegree () const
 If the basis is of polynomial or piecewise polynomial type, then this function returns the minimum polynomial degree.
 
virtual index_t nPieces () const
 Number of pieces in the domain of definition.
 
void numActive_into (const gsMatrix< T > &u, gsVector< index_t > &result) const
 The number of active basis functions at points u.
 
int numBreaks (int lvl, int k) const
 
virtual size_t numElements (boxSide const &s) const
 The number of elements on side s.
 
size_t numElements () const
 The number of elements.
 
int numKnots (int lvl, int k) const
 Returns the number of knots in direction k of level lvl.
 
index_t numLevels () const
 Returns the number of levels.
 
void only_insert_box (point const &k1, point const &k2, int lvl)
 Inserts a domain into the basis.
 
const gsBasis< T > & piece (const index_t k) const
 Returns the piece(s) of the function(s) at subdomain k.
 
virtual std::ostream & print (std::ostream &os) const =0
 Prints the object as a string.
 
void printBases (std::ostream &os=gsInfo) const
 Prints the spline-space hierarchy.
 
void printBasic (std::ostream &os=gsInfo) const
 Prints the spline-space hierarchy.
 
void printSpaces (std::ostream &os=gsInfo) const
 Prints the spline-space hierarchy.
 
void reduceContinuity (int const &i=1)
 Reduces spline continuity at interior knots by i.
 
virtual void refine (gsMatrix< T > const &boxes, int refExt)
 Refine the basis to levels and in the areas defined by boxes with an extension. More...
 
virtual void refine (gsMatrix< T > const &boxes)
 Refine the basis to levels and in the areas defined by boxes. More...
 
void refineBasisFunction (const index_t i)
 Refines the basis function with (hierarchical) index i.
 
virtual void refineElements (std::vector< index_t > const &boxes)
 Insert the given boxes into the quadtree. More...
 
void refineElements_withCoefs (gsMatrix< T > &coefs, std::vector< index_t > const &boxes)
 
void refineSide (const boxSide side, index_t lvl)
 Refines all the cells on the side side up to level lvl.
 
virtual void reverse ()
 Reverse the basis.
 
void setActiveToLvl (int level, std::vector< CMatrix > &x_matrix_lvl) const
 Creates characteristic matrices for basis where "level" is the maximum level i.e. ignoring higher level refinements.
 
void setDegree (short_t const &i)
 Set the degree of the basis (either elevate or reduce) in order to have degree equal to i wrt to each variable.
 
void setDegreePreservingMultiplicity (short_t const &i)
 Set the degree of the basis (either increase or decrecee) in order to have degree equal to i.
 
index_t size () const
 The number of basis functions in this basis.
 
virtual const gsBasissource () const
 
virtual gsBasissource ()
 
gsMatrix< T > support () const
 Returns the boundary basis for side s. More...
 
gsMatrix< T > support (const index_t &i) const
 Returns (a bounding box for) the support of the i-th basis function. More...
 
gsMatrix< T > supportInterval (index_t dir) const
 Returns an interval that contains the parameter values in direction dir. More...
 
virtual short_t targetDim () const
 Dimension of the target space. More...
 
virtual gsBasis::uPtr tensorize (const gsBasis &other) const
 Return a tensor basis of this and other.
 
tensorBasistensorLevel (index_t i) const
 Returns the tensor basis member of level i.
 
bool testPartitionOfUnity (const index_t npts=100, const T tol=1e-12) const
 Test the partition of unity. More...
 
virtual short_t totalDegree () const
 If the basis is of polynomial or piecewise polynomial type, then this function returns the total polynomial degree.
 
void transfer (const std::vector< gsSortedVector< index_t > > &old, gsSparseMatrix< T > &result)
 Returns transfer matrix between the hirarchical spline given by the characteristic matrix "old" and this.
 
const gsHDomain< d > & tree () const
 Returns a reference to m_tree.
 
gsHDomain< d > & tree ()
 Returns a reference to m_tree.
 
int treeSize () const
 The number of nodes in the tree representation.
 
virtual void uniformCoarsen (int numKnots=1)
 Coarsen the basis uniformly by removing groups of numKnots consecutive knots, each knot removed mul times. More...
 
void uniformCoarsen_withCoefs (gsMatrix< T > &coefs, int numKnots=1)
 Coarsen the basis uniformly. More...
 
virtual void uniformCoarsen_withTransfer (gsSparseMatrix< T, RowMajor > &transfer, int numKnots=1)
 Coarsen the basis uniformly and produce a sparse matrix which maps coarse coefficient vectors to refined ones. More...
 
virtual void uniformRefine (int numKnots=1, int mul=1, int dir=-1)
 Refine the basis uniformly by inserting numKnots new knots with multiplicity mul on each knot span.
 
void uniformRefine_withCoefs (gsMatrix< T > &coefs, int numKnots=1, int mul=1, int dir=-1)
 Refine the basis uniformly. More...
 
virtual void uniformRefine_withTransfer (gsSparseMatrix< T, RowMajor > &transfer, int numKnots=1, int mul=1)
 Refine the basis uniformly. More...
 
virtual void unrefineElements (std::vector< index_t > const &boxes)
 Clear the given boxes into the quadtree. More...
 
virtual const gsMatrix< T > & weights () const
 Only for compatibility reasons, with gsRationalBasis. It returns an empty matrix.
 
virtual gsMatrix< T > & weights ()
 Only for compatibility reasons, with gsRationalBasis. It returns an empty matrix.
 
virtual ~gsHTensorBasis ()
 Destructor.
 
Evaluation functions
gsMatrix< T > evalSingle (index_t i, const gsMatrix< T > &u) const
 Evaluate a single basis function i at points u.
 
gsMatrix< T > derivSingle (index_t i, const gsMatrix< T > &u) const
 Evaluate a single basis function i derivative at points u.
 
gsMatrix< T > deriv2Single (index_t i, const gsMatrix< T > &u) const
 Evaluate the second derivative of a single basis function i at points u.
 
gsVector< index_tnumActive (const gsMatrix< T > &u) const
 Number of active basis functions at an arbitrary parameter value. More...
 
virtual bool isActive (const index_t i, const gsVector< T > &u) const
 Returns true if there the point u with non-zero value or derivatives when evaluated at the basis function i.
 
virtual void activeCoefs_into (const gsVector< T > &u, const gsMatrix< T > &coefs, gsMatrix< T > &result) const
 Returns the matrix result of active coefficients at points u, each row being one coefficient. The order of the rows is the same as active_into and eval_into functions. More...
 
virtual void eval_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluates nonzero basis functions at point u into result. More...
 
virtual void evalSingle_into (index_t i, const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluate the i-th basis function at points u into result.
 
virtual void deriv_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluates the first partial derivatives of the nonzero basis function. More...
 
virtual void derivSingle_into (index_t i, const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluates the (partial) derivatives of the i-th basis function at points u into result. More...
 
virtual void deriv2_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluate the second derivatives of all active basis function at points u. More...
 
virtual void deriv2Single_into (index_t i, const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluate the (partial) derivatives of the i-th basis function at points u into result.
 
virtual void evalAllDers_into (const gsMatrix< T > &u, int n, std::vector< gsMatrix< T > > &result) const
 Evaluate the nonzero basis functions and their derivatives up to order n at points u into result. More...
 
virtual void evalAllDersSingle_into (index_t i, const gsMatrix< T > &u, int n, gsMatrix< T > &result) const
 Evaluate the basis function i and its derivatives up to order n at points u into result.
 
virtual void evalDerSingle_into (index_t i, const gsMatrix< T > &u, int n, gsMatrix< T > &result) const
 Evaluate the (partial) derivative(s) of order n the i-th basis function at points u into result.
 
virtual gsMatrix< T > laplacian (const gsMatrix< T > &u) const
 Compute the Laplacian of all nonzero basis functions at points u.
 

Static Public Attributes

static const short_t Dim
 Dimension of the parameter domain.
 

Protected Member Functions

void _diadicIndexToKnotIndex (const index_t level, gsVector< index_t, d > &diadicIndex) const
 Transfers the diadicIndex in the knot span in direction on level level to knot indices. More...
 
void _knotIndexToDiadicIndex (const index_t level, const index_t dir, index_t &knotIndex) const
 Transfers the knotIndex in the knot span in direction dir on level level to diadic indices. More...
 
void createMoreLevels (int numLevels) const
 Creates numLevels extra grids in the hierarchy.
 
void getBoxesAlongSlice (int dir, T par, std::vector< index_t > &boxes) const
 
void needLevel (int maxLevel) const
 Makes sure that there are numLevels grids computed in the hierarachy.
 
virtual void update_structure ()
 Updates the basis structure (eg. charact. matrices, etc), to be called after any modifications.
 

Protected Attributes

std::vector< tensorBasis * > m_bases
 The list of nested spaces. More...
 
hdomain_type m_tree
 The tree structure of the index space.
 
std::vector< std::vector
< std::vector< index_t > > > 
m_uIndices
 Store the indices of the element boundaries for each level (only if m_manualLevels==true)
 
std::vector< CMatrixm_xmatrix
 The characteristic matrices for each level. More...
 
std::vector< index_tm_xmatrix_offset
 Stores the offsets of active functions for all levels. More...
 

Private Member Functions

virtual gsSparseMatrix< T > coarsening (const std::vector< CMatrix > &old, const std::vector< CMatrix > &n, const gsSparseMatrix< T, RowMajor > &transfer) const =0
 returns a transfer matrix using the characteristic matrix of the old and new basis
 
std::vector< std::vector
< std::vector< index_t > > > 
domainBoundariesGeneric (std::vector< std::vector< std::vector< std::vector< index_t > > > > &indices, std::vector< std::vector< std::vector< std::vector< T > > > > &params, bool indicesFlag) const
 Implementation of the features common to domainBoundariesParams and domainBoundariesIndices. It takes both. More...
 
void functionOverlap (const point &boxLow, const point &boxUpp, const int level, point &actLow, point &actUpp)
 Returns the basis functions of level which have support on box, represented as an index box.
 
void insert_box (point const &k1, point const &k2, int lvl)
 Inserts a domain into the basis. More...
 

Geometry evaluation functions

These functions evaluate not the individual basis functions of the basis, but a geometry object which is represented by a coefficient matrix w.r.t. this basis object. For the format of the coefficient matrix, see gsGeometry.

These functions have default implementations which simply compute the basis function values and perform linear combination, but they may be overridden in derived classes if a higher-performance implementation is possible.

gsMatrix< T > evalFunc (const gsMatrix< T > &u, const gsMatrix< T > &coefs) const
 Evaluate the function described by coefs at points u. More...
 
virtual void evalFunc_into (const gsMatrix< T > &u, const gsMatrix< T > &coefs, gsMatrix< T > &result) const
 Evaluate the function described by coefs at points u, i.e., evaluates a linear combination of coefs x BasisFunctions, into result. More...
 
gsMatrix< T > derivFunc (const gsMatrix< T > &u, const gsMatrix< T > &coefs) const
 Evaluate the derivatives of the function described by coefs at points u. More...
 
virtual void derivFunc_into (const gsMatrix< T > &u, const gsMatrix< T > &coefs, gsMatrix< T > &result) const
 Evaluate the derivatives of the function described by coefs at points u. More...
 
virtual void jacobianFunc_into (const gsMatrix< T > &u, const gsMatrix< T > &coefs, gsMatrix< T > &result) const
 Evaluate the Jacobian of the function described by coefs at points u. Jacobian matrices are stacked in blocks.
 
gsMatrix< T > deriv2Func (const gsMatrix< T > &u, const gsMatrix< T > &coefs) const
 Evaluates the second derivatives of the function described by coefs at points u. More...
 
virtual void deriv2Func_into (const gsMatrix< T > &u, const gsMatrix< T > &coefs, gsMatrix< T > &result) const
 Evaluates the second derivatives of the function described by coefs at points u. More...
 
virtual void evalAllDersFunc_into (const gsMatrix< T > &u, const gsMatrix< T > &coefs, const unsigned n, std::vector< gsMatrix< T > > &result) const
 Evaluates all derivatives up to order n of the function described by coefs at points u. More...
 
static void linearCombination_into (const gsMatrix< T > &coefs, const gsMatrix< index_t > &actives, const gsMatrix< T > &values, gsMatrix< T > &result)
 Computes the linear combination coefs * values( actives ) More...
 

Constructor & Destructor Documentation

gsHTensorBasis ( gsTensorBSplineBasis< d, T > const &  tbasis,
gsMatrix< T > const &  boxes 
)
inline

gsHTensorBasis

Parameters
tbasis- tensor basis
boxes- matrix containing boxes - each 2x2 submatrix contains the lover left and upper right corner of the box
  • the level where the box should be inserted is one higher as the level where it is completely contained
gsHTensorBasis ( gsTensorBSplineBasis< d, T > const &  tbasis,
gsMatrix< T > const &  boxes,
const std::vector< index_t > &  levels 
)
inline

gsHTensorBasis

Parameters
tbasis
boxes- matrix containing boxes - eaxh 2x2 submatrix contains the lover left and upper right corner of the box
levels

Member Function Documentation

void _diadicIndexToKnotIndex ( const index_t  level,
gsVector< index_t, d > &  diadicIndex 
) const
protected

Transfers the diadicIndex in the knot span in direction on level level to knot indices.

Parameters
[in]levelThe level
diadicIndexThe diadic index
void _knotIndexToDiadicIndex ( const index_t  level,
const index_t  dir,
index_t knotIndex 
) const
protected

Transfers the knotIndex in the knot span in direction dir on level level to diadic indices.

Parameters
[in]levelThe level
[in]dirThe dir
knotIndexThe knot index
gsMatrix<index_t> active ( const gsMatrix< T > &  u) const
inlineinherited

Returns the indices of active (nonzero) functions at points u, as a list of indices.

See Also
active_into()
void active_into ( const gsMatrix< T > &  u,
gsMatrix< index_t > &  result 
) const
virtual

Returns the indices of active basis functions at points u, as a list of indices, in result. A function is said to be active in a point if this point lies in the closure of the function's support.

Parameters
[in]ugsMatrix containing evaluation points. Each column represents one evaluation point.
[out]resultFor every column i of u, a column containing the indices of the active basis functions at evaluation point u.col(i).

Reimplemented from gsBasis< T >.

Reimplemented in gsTHBSplineBasis< d, T >.

void activeBoundaryFunctionsOfLevel ( const unsigned  level,
const boxSide s,
std::vector< bool > &  actives 
) const

Fills the vector actives with booleans, that determine if a function of the given level is active. The functions on the boundary are ordered in ascending patchindex order.

Parameters
[in]level: level of the boundary functions
[in]s: boundary side
[out]actives: the result, true if its active, false if not
void activeCoefs_into ( const gsVector< T > &  u,
const gsMatrix< T > &  coefs,
gsMatrix< T > &  result 
) const
virtualinherited

Returns the matrix result of active coefficients at points u, each row being one coefficient. The order of the rows is the same as active_into and eval_into functions.

Parameters
[in]ugsVector containing an evaluation point.
[in]coefsgsMatrix is a coefficient matrix with as many rows as the size of the basis
[out]resultFor every column i of u, a column containing the indices of the active basis functions at evaluation point u.col(i).
gsMatrix<T> anchors ( ) const
inlineinherited

Returns the anchor points that represent the members of the basis. There is exactly one anchor point per basis function.

The exact definition of the anchor points depends on the particular basis. For instance, for a Bspline basis these are the Greville abscissae. In general, evaluating a function at the anchor points should provide enough information to interpolate that function using this basis.

virtual void anchors_into ( gsMatrix< T > &  result) const
inlinevirtual

Returns the anchors points that represent the members of the basis

Reimplemented from gsBasis< T >.

gsMatrix< index_t > boundaryOffset ( boxSide const &  s,
index_t  offset 
) const
virtual

Returns the indices of the basis functions that are nonzero at the domain boundary. If an offset is provided (the default is zero), it will return the indizes of the basis functions having this offset to the provided boxSide. Note that the offset cannot be bigger than the size of the basis in the direction orthogonal to boxSide.

Reimplemented from gsBasis< T >.

Reimplemented in gsTHBSplineBasis< d, T >.

gsSparseMatrix< T > collocationMatrix ( gsMatrix< T > const &  u) const
inlineinherited

Computes the collocation matrix w.r.t. points u.

The collocation matrix is a sparse matrix with u.cols rows and size() columns. The entry (i,j) is the value of basis function j at evaluation point i.

gsBasis< T >::uPtr componentBasis_withIndices ( boxComponent  b,
gsMatrix< index_t > &  indices,
bool  noBoundary = true 
) const
virtualinherited

Returns the basis that corresponds to the component.

Parameters
bThe component
indicesThe row vector where the indices are stored to
noBoundaryIf true, the transfer matrix does not include parts belonging to lower-order components (i.e., edges without corners or faces without corners and edges)
void compute ( const gsMatrix< T > &  in,
gsFuncData< T > &  out 
) const
virtualinherited

Computes function data.

This function evaluates the functions and their derivatives at the points in and writes them in the corresponding fields of out. Which field to write (and what to compute) is controlled by the out.flags (see also gsFuncData).

The input points in are expected to be compatible with the implementation/representation of the function, i.e. they should be points inside the domain of definitition of the function

Parameters
[in]in
[out]out

Reimplemented in gsGeometry< T >, and gsConstantFunction< T >.

void connectivity ( const gsMatrix< T > &  nodes,
gsMesh< T > &  mesh 
) const
virtual

Returns the connectivity structure of the basis The returned mesh has vertices the rows of matrix nodes

Reimplemented from gsBasis< T >.

void connectivityAtAnchors ( gsMesh< T > &  mesh) const
inlinevirtualinherited

Returns the connectivity structure of the basis The returned mesh has the anchor points as vertices

gsMatrix< T > deriv ( const gsMatrix< T > &  u) const
inherited

Evaluate the derivatives,.

See Also
deriv_into()
gsMatrix< T > deriv2 ( const gsMatrix< T > &  u) const
inherited

Evaluates the second derivatives of active (i.e., non-zero) basis at points u.

See documentation for deriv2_into() (the one without input parameter coefs) for details.

See Also
deriv2_into()
Parameters
[in]uEvaluation points in columns.
Returns
For every column of u, a column containing the second derivatives. See documentation for deriv2_into() (the one without input parameter coefs) for details.
void deriv2_into ( const gsMatrix< T > &  u,
gsMatrix< T > &  result 
) const
virtualinherited

Evaluate the second derivatives of all active basis function at points u.

Input parameter u is a gsMatrix of size d x N, where
d is the dimension of the parameter domain and
N is the number of evaluation points.
Each column of u corresponds to the coordinates of one evaluation point.

result is a gsMatrix of size (K * d) x N, where
K is the number of active basis functions at the evaluation point.
Each column of result corresponds to a column of u. It contains the "pure" and the mixed derivatives for each active basis function, "above" each other.

Example (bivariate): Let \(B_i(x,y)\), d = 2 be bivariate basis functions, and let the functions with indices 3,4,7, and 8 (K = 4) be active at an evaluation point u. Then, the corresponding column of result represents:
\( ( \partial_{xx}\, B_3(u), \partial_{yy}\, B_3(u), \partial_{xy}\, B_3(u), \partial_{xx}\, B_4(u), \partial_{yy}\, B_4(u), \partial_{xy}\, B_4(u), \partial_{xx}\, B_7(u), ... , \partial_{xy}\, B_8(u) )^T \)

Example (trivariate): Let \(B_i(x,y,z)\), d = 3 be trivariate basis functions, and let the functions with indices 3,4,7, and 8 be active at an evaluation point u. Then, the corresponding column of result represents:
\(( \partial_{xx}\, B_3(u), \partial_{yy}\, B_3(u), \partial_{zz}\, B_3(u), \partial_{xy}\, B_3(u), \partial_{xz}\, B_3(u), \partial_{yz}\, B_3(u), \partial_{xx}\, B_4(u), ... , \partial_{yz}\, B_8(u) )^T \)

Parameters
[in]uEvaluation points in columns (see above for format).
[in,out]resultFor every column of u, a column containing the second derivatives as described above.

See also deriv2() (the one without input parameter coefs).

Reimplemented from gsFunctionSet< T >.

Reimplemented in gsRationalBasis< SrcT >, gsTensorBSplineBasis< 1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsTHBSplineBasis< d, T >, gsLegendreBasis< T >, gsHBSplineBasis< d, T >, gsMonomialBasis< T >, and gsMvLegendreBasis< T >.

gsMatrix<T> deriv2Func ( const gsMatrix< T > &  u,
const gsMatrix< T > &  coefs 
) const
inlineinherited

Evaluates the second derivatives of the function described by coefs at points u.

See documentation for deriv2_into() (the one with input parameter coefs) for details.

Parameters
[in]uEvaluation points in columns.
[in]coefsCoefficient matrix describing the geometry in this basis.
Returns
For every column of u, a column containing the second derivatives. See documentation for deriv2_into() (the one with input parameter coefs) for details.
void deriv2Func_into ( const gsMatrix< T > &  u,
const gsMatrix< T > &  coefs,
gsMatrix< T > &  result 
) const
virtualinherited

Evaluates the second derivatives of the function described by coefs at points u.

...i.e., evaluates a linear combination of coefs * (2nd derivatives of basis functions), into result.

Evaluation points u are given as gsMatrix of size d x N, where
d is the dimension of the parameter domain and
N is the number of evaluation points.
Each column of u corresponds to the coordinates of one evaluation point.

The coefficients coefs are given as gsMatrix of size N x n, where
N is the number of points = number of basis functions and
n is the dimension of the physical domain.
Each row of coefs corresponds to the coordinates of one control point.

Let the function \( f: \mathbb R^3 \to \mathbb R^3\) be given by

\[ f = ( f_1, f_2, f_3)^T = \sum_{i=1}^N c_i B_i(x,y,z), \]

where \( B_i(x,y,z)\) are scalar basis functions and \(c_i\) are the corresponding (m-dimensional) control points. Then, for each column in u, the corresponding column in result represents

\[ ( \partial_{xx}\ f_1, \partial_{yy}\ f_1, \partial_{zz}\ f_1, \partial_{xy}\ f_1, \partial_{xz}\ f_1, \partial_{yz}\ f_1, \partial_{xx}\ f_2, \partial_{yy}\ f_2, \ldots , \partial_{xz}\ f_3, \partial_{yz}\ f_3)^T. \]

at the respective evaluation point.

Parameters
[in]uEvaluation points in columns (see above for format).
[in]coefsCoefficient matrix describing the geometry in this basis.
[out]resultFor every column of u, a column containing the second derivatives at the respective point in the format described above.

This function has a default implementation that may be overridden in derived classes for higher performance.
See also deriv2() (the one with input parameter coefs).

void deriv_into ( const gsMatrix< T > &  u,
gsMatrix< T > &  result 
) const
virtualinherited

Evaluates the first partial derivatives of the nonzero basis function.

Let
d denote the dimension of the parameter domain.
K denote the number of active (i.e., non-zero) basis functions (see active_into()).
N denote the number of evaluation points.
The N evaluation points u are given in a gsMatrix of size d x N. Each column of u represents one evaluation point.

The gsMatrix result contains the computed derivatives in the following form:
Column j of result corresponds to one evaluation point (specified by the j-th column of u). The column contains the gradients of all active functions "above" each other.
For example, for scalar basis functions \(B_i : (x,y,z)-> R\), a column represents
\((dx B_1, dy B_1, dz B_1, dx B_2, dy B_2, dz B_2, ... , dx B_n, dy B_N, dz B_N)^T\),
where the order the basis functions \(B_i\) is as returned by active() and active_into().

Parameters
[in]uEvaluation points given as gsMatrix of size d x N. See above for details.
[in,out]resultgsMatrix of size (K*d) x N. See above for details.
Todo:
Rename to _ grad_into

Reimplemented from gsFunctionSet< T >.

Reimplemented in gsRationalBasis< SrcT >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsTensorBSplineBasis< 1, T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsTHBSplineBasis< d, T >, gsHBSplineBasis< d, T >, gsLegendreBasis< T >, gsMvLegendreBasis< T >, gsMonomialBasis< T >, and gsConstantBasis< T >.

gsMatrix<T> derivFunc ( const gsMatrix< T > &  u,
const gsMatrix< T > &  coefs 
) const
inlineinherited

Evaluate the derivatives of the function described by coefs at points u.

Parameters
uevaluation points as N column vectors
coefscoefficient matrix describing the geometry in this basis, n columns
Returns
For every column of u, the result matrix will contain one Jacobian matrix of size d * n, such that the total size of the result is n x (d * n) x N
void derivFunc_into ( const gsMatrix< T > &  u,
const gsMatrix< T > &  coefs,
gsMatrix< T > &  result 
) const
virtualinherited

Evaluate the derivatives of the function described by coefs at points u.

Evaluates a linear combination of coefs*BasisFunctionDerivatives, into result.

This function has a default implementation that may be overridden in derived classes for higher performance.

Let the function \(f: \mathbb R^d \to \mathbb R^m \) be described by the coefficients coefs, i.e.,
each evaluation point is in \(\mathbb R^d\), and
each coefficient is a point in \(\mathbb R^m\).

The N evaluation points u are given in a gsMatrix of size d x N. Each column of u represents one evaluation point.

The K coefficients coefs are given as a gsMatrix of size K x m. Each row of coefs represents one coefficient in \(\mathbb R^m\).

The gsMatrix result contains the following data:
For every column of u, the corresponding column in the matrix result contains the gradients of the m components of the function above each other. Hence, the size of result is (d*m) x N.

Example 1:
Let \(f(s,t)\) be a bivariate scalar function, \(f:\mathbb R^2 \to \mathbb R\) (i.e., d=2, m=1), and let the evaluation point \( u_i\) be represented by the i-th column of u.
Then, result has the form

\[ \left( \begin{array}{cccc} \partial_s f(u_1) & \partial_s f(u_2) & \ldots & \partial_t f(u_{N}) \\ \partial_t f(u_1) & \partial_t f(u_2) & \ldots & \partial_t f(u_{N}) \end{array} \right) \]

Example 2:
Let \(f(s,t) = ( f_1(s,t), f_2(s,t), f_3(s,t) )\) represent a surface in space, \(f:\mathbb R^2 \to \mathbb R^3\) (i.e., d=2, m=3), and let the evaluation point \( u_i\) be represented by the i-th column of u.
Then, result has the form

\[ \left( \begin{array}{ccccccc} \partial_s f_1(u_1) & \partial_s f_1(u_2) & \ldots & \partial_s f_1(u_N) \\ \partial_t f_1(u_1) & \partial_t f_1(u_2) & \ldots & \partial_t f_1(u_N) \\ \partial_s f_2(u_1) & \partial_s f_2(u_2) & \ldots & \partial_s f_2(u_N) \\ \partial_t f_2(u_1) & \partial_t f_2(u_2) & \ldots & \partial_t f_2(u_N) \\ \partial_s f_3(u_1) & \partial_s f_3(u_2) & \ldots & \partial_s f_3(u_N) \\ \partial_t f_3(u_1) & \partial_t f_3(u_2) & \ldots & \partial_t f_3(u_N) \\ \end{array} \right) \]

Parameters
[in]uEvaluation points as d x N-matrix.
[in]coefsCoefficient matrix describing the geometry in this basis as K x m-matrix.
K should equal the size() of the basis, i.e., the number basis functions.
[in,out]resultgsMatrix of size d*m x N, see above for format.

where
d is the dimension of the parameter domain
m is the dimension of the physical domain
N is the number of evaluation points
K is the number of coefficients

void derivSingle_into ( index_t  i,
const gsMatrix< T > &  u,
gsMatrix< T > &  result 
) const
virtualinherited

Evaluates the (partial) derivatives of the i-th basis function at points u into result.

See deriv_into() for detailed documentation.

Todo:
rename grad_into

Reimplemented in gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsTensorBSplineBasis< 1, T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsTHBSplineBasis< d, T >, gsHBSplineBasis< d, T >, and gsMonomialBasis< T >.

gsDomain< T > * domain ( ) const
virtualinherited
std::vector< std::vector< std::vector< index_t > > > domainBoundariesGeneric ( std::vector< std::vector< std::vector< std::vector< index_t > > > > &  indices,
std::vector< std::vector< std::vector< std::vector< T > > > > &  params,
bool  indicesFlag 
) const
private

Implementation of the features common to domainBoundariesParams and domainBoundariesIndices. It takes both.

Parameters
indicesand
paramsbut fills in only one depending on
indicesFlag(if true, then it returns indices).
std::vector< std::vector< std::vector< index_t > > > domainBoundariesIndices ( std::vector< std::vector< std::vector< std::vector< index_t > > > > &  result) const

Gives polylines on the boundaries between different levels of the mesh.

Parameters
resultvariable where to write the polylines in the form < levels < polylines_in_one_level < one_polyline < one_segment (x1, y1, x2, y2) > > > > where <x1, y1, x2, y2 > are so that (x1, y1) <=LEX (x2, y2) and where x1, y1, x2 and y2 are indices of the knots with respect to m_maxInsLevel.
Returns
bounding boxes of the polylines in the form < levels < polylines_in_one_level < x_ll, y_ll, x_ur, y_ur > > >, where "ur" stands for "upper right" and "ll" for "lower left".
std::vector< std::vector< std::vector< index_t > > > domainBoundariesParams ( std::vector< std::vector< std::vector< std::vector< T > > > > &  result) const

Gives polylines on the boundaries between different levels of the mesh.

Parameters
resultvariable where to write the polylines in the form < levels < polylines_in_one_level < one_polyline < one_segment (x1, y1, x2, y2) > > > > , where <x1, y1, x2, y2 > are so that (x1, y1) <=LEX (x2, y2) and where x1, y1, x2 and y2 are parameters (knots).
Returns
bounding boxes of the polylines in the form < levels < polylines_in_one_level < x_ll, y_ll, x_ur, y_ur > > >, where "ur" stands for "upper right" and "ll" for "lower left".
gsMatrix< T > eval ( const gsMatrix< T > &  u) const
inherited

Evaluate the function,.

See Also
eval_into()
void eval_into ( const gsMatrix< T > &  u,
gsMatrix< T > &  result 
) const
virtualinherited

Evaluates nonzero basis functions at point u into result.

Let...
d denote the dimension of the parameter domain.
K denote the number of active (i.e., non-zero) basis functions (see active_into()). N denote the number of evaluation points.
The n evaluation points u are given in a gsMatrix of size d x N. Each column of u represents one evaluation point.

The gsMatrix result contains the computed function values in the following form:
Column j of result corresponds to one evaluation point (specified by the j-th column of u). The column contains the values of all active functions "above" each other.
For example, for scalar basis functions Bi : (x,y,z)-> R, a column represents
(B1, B2, ... , BN)^T,
where the order the basis functions Bi is as returned by active() and active_into().

Parameters
[in]uEvaluation points given as gsMatrix of size d x N. See above for details.
[in,out]resultgsMatrix of size K x N. See above for details.

Reimplemented from gsFunctionSet< T >.

Reimplemented in gsTHBSplineBasis< d, T >, gsRationalBasis< SrcT >, gsTensorBSplineBasis< 1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsHBSplineBasis< d, T >, gsLegendreBasis< T >, gsMvLegendreBasis< T >, gsConstantBasis< T >, and gsMonomialBasis< T >.

std::vector< gsMatrix< T > > evalAllDers ( const gsMatrix< T > &  u,
int  n 
) const
inherited

Evaluate all derivatives upto order n,.

See Also
evalAllDers_into
void evalAllDers_into ( const gsMatrix< T > &  u,
int  n,
std::vector< gsMatrix< T > > &  result 
) const
virtualinherited

Evaluate the nonzero basis functions and their derivatives up to order n at points u into result.

The derivatives (the 0-th derivative is the function value) are stored in a result. result is a std::vector, where result[i] is a gsMatrix which contains the i-th derivatives.

The entries in result[0], result[1], and result[2] are ordered as in eval_into(), deriv_into(), and deriv2_into(), respectively. For i > 2, the derivatives are stored in lexicographical order, e.g. for order i = 3 and dimension 2 the derivatives are stored as follows: \( \partial_{xxx}, \, \partial_{xxy}, \, \partial_{xyy}, \, \partial_{yyy}.\, \)

Parameters
[in]uEvaluation points, each column corresponds to one evaluation point.
[in]nAll derivatives up to order n are computed and stored in result.
[in,out]resultSee above for format.

Reimplemented from gsFunctionSet< T >.

Reimplemented in gsTensorBSplineBasis< 1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, and gsMappedSingleBasis< d, T >.

void evalAllDersFunc_into ( const gsMatrix< T > &  u,
const gsMatrix< T > &  coefs,
const unsigned  n,
std::vector< gsMatrix< T > > &  result 
) const
virtualinherited

Evaluates all derivatives up to order n of the function described by coefs at points u.

Evaluation points u are given as gsMatrix of size d x N, where
d is the dimension of the parameter domain and
N is the number of evaluation points.
Each column of u corresponds to the coordinates of one evaluation point.

The coefficients coefs are given as gsMatrix of size K x n, where
K is the number of (active) basis functions (=size()) and
n is the dimension of the physical domain.
Each row of coefs corresponds to the coordinates of one control point.

result is a std::vector, where the entry result[i] contains the gsMatrix corresponding to the i-th derivatives. The format of the respective entry is as in
evalFunc_into()
derivFunc_into()
deriv2Func_into()

Todo:
finish documentation
Parameters
[in]u
[in]coefs
[in]n
[out]result
gsMatrix<T> evalFunc ( const gsMatrix< T > &  u,
const gsMatrix< T > &  coefs 
) const
inlineinherited

Evaluate the function described by coefs at points u.

This function has a default implementation that may be overridden in derived classes for higher performance.

Parameters
uevaluation points as m column vectors
coefscoefficient matrix describing the geometry in this basis, n columns
Returns
a matrix of size n x m with one function value as a column vector per evaluation point
void evalFunc_into ( const gsMatrix< T > &  u,
const gsMatrix< T > &  coefs,
gsMatrix< T > &  result 
) const
virtualinherited

Evaluate the function described by coefs at points u, i.e., evaluates a linear combination of coefs x BasisFunctions, into result.

This function has a default implementation that may be overridden in derived classes for higher performance.

Parameters
uevaluation points as N column vectors
coefscoefficient matrix describing the geometry in this basis, n columns
[out]resulta matrix of size n x N with one function value as a column vector per evaluation point

Reimplemented in gsRationalBasis< SrcT >, gsTensorBSplineBasis< 1, T >, and gsMonomialBasis< T >.

void flatTensorIndexesToHierachicalIndexes ( gsSortedVector< int > &  indexes,
const int  level 
) const

transformes a sortedVector indexes of flat tensor index of the bspline basis of level to hierachical indexes in place. If a flat tensor index is not found, it will transform to -1.

Parameters
[in]indexesflat tensor indexes of the function in level
levelLevel of the basis.
index_t flatTensorIndexOf ( const index_t  i) const
inline

Returns the tensor index of the function indexed i (in continued indices).

Parameters
[in]iGlobal (continued) index of a basis function of the hierarchical basis.
Returns
The tensor index of this basis function with respect to the tensor-product basis of the corresponding level.
index_t flatTensorIndexOf ( const index_t  i,
const index_t  level 
) const
inline

Returns the tensor index of the function indexed i (in continued indices).

Parameters
[in]iGlobal (continued) index of a basis function of the hierarchical basis.
levelLevel of the i-th basis function.
Returns
The tensor index of this basis function with respect to the tensor-product basis of level.
int flatTensorIndexToHierachicalIndex ( index_t  index,
const int  level 
) const

takes a flat tensor index of the bspline basis of level and gives back the hierachical index. If a flat tensor index is not found, it will return -1.

Parameters
[in]indexflat tensor index of the function in level
levelLevel of the basis.
Returns
hierachical index, or -1 if it was not found
gsBasisFun< T > function ( index_t  i) const
inherited

Returns the i-th basis function as a gsFunction.

Note that the gsBasisFun object only holds a reference to the current basis, so it is invalidated when the basis is destroyed.

const std::vector<tensorBasis*>& getBases ( ) const
inline

Returns the tensor B-spline space of all levels.

Returns
void getBoxesAlongSlice ( int  dir,
par,
std::vector< index_t > &  boxes 
) const
protected

gets all the boxes along a slice in direction dir at parameter par. the boxes are given back in a std::vector<index_t> and are in the right format to be given to refineElements().

index_t getLevelAtIndex ( const point Pt) const
inline

Returns the level(s) at indexes in the parameter domain.

Parameters
[in]PtgsMatrix of size d x n, where
d is the dimension of the parameter domain and
n is the number of evaluation points.
Each column of Pts represents one evaluation point.
Returns
levels gsMatrix of size 1 x n.
levels(0,i) is the level of the point defined by the i-th column in Pts.
index_t getLevelAtPoint ( const gsMatrix< T > &  Pt) const
inline

Returns the level(s) at point(s) in the parameter domain.

Parameters
[in]PtgsMatrix of size d x n, where
d is the dimension of the parameter domain and
n is the number of evaluation points.
Each column of Pts represents one evaluation point.
Returns
levels gsMatrix of size 1 x n.
levels(0,i) is the level of the point defined by the i-th column in Pts.
void getLevelUniqueSpanAtPoints ( const gsMatrix< T > &  Pt,
gsVector< index_t > &  lvl,
gsMatrix< index_t > &  loIdx 
) const
inline

Returns the level(s) and knot span(s) at point(s) in the parameter domain.

Parameters
[in]PtgsMatrix of size d x n, where
d is the dimension of the parameter domain and
n is the number of evaluation points.
Each column of Pts represents one evaluation point.
[out]lvlgsVector of length n with the levels of the respective points.
[out]loIdxgsMatrix of size d x n.
Each column contains the lower corner of the knot span containing i-th point. The corner is given in unique knot span indices of level lvl[i].
void increaseMultiplicity ( index_t  lvl,
int  dir,
knotValue,
int  mult = 1 
)
virtual

Increases the multiplicity of a knot with the value knotValue in level lvl in direction dir by mult. If knotValue is not currently in the given knot vector its not added.

Parameters
[in]lvl: level
[in]dir: direction
[in]knotValue: value of the knot
[in]mult: multiplicity
void increaseMultiplicity ( index_t  lvl,
int  dir,
const std::vector< T > &  knotValue,
int  mult = 1 
)
virtual

Increases the multiplicity of several knots with the value knotValue in level lvl in direction dir by mult. If knotValue is not currently in the given knot vector its not added.

Parameters
[in]lvl: level
[in]dir: direction
[in]knotValue: value of the knot
[in]mult: multiplicity
void insert_box ( point const &  k1,
point const &  k2,
int  lvl 
)
inlineprivate

Inserts a domain into the basis.

private functions

void linearCombination_into ( const gsMatrix< T > &  coefs,
const gsMatrix< index_t > &  actives,
const gsMatrix< T > &  values,
gsMatrix< T > &  result 
)
staticinherited

Computes the linear combination coefs * values( actives )

Todo:
documentation
Parameters
[in]coefsgsMatrix of size K x m, where K should equal size() of the basis (i.e., the number of basis functions).
[in]activesgsMatrix of size numAct x numPts
[in]valuesgsMatrix of size stride*numAct x numPts
[out]resultgsMatrix of size stride x numPts
virtual memory::unique_ptr<gsBasis<T> > makeNonRational ( ) const
inlinevirtualinherited

Clone the source of this basis in case of rational basis, same as clone() otherwise

Reimplemented in gsRationalBasis< SrcT >, gsRationalBasis< gsBSplineTraits< d, T >::Basis >, and gsRationalBasis< gsBSplineBasis< T > >.

void matchWith ( const boundaryInterface bi,
const gsBasis< T > &  other,
gsMatrix< index_t > &  bndThis,
gsMatrix< index_t > &  bndOther,
index_t  offset 
) const
virtual

Computes the indices of DoFs that match on the interface bi. The interface is assumed to be a common face between this patch and other, with an offset offset. The output is two lists of indices bndThis and bndOther, with indices that match one-to-one on the boundary bi.

NOTE: bndThis will have offset but bndOther will NOT have an offset (hence offset 0)

Reimplemented from gsBasis< T >.

gsVector<index_t> numActive ( const gsMatrix< T > &  u) const
inlineinherited

Number of active basis functions at an arbitrary parameter value.

Usually, this is used for getting the active functions on one element, assuming that this number doesn't change for different parameters inside the element.

int numBreaks ( int  lvl,
int  k 
) const
inline

Returns the number of breaks (distinct knot values) in direction k of level lvl

void refine ( gsMatrix< T > const &  boxes,
int  refExt 
)
virtual

Refine the basis to levels and in the areas defined by boxes with an extension.

Parameters
[in]boxesgsMatrix of size d x n, where
n is the number of refinement boxes.
Every two consecutive columns specify the lower and upper corner of one refinement box (See also documentation of refine() for the format of box)
[in]refExtis an integer specifying how many cells should also be refined around the respective boxes.

Reimplemented from gsBasis< T >.

void refine ( gsMatrix< T > const &  boxes)
virtual

Refine the basis to levels and in the areas defined by boxes.

Parameters
[in]boxesgsMatrix of size d x n, where
n is the number of refinement boxes.
Every two consecutive columns specify the lower and upper corner of one refinement box (See also documentation of refine() for the format of box)
void refineElements ( std::vector< index_t > const &  boxes)
virtual

Insert the given boxes into the quadtree.

Each box is defined by 2d+1 indices, where d is the dimension of the parameter domain. The first index defines the level in which the box should be inserted, the next d indices the "coordinates" of the lower corner in the index space, and the last d indices the "coordinates" of the upper corner.

Example: Let d=3 and

\[ \mathsf{boxes} = [ L^1, \ell_x^1, \ell_y^1, \ell_z^1, u_x^1, u_y^1, u_z^1, L^2, \ell_x^2, \ell_y^2, \ell_z^2, u_x^2, u_y^2, u_z^2, L^3, \ell_x^3, \ell_y^3, \ldots ], \]

then, the first box will be inserted in level \(L^1\) and its lower and upper corner will have the indices \( (\ell_x^1, \ell_y^1, \ell_z^1)\) and \( (u_x^1, u_y^1, u_z^1) \) in the index space of level \(L^1\), respectively.

Parameters
boxesvector of size N (2d+1), where
N is the number of boxes,
d is the dimension of the parameter domain.
See description above for details on the format.

Reimplemented from gsBasis< T >.

void refineElements_withCoefs ( gsMatrix< T > &  coefs,
std::vector< index_t > const &  boxes 
)
virtual

Refine the basis and adjust the given matrix of coefficients accordingly.

Parameters
coefsis a matrix of coefficients as given, e.g., by gsTHBSpline<>::coefs();
boxesspecify where to refine; each 5-tuple gives the level of the box, then two indices (in the current level indexing) of the lower left corner and finally two indices of the upper right corner, see gsHTensorBasis::refineElements() for details.

Reimplemented from gsBasis< T >.

virtual const gsBasis& source ( ) const
inlinevirtualinherited

Applicable for rational bases: returns the underlying "source" (non-rational) basis

Reimplemented in gsRationalBasis< SrcT >, gsRationalBasis< gsBSplineTraits< d, T >::Basis >, and gsRationalBasis< gsBSplineBasis< T > >.

virtual gsBasis& source ( )
inlinevirtualinherited

Applicable for rational bases: returns the underlying "source" (non-rational) basis

Reimplemented in gsRationalBasis< SrcT >, gsRationalBasis< gsBSplineTraits< d, T >::Basis >, and gsRationalBasis< gsBSplineBasis< T > >.

gsMatrix< T > support ( ) const
virtual

Returns the boundary basis for side s.

Returns a bounding box for the basis' domain

Reimplemented from gsBasis< T >.

gsMatrix< T > support ( const index_t i) const
virtual

Returns (a bounding box for) the support of the i-th basis function.

Returns a dx2 matrix, containing the two diagonally extreme corners of a hypercube.

Reimplemented from gsBasis< T >.

gsMatrix< T > supportInterval ( index_t  dir) const
inherited

Returns an interval that contains the parameter values in direction dir.

Returns a 1x2 matrix, containing the two endpoints of the interval.

bool testPartitionOfUnity ( const index_t  npts = 100,
const T  tol = 1e-12 
) const

Test the partition of unity.

Parameters
[in]nptsThe number of points in each direction
Returns
True if the basis has the parition of unity property.
void uniformCoarsen ( int  numKnots = 1)
virtual

Coarsen the basis uniformly by removing groups of numKnots consecutive knots, each knot removed mul times.

This function is the oposite of gsBasis::uniformRefine

The execution of

basis->uniformRefine (nKnots, mul)
basis->uniformCoarsen(nKnots);

results in no overall change in "basis". However,

basis->uniformCoarsen(nKnots);
basis->uniformRefine (nKnots, mul)

is not guaranteed to keep "basis" unchanged.

See Also
gsBasis::uniformRefine

Reimplemented from gsBasis< T >.

void uniformCoarsen_withCoefs ( gsMatrix< T > &  coefs,
int  numKnots = 1 
)
virtual

Coarsen the basis uniformly.

The function simultainously updates the vector coefs, representing a function in the bases, such that its new version represents the same function.

This function is equivalent to

gsSparseMatrix<T,RowMajor> transfer;
basis->uniformCoarsen_withTransfer(transfer, numKnots);
coefs = transfer * coefs;
See Also
gsBasis::uniformRefine

Reimplemented from gsBasis< T >.

void uniformCoarsen_withTransfer ( gsSparseMatrix< T, RowMajor > &  transfer,
int  numKnots = 1 
)
virtualinherited

Coarsen the basis uniformly and produce a sparse matrix which maps coarse coefficient vectors to refined ones.

The function writes a sparse matrix into the variable transfer that indicates how the functions on the coarse grid are represented as linear combinations as fine grid functions

See Also
gsBasis::uniformCoarsen

Reimplemented in gsTensorBSplineBasis< 1, T >, gsTensorBasis< d, T >, and gsTensorBasis< 1, T >.

void uniformRefine_withCoefs ( gsMatrix< T > &  coefs,
int  numKnots = 1,
int  mul = 1,
int  dir = -1 
)
virtual

Refine the basis uniformly.

The function simultainously updates the vector coefs, representing a function in the bases, such that its new version represents the same function.

This function is equivalent to

gsSparseMatrix<T,RowMajor> transfer;
basis->uniformRefine_withTransfer(transfer, numKnots, mul);
coefs = transfer * coefs;
See Also
gsBasis::uniformRefine

Reimplemented from gsBasis< T >.

void uniformRefine_withTransfer ( gsSparseMatrix< T, RowMajor > &  transfer,
int  numKnots = 1,
int  mul = 1 
)
virtualinherited

Refine the basis uniformly.

The function writes a sparse matrix into the variable transfer that indicates how the functions on the coarse grid are represented as linear combinations as fine grid functions

See Also
gsBasis::uniformRefine

Reimplemented in gsTensorBSplineBasis< 1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, and gsRationalBasis< SrcT >.

void unrefineElements ( std::vector< index_t > const &  boxes)
virtual

Clear the given boxes into the quadtree.

Parameters
boxesSee refineElements
[in]refExtSee refineElements

Reimplemented from gsBasis< T >.

Member Data Documentation

std::vector<tensorBasis*> m_bases
mutableprotected

The list of nested spaces.

See documentation for the class for details on the underlying structure.

Recall that the hierarchical basis is built from a sequence of underlying bases \( B^0, B^1,\ldots, B^L\). These underlying bases are stored in gsHTensorBasis.m_bases, which is of type std::vector.
m_bases[k] stores the pointer to the (global) tensor-product basis \( B^k\).

std::vector< CMatrix > m_xmatrix
protected

The characteristic matrices for each level.

See documentation for the class for details on the underlying structure.

Characteristic matrices provide information on the relation between
the basis functions of this gsHTensorBasis \( H \) and
the tensor-product basis functions of the underlying tensor-product bases \( B^\ell \).

Let vk = m_xmatrix[k]. vk is a gsSortedVector. It contains a list of indices of the basis function of level k, i.e., of the basis functions which "are taken" from \(B^k\). These indices are stored as the global indices in \(B^k\).

std::vector<index_t> m_xmatrix_offset
protected

Stores the offsets of active functions for all levels.

See documentation for the class for details on the underlying structure. As mentioned there, the basis functions of the hierarchical basis \( H\) have a global numbering, where the functions from \(B^0\) come first, then those from \(B^1\), then \(B^2\), and so forth.

The entry m_xmatrix_offset[k] indicates the index from which the basis functions from level k (i.e., those taken from \( B^k \)) start.