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

Detailed Description

template<class T>
class gismo::gsMultiBasis< T >

Holds a set of patch-wise bases and their topology information.

Template Parameters
Tcoefficient type
+ Inheritance diagram for gsMultiBasis< T >:
+ Collaboration diagram for gsMultiBasis< T >:

Public Types

typedef BasisContainer::iterator iterator
 Type definitions.
 

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
 Indices of active (non-zero) function(s) for each point. More...
 
void addBasis (gsBasis< T > *g)
 Add a basis (ownership of the pointer is also acquired)
 
void addBasis (typename gsBasis< T >::uPtr g)
 Add a basis.
 
void addInterface (gsBasis< T > *g1, boxSide s1, gsBasis< T > *g2, boxSide s2)
 Add an interface joint between side s1 of geometry g1 side s2 of geometry g2. More...
 
void addPatchBoundary (gsBasis< T > *g, boxSide s)
 Add side s of patch g to the outer boundary of the domain.
 
const_reference at (size_t i) const
 Assess i-th parametric basis.
 
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.
 
const gsBasis< T > & basis (const size_t i) const
 Return the i-th basis block.
 
gsBasis< T > & basis (const size_t i)
 Return the i-th basis block.
 
const_iterator begin () const
 
iterator begin ()
 
void clear ()
 Clear (delete) all patches.
 
uPtr clone ()
 Clone methode. Produceds a deep copy inside a uPtr.
 
gsBasis< T >::uPtr componentBasis (patchComponent p) const
 Returns the basis that corresponds to the component.
 
gsBasis< T >::uPtr componentBasis_withIndices (patchComponent pc, const gsDofMapper &dm, gsMatrix< index_t > &indices, bool no_lower=true) const
 Returns the basis that corresponds to the component. More...
 
std::vector< typename gsBasis
< T >::uPtr
componentBasis_withIndices (const std::vector< patchComponent > &pc, const gsDofMapper &dm, gsMatrix< index_t > &indices, bool no_lower=true) const
 Returns the bases that correspond to the components. More...
 
virtual void compute (const gsMatrix< T > &in, gsFuncData< T > &out) const
 Computes function data. More...
 
short_t degree (size_t i=0, short_t comp=0) const
 Returns the polynomial degree of basis i in component j, if the basis is of polynomial or piecewise polynomial type.
 
void degreeDecrease (short_t const i=1, short_t const dir=-1)
 Increase the degree of every basis by the given amount. (keeping the multiplicity)
 
void degreeElevate (short_t const i=1, short_t const dir=-1)
 Elevate the degree of every basis by the given amount. (keeping the smoothness)
 
void degreeIncrease (short_t const i=1, short_t const dir=-1)
 Increase the degree of every basis by the given amount. (keeping the multiplicity)
 
void degreeReduce (short_t const i=1)
 Reduce the degree of the basis by the given amount.
 
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 void deriv2_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Second derivatives. More...
 
virtual void deriv_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 First derivatives. More...
 
short_t dim () const
 Dimension of the parameter domain (must match for all bases).
 
short_t domainDim () const
 Dimension of the (source) domain. More...
 
const_iterator end () const
 
iterator end ()
 
gsMatrix< T > eval (const gsMatrix< T > &u) const
 Evaluate the function,. 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.
 
int findBasisIndex (gsBasis< T > *g) const
 Search for the given basis and return its index.
 
const gsFunction< T > & function (const index_t k) const
 Helper which casts and returns the k-th piece of this function set as a gsFunction.
 
 gsMultiBasis ()
 Default empty constructor.
 
 gsMultiBasis (const gsMultiPatch< T > &mpatch, bool numeratorOnly=false)
 Create a multi-basis instance from a gsMultiPatch. More...
 
 gsMultiBasis (BasisContainer &bases, const gsBoxTopology &topology)
 Create from a vector of bases and topology.
 
 gsMultiBasis (const gsBasis< T > &geo, bool numeratorOnly=false)
 Create a single-basis instance.
 
 gsMultiBasis (BasisContainer &bases, const std::vector< patchSide > &boundary, const std::vector< boundaryInterface > &interfaces)
 Create from bases and boundary/interface information.
 
 gsMultiBasis (const gsMultiBasis &other)
 Copy constructor (makes deep copy)
 
