#include <QuadratureRuleContainer.h>
Public Member Functions | |
QuadratureRuleContainer (const QuadratureRuleAttributes &quadratureRuleAttributes, std::shared_ptr< const QuadratureRule > quadratureRule, std::shared_ptr< const basis::TriangulationBase > triangulation, const basis::CellMappingBase &cellMapping) | |
Constructor for assigning each cell with a common quadrature rule. More... | |
QuadratureRuleContainer (const QuadratureRuleAttributes &quadratureRuleAttributes, std::vector< std::shared_ptr< const QuadratureRule > > quadratureRuleVec, std::shared_ptr< const basis::TriangulationBase > triangulation, const basis::CellMappingBase &cellMapping) | |
Constructor for assigning each cell with different quadrature rule for each cell. More... | |
QuadratureRuleContainer (const QuadratureRuleAttributes &quadratureRuleAttributes, std::shared_ptr< const QuadratureRule > baseQuadratureRule, std::shared_ptr< const basis::TriangulationBase > triangulation, const basis::CellMappingBase &cellMapping, basis::ParentToChildCellsManagerBase &parentToChildCellsManager, std::vector< std::shared_ptr< const utils::ScalarSpatialFunctionReal > > functions, const std::vector< double > &absoluteTolerances, const std::vector< double > &relativeTolerances, const std::vector< double > &integralThresholds, const double smallestCellVolume=QuadratureRuleAdaptiveDefaults::SMALLEST_CELL_VOLUME, const unsigned int maxRecursion=QuadratureRuleAdaptiveDefaults::MAX_RECURSION) | |
Constructor for creating an adaptive quadrature rule in each cell based user-defined functions. More... | |
QuadratureRuleContainer (const QuadratureRuleAttributes &quadratureRuleAttributes, const size_type order1DMin, const size_type order1DMax, const size_type copies1DMax, std::shared_ptr< const basis::TriangulationBase > triangulation, const basis::CellMappingBase &cellMapping, std::vector< std::shared_ptr< const utils::ScalarSpatialFunctionReal > > functions, const std::vector< double > &absoluteTolerances, const std::vector< double > &relativeTolerances, const quadrature::QuadratureRuleContainer &quadratureRuleContainerReference, const utils::mpi::MPIComm &comm) | |
Constructor for creating a subdivided quadrature rule in each cell based user-defined functions and a reference quadrature. So one provides a minimum and maximum order (or number of 1D gauss points) and maimum number of copies of the orders one wants to go. Also one has a reference quadrature which can. More... | |
const QuadratureRuleAttributes & | getQuadratureRuleAttributes () const |
Returns the underlying QuadratureRuleAttributes. More... | |
size_type | nCells () const |
Returns the number of cells in the quadrature container. More... | |
const std::vector< dftefe::utils::Point > & | getRealPoints () const |
Function that returns a vector containing the real coordinates of the quad points in all cells. More... | |
std::vector< dftefe::utils::Point > | getCellRealPoints (const unsigned int cellId) const |
Function that returns a vector containing the real coordinates of the quad points in cell corresponding to the cellId. More... | |
const std::vector< dftefe::utils::Point > & | getCellParametricPoints (const unsigned int cellId) const |
Function that returns a vector containing the real coordinates of the quad points in cell corresponding to the cellId. More... | |
const std::vector< double > & | getCellQuadratureWeights (const unsigned int cellId) const |
Function that returns a vector containing the weight of the quad points in cell corresponding to the cellId. More... | |
const std::vector< double > & | getJxW () const |
Function that returns a vector containing the Jacobian times quadrature weight for all the quad points across all the cells in the triangulation. More... | |
std::vector< double > | getCellJxW (const unsigned int cellId) const |
Function that returns a vector containing the Jacobian times weight of the quad points in cell corresponding to the cellId. More... | |
const QuadratureRule & | getQuadratureRule (const unsigned int cellId) const |
Function that returns the handle to quadrature rule corresponding to the the cell Id. More... | |
size_type | nQuadraturePoints () const |
A function to return the total number of quadrature points in all the cells. More... | |
size_type | nCellQuadraturePoints (const unsigned int cellId) const |
A function to returns the number of quadrature points in cell corresponding to the the cell Id. More... | |
const std::vector< size_type > & | getCellQuadStartIds () const |
A function to return the starting index of the quadrature point of each cell. More... | |
size_type | getCellQuadStartId (const size_type cellId) const |
A function to return the starting index of the quadrature point of a given cell. More... | |
std::shared_ptr< const basis::TriangulationBase > | getTriangulation () const |
const basis::CellMappingBase & | getCellMapping () const |
Private Attributes | |
const QuadratureRuleAttributes & | d_quadratureRuleAttributes |
std::vector< std::shared_ptr< const QuadratureRule > > | d_quadratureRuleVec |
std::vector< size_type > | d_numCellQuadPoints |
std::vector< size_type > | d_cellQuadStartIds |
std::vector< dftefe::utils::Point > | d_realPoints |
std::vector< double > | d_JxW |
unsigned int | d_dim |
size_type | d_numQuadPoints |
size_type | d_numCells |
bool | d_storeJacobianInverse |
std::shared_ptr< const basis::TriangulationBase > | d_triangulation |
const basis::CellMappingBase & | d_cellMapping |
This class stores the quadrature points and corresponding JxW in each cell. This supports adaptive quadrature i.e each cell can have different quadrature rules. Further each cell can have arbitrary quadrature rule.
dftefe::quadrature::QuadratureRuleContainer::QuadratureRuleContainer | ( | const QuadratureRuleAttributes & | quadratureRuleAttributes, |
std::shared_ptr< const QuadratureRule > | quadratureRule, | ||
std::shared_ptr< const basis::TriangulationBase > | triangulation, | ||
const basis::CellMappingBase & | cellMapping | ||
) |
Constructor for assigning each cell with a common quadrature rule.
[in] | quadratureRule | The quadrature rule specifying the quad points in the parametric coordinates with its corresponding weights. |
[in] | triangulation | The triangulation that has information on the cell and its vertices |
[in] | cellMapping | cellMapping provides the the information on how the cell in real space is mapped to its parametric coordinates. This is required to calculate the JxW values at each quad point |
dftefe::quadrature::QuadratureRuleContainer::QuadratureRuleContainer | ( | const QuadratureRuleAttributes & | quadratureRuleAttributes, |
std::vector< std::shared_ptr< const QuadratureRule > > | quadratureRuleVec, | ||
std::shared_ptr< const basis::TriangulationBase > | triangulation, | ||
const basis::CellMappingBase & | cellMapping | ||
) |
Constructor for assigning each cell with different quadrature rule for each cell.
[in] | quadratureRulevec | vector of quadratureRule pointers specifying the quadrature rule for each cell |
[in] | triangulation | The triangulation that has information on the cell and its vertices |
[in] | cellMapping | cellMapping provides the the information on how the cell in real space is mapped to its parametric coordinates. This is required to calculate the JxW values at each quad point |
dftefe::quadrature::QuadratureRuleContainer::QuadratureRuleContainer | ( | const QuadratureRuleAttributes & | quadratureRuleAttributes, |
std::shared_ptr< const QuadratureRule > | baseQuadratureRule, | ||
std::shared_ptr< const basis::TriangulationBase > | triangulation, | ||
const basis::CellMappingBase & | cellMapping, | ||
basis::ParentToChildCellsManagerBase & | parentToChildCellsManager, | ||
std::vector< std::shared_ptr< const utils::ScalarSpatialFunctionReal > > | functions, | ||
const std::vector< double > & | absoluteTolerances, | ||
const std::vector< double > & | relativeTolerances, | ||
const std::vector< double > & | integralThresholds, | ||
const double | smallestCellVolume = QuadratureRuleAdaptiveDefaults::SMALLEST_CELL_VOLUME , |
||
const unsigned int | maxRecursion = QuadratureRuleAdaptiveDefaults::MAX_RECURSION |
||
) |
Constructor for creating an adaptive quadrature rule in each cell based user-defined functions.
[in] | baseQuadratureRule | The base quadrature rule to be used in constructing the adaptive quadrature rule |
[in] | triangulation | The triangulation that has information on the cell and its vertices |
[in] | cellMapping | cellMapping object that provides the the information on how the cell in real space is mapped to its parametric coordinates. This is required to calculate the JxW values at each quad point |
std::cout << "Function Eval in recursiveIntegrate function: " << timer["Function Eval"] / 1e6 << "\n" << std::flush; std::cout << "Child Cell Creation in recursiveIntegrate function: " << timer["Child Cell Creation"] / 1e6 << "\n" << std::flush; std::cout << "Cell Mapping to get jxw in recursiveIntegrate function: " << timer["Cell Mapping jxw"] / 1e6 << "\n" << std::flush; std::cout << "Cell Mapping to get real pts in recursiveIntegrate function: " << timer["Cell Mapping real"] / 1e6 << "\n" << std::flush;
dftefe::quadrature::QuadratureRuleContainer::QuadratureRuleContainer | ( | const QuadratureRuleAttributes & | quadratureRuleAttributes, |
const size_type | order1DMin, | ||
const size_type | order1DMax, | ||
const size_type | copies1DMax, | ||
std::shared_ptr< const basis::TriangulationBase > | triangulation, | ||
const basis::CellMappingBase & | cellMapping, | ||
std::vector< std::shared_ptr< const utils::ScalarSpatialFunctionReal > > | functions, | ||
const std::vector< double > & | absoluteTolerances, | ||
const std::vector< double > & | relativeTolerances, | ||
const quadrature::QuadratureRuleContainer & | quadratureRuleContainerReference, | ||
const utils::mpi::MPIComm & | comm | ||
) |
Constructor for creating a subdivided quadrature rule in each cell based user-defined functions and a reference quadrature. So one provides a minimum and maximum order (or number of 1D gauss points) and maimum number of copies of the orders one wants to go. Also one has a reference quadrature which can.
[in] | order1DMin | The minimum gauss 1D number of points from which to start |
[in] | order1DMax | The maximum gauss 1D number of points upto which to iterate |
[in] | copies1DMax | The maimum number of copies (starting from 1) to be done at each iteration of order. |
[in] | triangulation | The triangulation that has information on the cell and its vertices |
[in] | cellMapping | cellMapping object that provides the the information on how the cell in real space is mapped to its parametric coordinates. This is required to calculate the JxW values at each quad point |
std::vector< double > dftefe::quadrature::QuadratureRuleContainer::getCellJxW | ( | const unsigned int | cellId | ) | const |
Function that returns a vector containing the Jacobian times weight of the quad points in cell corresponding to the cellId.
[in] | cellId | the id to the cell |
const basis::CellMappingBase & dftefe::quadrature::QuadratureRuleContainer::getCellMapping | ( | ) | const |
const std::vector< dftefe::utils::Point > & dftefe::quadrature::QuadratureRuleContainer::getCellParametricPoints | ( | const unsigned int | cellId | ) | const |
Function that returns a vector containing the real coordinates of the quad points in cell corresponding to the cellId.
[in] | cellId | the id to the cell |
const std::vector< double > & dftefe::quadrature::QuadratureRuleContainer::getCellQuadratureWeights | ( | const unsigned int | cellId | ) | const |
Function that returns a vector containing the weight of the quad points in cell corresponding to the cellId.
[in] | cellId | the id to the cell |
size_type dftefe::quadrature::QuadratureRuleContainer::getCellQuadStartId | ( | const size_type | cellId | ) | const |
A function to return the starting index of the quadrature point of a given cell.
[in] | cellId | index of the cell |
const std::vector< size_type > & dftefe::quadrature::QuadratureRuleContainer::getCellQuadStartIds | ( | ) | const |
A function to return the starting index of the quadrature point of each cell.
std::vector< dftefe::utils::Point > dftefe::quadrature::QuadratureRuleContainer::getCellRealPoints | ( | const unsigned int | cellId | ) | const |
Function that returns a vector containing the real coordinates of the quad points in cell corresponding to the cellId.
[in] | cellId | the id to the cell |
const std::vector< double > & dftefe::quadrature::QuadratureRuleContainer::getJxW | ( | ) | const |
Function that returns a vector containing the Jacobian times quadrature weight for all the quad points across all the cells in the triangulation.
const QuadratureRule & dftefe::quadrature::QuadratureRuleContainer::getQuadratureRule | ( | const unsigned int | cellId | ) | const |
Function that returns the handle to quadrature rule corresponding to the the cell Id.
[in] | cellId | the id to the cell |
const QuadratureRuleAttributes & dftefe::quadrature::QuadratureRuleContainer::getQuadratureRuleAttributes | ( | ) | const |
Returns the underlying QuadratureRuleAttributes.
const std::vector< dftefe::utils::Point > & dftefe::quadrature::QuadratureRuleContainer::getRealPoints | ( | ) | const |
Function that returns a vector containing the real coordinates of the quad points in all cells.
std::shared_ptr< const basis::TriangulationBase > dftefe::quadrature::QuadratureRuleContainer::getTriangulation | ( | ) | const |
size_type dftefe::quadrature::QuadratureRuleContainer::nCellQuadraturePoints | ( | const unsigned int | cellId | ) | const |
A function to returns the number of quadrature points in cell corresponding to the the cell Id.
[in] | cellId | the id to the cell |
size_type dftefe::quadrature::QuadratureRuleContainer::nCells | ( | ) | const |
Returns the number of cells in the quadrature container.
size_type dftefe::quadrature::QuadratureRuleContainer::nQuadraturePoints | ( | ) | const |
A function to return the total number of quadrature points in all the cells.
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |