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

Detailed Description

template<class T>
class gismo::gsExprEvaluator< T >

Generic evaluator of isogeometric expressions.

The expressions may be scalar ot vector-valued. Computed quatities can be global or element.wise.

Note the following will not return correct result (return value is a reference): gsInfo << "Eval:"<< ev.eval(a) + ev.eval(b) <<"\n"; Use the evaluator once at a time.

+ Collaboration diagram for gsExprEvaluator< T >:

Public Member Functions

gsAsConstVector< T > allValues () const
 Returns a vector containing the last computed values per element.
 
gsAsConstMatrix< T > allValues (index_t nR, index_t nC) const
 Returns the last computed values per element, resized as a matrix.
 
void calcRoot (const index_t p)
 Calculates the p-th root of the lastly computed quantities (eg. integrals)
 
void calcSqrt ()
 Calculates the square root of the lastly computed quantities (eg. integrals)
 
const std::vector< T > & elementwise () const
 Returns an std::vector containing the last computed values per element.
 
template<class E >
void eval (const expr::_expr< E > &expr, gsGridIterator< T, mode, d > &git, const index_t patchInd=0)
 
template<class E >
gsAsConstMatrix< T > eval (const expr::_expr< E > &testExpr, const gsVector< T > &pt, const index_t patchInd=0)
 
element getElement ()
 Returns a handle to an isogeometric element.
 
geometryMap getMap (const gsMultiPatch< T > &mp)
 Registers mp as an isogeometric geometry map and return a handle to it.
 
geometryMap getMap (const gsFunctionSet< T > &gm)
 Registers g as an isogeometric geometry map and return a handle to it.
 
variable getVariable (const gsFunctionSet< T > &func, index_t dim=1)
 Registers func as a variable and returns a handle to it.
 
expr::gsComposition< T > getVariable (const gsFunctionSet< T > &func, geometryMap G)
 Registers func as a variable defined on G and returns a handle to it.
 
template<class E >
integral (const expr::_expr< E > &expr)
 Calculates the integral of the expression expr on the whole integration domain.
 
template<class E >
integralBdr (const expr::_expr< E > &expr)
 
template<class E >
integralBdr (const expr::_expr< E > &expr, const bContainer &bdrlist)
 
template<class E >
integralBdrBc (const bcRefList &BCs, const expr::_expr< E > &expr)
 
template<class E >
integralElWise (const expr::_expr< E > &expr)
 Calculates the integral of the expression expr on each element.
 
template<class E >
integralInterface (const expr::_expr< E > &expr)
 
template<class E >
integralInterface (const expr::_expr< E > &expr, const intContainer &iFaces)
 
template<class E >
max (const expr::_expr< E > &expr)
 
template<class E >
maxElWise (const expr::_expr< E > &expr)
 
template<class E >
maxInterface (const expr::_expr< E > &expr)
 
template<class E >
maxInterface (const expr::_expr< E > &expr, const intContainer &iFaces)
 
template<class E >
min (const expr::_expr< E > &expr)
 
template<class E >
minElWise (const expr::_expr< E > &expr)
 
template<class E >
minInterface (const expr::_expr< E > &expr)
 
template<class E >
minInterface (const expr::_expr< E > &expr, const intContainer &iFaces)
 
size_t nValues () const
 The number of lastly computed values.
 
void setIntegrationElements (const gsMultiBasis< T > &mesh)
 Sets the domain of integration. More...
 
template<class E >
void testEval (const expr::_expr< E > &expr, const gsVector< T > &pt, const index_t patchInd=0)
 
value () const
 Returns the last computed value.
 
template<class E >
void writeParaview (const expr::_expr< E > &expr, geometryMap G, std::string const &fn)
 Creates a paraview file named fn containing valies of the. More...
 
template<class E >
void writeParaview (const expr::_expr< E > &expr, std::string const &fn)
 Creates a paraview file named fn containing valies of the. More...
 

Private Member Functions

template<class E , class _op >
computeBdr_impl (const expr::_expr< E > &expr, const bContainer &bdrlist)
 

Member Function Documentation

T computeBdr_impl ( const expr::_expr< E > &  expr,
const bContainer &  bdrlist 
)
private

! not multipatch!

void eval ( const expr::_expr< E > &  expr,
gsGridIterator< T, mode, d > &  git,
const index_t  patchInd = 0 
)

Computes values of the expression expr at the grid points git of patch patchId

gsAsConstMatrix<T> eval ( const expr::_expr< E > &  testExpr,
const gsVector< T > &  pt,
const index_t  patchInd = 0 
)

Computes value of the expression testExpr at the point pt of patch patchId

T integralBdr ( const expr::_expr< E > &  expr)
inline

Calculates the integral of the expression expr on the boundary of the integration domain

T integralBdr ( const expr::_expr< E > &  expr,
const bContainer &  bdrlist 
)
inline

Calculates the integral of the expression expr on the boundaries contained in bdrlist

T integralBdrBc ( const bcRefList &  BCs,
const expr::_expr< E > &  expr 
)
inline

Calculates the integral of the expression expr on the boundaries contained in BCs taking into account boundary functions

T integralInterface ( const expr::_expr< E > &  expr)
inline

Calculates the integral of the expression expr on the interfaces of the (multi-basis) integration domain

T integralInterface ( const expr::_expr< E > &  expr,
const intContainer &  iFaces 
)
inline

Calculates the integral of the expression expr on the interfaces iFaces of the integration domain

T max ( const expr::_expr< E > &  expr)
inline

Calculates the maximum value of the expression expr by sampling over a finite number of points

T maxElWise ( const expr::_expr< E > &  expr)
inline

Calculates the maximum value of the expression expr by on each element by sampling over a finite number of points

T maxInterface ( const expr::_expr< E > &  expr)
inline

Calculates the maximum of the expression expr on the interfaces of the (multi-basis) integration domain

T maxInterface ( const expr::_expr< E > &  expr,
const intContainer &  iFaces 
)
inline

Calculates the maximum of the expression expr on the interfaces of the (multi-basis) integration domain

T min ( const expr::_expr< E > &  expr)
inline

Calculates the minimum value of the expression expr by sampling over a finite number of points

T minElWise ( const expr::_expr< E > &  expr)
inline

Calculates the minimum value of the expression expr on each element by sampling over a finite number of points

T minInterface ( const expr::_expr< E > &  expr)
inline

Calculates the minimum of the expression expr on the interfaces of the (multi-basis) integration domain

T minInterface ( const expr::_expr< E > &  expr,
const intContainer &  iFaces 
)
inline

Calculates the minimum of the expression expr on the interfaces of the (multi-basis) integration domain

void setIntegrationElements ( const gsMultiBasis< T > &  mesh)
inline

Sets the domain of integration.

Warning
Must be called before any computation is requested
void testEval ( const expr::_expr< E > &  expr,
const gsVector< T > &  pt,
const index_t  patchInd = 0 
)
inline

Computes value of the expression expr at the point pt of patch patchId, and displays the result

void writeParaview ( const expr::_expr< E > &  expr,
geometryMap  G,
std::string const &  fn 
)
inline

Creates a paraview file named fn containing valies of the.

Plotting properties are controlled by entries in the options

void writeParaview ( const expr::_expr< E > &  expr,
std::string const &  fn 
)
inline

Creates a paraview file named fn containing valies of the.

Plotting properties are controlled by entires in the options