void matchInterface (const boundaryInterface &bi, gsDofMapper &mapper) const
 Matches the degrees of freedom along an interface. More...
 
short_t maxCwiseDegree () const
 Maximum degree with respect to all variables.
 
short_t maxDegree (short_t k) const
 Maximum degree with respect to variable k.
 
short_t minCwiseDegree () const
 Minimum degree with respect to all variables.
 
short_t minDegree (short_t k) const
 Minimum degree with respect to variable k.
 
size_t nBases () const
 Number of patch-wise bases.
 
index_t nPieces () const
 Number of patch-wise bases.
 
gsMultiBasisoperator= (gsMultiBasis other)
 Assignment operator (uses copy-and-swap idiom)
 
const gsBasis< T > & piece (const index_t i) const
 Returns the piece(s) of the function(s) at subdomain k.
 
std::ostream & print (std::ostream &os) const
 Prints the object as a string.
 
void reduceContinuity (int const i=1)
 Reduce the continuity by i.
 
void refine (size_t k, gsMatrix< T > const &boxes, int refExt)
 Refine the are defined by boxes on patch k with extension refExt. More...
 
void refineElements (int k, std::vector< index_t > const &boxes)
 Refine the are defined by boxes on patch k. More...
 
bool repairInterface (const boundaryInterface &bi)
 Checks if the interface is fully matching, and if not, repairs it. More...
 
bool repairInterface2d (const boundaryInterface &bi)
 Checks if the 2D-interface is fully matching, and if not, repairs it. More...
 
template<short_t d>
bool repairInterfaceFindElements (const boundaryInterface &bi, std::vector< index_t > &refEltsFirst, std::vector< index_t > &refEltsSecond)
 Finds the elements that need to be refined in order to repair an interface. More...
 
void repairInterfaces (const std::vector< boundaryInterface > &bivec)
 Checks if the interfaces bivec are fully matching, and if not, repairs them, i.e., makes them fully matching. More...
 
void setDegree (short_t const &i)
 Set the degree of the basis.
 
int size (size_t i) const
 The number of basis functions in basis i.
 
index_t size () const
 size More...
 
void swap (gsMultiBasis &other)
 Swap with another gsMultiBasis.
 
short_t targetDim () const
 Dimension of the target space. More...
 
void tileParameters ()
 
size_t totalElements () const
 The total number of elements in all patches.
 
size_t totalSize () const
 The total number of basis functions in all bases.
 
void uniformCoarsen (int numKnots=1)
 Coarsen every basis uniformly. More...
 
void uniformCoarsen_withTransfer (gsSparseMatrix< T, RowMajor > &transfer, const gsBoundaryConditions< T > &boundaryConditions, const gsOptionList &assemblerOptions, int numKnots=1, index_t unk=0)
 Coarsen every basis uniformly. More...
 
void uniformRefine (int numKnots=1, int mul=1, int dir=-1)
 Refine every basis uniformly. More...
 
void uniformRefine_withTransfer (gsSparseMatrix< T, RowMajor > &transfer, const gsBoundaryConditions< T > &boundaryConditions, const gsOptionList &assemblerOptions, int numKnots=1, int mul=1, index_t unk=0)
 Refine every basis uniformly. More...
 
void uniformRefineComponent (int comp, int numKnots=1, int mul=1)
 Refine the component comp of every basis uniformly by inserting numKnots new knots on each knot span.
 
 ~gsMultiBasis ()
 Destructor.
 

Static Public Member Functions

static void combineTransferMatrices (const std::vector< gsSparseMatrix< T, RowMajor > > &localTransferMatrices, const gsDofMapper &coarseMapper, const gsDofMapper &fineMapper, gsSparseMatrix< T, RowMajor > &transferMatrix)
 This function takes local transfer matrices (per patch) and combines them using given DofMappers to a global transfer matrix. Simultanously, this function restricts the matrices to the free dofs, e.g., Dirichlet dofs might be eliminated. More...
 

