DOLFIN
DOLFIN C++ interface
|
This class defines quadrature rules for simplices. More...
#include <SimplexQuadrature.h>
Public Member Functions | |
SimplexQuadrature (std::size_t tdim, std::size_t order) | |
std::pair< std::vector< double >, std::vector< double > > | compute_quadrature_rule (const Cell &cell) const |
std::pair< std::vector< double >, std::vector< double > > | compute_quadrature_rule (const std::vector< Point > &coordinates, std::size_t gdim) const |
std::pair< std::vector< double >, std::vector< double > > | compute_quadrature_rule_interval (const std::vector< Point > &coordinates, std::size_t gdim) const |
std::pair< std::vector< double >, std::vector< double > > | compute_quadrature_rule_triangle (const std::vector< Point > &coordinates, std::size_t gdim) const |
std::pair< std::vector< double >, std::vector< double > > | compute_quadrature_rule_tetrahedron (const std::vector< Point > &coordinates, std::size_t gdim) const |
Static Public Member Functions | |
static std::vector< std::size_t > | compress (std::pair< std::vector< double >, std::vector< double >> &qr, std::size_t gdim, std::size_t quadrature_order) |
This class defines quadrature rules for simplices.
SimplexQuadrature::SimplexQuadrature | ( | std::size_t | tdim, |
std::size_t | order | ||
) |
Create SimplexQuadrature rules for reference simplex
Arguments tdim (std::size_t) The topological dimension of the simplex. order (std::size_t) The order of convergence of the quadrature rule.
|
static |
Compress a quadrature rule using algorithms from Compression of multivariate discrete measures and applications A. Sommariva, M. Vianello Numerical Functional Analysis and Optimization Volume 36, 2015 - Issue 9
Arguments qr (std::pair<std::vector<double>, std::vector<double>>) The quadrature rule to be compressed gdim (std::size_t) The geometric dimension quadrature_order (std::size_t) The order of the quadrature rule
Returns std::vector<std::size_t> The indices of the points that were kept (empty if no compression was made)
std::pair< std::vector< double >, std::vector< double > > SimplexQuadrature::compute_quadrature_rule | ( | const Cell & | cell | ) | const |
Compute quadrature rule for cell.
Arguments cell (Cell) The cell.
Returns std::pair<std::vector<double>, std::vector<double>> A flattened array of quadrature points and a corresponding array of quadrature weights.
std::pair< std::vector< double >, std::vector< double > > SimplexQuadrature::compute_quadrature_rule | ( | const std::vector< Point > & | coordinates, |
std::size_t | gdim | ||
) | const |
Compute quadrature rule for simplex.
Arguments coordinates (std::vector<Point>) Vertex coordinates for the simplex tdim (std::size_t) The topological dimension of the simplex. gdim (std::size_t) The geometric dimension.
Returns std::pair<std::vector<double>, std::vector<double>> A flattened array of quadrature points and a corresponding array of quadrature weights.
std::pair< std::vector< double >, std::vector< double > > SimplexQuadrature::compute_quadrature_rule_interval | ( | const std::vector< Point > & | coordinates, |
std::size_t | gdim | ||
) | const |
Compute quadrature rule for interval.
Arguments coordinates (std::vector<Point>) Vertex coordinates for the simplex gdim (std::size_t) The geometric dimension.
Returns std::pair<std::vector<double>, std::vector<double>> A flattened array of quadrature points and a corresponding array of quadrature weights.
std::pair< std::vector< double >, std::vector< double > > SimplexQuadrature::compute_quadrature_rule_tetrahedron | ( | const std::vector< Point > & | coordinates, |
std::size_t | gdim | ||
) | const |
Compute quadrature rule for tetrahedron.
Arguments coordinates (std::vector<Point>) Vertex coordinates for the simplex gdim (std::size_t) The geometric dimension.
Returns std::pair<std::vector<double>, std::vector<double>> A flattened array of quadrature points and a corresponding array of quadrature weights.
std::pair< std::vector< double >, std::vector< double > > SimplexQuadrature::compute_quadrature_rule_triangle | ( | const std::vector< Point > & | coordinates, |
std::size_t | gdim | ||
) | const |
Compute quadrature rule for triangle.
Arguments coordinates (std::vector<Point>) Vertex coordinates for the simplex gdim (std::size_t) The geometric dimension.
Returns std::pair<std::vector<double>, std::vector<double>> A flattened array of quadrature points and a corresponding array of quadrature weights.