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

Detailed Description

template<class SrcT>
class gismo::gsRationalBasis< SrcT >

Class that creates a rational counterpart for a given basis.

A rational basis holds an inner (referred to as "source") basis of the given source type, and a matrix of coefficients defining a weight function in terms of the source basis.

If w_i is the i-th weight coefficient and b_i the i-th basis function of the source basis, then the i-th basis function of the resulting rational basis is given by r_i(u) = b_i(u) w_i / (sum(b_j(u) w_j, j = 1, ..., size))

If the weights are all equal to one (or all equal to a constant), the rational basis is identical (or a scalar multiple) to the source basis.

Example: the rational version of a B-spline basis is a NURBS basis.

Template Parameters
SrcTthe basis of which to create a rational version
+ Inheritance diagram for gsRationalBasis< SrcT >:
+ Collaboration diagram for gsRationalBasis< SrcT >:

Public Types

typedef SrcT SourceBasis
 Associated source basis type.
 

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...
 
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...
 
virtual void active_into (const gsMatrix< T > &u, gsMatrix< index_t > &result) const
 Indices of active (non-zero) function(s) for each point. More...
 
gsMatrix< index_tallBoundary () const
 Returns the indices of the basis functions that are nonzero at the domain boundary.
 
void anchor_into (index_t i, gsMatrix< T > &result) const
 Returns the anchor point for member i of the basis.
 
gsMatrix< SrcT::Scalar_t > anchors () const
 Returns the anchor points that represent the members of the basis. There is exactly one anchor point per basis function. More...
 
void anchors_into (gsMatrix< T > &result) const
 Returns the anchor points that represent the members of the basis in result. There is exactly one anchor point per basis function. More...
 
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< SrcT::Scalar_t >::uPtr boundaryBasis (boxSide const &s)
 Returns the boundary basis for side s.
 
gsMatrix< index_tboundaryOffset (boxSide const &s, index_t offset) const
 
bool check () const
 Check the rational basis for consistency.
 
uPtr clone ()
 Clone methode. Produceds a deep copy inside a uPtr.
 
gsSparseMatrix< SrcT::Scalar_t > collocationMatrix (gsMatrix< SrcT::Scalar_t > const &u) const
 Computes the collocation matrix w.r.t. points u. More...
 
virtual const gsBasis< T > & component (short_t i) const
 For a tensor product basis, return the (const) 1-d basis for the i-th parameter component.
 
virtual gsBasis< SrcT::Scalar_t > & component (short_t i)
 For a tensor product basis, return the 1-d basis for the i-th parameter component.
 
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...
 
void connectivity (const gsMatrix< T > &nodes, gsMesh< T > &mesh) const
 
virtual void connectivityAtAnchors (gsMesh< SrcT::Scalar_t > &mesh) const
 
virtual gsBasis::uPtr create () const
 Create an empty basis of the derived type and return a pointer to it.
 
short_t degree (short_t i=0) const
 Degree with respect to the i-th variable. If the basis is a tensor product of (piecewise) polynomial bases, then this function returns the polynomial degree of the i-th component.
 
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.
 
void degreeElevate (short_t const &i=1, short_t const dir=-1)
 Elevate the degree of the basis by the given amount, preserve smoothness.
 
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.
 
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...
 
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 deriv2_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Second derivatives. More...
 
void deriv_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluates the first partial derivatives of the nonzero basis function. More...
 
virtual void deriv_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 First derivatives. More...
 
virtual std::string detail () const
 Prints the object as a string with extended details.
 
gsDomain< T > * domain () const
 
short_t domainDim () const
 Dimension of the (source) domain. More...
 
size_t elementIndex (const gsVector< T > &u) const
 See gsBasis for a description.
 
virtual gsMatrix< SrcT::Scalar_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...
 
void eval_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluates nonzero basis functions at point u into result. More...
 
virtual void eval_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluates the function(s). More...
 
std::vector< gsMatrix< T > > evalAllDers (const gsMatrix< T > &u, int n) const
 Evaluate all derivatives upto order n,. More...
 