Constructor & Destructor Documentation

gsMultiBasis ( const gsMultiPatch< T > &  mpatch,
bool  numeratorOnly = false 
)
explicit

Create a multi-basis instance from a gsMultiPatch.

Parameters
mpatchused gsMultiPatch
numeratorOnlyIf true, and the bases are derived from gsRationalBasis, then only the source bases (numerators) are returned
Note
In the case of NURBS, the numerator possess the same approximation power, while the evaluation of values and partial derivatives are much less expensive

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
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 >.

void addInterface ( gsBasis< T > *  g1,
boxSide  s1,
gsBasis< T > *  g2,
boxSide  s2 
)

Add an interface joint between side s1 of geometry g1 side s2 of geometry g2.

Todo:
add orientation information
const_iterator begin ( ) const
inline

Get a const-iterator to the patches

Returns
an iterator to the beginning of the patches
iterator begin ( )
inline

Get an iterator to the beginning of the patches

Returns
an iterator to the beginning of the patches
void combineTransferMatrices ( const std::vector< gsSparseMatrix< T, RowMajor > > &  localTransferMatrices,
const gsDofMapper coarseMapper,
const gsDofMapper fineMapper,
gsSparseMatrix< T, RowMajor > &  transferMatrix 
)
static

This function takes local transfer matrices (per patch) and combines them using given DofMappers to a global transfer matrix. Simultanously, this function restricts the matrices to the free dofs, e.g., Dirichlet dofs might be eliminated.

Parameters
[in]localTransferMatricesThe local and full (also non-free dofs) transfer matrices per patch
[in]coarseMapperThe DofMapper on the coarse grid
[in]fineMapperThe DofMapper on the fine grid
[out]transferMatrixThe combined transfer matrix restricted to the free dofs
gsBasis< T >::uPtr componentBasis_withIndices ( patchComponent  pc,
const gsDofMapper dm,
gsMatrix< index_t > &  indices,
bool  no_lower = true 
) const

Returns the basis that corresponds to the component.

Parameters
pcThe component
dmThe dof mapper to be used
indicesThe row vector where the indices are stored to
no_lowerIf true, the transfer matrix does not include parts belonging to lower-order components (i.e., edges without corners or faces without corners and edges)
std::vector< typename gsBasis< T >::uPtr > componentBasis_withIndices ( const std::vector< patchComponent > &  pc,
const gsDofMapper dm,
gsMatrix< index_t > &  indices,
bool  no_lower = true 
) const

Returns the bases that correspond to the components.

Parameters
pcThe components
dmThe dof mapper to be used
indicesThe row vector where the indices are stored to
no_lowerIf 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 >.

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

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 >.

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 >.

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 >.

const_iterator end ( ) const
inline

Get a const iterator to the end of the patches

Returns
an iterator to the end of the patches
iterator end ( )
inline

Get an iterator to the end of the patches

Returns
an iterator to the end of the patches
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 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
void matchInterface ( const boundaryInterface bi,
gsDofMapper mapper 
) const

Matches the degrees of freedom along an interface.

The boundaryInterface specifying the interface is passed as argument.
The degrees of freedom (DOFs) along the interface from both patches are matched in the sense that they are then treated as one global DOF by the mapper.

Todo:
Check if this description is accurate.
Remarks
This function assumes tensor-product-structure with matching mesh along the interface! If the gsMultiBasis contains bases of class gsHTensorBasis (or derived), it calls the function matchInterfaceHTensor().
Parameters
bispecifies the interface to be matched
mapperthe gsDofMapper which should know that these interface-DOFs are matched.
void refine ( size_t  k,
gsMatrix< T > const &  boxes,
int  refExt 
)
inline

Refine the are defined by boxes on patch k with extension refExt.

See gsHTensorBasis::refineWithExtension() for further documentation.

void refineElements ( int  k,
std::vector< index_t > const &  boxes 
)
inline

Refine the are defined by boxes on patch k.

See gsHTensorBasis::refineElements() for further documentation.

bool repairInterface ( const boundaryInterface bi)

Checks if the interface is fully matching, and if not, repairs it.

Remarks
Designed for gsHTensorBasis and derived bases. Assumes that the respective meshes on all levels of the gsHTensorBasis are fully matching.
Returns
true, if something was repaired, i.e., if the mesh on the interface was changed.
bool repairInterface2d ( const boundaryInterface bi)

Checks if the 2D-interface is fully matching, and if not, repairs it.

Same as repairInterface(), but only for 2D and a bit more efficient.

Returns
true, if something was repaired, i.e., if the mesh on the interface was changed.
bool repairInterfaceFindElements ( const boundaryInterface bi,
std::vector< index_t > &  refEltsFirst,
std::vector< index_t > &  refEltsSecond 
)

Finds the elements that need to be refined in order to repair an interface.

This function compares the hierarchical meshes on both patches associated with the boundaryInterface bi. The elements that need to be refined on bi.first() and bi.second() are stored in refEltsFirst and refEltsSecond, respectively.

Subsequent calls of the functions
m_bases[ bi.first().patch ]->refineElements( refEltsFirst )
m_bases[ bi.second().patch ]->refineElements( refEltsSecond )
will repair the interface in the sense that the resulting meshes are fully matching (this is done in repairInterface()).

Parameters
[in]bibondaryInterface to be checked.
[out]refEltsFirstContains elements (and levels) specifying needed refinement on patch bi.first().
[out]refEltsSecondContains elements (and levels) specifying needed refinement on patch bi.second().
Returns
True, if anything needs to be refined (i.e., if refEltsFirst.size() > 0 or refEltsSecond.size() > 0).
False if the patches are already fully matching at interface bi.

Is called by repairInterface(), templated over dimension.

void repairInterfaces ( const std::vector< boundaryInterface > &  bivec)
inline

Checks if the interfaces bivec are fully matching, and if not, repairs them, i.e., makes them fully matching.

Remarks
Designed for gsHTensorBasis and derived bases. Assumes that the meshes on all levels of the gsHTensorBasis are fully matching.

Calls repairInterface() for each boundaryInterface in bivec.

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 >.

short_t targetDim ( ) const
inlinevirtual

Dimension of the target space.

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

Reimplemented from gsFunctionSet< T >.

void tileParameters ( )

Tile the parameter domains of the pieces according to the topology

void uniformCoarsen ( int  numKnots = 1)
inline

Coarsen every basis uniformly.

This calls gsBasis::uniformCoarsen(numKnots) for all patches

void uniformCoarsen_withTransfer ( gsSparseMatrix< T, RowMajor > &  transfer,
const gsBoundaryConditions< T > &  boundaryConditions,
const gsOptionList assemblerOptions,
int  numKnots = 1,
index_t  unk = 0 
)

Coarsen every 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.

For computing the transfer matrix (but not for refinement), the boundaryConditions and the assemblerOptions have to be provided. By deault, the boundary conditions for unknown 0 are chosen. Use the parameter unk to choose another one.

See Also
gsMultiBasis::uniformCoarsen
void uniformRefine ( int  numKnots = 1,
int  mul = 1,
int  dir = -1 
)
inline

Refine every basis uniformly.

This calls gsBasis::uniformRefine(numKnots,mul) for all patches

void uniformRefine_withTransfer ( gsSparseMatrix< T, RowMajor > &  transfer,
const gsBoundaryConditions< T > &  boundaryConditions,
const gsOptionList assemblerOptions,
int  numKnots = 1,
int  mul = 1,
index_t  unk = 0 
)

Refine every 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.

For computing the transfer matrix (but not for refinement), the boundaryConditions and the assemblerOptions have to be provided. By deault, the boundary conditions for unknown 0 are chosen. Use the parameter unk to choose another one.

See Also
gsMultiBasis::uniformRefine