G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsQuadrature Struct Reference

Detailed Description

Helper class for obtaining a quadrature rule.

Public Types

enum  rule {
  GaussLegendre ,
  GaussLobatto ,
  PatchRule
}
 Quadrature rule types. More...
 

Static Public Member Functions

template<class T >
static gsQuadRule< T > get (const gsBasis< T > &basis, const gsOptionList &options, short_t fixDir=-1)
 Constructs a quadrature rule based on input options.
 
template<class T >
static gsQuadRule< T > get (index_t qu, gsVector< index_t > const &numNodes, unsigned digits=0)
 Constructs a quadrature rule based on input options.
 
template<class T >
static gsMatrix< T > getAllNodes (const gsBasis< T > &basis, const gsGeometry< T > &geom, const gsOptionList &options, const patchSide side)
 Get all quadrature nodes for a specified side of a basis and evaluates them using a geometry.
 
template<class T >
static gsMatrix< T > getAllNodes (const gsBasis< T > &basis, const gsGeometry< T > &geom, const gsOptionList &options, const std::vector< patchSide > sides)
 Collects and evaluates all quadrature nodes for multiple sides of a given basis.
 
template<class T >
static gsMatrix< T > getAllNodes (const gsBasis< T > &basis, const gsOptionList &options)
 Retrieves all quadrature nodes for the given basis.
 
template<class T >
static gsMatrix< T > getAllNodes (const gsBasis< T > &basis, const gsOptionList &options, const patchSide side)
 Get all quadrature nodes for a specified side of a given basis.
 
template<class T >
static gsMatrix< T > getAllNodes (const gsBasis< T > &basis, const gsOptionList &options, const std::vector< patchSide > sides)
 Retrieves all quadrature nodes for multiple sides of a given basis.
 
template<class T >
static gsMatrix< T > getAllNodes (const gsMultiBasis< T > &bases, const gsMultiPatch< T > &mp, const gsOptionList &options, const std::vector< patchSide > sides)
 Gets all quadrature nodes for several sides of a multi-basis for multi-patch geometry.
 
template<class T >
static std::vector< gsMatrix< T > > getAllNodes (const gsMultiBasis< T > &bases, const gsOptionList &options, const std::vector< patchSide > sides)
 Collects all quadrature nodes for a multi-basis.
 
template<class T >
static gsQuadRule< T >::uPtr getPtr (const gsBasis< T > &basis, const gsOptionList &options, short_t fixDir=-1)
 Constructs a quadrature rule based on input options.
 
template<class T >
static gsQuadRule< T > getUnivariate (index_t qu, index_t numNodes, unsigned digits=0)
 Constructs a quadrature rule based on input options.
 
template<class T >
static gsVector< index_tnumNodes (const gsBasis< T > &basis, const Real quA, const index_t quB, short_t fixDir=-1)
 

Member Enumeration Documentation

◆ rule

enum rule

Quadrature rule types.

Enumerator
GaussLegendre 

Gauss-Legendre quadrature.

GaussLobatto 

Gauss-Lobatto quadrature.

PatchRule 

Patch-wise quadrature rule (Johannessen 2017)

Member Function Documentation

◆ getAllNodes() [1/3]

template<class T >
static gsMatrix< T > getAllNodes ( const gsBasis< T > &  basis,
const gsOptionList options 
)
inlinestatic

Retrieves all quadrature nodes for the given basis.

This function computes and returns all quadrature nodes for a given basis, using the provided options to determine the quadrature rules.

Template Parameters
TReal type.
Parameters
[in]basisThe basis for which quadrature nodes are computed.
[in]optionsOptions specifying the quadrature rule.
Returns
A matrix where each column represents a quadrature node in the parametric domain.

◆ getAllNodes() [2/3]

template<class T >
static gsMatrix< T > getAllNodes ( const gsBasis< T > &  basis,
const gsOptionList options,
const patchSide  side 
)
inlinestatic

Get all quadrature nodes for a specified side of a given basis.

Template Parameters
TReal type.
Parameters
[in]basisThe basis for which the quadrature nodes are to be collected.
[in]optionsQuadrature rule.
[in]sideThe side of the basis.
Returns
result A matrix of quadrature nodes, where each column corresponds to a quadrature node.

◆ getAllNodes() [3/3]

template<class T >
static gsMatrix< T > getAllNodes ( const gsMultiBasis< T > &  bases,
const gsMultiPatch< T > &  mp,
const gsOptionList options,
const std::vector< patchSide sides 
)
inlinestatic

Gets all quadrature nodes for several sides of a multi-basis for multi-patch geometry.

Template Parameters
TData type for computations.
Parameters
[in]basesThe multi-basis for which the quadrature nodes are to be collected.
[in]optionsOptions specifying the quadrature rule and other settings.
[in]sidesA vector of sides corresponding to the patches in the multi-basis.
[in]mpA multi-patch geometry.
Returns
A vector of matrices, where each matrix contains quadrature nodes for a specific patch.

◆ numNodes()

template<class T >
static gsVector< index_t > numNodes ( const gsBasis< T > &  basis,
const Real  quA,
const index_t  quB,
short_t  fixDir = -1 
)
inlinestatic

Computes and integer quA*deg_i + quB where deg_i is the degree of basis