virtual void evalAllDers_into (const gsMatrix< T > &u, int n, std::vector< gsMatrix< T > > &result) const
 Evaluate the nonzero functions and their derivatives up to order n. If n is -1 then no computation is performed.
 
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...
 
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.
 
gsBasisFun< SrcT::Scalar_t > function (index_t i) const
 Returns the i-th basis function as a gsFunction. More...
 
virtual SrcT::Scalar_t getMaxCellLength () const
 Get the maximum mesh size, as expected for approximation error estimates.
 
virtual SrcT::Scalar_t getMinCellLength () const
 Get the minimum mesh size, as expected for inverse inequalities.
 
 gsRationalBasis ()
 Default empty constructor.
 
 gsRationalBasis (SrcT *basis)
 Construct a rational counterpart of basis.
 
 gsRationalBasis (const SrcT &basis)
 Construct a rational counterpart of basis.
 
 gsRationalBasis (SrcT *basis, gsMatrix< T > w)
 Construct a rational counterpart of basis.
 
 gsRationalBasis (const gsRationalBasis &o)
 Copy Constructor.
 
virtual memory::unique_ptr
< gsGeometry< SrcT::Scalar_t > > 
interpolateAtAnchors (gsMatrix< SrcT::Scalar_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
< SrcT::Scalar_t > > 
interpolateData (gsMatrix< SrcT::Scalar_t > const &vals, gsMatrix< SrcT::Scalar_t > const &pts) const
 Applies interpolation given the parameter values pts and values vals.
 
virtual bool isRational () const
 Returns true, since by definition a gsRationalBasis is rational.
 
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< SrcT::Scalar_t > > 
makeGeometry (gsMatrix< SrcT::Scalar_t > coefs) const =0
 Create a gsGeometry of proper type for this basis with the given coefficient matrix.
 
memory::unique_ptr< gsBasis< T > > makeNonRational () const
 
virtual 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.
 
virtual void matchWith (const boundaryInterface &bi, const gsBasis< SrcT::Scalar_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.
 
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.
 
size_t numElements () const
 The number of elements.
 
size_t numElements (boxSide const &s) const
 The number of elements on side s.
 
gsRationalBasisoperator= (const gsRationalBasis &o)
 Assignment operator.
 
const gsBasis< SrcT::Scalar_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.
 
gsMatrix< T > projectiveCoefs (const gsMatrix< T > &coefs) const
 
virtual void reduceContinuity (int const &i=1)
 Reduces the continuity of the basis along element boundaries.
 
virtual void refine (gsMatrix< SrcT::Scalar_t > const &boxes, int refExt=0)
 Refine the basis on the area defined by the matrix boxes. More...
 
void refineElements (std::vector< index_t > const &boxes)
 Refines specified areas or boxes, depending on underlying basis. More...
 
void refineElements_withCoefs (gsMatrix< T > &coefs, std::vector< index_t > const &boxes)
 Refines specified areas or boxes, depending on underlying basis. More...
 
virtual void reverse ()
 Reverse the basis.
 
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.
 
void setWeights (gsMatrix< T > const &w)
 Set weights.
 
index_t size () const
 size More...
 
const SrcT & source () const
 Returns the source basis of the rational basis.
 
SrcT & source ()
 
gsMatrix< T > support () const
 Returns (a bounding box for) the domain of the whole basis. More...
 
gsMatrix< T > support (const index_t &i) const
 Returns (a bounding box for) the support of the i-th basis function. More...
 
gsMatrix< SrcT::Scalar_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.
 
short_t totalDegree () const
 If the basis is of polynomial or piecewise polynomial type, then this function returns the total polynomial degree.
 
virtual void uniformCoarsen (int numKnots=1)
 Coarsen the basis uniformly by removing groups of numKnots consecutive knots, each knot removed mul times. More...
 
virtual void uniformCoarsen_withCoefs (gsMatrix< SrcT::Scalar_t > &coefs, int numKnots=1)
 Coarsen the basis uniformly. More...
 
virtual void uniformCoarsen_withTransfer (gsSparseMatrix< SrcT::Scalar_t, RowMajor > &transfer, int numKnots=1)
 Coarsen the basis uniformly and produce a sparse matrix which maps coarse coefficient vectors to refined ones. More...
 
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...
 
void uniformRefine_withTransfer (gsSparseMatrix< T, RowMajor > &transfer, int numKnots=1, int mul=1)
 Refine the basis uniformly. More...
 
T & weight (int i)
 Access to i-th weight.
 
const T & weight (int i) const
 Const access to i-th weight.
 
const gsMatrix< T > & weights () const
 Returns the weights of the rational basis.
 
gsMatrix< T > & weights ()
 Returns the weights of the rational basis.
 
Evaluation functions
gsMatrix< SrcT::Scalar_t > evalSingle (index_t i, const gsMatrix< SrcT::Scalar_t > &u) const
 Evaluate a single basis function i at points u.
 
gsMatrix< SrcT::Scalar_t > derivSingle (index_t i, const gsMatrix< SrcT::Scalar_t > &u) const
 Evaluate a single basis function i derivative at points u.
 
gsMatrix< SrcT::Scalar_t > deriv2Single (index_t i, const gsMatrix< SrcT::Scalar_t > &u) const
 Evaluate the second derivative of a single basis function i at points u.
 
gsVector< index_tnumActive (const gsMatrix< SrcT::Scalar_t > &u) const
 Number of active basis functions at an arbitrary parameter value. More...
 
virtual void numActive_into (const gsMatrix< SrcT::Scalar_t > &u, gsVector< index_t > &result) const
 Returns the number of active (nonzero) basis functions at points u in result.
 
virtual bool isActive (const index_t i, const gsVector< SrcT::Scalar_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< SrcT::Scalar_t > &u, const gsMatrix< SrcT::Scalar_t > &coefs, gsMatrix< SrcT::Scalar_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 derivSingle_into (index_t i, const gsMatrix< SrcT::Scalar_t > &u, gsMatrix< SrcT::Scalar_t > &result) const
 Evaluates the (partial) derivatives of the i-th basis function at points u into result. More...
 
virtual void deriv2Single_into (index_t i, const gsMatrix< SrcT::Scalar_t > &u, gsMatrix< SrcT::Scalar_t > &result) const
 Evaluate the (partial) derivatives of the i-th basis function at points u into result.
 
virtual void evalAllDers_into (const gsMatrix< SrcT::Scalar_t > &u, int n, std::vector< gsMatrix< SrcT::Scalar_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< SrcT::Scalar_t > &u, int n, gsMatrix< SrcT::Scalar_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< SrcT::Scalar_t > &u, int n, gsMatrix< SrcT::Scalar_t > &result) const
 Evaluate the (partial) derivative(s) of order n the i-th basis function at points u into result.
 
virtual gsMatrix< SrcT::Scalar_t > laplacian (const gsMatrix< SrcT::Scalar_t > &u) const
 Compute the Laplacian of all nonzero basis functions at points u.
 

Static Public Member Functions

static gsMatrix< T > projectiveCoefs (const gsMatrix< T > &coefs, const gsMatrix< T > &weights)
 
static void setFromProjectiveCoefs (const gsMatrix< T > &pr_coefs, gsMatrix< T > &coefs, gsMatrix< T > &weights)
 

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< SrcT::Scalar_t > evalFunc (const gsMatrix< SrcT::Scalar_t > &u, const gsMatrix< SrcT::Scalar_t > &coefs) const
 Evaluate the function described by coefs at points u. More...
 
gsMatrix< SrcT::Scalar_t > derivFunc (const gsMatrix< SrcT::Scalar_t > &u, const gsMatrix< SrcT::Scalar_t > &coefs) const
 Evaluate the derivatives of the function described by coefs at points u. More...
 
virtual void derivFunc_into (const gsMatrix< SrcT::Scalar_t > &u, const gsMatrix< SrcT::Scalar_t > &coefs, gsMatrix< SrcT::Scalar_t > &result) const
 Evaluate the derivatives of the function described by coefs at points u. More...
 
virtual void jacobianFunc_into (const gsMatrix< SrcT::Scalar_t > &u, const gsMatrix< SrcT::Scalar_t > &coefs, gsMatrix< SrcT::Scalar_t > &result) const
 Evaluate the Jacobian of the function described by coefs at points u. Jacobian matrices are stacked in blocks.
 
gsMatrix< SrcT::Scalar_t > deriv2Func (const gsMatrix< SrcT::Scalar_t > &u, const gsMatrix< SrcT::Scalar_t > &coefs) const
 Evaluates the second derivatives of the function described by coefs at points u. More...
 
virtual void deriv2Func_into (const gsMatrix< SrcT::Scalar_t > &u, const gsMatrix< SrcT::Scalar_t > &coefs, gsMatrix< SrcT::Scalar_t > &result) const
 Evaluates the second derivatives of the function described by coefs at points u. More...
 
virtual void evalAllDersFunc_into (const gsMatrix< SrcT::Scalar_t > &u, const gsMatrix< SrcT::Scalar_t > &coefs, const unsigned n, std::vector< gsMatrix< SrcT::Scalar_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< SrcT::Scalar_t > &coefs, const gsMatrix< index_t > &actives, const gsMatrix< SrcT::Scalar_t > &values, gsMatrix< SrcT::Scalar_t > &result)
 Computes the linear combination coefs * values( actives ) More...
 

Member Function Documentation

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
inlinevirtual

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< SrcT::Scalar_t >.

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

Indices of active (non-zero) function(s) for each point.

The columns are sorted in increasing order, if on a point there are less active then the number of rows in the result matrix (some other point has more actives) then the rest of the column is filled with 0s.

Parameters
u
result

Reimplemented in gsHTensorBasis< d, T >, gsBasis< T >, gsBasis< Scalar >, gsBasis< real_t >, gsTensorBSplineBasis< 1, T >, gsTensorBSplineBasis< d, T >, gsTensorBSplineBasis< domainDim+1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsTHBSplineBasis< d, T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsFunction< T >, gsConstantBasis< T >, gsLegendreBasis< T >, gsMonomialBasis< T >, and gsMvLegendreBasis< T >.

virtual void activeCoefs_into ( const gsVector< SrcT::Scalar_t > &  u,
const gsMatrix< SrcT::Scalar_t > &  coefs,
gsMatrix< SrcT::Scalar_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<SrcT::Scalar_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.

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

Returns the anchor points that represent the members of the basis in result. 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.

Reimplemented from gsBasis< SrcT::Scalar_t >.

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

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< SrcT::Scalar_t >.

gsSparseMatrix<SrcT::Scalar_t > collocationMatrix ( gsMatrix< SrcT::Scalar_t > const &  u) const
inherited

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.

virtual 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
inlinevirtual

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

Reimplemented from gsBasis< SrcT::Scalar_t >.

virtual void connectivityAtAnchors ( gsMesh< SrcT::Scalar_t > &  mesh) const
virtualinherited

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
virtual

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 gsBasis< SrcT::Scalar_t >.

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

Second derivatives.

For scalar valued functions \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}\) format is:

\[ \left[ \begin{array}{ccccc} \partial_{1}\partial_{1}f_1(p_1) & \partial_{1}\partial_{1}f_1(p_2) & \ldots & \partial_{1}\partial_{1}f_1(p_N)\\ \partial_{2}\partial_{2}f_1(p_1) & \partial_{2}\partial_{2}f_1(p_2) & \ldots & \partial_{2}\partial_{2}f_1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}\partial_{k}f_1(p_1) & \partial_{k}\partial_{k}f_1(p_2) & \ldots & \partial_{k}\partial_{k}f_1(p_N)\\ \partial_{1}\partial_{2}f_1(p_1) & \partial_{1}\partial_{2}f_1(p_2) & \ldots & \partial_{1}\partial_{2}f_1(p_N)\\ \partial_{1}\partial_{3}f_1(p_1) & \partial_{1}\partial_{3}f_1(p_2) & \ldots & \partial_{1}\partial_{3}f_1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{1}\partial_{k}f_1(p_1) & \partial_{1}\partial_{k}f_1(p_2) & \ldots & \partial_{1}\partial_{k}f_1(p_N)\\ \partial_{2}\partial_{3}f_1(p_1) & \partial_{2}\partial_{3}f_1(p_2) & \ldots & \partial_{2}\partial_{3}f_1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{2}\partial_{k}f_1(p_1) & \partial_{2}\partial_{k}f_1(p_2) & \ldots & \partial_{2}\partial_{k}f_1(p_N)\\ \partial_{3}\partial_{4}f_1(p_1) & \partial_{3}\partial_{4}f_1(p_2) & \ldots & \partial_{3}\partial_{4}f_1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k-1}\partial_{k}f_1(p_1) & \partial_{k-1}\partial_{k}f_1(p_2) & \ldots & \partial_{k-1}\partial_{k}f_1(p_N)\\ \partial_{1}\partial_{1}f_2(p_1) & \partial_{1}\partial_{1}f_2(p_2) & \ldots & \partial_{1}\partial_{1}f_2(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k-1}\partial_{k}f_S(p_1) & \partial_{k-1}\partial_{k}f_S(p_2) & \ldots & \partial_{k-1}\partial_{k}f_S(p_N)\\ \end{array} \right] \]

