3#ifndef DUNE_MMESH_INTERFACE_GEOMETRY_HH
4#define DUNE_MMESH_INTERFACE_GEOMETRY_HH
12#include <dune/common/fmatrix.hh>
13#include <dune/common/typetraits.hh>
14#include <dune/geometry/affinegeometry.hh>
15#include <dune/grid/common/geometry.hh>
26template <
int mydim,
int coorddim,
class Gr
idImp>
31template <
int md,
class Gr
idImp>
33 :
public AffineGeometry<typename GridImp::ctype, md, 2> {
35 static constexpr int mydim = md;
36 static constexpr int coorddim = 2;
37 typedef AffineGeometry<typename GridImp::ctype, mydim, coorddim> BaseType;
38 typedef FieldVector<typename GridImp::ctype, coorddim> FVector;
40 enum { dimension = GridImp::dimension };
41 enum { dimensionworld = GridImp::dimensionworld };
42 enum { coorddimension = coorddim };
43 enum { mydimension = mydim };
47 const typename GridImp::template MMeshInterfaceEntity<0>& hostEntity)
48 : BaseType(GeometryTypes::simplex(mydim), getVertices(hostEntity)) {}
52 : BaseType(GeometryTypes::simplex(mydim), points) {}
56 const typename GridImp::template MMeshInterfaceEntity<1>& hostEntity)
58 GeometryTypes::simplex(mydim),
59 std::array<FVector, 1>({makeFieldVector(hostEntity->point())})) {}
64 : BaseType(GeometryTypes::simplex(1),
65 std::array<FVector, 2>(
66 {FVector({i == 2 ? 1.0 : 0.0, 0.0}),
67 FVector({i == 0 ? 1.0 : 0.0, i > 0 ? 1.0 : 0.0})})) {}
73 static inline std::array<FVector, 2> getVertices(
74 const typename GridImp::template MMeshInterfaceEntity<0>& hostEntity) {
75 const auto& cgalIdx = MMeshInterfaceImpl::computeCGALIndices<
76 typename GridImp::template MMeshInterfaceEntity<0>, 1>(hostEntity);
78 std::array<FVector, 2> vertices(
79 {makeFieldVector(hostEntity.first->vertex(cgalIdx[0])->point()),
80 makeFieldVector(hostEntity.first->vertex(cgalIdx[1])->point())});
88template <
int md,
class Gr
idImp>
90 :
public AffineGeometry<typename GridImp::ctype, md, 3> {
92 static constexpr int mydim = md;
93 static constexpr int coorddim = 3;
94 typedef AffineGeometry<typename GridImp::ctype, mydim, coorddim> BaseType;
95 typedef typename GridImp::ctype ctype;
96 typedef FieldVector<ctype, coorddim> FVector;
98 enum { dimension = GridImp::dimension };
99 enum { dimensionworld = GridImp::dimensionworld };
100 enum { coorddimension = coorddim };
101 enum { mydimension = mydim };
105 const typename GridImp::template MMeshInterfaceEntity<0>& hostEntity)
106 : BaseType(GeometryTypes::simplex(mydim), getVertices<0>(hostEntity)) {}
110 : BaseType(GeometryTypes::simplex(mydim), points) {}
114 const typename GridImp::template MMeshInterfaceEntity<1>& hostEntity)
115 : BaseType(GeometryTypes::simplex(mydim), getVertices<1>(hostEntity)) {}
119 const typename GridImp::template MMeshInterfaceEntity<2>& hostEntity)
121 GeometryTypes::simplex(mydim),
122 std::array<FVector, 1>({makeFieldVector(hostEntity->point())})) {}
128 template <
int codim,
typename Enable = std::enable_if_t<codim == 0> >
129 static inline std::array<FVector, 3> getVertices(
130 const typename GridImp::template MMeshInterfaceEntity<0>& hostEntity) {
131 const auto& cgalIdx = MMeshInterfaceImpl::computeCGALIndices<
132 typename GridImp::template MMeshInterfaceEntity<0>, 2>(hostEntity);
134 std::array<FVector, 3> vertices;
136 const auto& cell = hostEntity.first;
139 for (
int i = 0; i < 3; ++i)
140 vertices[i] = makeFieldVector(cell->vertex(cgalIdx[i])->point());
145 template <
int codim,
typename Enable = std::enable_if_t<codim == 1> >
146 static inline std::array<FVector, 2> getVertices(
147 const typename GridImp::template MMeshInterfaceEntity<1>& hostEntity) {
148 const auto& cell = hostEntity.first;
149 const auto& i = hostEntity.second;
150 const auto& j = hostEntity.third;
152 std::array<FVector, 2> vertices;
153 vertices[0] = makeFieldVector(cell->vertex(i)->point());
154 vertices[1] = makeFieldVector(cell->vertex(j)->point());
156 if (cell->vertex(i)->info().id > cell->vertex(j)->info().id)
157 std::swap(vertices[0], vertices[1]);
164template <
int md,
class Gr
idImp>
166 :
public AffineGeometry<typename GridImp::ctype, md, 1> {
168 static constexpr int mydim = md;
169 static constexpr int coorddim = 1;
170 typedef AffineGeometry<typename GridImp::ctype, mydim, coorddim> BaseType;
171 typedef FieldVector<typename GridImp::ctype, coorddim> FVector;
173 enum { dimension = GridImp::dimension };
174 enum { dimensionworld = GridImp::dimensionworld };
175 enum { coorddimension = coorddim };
176 enum { mydimension = mydim };
180 : BaseType(GeometryTypes::simplex(mydim),
std::array<FVector, 1>({i})) {}
184 : BaseType(GeometryTypes::simplex(mydim), points) {}
MMeshInterfaceGridGeometry(const std::array< FVector, 2 > &points)
Constructor from vertex array.
Definition: geometry.hh:183
MMeshInterfaceGridGeometry(int i)
Constructor for local geometry of intersection from intersection index.
Definition: geometry.hh:179
const FVector circumcenter() const
Obtain the circumcenter.
Definition: geometry.hh:70
MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity< 1 > &hostEntity)
Constructor from host geometry with codim 1.
Definition: geometry.hh:55
MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity< 0 > &hostEntity)
Constructor from host geometry with codim 0.
Definition: geometry.hh:46
MMeshInterfaceGridGeometry(int i)
Definition: geometry.hh:63
MMeshInterfaceGridGeometry(const std::array< FVector, 2 > &points)
Constructor from vertex array.
Definition: geometry.hh:51
MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity< 0 > &hostEntity)
Constructor from host geometry with codim 0.
Definition: geometry.hh:104
MMeshInterfaceGridGeometry(const std::array< FVector, 3 > &points)
Constructor from vertex list.
Definition: geometry.hh:109
MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity< 2 > &hostEntity)
Constructor from host geometry with codim 2.
Definition: geometry.hh:118
const FVector circumcenter() const
Obtain the circumcenter.
Definition: geometry.hh:125
MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity< 1 > &hostEntity)
Constructor from host geometry with codim 1.
Definition: geometry.hh:113
Geometry.
Definition: geometry.hh:27
Some common helper methods.
Helpers for conversion from CGAL::Point_x to DUNE::FieldVector.
auto computeCircumcenter(const Geometry &geo)
Compute circumcenter.
Definition: pointfieldvector.hh:78