DFT-EFE
 
Loading...
Searching...
No Matches
TriangulationCellBase.h
Go to the documentation of this file.
1#ifndef dftefeTriangulationCellBase_h
2#define dftefeTriangulationCellBase_h
3
4#include <utils/TypeConfig.h>
5#include <utils/Point.h>
6#include "CellMappingBase.h"
7
8#include <memory>
9namespace dftefe
10{
11 namespace basis
12 {
20 {
21 public:
22 virtual ~TriangulationCellBase() = default;
23 virtual void
24 getVertices(std::vector<utils::Point> &points) const = 0;
25
26 virtual void
27 getVertex(size_type i, utils::Point &point) const = 0;
28
29 virtual size_type
30 getId() const = 0;
31
32 virtual bool
33 isPointInside(const utils::Point &point) const = 0;
34
35 virtual bool
36 isAtBoundary(const unsigned int i) const = 0;
37
38 virtual bool
39 isAtBoundary() const = 0;
40
41 virtual bool
42 hasPeriodicNeighbor(const unsigned int i) const = 0;
43
44 virtual size_type
45 getDim() const = 0;
46
47 virtual double
48 diameter() const = 0;
49
50 virtual void
51 center(dftefe::utils::Point &centerPoint) const = 0;
52
53 virtual void
55
56 virtual void
58
59 virtual double
61
62 virtual double
63 distanceToUnitCell(dftefe::utils::Point &parametricPoint) const = 0;
64
65 virtual void
67 const CellMappingBase &cellMapping,
68 utils::Point & parametricPoint) const = 0;
69
70 virtual void
71 getRealPoint(const utils::Point & parametricPoint,
72 const CellMappingBase &cellMapping,
73 utils::Point & realPoint) const = 0;
74
75 }; // end of class TriaCellBase
76 } // end of namespace basis
77
78} // end of namespace dftefe
79#endif // dftefeTriaCellBase_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
virtual void getVertex(size_type i, utils::Point &point) const =0
virtual void getParametricPoint(const utils::Point &realPoint, const CellMappingBase &cellMapping, utils::Point &parametricPoint) const =0
virtual bool isPointInside(const utils::Point &point) const =0
virtual void getVertices(std::vector< utils::Point > &points) const =0
virtual bool isAtBoundary() const =0
virtual double diameter() const =0
virtual size_type getDim() const =0
virtual bool hasPeriodicNeighbor(const unsigned int i) const =0
virtual void center(dftefe::utils::Point &centerPoint) const =0
virtual ~TriangulationCellBase()=default
virtual double minimumVertexDistance() const =0
virtual bool isAtBoundary(const unsigned int i) const =0
virtual size_type getId() const =0
virtual void getRealPoint(const utils::Point &parametricPoint, const CellMappingBase &cellMapping, utils::Point &realPoint) const =0
virtual double distanceToUnitCell(dftefe::utils::Point &parametricPoint) const =0
Definition: PointImpl.h:13
dealii includes
Definition: AtomFieldDataSpherical.cpp:31
unsigned int size_type
Definition: TypeConfig.h:8