chrono::ChQuadrature Class Reference

Description

Class to perform Gauss-Legendre quadrature, in 1D, 2D, 3D.

It integrates a function on a nD domain using the Gauss quadrature; this is mostly useful in the case that the integrand is polynomial, because the result is exact if the order of quadrature is greater or equal to the degree of the polynomial. Integrate() functions are static; no need to allocate an instance of this class.

#include <ChQuadrature.h>

Static Public Member Functions

template<class T >
static void Integrate1D (T &result, ChIntegrable1D< T > &integrand, const double a, const double b, const int order)
 Integrate the integrand T = f(x) over the 1D interval [xA, xB], with desired order of quadrature. More...
 
template<class T >
static void Integrate2D (T &result, ChIntegrable2D< T > &integrand, const double Xa, const double Xb, const double Ya, const double Yb, const int order)
 Integrate the integrand T = f(x,y) over the 2D interval [xA, xB][yA, yB], with desired order of quadrature. More...
 
template<class T >
static void Integrate3D (T &result, ChIntegrable3D< T > &integrand, const double Xa, const double Xb, const double Ya, const double Yb, const double Za, const double Zb, const int order)
 Integrate the integrand T = f(x,y,z) over the 3D interval [xA, xB][yA, yB][zA, zB], with desired order of quadrature. More...
 
template<class T >
static void Integrate2Dtriangle (T &result, ChIntegrable2D< T > &integrand, const int order)
 Special case of 2D integration: integrate the integrand T = f(u,v) over a triangle, with desired order of quadrature. More...
 
template<class T >
static void Integrate3Dtetrahedron (T &result, ChIntegrable3D< T > &integrand, const int order)
 Special case of 3D integration: integrate the integrand T = f(u,v,w) over a tetrahedron, with desired order of quadrature. More...
 
static ChQuadratureTablesGetStaticTables ()
 Access a statically-allocated set of tables, from 0 to a 10th order, with precomputed tables.
 
static ChQuadratureTablesTriangleGetStaticTablesTriangle ()
 Access a statically-allocated set of tables for tetrahedron quadrature, with 5 precomputed tables. More...
 
static ChQuadratureTablesTetrahedronGetStaticTablesTetrahedron ()
 Access a statically-allocated set of tables for tetrahedron quadrature, with 5 precomputed tables. More...
 

Member Function Documentation

◆ GetStaticTablesTetrahedron()

ChQuadratureTablesTetrahedron * chrono::ChQuadrature::GetStaticTablesTetrahedron ( )
static

Access a statically-allocated set of tables for tetrahedron quadrature, with 5 precomputed tables.

Use Dunavant theory.

◆ GetStaticTablesTriangle()

ChQuadratureTablesTriangle * chrono::ChQuadrature::GetStaticTablesTriangle ( )
static

Access a statically-allocated set of tables for tetrahedron quadrature, with 5 precomputed tables.


◆ Integrate1D()

template<class T >
static void chrono::ChQuadrature::Integrate1D ( T &  result,
ChIntegrable1D< T > &  integrand,
const double  a,
const double  b,
const int  order 
)
inlinestatic

Integrate the integrand T = f(x) over the 1D interval [xA, xB], with desired order of quadrature.

Best if integrand is polynomial. For order in 1..10, precomputed polynomial coefficients are used for max speed.

Parameters
resultresult is returned here
integrandthis is the integrand
amin limit for x domain
bmin limit for x domain
orderorder of integration

◆ Integrate2D()

template<class T >
static void chrono::ChQuadrature::Integrate2D ( T &  result,
ChIntegrable2D< T > &  integrand,
const double  Xa,
const double  Xb,
const double  Ya,
const double  Yb,
const int  order 
)
inlinestatic

Integrate the integrand T = f(x,y) over the 2D interval [xA, xB][yA, yB], with desired order of quadrature.

Best if integrand is polynomial. For order in 1..10, precomputed polynomial coefficients are used for max speed.

Parameters
resultresult is returned here
integrandthis is the integrand
Xamin limit for x domain
Xbmin limit for x domain
Yamin limit for y domain
Ybmin limit for y domain
orderorder of integration

◆ Integrate2Dtriangle()

template<class T >
static void chrono::ChQuadrature::Integrate2Dtriangle ( T &  result,
ChIntegrable2D< T > &  integrand,
const int  order 
)
inlinestatic

Special case of 2D integration: integrate the integrand T = f(u,v) over a triangle, with desired order of quadrature.

Best if integrand is polynomial. Two triangle coordinates are assumed to be 'area' coordinates u,v in [0...1]. The third is assumed 1-u-v.
For order in 1..5. Use precomputed polynomial coefficients for max speed.

Parameters
resultresult is returned here
integrandthis is the integrand
orderorder of integration

◆ Integrate3D()

template<class T >
static void chrono::ChQuadrature::Integrate3D ( T &  result,
ChIntegrable3D< T > &  integrand,
const double  Xa,
const double  Xb,
const double  Ya,
const double  Yb,
const double  Za,
const double  Zb,
const int  order 
)
inlinestatic

Integrate the integrand T = f(x,y,z) over the 3D interval [xA, xB][yA, yB][zA, zB], with desired order of quadrature.

Best if integrand is polynomial. For order in 1..10, precomputed polynomial coefficients are used for max speed.

Parameters
resultresult is returned here
integrandthis is the integrand
Xamin limit for x domain
Xbmin limit for x domain
Yamin limit for y domain
Ybmin limit for y domain
Zamin limit for z domain
Zbmin limit for z domain
orderorder of integration

◆ Integrate3Dtetrahedron()

template<class T >
static void chrono::ChQuadrature::Integrate3Dtetrahedron ( T &  result,
ChIntegrable3D< T > &  integrand,
const int  order 
)
inlinestatic

Special case of 3D integration: integrate the integrand T = f(u,v,w) over a tetrahedron, with desired order of quadrature.

Best if integrand is polynomial. Three tetrahedron coordinates are assumed to be 'volume' coordinates u,v,w in [0...1]. The 4th is assumed 1-u-v-w.
For order in 1..5. Use precomputed polynomial coefficients for max speed.

Parameters
resultresult is returned here
integrandthis is the integrand
orderorder of integration

The documentation for this class was generated from the following files:
  • /builds/uwsbel/chrono/src/chrono/core/ChQuadrature.h
  • /builds/uwsbel/chrono/src/chrono/core/ChQuadrature.cpp