DFT-EFE
 
Loading...
Searching...
No Matches
LinearCellMappingDealii.h
Go to the documentation of this file.
1#ifndef dftefeLinearCellMappingDealii_h
2#define dftefeLinearCellMappingDealii_h
3
4#include "CellMappingBase.h"
6#include <vector>
7#include <utils/Point.h>
8#include <deal.II/fe/mapping.h>
9#include <deal.II/fe/fe_q.h>
10#include <deal.II/fe/mapping_q1.h>
11namespace dftefe
12{
13 namespace basis
14 {
15 template <unsigned int dim>
17 {
18 public:
21
22 void
23 getJxW(const TriangulationCellBase & triaCellBase,
24 const std::vector<dftefe::utils::Point> &paramPoints,
25 const std::vector<double> & weights,
26 std::vector<double> & valuesJxW) const override;
27 void
29 const TriangulationCellBase &triaCellBase,
30 dftefe::utils::Point & parametricPoint,
31 bool &isPointInside) const override;
32
33 void
34 getParametricPoints(const std::vector<dftefe::utils::Point> &realPoints,
35 const TriangulationCellBase & triaCellBase,
36 std::vector<utils::Point> &parametricPoints,
37 std::vector<bool> &arePointsInside) const override;
38
39 void
40 getRealPoint(const dftefe::utils::Point & parametricPoint,
41 const TriangulationCellBase &triaCellBase,
42 dftefe::utils::Point & realPoint) const override;
43
44 void
46 const std::vector<dftefe::utils::Point> &parametricPoints,
47 const TriangulationCellBase & triaCellBase,
48 std::vector<dftefe::utils::Point> & realPoints) const override;
49
50 private:
51 dealii::MappingQ1<dim> d_mappingDealii;
52 dealii::FE_Q<dim> d_fe;
53
54
55 }; // end of class LinearCellMappingDealii
56
57 } // end of namespace basis
58} // end of namespace dftefe
59
61#endif /* dftefeLinearCellMappingDealii_h */
An abstract class to map a real point to parametric point and vice-versa.
Definition: CellMappingBase.h:27
Definition: LinearCellMappingDealii.h:17
void getJxW(const TriangulationCellBase &triaCellBase, const std::vector< dftefe::utils::Point > &paramPoints, const std::vector< double > &weights, std::vector< double > &valuesJxW) const override
Definition: LinearCellMappingDealii.t.cpp:27
dealii::FE_Q< dim > d_fe
Definition: LinearCellMappingDealii.h:52
void getParametricPoints(const std::vector< dftefe::utils::Point > &realPoints, const TriangulationCellBase &triaCellBase, std::vector< utils::Point > &parametricPoints, std::vector< bool > &arePointsInside) const override
Definition: LinearCellMappingDealii.t.cpp:86
void getParametricPoint(const dftefe::utils::Point &realPoint, const TriangulationCellBase &triaCellBase, dftefe::utils::Point &parametricPoint, bool &isPointInside) const override
Definition: LinearCellMappingDealii.t.cpp:55
void getRealPoints(const std::vector< dftefe::utils::Point > &parametricPoints, const TriangulationCellBase &triaCellBase, std::vector< dftefe::utils::Point > &realPoints) const override
Definition: LinearCellMappingDealii.t.cpp:138
LinearCellMappingDealii()
Definition: LinearCellMappingDealii.t.cpp:14
void getRealPoint(const dftefe::utils::Point &parametricPoint, const TriangulationCellBase &triaCellBase, dftefe::utils::Point &realPoint) const override
Definition: LinearCellMappingDealii.t.cpp:118
~LinearCellMappingDealii()
Definition: LinearCellMappingDealii.t.cpp:21
dealii::MappingQ1< dim > d_mappingDealii
Definition: LinearCellMappingDealii.h:51
An abstract class for an geometric cell. This is done to prevent the template (as required by deal....
Definition: TriangulationCellBase.h:20
Definition: PointImpl.h:13
dealii includes
Definition: AtomFieldDataSpherical.cpp:31