For vector valued functions function \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}^{m}\) the format is:

\[ \left[ \begin{array}{ccccc} \partial_{1}\partial_{1}f_1^{(1)}(p_1) & \partial_{1}\partial_{1}f_1^{(1)}(p_2) & \ldots & \partial_{1}\partial_{1}f_1^{(1)}(p_N)\\ \partial_{2}\partial_{2}f_1^{(1)}(p_1) & \partial_{2}\partial_{2}f_1^{(1)}(p_2) & \ldots & \partial_{2}\partial_{2}f_1^{(1)}(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}\partial_{k}f_1^{(1)}(p_1) & \partial_{k}\partial_{k}f_1^{(1)}(p_2) & \ldots & \partial_{k}\partial_{k}f_1^{(1)}(p_N)\\ \partial_{1}\partial_{2}f_1^{(1)}(p_1) & \partial_{1}\partial_{2}f_1^{(1)}(p_2) & \ldots & \partial_{1}\partial_{2}f_1^{(1)}(p_N)\\ \partial_{1}\partial_{3}f_1^{(1)}(p_1) & \partial_{1}\partial_{3}f_1^{(1)}(p_2) & \ldots & \partial_{1}\partial_{3}f_1^{(1)}(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k-1}\partial_{k}f_1{(1)}(p_1) & \partial_{k-1}\partial_{k}f_1^{(1)}(p_2) & \ldots & \partial_{k-1}\partial_{k}f_1^{(1)}(p_N)\\ \partial_{1}\partial_{1}f_1^{(2)}(p_1) & \partial_{1}\partial_{1}f_1^{(2)}(p_2) & \ldots & \partial_{1}\partial_{1}f_1^{(2)}(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k-1}\partial_{k}f_1^{(m)}(p_1) & \partial_{k-1}\partial_{k}f_1^{(m)}(p_2) & \ldots & \partial_{k-1}\partial_{k}f_1^{(m)}(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k-1}\partial_{k}f_S^{(m)}(p_1) & \partial_{k-1}\partial_{k}f_S^{(m)}(p_2) & \ldots & \partial_{k-1}\partial_{k}f_S^{(m)}(p_N)\\ \end{array} \right] \]

