DFT-EFE
 
Loading...
Searching...
No Matches
TriangulationCellDealii.h
Go to the documentation of this file.
1#ifndef dftefeTriangulationCellDealii_h
2#define dftefeTriangulationCellDealii_h
3
4#include <utils/Point.h>
5#include <utils/TypeConfig.h>
6#include "CellMappingBase.h"
8
9#include <deal.II/grid/tria.h>
10#include <deal.II/grid/tria_accessor.h>
11#include <deal.II/grid/tria_iterator.h>
12
13#include <memory>
14namespace dftefe
15{
16 namespace basis
17 {
21 template <unsigned int dim>
23 {
24 /*
25 * typedefs
26 */
27 public:
29 typename dealii::Triangulation<dim>::active_cell_iterator;
30
31 public:
34
35 void
36 getVertices(std::vector<utils::Point> &points) const override;
37
38 void
39 getVertex(size_type i, utils::Point &point) const override;
40
42 getId() const override;
43
44 bool
45 isPointInside(const utils::Point &point) const override;
46
47 bool
48 isAtBoundary(const unsigned int i) const override;
49
50 bool
51 isAtBoundary() const override;
52
53 unsigned int
54 getDim() const override;
55
56 double
57 diameter() const override;
58
59 void
60 center(dftefe::utils::Point &centerPoint) const override;
61
62 void
63 setRefineFlag() override;
64 /*
65 * \todo
66 * TODO : Should implement the cellMapping before implementation
67 */
68
69 void
70 clearRefineFlag() override;
71
72 double
73 minimumVertexDistance() const override;
74
75 double
76 distanceToUnitCell(dftefe::utils::Point &parametricPoint) const override;
77
78 void
80 const CellMappingBase & cellMapping,
81 dftefe::utils::Point &parametricPoint) const override;
82
83 /*
84 * \todo
85 * TODO : Should implement the cellMapping before implementation
86 */
87 void
88 getRealPoint(const utils::Point & parametricPoint,
89 const CellMappingBase &cellMapping,
90 utils::Point & realPoint) const override;
91
94
95 private:
97
98 }; // end of class TriaCellDealii
99 } // end of namespace basis
100
101} // end of namespace dftefe
103#endif // dftefeTriaCellDealii_h
An abstract class to map a real point to parametric point and vice-versa.
Definition: CellMappingBase.h:27
An abstract class for an geometric cell. This is done to prevent the template (as required by deal....
Definition: TriangulationCellBase.h:20
An interface to deal.ii geometric cell.
Definition: TriangulationCellDealii.h:23
void getParametricPoint(const dftefe::utils::Point &realPoint, const CellMappingBase &cellMapping, dftefe::utils::Point &parametricPoint) const override
Definition: TriangulationCellDealii.t.cpp:135
double minimumVertexDistance() const override
Definition: TriangulationCellDealii.t.cpp:118
void getVertices(std::vector< utils::Point > &points) const override
Definition: TriangulationCellDealii.t.cpp:22
void clearRefineFlag() override
Definition: TriangulationCellDealii.t.cpp:111
double distanceToUnitCell(dftefe::utils::Point &parametricPoint) const override
Definition: TriangulationCellDealii.t.cpp:125
DealiiTriangulationCellIterator & getCellIterator()
Definition: TriangulationCellDealii.t.cpp:161
bool isAtBoundary() const override
Definition: TriangulationCellDealii.t.cpp:74
double diameter() const override
Definition: TriangulationCellDealii.t.cpp:87
void getRealPoint(const utils::Point &parametricPoint, const CellMappingBase &cellMapping, utils::Point &realPoint) const override
Definition: TriangulationCellDealii.t.cpp:149
void setRefineFlag() override
Definition: TriangulationCellDealii.t.cpp:104
DealiiTriangulationCellIterator d_cellItr
Definition: TriangulationCellDealii.h:96
bool isPointInside(const utils::Point &point) const override
Definition: TriangulationCellDealii.t.cpp:57
void center(dftefe::utils::Point &centerPoint) const override
Definition: TriangulationCellDealii.t.cpp:94
size_type getId() const override
Definition: TriangulationCellDealii.t.cpp:48
typename dealii::Triangulation< dim >::active_cell_iterator DealiiTriangulationCellIterator
Definition: TriangulationCellDealii.h:29
~TriangulationCellDealii()
Definition: TriangulationCellDealii.t.cpp:17
void getVertex(size_type i, utils::Point &point) const override
Definition: TriangulationCellDealii.t.cpp:40
unsigned int getDim() const override
Definition: TriangulationCellDealii.t.cpp:81
Definition: PointImpl.h:13
dealii includes
Definition: AtomFieldDataSpherical.cpp:31
unsigned int size_type
Definition: TypeConfig.h:8