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

Detailed Description

template<class T>
class gismo::gsBasisFun< T >

Represents an individual function in a function set, or a certain component of a vector-valued function.

It is constructed using a function set and the index of a single function in that basis. Note that it does not own the underlying function set (eg. basis). The lifetime of the parent object should be at least the lifetime of gsBasisFun.

+ Inheritance diagram for gsBasisFun< T >:
+ Collaboration diagram for gsBasisFun< T >:

Public Types

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

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
 Indices of active (non-zero) function(s) for each point. 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.
 
uPtr clone ()
 Clone methode. Produceds a deep copy inside a uPtr.
 
virtual void compute (const gsMatrix< T > &in, gsFuncData< T > &out) const
 Computes function data. More...
 
virtual void computeMap (gsMapData< T > &InOut) const
 Computes map function data. More...
 
gsFuncCoordinate< T > coord (const index_t c) const
 Returns the scalar function giving the i-th coordinate of this function.
 
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) functions at points u. More...
 
void deriv_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluate derivatives of the function \(f:\mathbb{R}^n\rightarrow\mathbb{R}^m\) at points u into result. More...
 
virtual T distanceL2 (gsFunction< T > const &) const
 Computes the L2-distance between this function and the field and a function func.
 
short_t domainDim () const
 Dimension of the (source) domain. More...
 
gsMatrix< T > eval (const gsMatrix< T > &u) const
 Evaluate the function,. More...
 
void eval_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluate the function at points u into result. More...
 
std::vector< gsMatrix< T > > evalAllDers (const gsMatrix< T > &u, int n, bool sameElement=false) const
 Evaluate all derivatives upto order n,. More...
 
virtual void evalAllDers_into (const gsMatrix< T > &u, int n, std::vector< gsMatrix< T > > &result, bool sameElement=false) const
 Evaluate the nonzero functions and their derivatives up to order n at points u into result. More...
 
void first ()
 Point to the first basis function of the basis.
 
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.
 
 gsBasisFun (const gsBasis< T > &basis, unsigned const i)
 Construct a basis function by a pointer to a basis and an index i.
 
virtual void invertPoints (const gsMatrix< T > &points, gsMatrix< T > &result, const T accuracy=1e-6, const bool useInitialPoint=false) const
 
int newtonRaphson (const gsVector< T > &value, gsVector< T > &arg, bool withSupport=true, const T accuracy=1e-6, int max_loop=100, T damping_factor=1) const
 
bool next ()
 
virtual index_t nPieces () const
 Number of pieces in the domain of definition.
 
virtual gsMatrix< T > parameterCenter () const
 Returns a "central" point inside inside the parameter domain.
 
gsMatrix< T > parameterCenter (const boxCorner &bc) const
 Get coordinates of the boxCorner bc in the parameter domain.
 
gsMatrix< T > parameterCenter (const boxSide &bs) const
 Get coordinates of the midpoint of the boxSide bs in the parameter domain.
 
virtual const gsBasisFunpiece (const index_t) const
 Returns the piece(s) of the function(s) at subdomain k.
 
virtual std::ostream & print (std::ostream &os) const
 Prints the object as a string.
 
void recoverPoints (gsMatrix< T > &xyz, gsMatrix< T > &uv, index_t k, const T accuracy=1e-6) const
 
void setFunction (unsigned const &i)
 
index_t size () const
 size More...
 
short_t targetDim () const
 Dimension of the target space. More...
 
bool valid ()
 Return false if the the function iteration is invalidated.
 
Evaluation functions

These functions allow one to evaluate the function as well as its derivatives at one or more points in the parameter space. See also Evaluation members.

virtual void eval_component_into (const gsMatrix< T > &u, const index_t comp, gsMatrix< T > &result) const
 Evaluate the function for component comp in the target dimension at points u into result.
 
virtual void jacobian_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Computes for each point u a block of result containing the Jacobian matrix.
 
void div_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Computes for each point u a block of result containing the divergence matrix.
 
gsMatrix< T > jacobian (const gsMatrix< T > &u) const
 
virtual void deriv2_into (const gsMatrix< T > &u, gsMatrix< T > &result) const
 Evaluate second derivatives of the function at points u into result. More...
 
virtual void hessian_into (const gsMatrix< T > &u, gsMatrix< T > &result, index_t coord=0) const
 
virtual gsMatrix< T > hessian (const gsMatrix< T > &u, index_t coord=0) const
 
virtual gsMatrix< T > laplacian (const gsMatrix< T > &u) const
 Evaluate the Laplacian at points u. More...
 

Private Member Functions

 gsBasisFun ()
 Default empty constructor.
 

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
inlinevirtualinherited

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 from gsFunctionSet< T >.

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 computeMap ( gsMapData< T > &  InOut) const
virtualinherited

Computes map function data.

This function evaluates the functions and their derivatives at the points InOut.points and writes them in the corresponding fields of InOut. Which field to write (and what to compute) is controlled by the InOut.flags (see also gsMapData). This is intended for parametrizations only and it works on functions sets of cardinality 1 only.

Parameters
[in,out]InOut
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) functions 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 second derivatives of the function at points u into result.

Let n be the dimension of the source space ( n = domainDim() ).
Let m be the dimension of the image/target space ( m = targetDim() ).
Let N denote the number of evaluation points.