where \( f^{(i)}_j\) is the \(i\)-th component of function \( f_j\) of the set.

Parameters
u
result

Reimplemented in gsBasis< T >, gsBasis< Scalar >, gsBasis< real_t >, gsTensorBSplineBasis< 1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsGeometry< T >, gsFunctionExpr< T >, gsFunction< T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsConstantFunction< T >, gsTHBSplineBasis< d, T >, gsPreCICEFunction< T >, gsLegendreBasis< T >, gsAffineFunction< T >, gsHBSplineBasis< d, T >, gsMappedSingleSpline< d, T >, gsMonomialBasis< T >, gsMvLegendreBasis< T >, gsGeometryTransform< T >, and gsFuncCoordinate< T >.

gsMatrix<SrcT::Scalar_t > deriv2Func ( const gsMatrix< SrcT::Scalar_t > &  u,
const gsMatrix< SrcT::Scalar_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.
virtual void deriv2Func_into ( const gsMatrix< SrcT::Scalar_t > &  u,
const gsMatrix< SrcT::Scalar_t > &  coefs,
gsMatrix< SrcT::Scalar_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
virtual

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 gsBasis< SrcT::Scalar_t >.

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

First derivatives.

For scalar valued functions \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}\) format is:

\[ \left[ \begin{array}{ccccc} \partial_{1}f_1(p_1) & \partial_{1}f_1(p_2) & \ldots & \partial_{1}f_1(p_N)\\ \partial_{2}f_1(p_1) & \partial_{2}f_1(p_2) & \ldots & \partial_{2}f_1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}f_1(p_1) & \partial_{k}f_1(p_2) & \ldots & \partial_{k}f_1(p_N)\\ \partial_{1}f_2(p_1) & \partial_{1}f_2(p_2) & \ldots & \partial_{1}f_2(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}f_S(p_1) & \partial_{k}f_S(p_2) & \ldots & \partial_{k}f_S(p_N)\\ \end{array} \right] \]

