Loading [MathJax]/extensions/tex2jax.js

DUNE-GRID-GLUE (2.10)

projection.hh
1// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LGPL-3.0-or-later OR LicenseRef-GPL-2.0-only-with-dune-grid-glue-exception
3#ifndef DUNE_GRIDGLUE_COMMON_PROJECTIONHELPER2_HH
4#define DUNE_GRIDGLUE_COMMON_PROJECTIONHELPER2_HH
5
6#include <array>
7#include <bitset>
8#include <tuple>
9
10namespace Dune {
11namespace GridGlue {
12
19template<typename Coordinate>
21{
22public:
29 {
33 std::array<unsigned, 2> edge;
34
41 std::array<Coordinate, 2> local;
42 };
43
47 constexpr static unsigned dim = Coordinate::dimension;
48
54 constexpr static unsigned maxEdgeIntersections = dim == 3 ? 9 : 0;
55
56 static_assert(dim == 2 || dim == 3, "Projection only implemented for dim=2 or dim=3");
57
61 typedef typename Coordinate::field_type Field;
62
70 typedef std::array<Coordinate, dim> Images;
71
79
80private:
84 const Field m_overlap;
85
94 const Field m_max_normal_product;
95
101 Field m_epsilon = Field(1e-12);
102
104 std::tuple<Images, Preimages> m_images;
105
107 std::tuple<std::bitset<dim>, std::bitset<dim> > m_success;
108
110 unsigned m_number_of_edge_intersections;
111
113 std::array<EdgeIntersection, maxEdgeIntersections> m_edge_intersections;
114
126 bool m_projection_valid;
127
133 template<typename Corners, typename Normals>
134 void doProjection(const std::tuple<Corners&, Corners&>& corners, const std::tuple<Normals&, Normals&>& normals);
135
144 template<typename Corners, typename Normals>
145 void doInverseProjection(const std::tuple<Corners&, Corners&>& corners, const std::tuple<Normals&, Normals&>& normals);
146
155 template<typename Corners, typename Normals>
156 void doEdgeIntersection(const std::tuple<Corners&, Corners&>& corners, const std::tuple<Normals&, Normals&>& normals);
157
183 template<typename Corners, typename Normals>
184 inline bool projectionFeasible(const Coordinate& x, const Coordinate& nx, const Coordinate& px, const Corners& corners, const Normals& normals) const;
185
186public:
191 Projection(const Field overlap = Field(0), const Field max_normal_product = Field(-0.1));
192
198 void epsilon(const Field epsilon);
199
212 template<typename Corners, typename Normals>
213 void project(const std::tuple<Corners&, Corners&>& corners, const std::tuple<Normals&, Normals&>& normals);
214
235 const std::tuple<Images, Preimages>& images() const
236 { return m_images; }
237
252 const std::tuple<std::bitset<dim>, std::bitset<dim> >& success() const
253 { return m_success; }
254
263 { return m_number_of_edge_intersections; }
264
273 const std::array<EdgeIntersection, maxEdgeIntersections>& edgeIntersections() const
274 { return m_edge_intersections; }
275};
276
277} /* namespace GridGlue */
278} /* namespace Dune */
279
280#include "projection_impl.hh"
281
282#endif
Projection of a line (triangle) on another line (triangle).
Definition: projection.hh:21
Coordinate::field_type Field
Scalar type.
Definition: projection.hh:61
Projection(const Field overlap=Field(0), const Field max_normal_product=Field(-0.1))
Definition: projection_impl.hh:130
const std::tuple< std::bitset< dim >, std::bitset< dim > > & success() const
Indicate whether projection (inverse projection) is valid for each corner or not.
Definition: projection.hh:252
static constexpr unsigned maxEdgeIntersections
maximum number of edge-edge intersections
Definition: projection.hh:54
Images Preimages
Definition: projection.hh:78
std::array< Coordinate, dim > Images
List of corner images.
Definition: projection.hh:70
void epsilon(const Field epsilon)
Set epsilon used for floating-point comparisons.
Definition: projection_impl.hh:140
void project(const std::tuple< Corners &, Corners & > &corners, const std::tuple< Normals &, Normals & > &normals)
Do the actual projection.
Definition: projection_impl.hh:469
static constexpr unsigned dim
dimension of coordinates
Definition: projection.hh:47
unsigned numberOfEdgeIntersections() const
Number of edge intersections.
Definition: projection.hh:262
const std::tuple< Images, Preimages > & images() const
Images and preimages of corners.
Definition: projection.hh:235
const std::array< EdgeIntersection, maxEdgeIntersections > & edgeIntersections() const
Edge-edge intersections.
Definition: projection.hh:273
Intersection between two edges of a triangle.
Definition: projection.hh:29
std::array< Coordinate, 2 > local
Local coordinates of intersection and distance along normals.
Definition: projection.hh:41
std::array< unsigned, 2 > edge
Edge numbers in image and preimage triangle.
Definition: projection.hh:33
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 15, 23:04, 2025)