chrono::ChQuadrature Class Reference

Description

Gauss-Legendre quadrature in 1D, 2D, or 3D.

This class provides methods to integrate a function (the integrand) on a domain using the Gauss quadrature. This method is most useful for polynomial integrands, in which case the result is exact if the order of quadrature is greater or equal to the degree of the polynomial.

#include <ChQuadrature.h>

Static Public Member Functions

template<class T >
static void Integrate1D (T &result, ChIntegrand1D< T > &integrand, const double x_min, const double x_max, const int order)
 Integrate the integrand T = f(x) over the 1D interval [xA, xB], with specified order of quadrature. More...
 
template<class T >
static void Integrate2D (T &result, ChIntegrand2D< T > &integrand, const double x_min, const double x_max, const double y_min, const double y_max, 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, ChIntegrand3D< T > &integrand, const double x_min, const double x_max, const double y_min, const double y_max, const double z_min, const double z_max, 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, ChIntegrand2D< T > &integrand, const int order)
 Integrate the 2D integrand T = f(u,v) over a triangle, with desired order of quadrature. More...
 
template<class T >
static void Integrate3Dtetrahedron (T &result, ChIntegrand3D< T > &integrand, const int order)
 Integrate the 3D 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 10th order, with precomputed tables.
 
static ChQuadratureTablesTriangleGetStaticTablesTriangle ()
 Access a statically-allocated set of tables for tetrahedron quadrature, with 5 precomputed tables.
 
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.

◆ Integrate1D()

template<class T >
static void chrono::ChQuadrature::Integrate1D ( T &  result,
ChIntegrand1D< T > &  integrand,
const double  x_min,
const double  x_max,
const int  order 
)
inlinestatic

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

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

Parameters
resultresult is returned here
integrandthis is the integrand
x_minmin limit for x domain
x_maxmin limit for x domain
orderorder of integration

◆ Integrate2D()

template<class T >
static void chrono::ChQuadrature::Integrate2D ( T &  result,
ChIntegrand2D< T > &  integrand,
const double  x_min,
const double  x_max,
const double  y_min,
const double  y_max,
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.

Parameters
resultresult is returned here
integrandthis is the integrand
x_minmin limit for x domain
x_maxmin limit for x domain
y_minmin limit for y domain
y_maxmin limit for y domain
orderorder of integration

◆ Integrate2Dtriangle()

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

Integrate the 2D 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 and v in [0...1], with the 3rd assumed to be 1-u-v. For order between 1 and 5, use precomputed polynomial coefficients.

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

◆ Integrate3D()

template<class T >
static void chrono::ChQuadrature::Integrate3D ( T &  result,
ChIntegrand3D< T > &  integrand,
const double  x_min,
const double  x_max,
const double  y_min,
const double  y_max,
const double  z_min,
const double  z_max,
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.

Parameters
resultresult is returned here
integrandthis is the integrand
x_minmin limit for x domain
x_maxmin limit for x domain
y_minmin limit for y domain
y_maxmin limit for y domain
z_minmin limit for z domain
z_maxmin limit for z domain
orderorder of integration

◆ Integrate3Dtetrahedron()

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

Integrate the 3D 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, and w in [0...1], with the 4th assumed to be1-u-v-w. For order between 1 and 5 use precomputed polynomial coefficients.

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