For vector valued functions function \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}^{m}\) the format is:

\[ \left[ \begin{array}{ccccc} \partial_{1}f_1^1(p_1) & \partial_{1}f_1^1(p_2) & \ldots & \partial_{1}f_1^1(p_N)\\ \partial_{2}f_1^1(p_1) & \partial_{1}f_2^1(p_2) & \ldots & \partial_{1}f_2^1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}f_1^1(p_1) & \partial_{k}f_1^1(p_2) & \ldots & \partial_{k}f_1^1(p_N)\\ \partial_{1}f_1^2(p_1) & \partial_{1}f_1^2(p_2) & \ldots & \partial_{1}f_1^2(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}f_1^2(p_1) & \partial_{k}f_1^2(p_2) & \ldots & \partial_{k}f_1^2(p_N)\\ \partial_{1}f_2^1(p_1) & \partial_{1}f_2^1(p_2) & \ldots & \partial_{1}f_2^1(p_N)\\ \vdots & \vdots & & \vdots\\ \partial_{k}f_S^{(m)}(p_1) & \partial_{k}f_S^{(m)}(p_2) & \ldots & \partial_{k}f_S^{(m)}(p_N) \end{array} \right] \]

where \(f^{(i)}_j\) is the \(i\)-th component of function \(f_j\) of the set.

