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

Detailed Description

template<class T>
class gismo::gsOverIntegrateRule< T >

Class that defines a mixed quadrature rule with different rules for the interior and the boundaries.

This class is defined using two other quadrature rules (only works for gsGaussRule or gsLobattoRule) for the interior and for the boundary. Depending on the location of the considered element, it will either use the interior rule or the boundary rule. In this way, one can for example use full (exact) integration on the boundary points and reduced integration in the interior.

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

Public Member Functions

index_t dim () const
 Dimension of the rule.
 
 gsOverIntegrateRule ()
 Default empty constructor.
 
 gsOverIntegrateRule (const gsBasis< T > &basis, const std::vector< gsQuadRule< T > > &quadInterior, const std::vector< gsQuadRule< T > > &quadBoundary)
 Constructor. More...
 
void mapTo (const gsVector< T > &lower, const gsVector< T > &upper, gsMatrix< T > &nodes, gsVector< T > &weights) const
 Maps the points in the d-dimensional cube with points lower and upper. More...
 
virtual void mapTo (T startVal, T endVal, gsMatrix< T > &nodes, gsVector< T > &weights) const
 Maps a univariate quadrature rule (i.e., points and weights) from the reference interval to an arbitrary interval.
 
void mapToAll (const std::vector< T > &breaks, gsMatrix< T > &nodes, gsVector< T > &weights) const
 Maps a univariate quadrature rule (i.e., points and weights) from the reference interval to a number of intervals, implied by the breaks.
 
index_t numNodes () const
 Number of nodes in the currently kept rule.
 
const gsMatrix< T > & referenceNodes () const
 Returns reference nodes for the currently kept rule. More...
 
const gsVector< T > & referenceWeights () const
 Returns reference weights for the currently kept rule. More...
 
virtual void setNodes (gsVector< index_t > const &, unsigned=0)
 Initialize quadrature rule with numNodes number of quadrature points per integration variable. More...
 
void setNodes (index_t numNodes, unsigned digits=0)
 Initialize a univariate quadrature rule with numNodes quadrature points. More...
 

Static Public Member Functions

static uPtr make (const gsBasis< T > &basis, const std::vector< gsQuadRule< T > > &quadInterior, const std::vector< gsQuadRule< T > > &quadBoundary)
 Construct a smart-pointer to the rule. More...
 

Protected Member Functions

void computeTensorProductRule (const std::vector< gsVector< T > > &nodes, const std::vector< gsVector< T > > &weights)
 Computes the tensor product rule from coordinate-wise 1D nodes and weights.
 

Protected Attributes

gsMatrix< T > m_nodes
 Reference quadrature nodes (on the interval [-1,1]).
 
gsVector< T > m_weights
 Reference quadrature weights (corresponding to interval [-1,1]).
 

Constructor & Destructor Documentation

gsOverIntegrateRule ( const gsBasis< T > &  basis,
const std::vector< gsQuadRule< T > > &  quadInterior,
const std::vector< gsQuadRule< T > > &  quadBoundary 
)
inline

Constructor.

Parameters
[in]basisThe basis
[in]quadInteriorvector with the univariate INTERIOR quadrature rules for each dimension of the basis.
[in]quadBoundaryvector with the univariate BOUNDARY quadrature rules for each dimension of the basis.
Note
: Only works for QuadRules that do not re-implement mapTo (slicing occurs)

Member Function Documentation

static uPtr make ( const gsBasis< T > &  basis,
const std::vector< gsQuadRule< T > > &  quadInterior,
const std::vector< gsQuadRule< T > > &  quadBoundary 
)
inlinestatic

Construct a smart-pointer to the rule.

Parameters
[in]basisThe basis
[in]quadInteriorThe rule used for the interior
[in]quadBoundaryThe rule used for the boundary
Note
: Only works for QuadRules that do not re-implement mapTo (slicing occurs)
void mapTo ( const gsVector< T > &  lower,
const gsVector< T > &  upper,
gsMatrix< T > &  nodes,
gsVector< T > &  weights 
) const
inlinevirtual

Maps the points in the d-dimensional cube with points lower and upper.

Parameters
[in]lowerThe lower corner
[in]upperThe upper corner
nodesQuadrature points
weightsQuadrature weights

Reimplemented from gsQuadRule< T >.

const gsMatrix<T>& referenceNodes ( ) const
inlineinherited

Returns reference nodes for the currently kept rule.

Returns
v Vector of length d = dim(), where each entry v_i of v is again a vector. v_i contains the reference quadrature points for the i-th coordinate direction.
const gsVector<T>& referenceWeights ( ) const
inlineinherited

Returns reference weights for the currently kept rule.

Returns
v Vector of length d = dim(), where each entry v_i of v is again a vector. v_i contains the reference quadrature weights for the i-th coordinate direction.
virtual void setNodes ( gsVector< index_t > const &  ,
unsigned  = 0 
)
inlinevirtualinherited

Initialize quadrature rule with numNodes number of quadrature points per integration variable.

The call of this function initializes the quadrature rule for further use, i.e., the quadrature points and weights on the reference cube are set up and stored. The dimension d of the reference cube is specified by the length of the vector numNodes.

Example: numNodes = [2,5] corresponds to quadrature in 2D where two quadrature points are used for the first coordinate and five quadrature points for the second coordinate. In this case, 2*5 = 10 quadrature nodes and weights are set up.


Parameters
numNodes vector of length digits containing numbers of nodes to be used (per integration variable).
digits accuracy of nodes and weights. If 0, the quadrature rule will use precomputed tables if possible. If greater then 0, it will force to compute it everytime.

Reimplemented in gsGaussRule< T >, gsNewtonCotesRule< T >, and gsLobattoRule< T >.

void setNodes ( index_t  numNodes,
unsigned  digits = 0 
)
inlineinherited

Initialize a univariate quadrature rule with numNodes quadrature points.

Parameters
numNodesquadrature point
digitsaccuracy of nodes and weights. If 0, the quadrature rule will use precomputed tables if possible. If greater then 0, it will force to compute it everytime.