DFT-EFE
 
Loading...
Searching...
No Matches
QuadratureRuleContainer.h
Go to the documentation of this file.
1#ifndef dftefeQuadratureRuleContainer_h
2#define dftefeQuadratureRuleContainer_h
3
6#include <utils/TypeConfig.h>
7#include <utils/Point.h>
12#include <memory>
13#include <utils/MPITypes.h>
14#include <utils/MPIWrapper.h>
15#include <quadrature/Defaults.h>
16namespace dftefe
17{
18 namespace quadrature
19 {
26 {
27 public:
39 const QuadratureRuleAttributes & quadratureRuleAttributes,
40 std::shared_ptr<const QuadratureRule> quadratureRule,
41 std::shared_ptr<const basis::TriangulationBase> triangulation,
42 const basis::CellMappingBase & cellMapping);
43
55 const QuadratureRuleAttributes &quadratureRuleAttributes,
56 std::vector<std::shared_ptr<const QuadratureRule>> quadratureRuleVec,
57 std::shared_ptr<const basis::TriangulationBase> triangulation,
58 const basis::CellMappingBase & cellMapping);
59
73 const QuadratureRuleAttributes & quadratureRuleAttributes,
74 std::shared_ptr<const QuadratureRule> baseQuadratureRule,
75 std::shared_ptr<const basis::TriangulationBase> triangulation,
76 const basis::CellMappingBase & cellMapping,
77 basis::ParentToChildCellsManagerBase &parentToChildCellsManager,
78 std::vector<std::shared_ptr<const utils::ScalarSpatialFunctionReal>>
79 functions,
80 const std::vector<double> &absoluteTolerances,
81 const std::vector<double> &relativeTolerances,
82 const std::vector<double> &integralThresholds,
83 const double smallestCellVolume =
85 const unsigned int maxRecursion =
87
118 const QuadratureRuleAttributes &quadratureRuleAttributes,
119 const size_type order1DMin,
120 const size_type order1DMax,
121 const size_type copies1DMax,
122 std::shared_ptr<const basis::TriangulationBase> triangulation,
123 const basis::CellMappingBase & cellMapping,
124 std::vector<std::shared_ptr<const utils::ScalarSpatialFunctionReal>>
125 functions,
126 const std::vector<double> &absoluteTolerances,
127 const std::vector<double> &relativeTolerances,
129 & quadratureRuleContainerReference,
130 const utils::mpi::MPIComm &comm);
131
133
140
146 nCells() const;
147
153 const std::vector<dftefe::utils::Point> &
154 getRealPoints() const;
155
163 std::vector<dftefe::utils::Point>
164 getCellRealPoints(const unsigned int cellId) const;
165
173 const std::vector<dftefe::utils::Point> &
174 getCellParametricPoints(const unsigned int cellId) const;
175
183 const std::vector<double> &
184 getCellQuadratureWeights(const unsigned int cellId) const;
185
192 const std::vector<double> &
193 getJxW() const;
194
202 std::vector<double>
203 getCellJxW(const unsigned int cellId) const;
204
212 const QuadratureRule &
213 getQuadratureRule(const unsigned int cellId) const;
214
221 nQuadraturePoints() const;
222
231 nCellQuadraturePoints(const unsigned int cellId) const;
232
238 const std::vector<size_type> &
239 getCellQuadStartIds() const;
240
248 getCellQuadStartId(const size_type cellId) const;
249
250 std::shared_ptr<const basis::TriangulationBase>
251 getTriangulation() const;
252
254 getCellMapping() const;
255
256 private:
258 std::vector<std::shared_ptr<const QuadratureRule>> d_quadratureRuleVec;
259 std::vector<size_type> d_numCellQuadPoints;
260 std::vector<size_type> d_cellQuadStartIds;
261 std::vector<dftefe::utils::Point> d_realPoints;
262 std::vector<double> d_JxW;
263 unsigned int d_dim;
267 std::shared_ptr<const basis::TriangulationBase> d_triangulation;
269 };
270 } // end of namespace quadrature
271
272} // end of namespace dftefe
273
274#endif // dftefeQuadratureRuleContainer_h
An abstract class to map a real point to parametric point and vice-versa.
Definition: CellMappingBase.h:27
Definition: ParentToChildCellsManagerBase.h:13
static const double SMALLEST_CELL_VOLUME
Setting all the QuadratureRuleAdaptiveDefaults.
Definition: Defaults.h:42
static const unsigned int MAX_RECURSION
Definition: Defaults.h:47
Definition: QuadratureAttributes.h:135
Definition: QuadratureRuleContainer.h:26
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 ...
Definition: QuadratureRuleContainer.cpp:690
const std::vector< double > & getJxW() const
Function that returns a vector containing the Jacobian times quadrature weight for all the quad point...
Definition: QuadratureRuleContainer.cpp:697
bool d_storeJacobianInverse
Definition: QuadratureRuleContainer.h:266
const basis::CellMappingBase & getCellMapping() const
Definition: QuadratureRuleContainer.cpp:754
const std::vector< size_type > & getCellQuadStartIds() const
A function to return the starting index of the quadrature point of each cell.
Definition: QuadratureRuleContainer.cpp:736
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.
Definition: QuadratureRuleContainer.cpp:729
std::vector< std::shared_ptr< const QuadratureRule > > d_quadratureRuleVec
Definition: QuadratureRuleContainer.h:258
const QuadratureRuleAttributes & getQuadratureRuleAttributes() const
Returns the underlying QuadratureRuleAttributes.
Definition: QuadratureRuleContainer.cpp:648
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 correspondi...
Definition: QuadratureRuleContainer.cpp:666
std::shared_ptr< const basis::TriangulationBase > d_triangulation
Definition: QuadratureRuleContainer.h:267
size_type nCells() const
Returns the number of cells in the quadrature container.
Definition: QuadratureRuleContainer.cpp:654
const basis::CellMappingBase & d_cellMapping
Definition: QuadratureRuleContainer.h:268
size_type nQuadraturePoints() const
A function to return the total number of quadrature points in all the cells.
Definition: QuadratureRuleContainer.cpp:723
const std::vector< dftefe::utils::Point > & getRealPoints() const
Function that returns a vector containing the real coordinates of the quad points in all cells.
Definition: QuadratureRuleContainer.cpp:660
std::vector< size_type > d_cellQuadStartIds
Definition: QuadratureRuleContainer.h:260
std::vector< double > d_JxW
Definition: QuadratureRuleContainer.h:262
std::vector< dftefe::utils::Point > d_realPoints
Definition: QuadratureRuleContainer.h:261
size_type d_numCells
Definition: QuadratureRuleContainer.h:265
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 correspondi...
Definition: QuadratureRuleContainer.cpp:683
std::vector< size_type > d_numCellQuadPoints
Definition: QuadratureRuleContainer.h:259
unsigned int d_dim
Definition: QuadratureRuleContainer.h:263
size_type getCellQuadStartId(const size_type cellId) const
A function to return the starting index of the quadrature point of a given cell.
Definition: QuadratureRuleContainer.cpp:742
const QuadratureRule & getQuadratureRule(const unsigned int cellId) const
Function that returns the handle to quadrature rule corresponding to the the cell Id.
Definition: QuadratureRuleContainer.cpp:717
const QuadratureRuleAttributes & d_quadratureRuleAttributes
Definition: QuadratureRuleContainer.h:257
size_type d_numQuadPoints
Definition: QuadratureRuleContainer.h:264
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 corres...
Definition: QuadratureRuleContainer.cpp:703
std::shared_ptr< const basis::TriangulationBase > getTriangulation() const
Definition: QuadratureRuleContainer.cpp:748
Definition: QuadratureRule.h:18
int MPIComm
Definition: MPITypes.h:83
dealii includes
Definition: AtomFieldDataSpherical.cpp:31
unsigned int size_type
Definition: TypeConfig.h:8