Parameters
u
result

Reimplemented in gsBasis< T >, gsBasis< Scalar >, gsBasis< real_t >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsTensorBSplineBasis< 1, T >, gsGeometry< T >, gsFunctionExpr< T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsFunction< T >, gsConstantFunction< T >, gsTHBSplineBasis< d, T >, gsPreCICEFunction< T >, gsAffineFunction< T >, gsHBSplineBasis< d, T >, gsLegendreBasis< T >, gsMappedSingleSpline< d, T >, gsMvLegendreBasis< T >, gsMonomialBasis< T >, gsConstantBasis< T >, gsGeometrySlice< T >, gsSquaredDistance< T >, gsBasisFun< T >, gsGeometryTransform< T >, and gsFuncCoordinate< T >.

gsMatrix<SrcT::Scalar_t > derivFunc ( const gsMatrix< SrcT::Scalar_t > &  u,
const gsMatrix< SrcT::Scalar_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
virtual void derivFunc_into ( const gsMatrix< SrcT::Scalar_t > &  u,
const gsMatrix< SrcT::Scalar_t > &  coefs,
gsMatrix< SrcT::Scalar_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

virtual void derivSingle_into ( index_t  i,
const gsMatrix< SrcT::Scalar_t > &  u,
gsMatrix< SrcT::Scalar_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
gsDomain<T>* domain ( ) const
inlinevirtual

Return the gsDomain which represents the parameter domain of this basis. Currently unused.

Reimplemented from gsBasis< SrcT::Scalar_t >.

short_t domainDim ( ) const
inlinevirtual

Dimension of the (source) domain.

Returns
For \(f:\mathbb{R}^n\rightarrow\mathbb{R}^m\), returns \(n\).

Implements gsFunctionSet< T >.

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
virtual

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 gsBasis< SrcT::Scalar_t >.

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

Evaluates the function(s).

For scalar valued functions \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}\) format is:

\[ \left[ \begin{array}{ccccc} f_1(p_1) & f_1(p_2) & \ldots & f_1(p_N)\\ f_2(p_1) & f_2(p_2) & \ldots & f_2(p_N)\\ \vdots & \vdots & & \vdots\\ f_S(p_1) & f_S(p_2) & \ldots & f_S(p_N) \end{array} \right] \]

For vector valued functions function \(f_1, \ldots, f_S\) from \(\mathbb{R}^n\rightarrow\mathbb{R}^m\) the format is:

\[ \left[ \begin{array}{ccccc} f_1^1(p_1) & f_1^{(1)}(p_2) & \ldots & f_1^{(1)}(p_N)\\ f_1^2(p_1) & f_1^{(2)}(p_2) & \ldots & f_1^{(2)}(p_N)\\ \vdots & \vdots & & \vdots\\ f_1^{(m)}(p_1) & f_1^{(m)}(p_2) & \ldots & f_1^{(m)}(p_N)\\ f_2^{(1)}(p_1) & f_2^{(1)}(p_2) & \ldots & f_2^{(1)}(p_N)\\ \vdots & \vdots & & \vdots\\ f_S^{(m)}(p_1) & f_S^{(m)}(p_2) & \ldots & f_S^{(m)}(p_N) \end{array} \right] \]

where \(f^{(i)}_j\) is the \(i\)-th component of function \(f_j\) of the set.

Parameters
u
result

Reimplemented in gsBasis< T >, gsBasis< Scalar >, gsBasis< real_t >, gsTHBSplineBasis< d, T >, gsPatchIdField< T >, gsMaterialMatrixEvalSingle< T, out >, gsParamField< T >, gsTensorBSplineBasis< 1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsShellStressFunction< T >, gsMaterialMatrixIntegrateSingle< T, out >, gsNormalField< T >, gsFsiLoad< T >, gsFunctionExpr< T >, gsGeometry< T >, gsLagrangeBasis< T >, gsMappedSingleBasis< d, T >, gsPiecewiseFunction< T >, gsConstantFunction< T >, gsJacDetField< T >, gsDetFunction< T >, gsRemapInterface< T >, gsFunction< T >, gsGradientField< T >, gsCPPInterface< T >, gsHBSplineBasis< d, T >, gsAffineFunction< T >, gsPreCICEFunction< T >, gsLegendreBasis< T >, gsMvLegendreBasis< T >, gsMappedSingleSpline< d, T >, gsConstantBasis< T >, gsMonomialBasis< T >, gsCauchyStressFunction< T >, gsBasisFun< T >, gsGeometrySlice< T >, gsGeometryTransform< T >, gsFuncCoordinate< T >, gsAbsError< T >, gsElementErrorPlotter< T >, and gsSquaredDistance< T >.

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

Evaluate all derivatives upto order n,.

See Also
evalAllDers_into
virtual void evalAllDers_into ( const gsMatrix< SrcT::Scalar_t > &  u,
int  n,
std::vector< gsMatrix< SrcT::Scalar_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.
virtual void evalAllDersFunc_into ( const gsMatrix< SrcT::Scalar_t > &  u,
const gsMatrix< SrcT::Scalar_t > &  coefs,
const unsigned  n,
std::vector< gsMatrix< SrcT::Scalar_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<SrcT::Scalar_t > evalFunc ( const gsMatrix< SrcT::Scalar_t > &  u,
const gsMatrix< SrcT::Scalar_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
virtual

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 from gsBasis< SrcT::Scalar_t >.

gsBasisFun<SrcT::Scalar_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.

static void linearCombination_into ( const gsMatrix< SrcT::Scalar_t > &  coefs,
const gsMatrix< index_t > &  actives,
const gsMatrix< SrcT::Scalar_t > &  values,
gsMatrix< SrcT::Scalar_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
memory::unique_ptr<gsBasis<T> > makeNonRational ( ) const
inlinevirtual

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

Reimplemented from gsBasis< SrcT::Scalar_t >.

virtual void matchWith ( const boundaryInterface bi,
const gsBasis< SrcT::Scalar_t > &  other,
gsMatrix< index_t > &  bndThis,
gsMatrix< index_t > &  bndOther,
index_t  offset 
) const
virtualinherited

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)

gsVector<index_t> numActive ( const gsMatrix< SrcT::Scalar_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.

gsMatrix<T> projectiveCoefs ( const gsMatrix< T > &  coefs) const
inline

Returns a matrix of projective coefficients. The input coefs are affine coefficients for this basis

static gsMatrix<T> projectiveCoefs ( const gsMatrix< T > &  coefs,
const gsMatrix< T > &  weights 
)
inlinestatic

Returns a matrix of projective coefficients. The input coefs are affine coefficients and weights

virtual void refine ( gsMatrix< SrcT::Scalar_t > const &  boxes,
int  refExt = 0 
)
virtualinherited

Refine the basis on the area defined by the matrix boxes.

boxes is a d x n-matrix (n even), where d is the dimension of the parameter domain.
n must be even, and every 2 successive columns in the matrix define a box in the parameter domain (the first column represents the coordinates of the lower corner, the second column the coordinates of the upper corner).

Example: The input of the matrix

\[ \left[\begin{array}{cccc} 0 & 0.2 & 0.8 & 1 \\ 0.4 & 0.6 & 0.2 & 0.4 \end{array} \right] \]

results in refinement of the two boxes \([0,0.2]\times[0.4,0.6]\) and \([0.8,1]\times[0.2,0.4]\).

Parameters
[in]boxesgsMatrix of size d x n, see above for description of size and meaning.
[in]refExtExtension to be applied to the refinement boxes
void refineElements ( std::vector< index_t > const &  boxes)
inlinevirtual

Refines specified areas or boxes, depending on underlying basis.

Parameters
boxesSee the function gsBasis::refineElements() of the underlying basis for syntax.

Reimplemented from gsBasis< SrcT::Scalar_t >.

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

Refines specified areas or boxes, depending on underlying basis.

Parameters
coefsCoefficients, given as gsMatrix of size \( n \times d\), where \(n\) is the number of basis functions and \(d\) is the target dimension.
boxesSee the function gsBasis::refineElements() of the underlying basis for syntax.

Reimplemented from gsBasis< SrcT::Scalar_t >.

static void setFromProjectiveCoefs ( const gsMatrix< T > &  pr_coefs,
gsMatrix< T > &  coefs,
gsMatrix< T > &  weights 
)
inlinestatic

Sets the weights and the coefs to be the affine coefficients corresponding to the projective coefficients pr_coefs

index_t size ( ) const
inlinevirtual

size

Warning
gsFunction and gsGeometry have size() == 1. This should not be confused with the size eg. of gsGeometry::basis(), which is the number of basis functions in the basis
Returns
the size of the function set: the total number of functions

Reimplemented from gsFunctionSet< T >.

SrcT& source ( )
inlinevirtual

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

Reimplemented from gsBasis< SrcT::Scalar_t >.

gsMatrix<T> support ( ) const
inlinevirtual

Returns (a bounding box for) the domain of the whole basis.

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

Reimplemented from gsBasis< SrcT::Scalar_t >.

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

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< SrcT::Scalar_t >.

gsMatrix<SrcT::Scalar_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.

virtual void uniformCoarsen ( int  numKnots = 1)
virtualinherited

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
virtual void uniformCoarsen_withCoefs ( gsMatrix< SrcT::Scalar_t > &  coefs,
int  numKnots = 1 
)
virtualinherited

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
virtual void uniformCoarsen_withTransfer ( gsSparseMatrix< SrcT::Scalar_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
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< SrcT::Scalar_t >.

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

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 from gsBasis< SrcT::Scalar_t >.