Parameters
[in]ugsMatrix of size n x N, where each column of u represents one evaluation point.
[out]resultgsMatrix of size (S*m) x N, where S=n*(n+1)/2.
Each column in result corresponds to one point (i.e., one column in u)
and contains the following values (for n=3, m=3):
\( (\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)},\ldots,\partial_{yz} f^{(3)} )^T\)
Warning
By default uses central finite differences with h=0.00001! One must override this function in derived classes to get proper results.

Reimplemented from gsFunctionSet< T >.

Reimplemented in gsGeometry< T >, gsFunctionExpr< T >, gsConstantFunction< T >, gsPreCICEFunction< T >, gsAffineFunction< T >, gsMappedSingleSpline< d, T >, gsGeometryTransform< T >, and gsFuncCoordinate< T >.

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

Evaluate derivatives of the function \(f:\mathbb{R}^n\rightarrow\mathbb{R}^m\) at points u into result.

Let n be the dimension of the source space ( n = domainDim() ).
Let m be the dimension of the image/target space ( m = targetDim() ).
Let N denote the number of evaluation points.

Let \( f:\mathbb R^2 \rightarrow \mathbb R^3 \), i.e., \( f(x,y) = ( f^{(1)}(x,y), f^{(2)}(x,y), f^{(3)}(x,y) )^T\),
and let \( u = ( u_1, \ldots, u_N) = ( (x_1,y_1)^T, \ldots, (x_N, y_N)^T )\).
Then, result is of the form

\[ \left[ \begin{array}{cccc} \partial_x f^{(1)}(u_1) & \partial_x f^{(1)}(u_2) & \ldots & \partial_x f^{(1)}(u_N) \\ \partial_y f^{(1)}(u_1) & \partial_y f^{(1)}(u_2) & \ldots & \partial_y f^{(1)}(u_N) \\ \partial_x f^{(2)}(u_1) & \partial_x f^{(2)}(u_2) & \ldots & \partial_x f^{(2)}(u_N) \\ \partial_y f^{(2)}(u_1) & \partial_y f^{(2)}(u_2) & \ldots & \partial_x f^{(2)}(u_N) \\ \partial_x f^{(3)}(u_1) & \partial_x f^{(3)}(u_2) & \ldots & \partial_x f^{(3)}(u_N)\\ \partial_y f^{(3)}(u_1) & \partial_y f^{(3)}(u_2) & \ldots & \partial_y f^{(3)}(u_N) \end{array} \right] \]

Parameters
[in]ugsMatrix of size n x N, where each column of u represents one evaluation point.
[out]resultgsMatrix of size (n * m) x N. Each row of result corresponds to one component in the target space and contains the gradients for each evaluation point, as row vectors, one after the other (see above for details on the format).
Warning
By default, gsFunction uses central finite differences with h=0.00001! One must override this function in derived classes to get proper results.

Reimplemented from gsFunction< 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

Evaluate the function at points u into result.

Let n be the dimension of the source space ( n = domainDim() ).
Let m be the dimension of the image/target space ( m = targetDim() ).
Let N denote the number of evaluation points.

Parameters
[in]ugsMatrix of size n x N, where each column of u represents one evaluation point.
[out]resultgsMatrix of size m x N, where each column of u represents the result of the function at the respective valuation point.

Implements gsFunction< T >.

std::vector< gsMatrix< T > > evalAllDers ( const gsMatrix< T > &  u,
int  n,
bool  sameElement = false 
) 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,
bool  sameElement = false 
) const
virtualinherited

Evaluate the nonzero 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 in gsTensorBSplineBasis< 1, T >, gsTensorBasis< d, T >, gsTensorBasis< 1, T >, gsGeometry< T >, gsMappedSingleBasis< d, T >, gsConstantFunction< T >, gsTHBSplineBasis< d, T >, gsMappedSingleSpline< d, T >, and gsSquaredDistance< T >.

virtual gsMatrix<T> hessian ( const gsMatrix< T > &  u,
index_t  coord = 0 
) const
inlinevirtualinherited

Evaluates the Hessian (matrix of second partial derivatives) of coordinate coord at points u.

void invertPoints ( const gsMatrix< T > &  points,
gsMatrix< T > &  result,
const T  accuracy = 1e-6,
const bool  useInitialPoint = false 
) const
virtualinherited

Takes the physical points and computes the corresponding parameter values. If the point cannot be inverted (eg. is not part of the geometry) the corresponding parameter values will be undefined

gsMatrix< T > laplacian ( const gsMatrix< T > &  u) const
virtualinherited

Evaluate the Laplacian at points u.

By default uses central finite differences with h=0.00001

Reimplemented in gsFunctionExpr< T >.

int newtonRaphson ( const gsVector< T > &  value,
gsVector< T > &  arg,
bool  withSupport = true,
const T  accuracy = 1e-6,
int  max_loop = 100,
damping_factor = 1 
) const
inherited

Newton-Raphson method to find a solution of the equation f(arg) = value with starting vector arg. If the point cannot be inverted the corresponding parameter values will be undefined

bool next ( )

Points to the first basis function of the basis or return false if this is the last basis function

void recoverPoints ( gsMatrix< T > &  xyz,
gsMatrix< T > &  uv,
index_t  k,
const T  accuracy = 1e-6 
) const
inherited

Recovers a point on the (geometry) together with its parameters uv, assuming that the k-th coordinate of the point xyz is not known (and has a random value as input argument).

void setFunction ( unsigned const &  i)

The gsBasisFun points to the i-th basis function of m_basis after calling this setter.

index_t size ( ) const
inlinevirtualinherited

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

Reimplemented in gsPiecewiseFunction< 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 >.