Loading [MathJax]/extensions/tex2jax.js

dune-mmesh (1.4)

geometry.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_MMESH_INTERFACE_GEOMETRY_HH
4#define DUNE_MMESH_INTERFACE_GEOMETRY_HH
5
11// Dune includes
12#include <dune/common/fmatrix.hh>
13#include <dune/common/typetraits.hh>
14#include <dune/grid/common/geometry.hh>
15#include <dune/geometry/affinegeometry.hh>
16
17// local includes
19
20namespace Dune
21{
22
26 template<int mydim, int coorddim, class GridImp>
28
30
31 template<int md, class GridImp>
32 class MMeshInterfaceGridGeometry<md, 2, GridImp> :
33 public AffineGeometry <typename GridImp::ctype, md, 2>
34 {
35 public:
36 static constexpr int mydim = md;
37 static constexpr int coorddim = 2;
38 typedef AffineGeometry <typename GridImp::ctype, mydim, coorddim> BaseType;
39 typedef FieldVector<typename GridImp::ctype, coorddim> FVector;
40
41 enum { dimension = GridImp::dimension };
42 enum { dimensionworld = GridImp::dimensionworld };
43 enum { coorddimension = coorddim };
44 enum { mydimension = mydim };
45
47 MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity<0>& hostEntity)
48 : BaseType( GeometryTypes::simplex(mydim), getVertices(hostEntity) )
49 {}
50
52 MMeshInterfaceGridGeometry(const std::array<FVector, 2>& points)
53 : BaseType( GeometryTypes::simplex(mydim), points )
54 {}
55
57 MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity<1>& hostEntity)
58 : BaseType( GeometryTypes::simplex(mydim), std::array<FVector, 1>( { makeFieldVector( hostEntity->point() ) } ) )
59 {}
60
63 : BaseType( GeometryTypes::simplex(1),
64 std::array<FVector, 2>(
65 {
66 FVector ( { i==2 ? 1.0 : 0.0, 0.0 } ),
67 FVector ( { i==0 ? 1.0 : 0.0, i>0 ? 1.0 : 0.0 } )
68 } ) )
69 {}
70
72 const FVector circumcenter () const
73 {
74 return computeCircumcenter(*this);
75 }
76
77 private:
78 static inline std::array<FVector, 2> getVertices ( const typename GridImp::template MMeshInterfaceEntity<0>& hostEntity )
79 {
80 const auto& cgalIdx = MMeshInterfaceImpl::computeCGALIndices<typename GridImp::template MMeshInterfaceEntity<0>, 1> ( hostEntity );
81
82 std::array<FVector, 2> vertices ( {
83 makeFieldVector( hostEntity.first->vertex( cgalIdx[0] )->point() ),
84 makeFieldVector( hostEntity.first->vertex( cgalIdx[1] )->point() )
85 } );
86
87 return vertices;
88 }
89 };
90
92
93 template<int md, class GridImp>
94 class MMeshInterfaceGridGeometry<md, 3, GridImp> :
95 public AffineGeometry <typename GridImp::ctype, md, 3>
96 {
97 public:
98 static constexpr int mydim = md;
99 static constexpr int coorddim = 3;
100 typedef AffineGeometry <typename GridImp::ctype, mydim, coorddim> BaseType;
101 typedef typename GridImp::ctype ctype;
102 typedef FieldVector<ctype, coorddim> FVector;
103
104 enum { dimension = GridImp::dimension };
105 enum { dimensionworld = GridImp::dimensionworld };
106 enum { coorddimension = coorddim };
107 enum { mydimension = mydim };
108
110 MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity<0>& hostEntity)
111 : BaseType( GeometryTypes::simplex(mydim), getVertices<0>(hostEntity) )
112 {}
113
115 MMeshInterfaceGridGeometry(const std::array<FVector, 3>& points)
116 : BaseType( GeometryTypes::simplex(mydim), points )
117 {}
118
120 MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity<1>& hostEntity)
121 : BaseType( GeometryTypes::simplex(mydim), getVertices<1>( hostEntity ) )
122 {}
123
125 MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity<2>& hostEntity)
126 : BaseType( GeometryTypes::simplex(mydim), std::array<FVector, 1>( { makeFieldVector( hostEntity->point() ) } ) )
127 {}
128
130 const FVector circumcenter () const
131 {
132 return computeCircumcenter(*this);
133 }
134
135 private:
136
137 template< int codim, typename Enable = std::enable_if_t< codim == 0 > >
138 static inline std::array<FVector, 3> getVertices ( const typename GridImp::template MMeshInterfaceEntity<0>& hostEntity )
139 {
140 const auto& cgalIdx = MMeshInterfaceImpl::computeCGALIndices<typename GridImp::template MMeshInterfaceEntity<0>, 2> ( hostEntity );
141
142 std::array<FVector, 3> vertices;
143
144 const auto& cell = hostEntity.first;
145
146 // use the CGAL index convention to obtain the vertices
147 for ( int i = 0; i < 3; ++i )
148 vertices[i] = makeFieldVector( cell->vertex( cgalIdx[i] )->point() );
149
150 return vertices;
151 }
152
153 template< int codim, typename Enable = std::enable_if_t< codim == 1 > >
154 static inline std::array<FVector, 2> getVertices ( const typename GridImp::template MMeshInterfaceEntity<1>& hostEntity )
155 {
156 const auto& cell = hostEntity.first;
157 const auto& i = hostEntity.second;
158 const auto& j = hostEntity.third;
159
160 std::array<FVector, 2> vertices;
161 vertices[0] = makeFieldVector( cell->vertex( i )->point() );
162 vertices[1] = makeFieldVector( cell->vertex( j )->point() );
163
164 if ( cell->vertex( i )->info().id > cell->vertex( j )->info().id )
165 std::swap( vertices[0], vertices[1] );
166
167 return vertices;
168 }
169 };
170
172 template<int md, class GridImp>
173 class MMeshInterfaceGridGeometry<md, 1, GridImp> :
174 public AffineGeometry <typename GridImp::ctype, md, 1>
175 {
176 public:
177 static constexpr int mydim = md;
178 static constexpr int coorddim = 1;
179 typedef AffineGeometry <typename GridImp::ctype, mydim, coorddim> BaseType;
180 typedef FieldVector<typename GridImp::ctype, coorddim> FVector;
181
182 enum { dimension = GridImp::dimension };
183 enum { dimensionworld = GridImp::dimensionworld };
184 enum { coorddimension = coorddim };
185 enum { mydimension = mydim };
186
189 : BaseType( GeometryTypes::simplex(mydim), std::array<FVector, 1>( { i } ) )
190 {}
191
193 MMeshInterfaceGridGeometry(const std::array<FVector, 2>& points)
194 : BaseType( GeometryTypes::simplex(mydim), points )
195 {}
196
197 };
198
199} // namespace Dune
200
201#endif
MMeshInterfaceGridGeometry(const std::array< FVector, 2 > &points)
Constructor from vertex array.
Definition: geometry.hh:193
MMeshInterfaceGridGeometry(int i)
Constructor for local geometry of intersection from intersection index.
Definition: geometry.hh:188
const FVector circumcenter() const
Obtain the circumcenter.
Definition: geometry.hh:72
MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity< 1 > &hostEntity)
Constructor from host geometry with codim 1.
Definition: geometry.hh:57
MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity< 0 > &hostEntity)
Constructor from host geometry with codim 0.
Definition: geometry.hh:47
MMeshInterfaceGridGeometry(int i)
Constructor for local geometry of intersection from intersection index for 3D.
Definition: geometry.hh:62
MMeshInterfaceGridGeometry(const std::array< FVector, 2 > &points)
Constructor from vertex array.
Definition: geometry.hh:52
MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity< 0 > &hostEntity)
Constructor from host geometry with codim 0.
Definition: geometry.hh:110
MMeshInterfaceGridGeometry(const std::array< FVector, 3 > &points)
Constructor from vertex list.
Definition: geometry.hh:115
MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity< 2 > &hostEntity)
Constructor from host geometry with codim 2.
Definition: geometry.hh:125
const FVector circumcenter() const
Obtain the circumcenter.
Definition: geometry.hh:130
MMeshInterfaceGridGeometry(const typename GridImp::template MMeshInterfaceEntity< 1 > &hostEntity)
Constructor from host geometry with codim 1.
Definition: geometry.hh:120
Geometry.
Definition: geometry.hh:27
STL namespace.
Helpers for conversion from CGAL::Point_x to DUNE::FieldVector.
auto computeCircumcenter(const Geometry &geo)
Compute circumcenter.
Definition: pointfieldvector.hh:82
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 6, 22:49, 2025)