DUNE-GRID-GLUE (2.10)

intersection.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// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LGPL-3.0-or-later OR LicenseRef-GPL-2.0-only-with-dune-grid-glue-exception
11#ifndef DUNE_GRIDGLUE_ADAPTER_INTERSECTION_HH
12#define DUNE_GRIDGLUE_ADAPTER_INTERSECTION_HH
13
14#include <algorithm>
15#include <optional>
16#include <tuple>
17
18#include <dune/common/deprecated.hh>
19#include <dune/common/version.hh>
20#include <dune/geometry/affinegeometry.hh>
21#include <dune/geometry/referenceelements.hh>
23
24#define ONLY_SIMPLEX_INTERSECTIONS
25
26namespace Dune {
27 namespace GridGlue {
28
29 // forward declaration
30 template<typename P0, typename P1>
31 class IntersectionIndexSet;
32
36 template<typename P0, typename P1>
38 {
39 public:
40 typedef ::Dune::GridGlue::GridGlue<P0, P1> GridGlue;
41
42 typedef typename GridGlue::IndexType IndexType;
43
45 static constexpr int coorddim = GridGlue::dimworld;
46
47 private:
48 // intermediate quantities
49 template<int side>
50 static constexpr int dim() { return GridGlue::template GridView<side>::Grid::dimension - GridGlue::template GridPatch<side>::codim; }
51
52 public:
54 static constexpr int mydim = dim<0>() < dim<1>() ? dim<0>() : dim<1>();
55
56 template<int side>
57 using GridLocalGeometry = AffineGeometry<
58 typename GridGlue::template GridView<side>::ctype, mydim, GridGlue::template GridView<side>::dimension>;
59
60 using Grid0LocalGeometry [[deprecated("please use GridLocalGeometry<0> instead")]] = GridLocalGeometry<0>;
61 using Grid1LocalGeometry [[deprecated("please use GridLocalGeometry<1> instead")]] = GridLocalGeometry<1>;
62
63 template<int side>
64 using GridGeometry = AffineGeometry<
65 typename GridGlue::template GridView<side>::ctype, mydim, GridGlue::template GridView<side>::dimensionworld>;
66
67 using Grid0Geometry [[deprecated("please use GridGeometry<0> instead")]] = GridGeometry<0>;
68 using Grid1Geometry [[deprecated("please use GridGeometry<1> instead")]] = GridGeometry<1>;
69
70 template<int side>
71 using GridIndexType = typename GridGlue::template GridView<side>::IndexSet::IndexType;
72
73 using Grid0IndexType [[deprecated("please use GridIndexType<0> instead")]] = GridIndexType<0>;
74 using Grid1IndexType [[deprecated("please use GridIndexType<1> instead")]] = GridIndexType<1>;
75
77 IntersectionData(const GridGlue& glue, unsigned int mergeindex, unsigned int offset, bool grid0local, bool grid1local);
78
80 IntersectionData() = default;
81
82 /* Accessor Functions */
83
84 template<int side>
85 const GridLocalGeometry<side>& localGeometry(unsigned int parentId = 0) const
86 { return *std::get<side>(sideData_).gridlocalgeom[parentId]; }
87
88 template<int side>
89 const GridGeometry<side>& geometry() const
90 { return *std::get<side>(sideData_).gridgeom; }
91
92 template<int side>
93 bool local() const
94 { return std::get<side>(sideData_).gridlocal; }
95
96 template<int side>
97 IndexType index(unsigned int parentId = 0) const
98 { return std::get<side>(sideData_).gridindices[parentId]; }
99
100 template<int side>
101 IndexType parents() const
102 { return std::get<side>(sideData_).gridindices.size(); }
103
104 private:
105 template<int side>
106 void initializeGeometry(const GridGlue& glue, unsigned mergeindex);
107
108 /* M E M B E R V A R I A B L E S */
109
110 public:
112 IndexType index_;
113
114 private:
115 template<int side>
116 struct SideData {
118 bool gridlocal = false;
119
121 std::vector< GridIndexType<side> > gridindices;
122
124 std::vector< std::optional< GridLocalGeometry<side> > > gridlocalgeom;
125
133 std::optional< GridGeometry<side> > gridgeom;
134 };
135
136 std::tuple< SideData<0>, SideData<1> > sideData_;
137 };
138
139 template<typename P0, typename P1>
140 template<int side>
141 void IntersectionData<P0, P1>::initializeGeometry(const GridGlue& glue, unsigned mergeindex)
142 {
143 auto& data = std::get<side>(sideData_);
144
145 const unsigned n_parents = glue.merger_->template parents<side>(mergeindex);
146
147 // init containers
148 data.gridindices.resize(n_parents);
149 data.gridlocalgeom.resize(n_parents);
150
151 // default values
152 data.gridindices[0] = 0;
153
154 static constexpr int nSimplexCorners = mydim + 1;
155 using ctype = typename GridGlue::ctype;
156
157 // initialize the local and the global geometries of grid `side`
158
159 // compute the coordinates of the subface's corners in codim 0 entity local coordinates
160 static constexpr int elementdim = GridGlue::template GridView<side>::template Codim<0>::Geometry::mydimension;
161
162 // coordinates within the subentity that contains the remote intersection
163 std::array<Dune::FieldVector< ctype, dim<side>() >, nSimplexCorners> corners_subEntity_local;
164
165 for (unsigned int par = 0; par < n_parents; ++par) {
166 for (int i = 0; i < nSimplexCorners; ++i)
167 corners_subEntity_local[i] = glue.merger_->template parentLocal<side>(mergeindex, i, par);
168
169 // Coordinates of the remote intersection corners wrt the element coordinate system
170 std::array<Dune::FieldVector<ctype, elementdim>, nSimplexCorners> corners_element_local;
171
172 if (data.gridlocal) {
173 data.gridindices[par] = glue.merger_->template parent<side>(mergeindex,par);
174
175 typename GridGlue::template GridPatch<side>::LocalGeometry
176 gridLocalGeometry = glue.template patch<side>().geometryLocal(data.gridindices[par]);
177 for (std::size_t i=0; i<corners_subEntity_local.size(); i++) {
178 corners_element_local[i] = gridLocalGeometry.global(corners_subEntity_local[i]);
179 }
180
181 // set the corners of the local geometry
182#ifdef ONLY_SIMPLEX_INTERSECTIONS
183 const Dune::GeometryType type = Dune::GeometryTypes::simplex(mydim);
184#else
185#error Not Implemented
186#endif
187 data.gridlocalgeom[par].emplace(type, corners_element_local);
188
189 // Add world geometry only for 0th parent
190 if (par == 0) {
191 typename GridGlue::template GridPatch<side>::Geometry
192 gridWorldGeometry = glue.template patch<side>().geometry(data.gridindices[par]);
193
194 // world coordinates of the remote intersection corners
195 std::array<Dune::FieldVector<ctype, GridGlue::template GridView<side>::dimensionworld>, nSimplexCorners> corners_global;
196
197 for (std::size_t i=0; i<corners_subEntity_local.size(); i++) {
198 corners_global[i] = gridWorldGeometry.global(corners_subEntity_local[i]);
199 }
200
201 data.gridgeom.emplace(type, corners_global);
202 }
203 }
204 }
205 }
206
208 template<typename P0, typename P1>
209 IntersectionData<P0, P1>::IntersectionData(const GridGlue& glue, unsigned int mergeindex, unsigned int offset,
210 bool grid0local, bool grid1local)
211 : index_(mergeindex+offset)
212 {
213 // if an invalid index is given do not proceed!
214 // (happens when the parent GridGlue initializes the "end"-Intersection)
215 assert (0 <= mergeindex || mergeindex < glue.index__sz);
216
217 std::get<0>(sideData_).gridlocal = grid0local;
218 std::get<1>(sideData_).gridlocal = grid1local;
219
220 initializeGeometry<0>(glue, mergeindex);
221 initializeGeometry<1>(glue, mergeindex);
222 }
223
228 template<typename P0, typename P1, int inside, int outside>
230 {
233
234 using InsideGridView = typename GridGlue::template GridView<inside>;
235 using OutsideGridView = typename GridGlue::template GridView<outside>;
236
237 using InsideLocalGeometry = typename IntersectionData::template GridLocalGeometry<inside>;
238 using OutsideLocalGeometry = typename IntersectionData::template GridLocalGeometry<outside>;
239
240 using Geometry = typename IntersectionData::template GridGeometry<inside>;
241 using OutsideGeometry = typename IntersectionData::template GridGeometry<outside>;
242
243 static constexpr auto coorddim = IntersectionData::coorddim;
244 static constexpr auto mydim = IntersectionData::mydim;
245 static constexpr int insidePatch = inside;
246 static constexpr int outsidePatch = outside;
247
248 using ctype = typename GridGlue::ctype;
249 using LocalCoordinate = Dune::FieldVector<ctype, mydim>;
250 using GlobalCoordinate = Dune::FieldVector<ctype, coorddim>;
251 };
252
255 template<typename P0, typename P1, int I, int O>
257 {
258
259 public:
260
262
263 typedef typename Traits::GridGlue GridGlue;
265
266
267 typedef typename Traits::InsideGridView InsideGridView;
268 typedef typename Traits::InsideLocalGeometry InsideLocalGeometry;
269
270 typedef typename Traits::OutsideGridView OutsideGridView;
271 typedef typename Traits::OutsideLocalGeometry OutsideLocalGeometry;
272 typedef typename Traits::OutsideGeometry OutsideGeometry;
273
274 typedef typename Traits::Geometry Geometry;
275 typedef typename Traits::ctype ctype;
276
277 typedef typename InsideGridView::Traits::template Codim<0>::Entity InsideEntity;
278 typedef typename OutsideGridView::Traits::template Codim<0>::Entity OutsideEntity;
279
280 typedef typename Traits::LocalCoordinate LocalCoordinate;
281 typedef typename Traits::GlobalCoordinate GlobalCoordinate;
282
284 static constexpr auto coorddim = Traits::coorddim;
285
287 static constexpr auto mydim = Traits::mydim;
288
290 static constexpr int insidePatch = Traits::insidePatch;
291
293 static constexpr int outsidePatch = Traits::outsidePatch;
294
295 // typedef unsigned int IndexType;
296
297 private:
301 static constexpr int codimensionWorld = coorddim - mydim;
302
303 public:
304 /* C O N S T R U C T O R S */
305
307 Intersection(const GridGlue* glue, const IntersectionData* i) :
308 glue_(glue), i_(i) {}
309
310 /* F U N C T I O N A L I T Y */
311
314 InsideEntity
315 inside(unsigned int parentId = 0) const
316 {
317 assert(self());
318 return glue_->template patch<I>().element(i_->template index<I>(parentId));
319 }
320
323 OutsideEntity
324 outside(unsigned int parentId = 0) const
325 {
326 assert(neighbor());
327 return glue_->template patch<O>().element(i_->template index<O>(parentId));
328 }
329
331 bool conforming() const
332 {
333 throw Dune::NotImplemented();
334 }
335
338 const InsideLocalGeometry& geometryInInside(unsigned int parentId = 0) const
339 {
340 return i_->template localGeometry<I>(parentId);
341 }
342
345 const OutsideLocalGeometry& geometryInOutside(unsigned int parentId = 0) const
346 {
347 return i_->template localGeometry<O>(parentId);
348 }
349
356 const Geometry& geometry() const
357 {
358 return i_->template geometry<I>();
359 }
360
367 const OutsideGeometry& geometryOutside() const // DUNE_DEPRECATED
368 {
369 return i_->template geometry<O>();
370 }
371
373 Dune::GeometryType type() const
374 {
375 #ifdef ONLY_SIMPLEX_INTERSECTIONS
376 return Dune::GeometryTypes::simplex(mydim);
377 #else
378 #error Not Implemented
379 #endif
380 }
381
382
384 bool self() const
385 {
386 return i_->template local<I>();
387 }
388
390 size_t neighbor(unsigned int g = 0) const
391 {
392 if (g == 0 && i_->template local<O>()) {
393 return i_->template parents<O>();
394 } else if (g == 1 && i_->template local<I>()) {
395 return i_->template parents<I>();
396 }
397 return 0;
398 }
399
401 int indexInInside(unsigned int parentId = 0) const
402 {
403 assert(self());
404 return glue_->template patch<I>().indexInInside(i_->template index<I>(parentId));
405 }
406
408 int indexInOutside(unsigned int parentId = 0) const
409 {
410 assert(neighbor());
411 return glue_->template patch<O>().indexInInside(i_->template index<O>(parentId));
412 }
413
418 GlobalCoordinate outerNormal(const LocalCoordinate &local) const
419 {
420 GlobalCoordinate normal;
421
422 if (codimensionWorld == 0)
423 DUNE_THROW(Dune::Exception, "There is no normal vector to a full-dimensional intersection");
424 else if (codimensionWorld == 1) {
425 /* TODO: Implement the general n-ary cross product here */
426 const auto jacobianTransposed = geometry().jacobianTransposed(local);
427 if (mydim==1) {
428 normal[0] = - jacobianTransposed[0][1];
429 normal[1] = jacobianTransposed[0][0];
430 } else if (mydim==2) {
431 normal[0] = (jacobianTransposed[0][1] * jacobianTransposed[1][2] - jacobianTransposed[0][2] * jacobianTransposed[1][1]);
432 normal[1] = - (jacobianTransposed[0][0] * jacobianTransposed[1][2] - jacobianTransposed[0][2] * jacobianTransposed[1][0]);
433 normal[2] = (jacobianTransposed[0][0] * jacobianTransposed[1][1] - jacobianTransposed[0][1] * jacobianTransposed[1][0]);
434 } else
435 DUNE_THROW(Dune::NotImplemented, "Remote intersections don't implement the 'outerNormal' method for " << mydim << "-dimensional intersections yet");
436 } else
437 DUNE_THROW(Dune::NotImplemented, "Remote intersections don't implement the 'outerNormal' method for intersections with codim >= 2 yet");
438
439 return normal;
440 }
441
446 GlobalCoordinate unitOuterNormal(const LocalCoordinate &local) const
447 {
448 GlobalCoordinate normal = outerNormal(local);
449 normal /= normal.two_norm();
450 return normal;
451 }
452
457 GlobalCoordinate integrationOuterNormal(const LocalCoordinate &local) const
458 {
459 return (unitOuterNormal(local) *= geometry().integrationElement(local));
460 }
461
466 GlobalCoordinate centerUnitOuterNormal () const
467 {
468 return unitOuterNormal(ReferenceElements<ctype,mydim>::general(type()).position(0,0));
469 }
470
475 {
476 return Intersection<P0,P1,O,I>(glue_,i_);
477 }
478
479#ifdef QUICKHACK_INDEX
480 typedef typename GridGlue::IndexType IndexType;
481
482 IndexType index() const
483 {
484 return i_->index_;
485 }
486
487#endif
488
489 private:
490
491 friend class IntersectionIndexSet<P0,P1>;
492
493 /* M E M B E R V A R I A B L E S */
494
496 const GridGlue* glue_;
497
499 const IntersectionData* i_;
500 };
501
502
503 } // end namespace GridGlue
504} // end namespace Dune
505
506#endif // DUNE_GRIDGLUE_ADAPTER_INTERSECTION_HH
sequential adapter to couple two grids at specified close together boundaries
Definition: gridglue.hh:67
unsigned int IndexType
Definition: gridglue.hh:147
PromotionTraits< typenameGridView< 0 >::ctype, typenameGridView< 1 >::ctype >::PromotedType ctype
The type used for coordinates.
Definition: gridglue.hh:171
static constexpr int dimworld
export the world dimension This is the maximum of the extractors' world dimensions.
Definition: gridglue.hh:166
storage class for Dune::GridGlue::Intersection related data
Definition: intersection.hh:38
static constexpr int coorddim
Dimension of the world space of the intersection.
Definition: intersection.hh:45
static constexpr int mydim
Dimension of the intersection.
Definition: intersection.hh:54
IndexType index_
index of this intersection after GridGlue interface
Definition: intersection.hh:112
IntersectionData()=default
Default Constructor.
The intersection of two entities of the two patches of a GridGlue.
Definition: intersection.hh:257
Intersection< P0, P1, O, I > flip() const
Return a copy of the intersection with inside and outside switched.
Definition: intersection.hh:474
bool conforming() const
Return true if intersection is conforming.
Definition: intersection.hh:331
int indexInOutside(unsigned int parentId=0) const
Local number of codim 1 entity in outside() Entity where intersection is contained in.
Definition: intersection.hh:408
int indexInInside(unsigned int parentId=0) const
Local number of codim 1 entity in the inside() Entity where intersection is contained in.
Definition: intersection.hh:401
bool self() const
For parallel computations: Return true if inside() entity exists locally.
Definition: intersection.hh:384
static constexpr auto mydim
dimension of the intersection
Definition: intersection.hh:287
static constexpr int outsidePatch
outside patch
Definition: intersection.hh:293
Dune::GeometryType type() const
Type of reference element for this intersection.
Definition: intersection.hh:373
InsideEntity inside(unsigned int parentId=0) const
Return element on the inside of this intersection.
Definition: intersection.hh:315
const OutsideLocalGeometry & geometryInOutside(unsigned int parentId=0) const
Geometric information about this intersection in local coordinates of the outside() element.
Definition: intersection.hh:345
GlobalCoordinate unitOuterNormal(const LocalCoordinate &local) const
Return a unit outer normal.
Definition: intersection.hh:446
size_t neighbor(unsigned int g=0) const
Return number of embeddings into local grid0 (grid1) entities.
Definition: intersection.hh:390
Intersection(const GridGlue *glue, const IntersectionData *i)
Constructor for a given Dataset.
Definition: intersection.hh:307
static constexpr int insidePatch
inside patch
Definition: intersection.hh:290
OutsideEntity outside(unsigned int parentId=0) const
Return element on the outside of this intersection.
Definition: intersection.hh:324
const Geometry & geometry() const
Geometric information about this intersection as part of the inside grid.
Definition: intersection.hh:356
GlobalCoordinate outerNormal(const LocalCoordinate &local) const
Return an outer normal (length not necessarily 1)
Definition: intersection.hh:418
const OutsideGeometry & geometryOutside() const
Geometric information about this intersection as part of the outside grid.
Definition: intersection.hh:367
GlobalCoordinate integrationOuterNormal(const LocalCoordinate &local) const
Return an outer normal with the length of the integration element.
Definition: intersection.hh:457
const InsideLocalGeometry & geometryInInside(unsigned int parentId=0) const
Geometric information about this intersection in local coordinates of the inside() element.
Definition: intersection.hh:338
GlobalCoordinate centerUnitOuterNormal() const
Unit outer normal at the center of the intersection.
Definition: intersection.hh:466
static constexpr auto coorddim
dimension of the world space of the intersection
Definition: intersection.hh:284
Central component of the module implementing the coupling of two grids.
Definition: intersection.hh:230
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 3, 22:46, 2025)