DFT-EFE
 
Loading...
Searching...
No Matches
dftefe::quadrature::QuadratureRuleContainer Class Reference

#include <QuadratureRuleContainer.h>

Collaboration diagram for dftefe::quadrature::QuadratureRuleContainer:

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 QuadratureRuleAttributesgetQuadratureRuleAttributes () 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::PointgetCellRealPoints (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 QuadratureRulegetQuadratureRule (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::TriangulationBasegetTriangulation () const
 
const basis::CellMappingBasegetCellMapping () const
 

Private Attributes

const QuadratureRuleAttributesd_quadratureRuleAttributes
 
std::vector< std::shared_ptr< const QuadratureRule > > d_quadratureRuleVec
 
std::vector< size_typed_numCellQuadPoints
 
std::vector< size_typed_cellQuadStartIds
 
std::vector< dftefe::utils::Pointd_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::TriangulationBased_triangulation
 
const basis::CellMappingBased_cellMapping
 

Detailed Description

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.

Constructor & Destructor Documentation

◆ QuadratureRuleContainer() [1/4]

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.

Parameters
[in]quadratureRuleThe quadrature rule specifying the quad points in the parametric coordinates with its corresponding weights.
[in]triangulationThe triangulation that has information on the cell and its vertices
[in]cellMappingcellMapping 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
Here is the call graph for this function:

◆ QuadratureRuleContainer() [2/4]

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.

Parameters
[in]quadratureRulevecvector of quadratureRule pointers specifying the quadrature rule for each cell
[in]triangulationThe triangulation that has information on the cell and its vertices
[in]cellMappingcellMapping 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
Here is the call graph for this function:

◆ QuadratureRuleContainer() [3/4]

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.

Parameters
[in]baseQuadratureRuleThe base quadrature rule to be used in constructing the adaptive quadrature rule
[in]triangulationThe triangulation that has information on the cell and its vertices
[in]cellMappingcellMapping 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;

Here is the call graph for this function:

◆ QuadratureRuleContainer() [4/4]

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.

  1. either be constructed from highly refined tensor structures of GAUSS or GLL. from spatally variable quadrature rules like GAUSS_VARIABLE, GLL_VARIABLE, or ADAPTIVE. Care should be taken that for ADAPTIVE it is assumed that the quadrature is optimal and hence only cells with high quadrature density are traversed. Whereas for all others all the cells are traversed. Then the {order, copy} pair with minimum quad points per cell is chosen. Then again the maximum of this points over all processors are taken. In one statement the condition can be written as sup_{processors} inf_{all pairs in a processor} (QuadPoints in Cells statisfing the given tolerances w.r.t. the reference quadrature)
    Parameters
    [in]order1DMinThe minimum gauss 1D number of points from which to start
    [in]order1DMaxThe maximum gauss 1D number of points upto which to iterate
    [in]copies1DMaxThe maimum number of copies (starting from 1) to be done at each iteration of order.
    [in]triangulationThe triangulation that has information on the cell and its vertices
    [in]cellMappingcellMapping 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
Here is the call graph for this function:

Member Function Documentation

◆ getCellJxW()

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.

Parameters
[in]cellIdthe id to the cell
Returns
a vector (double) of Jacobian times weight
Here is the caller graph for this function:

◆ getCellMapping()

const basis::CellMappingBase & dftefe::quadrature::QuadratureRuleContainer::getCellMapping ( ) const
Here is the caller graph for this function:

◆ getCellParametricPoints()

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.

Parameters
[in]cellIdthe id to the cell
Returns
a vector of dftefe::utils::Point
Here is the caller graph for this function:

◆ getCellQuadratureWeights()

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.

Parameters
[in]cellIdthe id to the cell
Returns
a vector of weights double
Here is the caller graph for this function:

◆ getCellQuadStartId()

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.

Parameters
[in]cellIdindex of the cell
Returns
the starting index of the quadrature point of the cell

◆ getCellQuadStartIds()

const std::vector< size_type > & dftefe::quadrature::QuadratureRuleContainer::getCellQuadStartIds ( ) const

A function to return the starting index of the quadrature point of each cell.

Returns
vector storing the starting index of the quadrature point of each cell

◆ getCellRealPoints()

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.

Parameters
[in]cellIdthe id to the cell
Returns
a vector of dftefe::utils::Point

◆ getJxW()

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.

Returns
vectors of Jacobian times quadrature weight
Here is the caller graph for this function:

◆ getQuadratureRule()

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.

Parameters
[in]cellIdthe id to the cell
Returns
Const reference to QuadratureRule

◆ getQuadratureRuleAttributes()

const QuadratureRuleAttributes & dftefe::quadrature::QuadratureRuleContainer::getQuadratureRuleAttributes ( ) const

Returns the underlying QuadratureRuleAttributes.

Returns
const reference to the QuadratureAttributes
Here is the caller graph for this function:

◆ getRealPoints()

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.

Returns
a vector of dftefe::utils::Point
Here is the caller graph for this function:

◆ getTriangulation()

std::shared_ptr< const basis::TriangulationBase > dftefe::quadrature::QuadratureRuleContainer::getTriangulation ( ) const
Here is the caller graph for this function:

◆ nCellQuadraturePoints()

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.

Parameters
[in]cellIdthe id to the cell
Returns
number of quadrature points
Here is the caller graph for this function:

◆ nCells()

size_type dftefe::quadrature::QuadratureRuleContainer::nCells ( ) const

Returns the number of cells in the quadrature container.

Returns
number of cells in the quadrature container

◆ nQuadraturePoints()

size_type dftefe::quadrature::QuadratureRuleContainer::nQuadraturePoints ( ) const

A function to return the total number of quadrature points in all the cells.

Returns
the number of quadrature points in all the cells
Here is the caller graph for this function:

Member Data Documentation

◆ d_cellMapping

const basis::CellMappingBase& dftefe::quadrature::QuadratureRuleContainer::d_cellMapping
private

◆ d_cellQuadStartIds

std::vector<size_type> dftefe::quadrature::QuadratureRuleContainer::d_cellQuadStartIds
private

◆ d_dim

unsigned int dftefe::quadrature::QuadratureRuleContainer::d_dim
private

◆ d_JxW

std::vector<double> dftefe::quadrature::QuadratureRuleContainer::d_JxW
private

◆ d_numCellQuadPoints

std::vector<size_type> dftefe::quadrature::QuadratureRuleContainer::d_numCellQuadPoints
private

◆ d_numCells

size_type dftefe::quadrature::QuadratureRuleContainer::d_numCells
private

◆ d_numQuadPoints

size_type dftefe::quadrature::QuadratureRuleContainer::d_numQuadPoints
private

◆ d_quadratureRuleAttributes

const QuadratureRuleAttributes& dftefe::quadrature::QuadratureRuleContainer::d_quadratureRuleAttributes
private

◆ d_quadratureRuleVec

std::vector<std::shared_ptr<const QuadratureRule> > dftefe::quadrature::QuadratureRuleContainer::d_quadratureRuleVec
private

◆ d_realPoints

std::vector<dftefe::utils::Point> dftefe::quadrature::QuadratureRuleContainer::d_realPoints
private

◆ d_storeJacobianInverse

bool dftefe::quadrature::QuadratureRuleContainer::d_storeJacobianInverse
private

◆ d_triangulation

std::shared_ptr<const basis::TriangulationBase> dftefe::quadrature::QuadratureRuleContainer::d_triangulation
private

The documentation for this class was generated from the following files: