3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
localintersectionindexmapper.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_LOCAL_INTERSECTION_INDEX_MAPPER_HH
25#define DUMUX_DISCRETIZATION_LOCAL_INTERSECTION_INDEX_MAPPER_HH
26
27#include <array>
28#include <numeric>
29#include <dune/common/float_cmp.hh>
31#include <dumux/common/math.hh>
35
36namespace Dumux {
37
38namespace Detail {
39
40template<class GridView, bool consistentlyOrientedGrid>
42
49template<class GridView>
51{
52 using SmallLocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex;
53 using Element = typename GridView::template Codim<0>::Entity;
54 static constexpr auto numElementFaces = GridView::Grid::dimension * 2;
55public:
56 void update(const GridView& gv, const Element& element)
57 {
58 static const bool makeConsistentlyOriented = getParam<bool>("Grid.MakeConsistentlyOriented", true);
59 if (!makeConsistentlyOriented)
60 {
61 std::iota(realToRefMap_.begin(), realToRefMap_.end(), 0);
62 refToRealMap_ = realToRefMap_;
63 return;
64 }
65
66 for (const auto& is : intersections(gv, element))
67 {
68 const auto& otherOuterNormal = is.centerUnitOuterNormal();
69 const auto idx = normalAxis(otherOuterNormal);
70 const int positveOrientation = !std::signbit(otherOuterNormal[idx]);
71 const auto refIdx = idx * 2 + positveOrientation;
72 const auto realIdx = is.indexInInside();
73 realToRefMap_[realIdx] = refIdx;
74 refToRealMap_[refIdx] = realIdx;
75 }
76 }
77
79 SmallLocalIndexType realToRefIdx(const SmallLocalIndexType localIsIdx) const
80 { return realToRefMap_[localIsIdx]; }
81
83 SmallLocalIndexType refToRealIdx(const SmallLocalIndexType localIsIdx) const
84 { return refToRealMap_[localIsIdx]; }
85
86private:
87 std::array<SmallLocalIndexType, numElementFaces> realToRefMap_ = {};
88 std::array<SmallLocalIndexType, numElementFaces> refToRealMap_ = {};
89};
90
99template<class GridView>
101{
102 using SmallLocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex;
103 using Element = typename GridView::template Codim<0>::Entity;
104public:
105
107 void update(const GridView&, const Element&) {}
108
110 SmallLocalIndexType realToRefIdx(const SmallLocalIndexType localIsIdx) const
111 { return localIsIdx; }
112
114 SmallLocalIndexType refToRealIdx(const SmallLocalIndexType localIsIdx) const
115 { return localIsIdx; }
116};
117
118} // end namespace Detail
119
127template<class GridView>
130
131} // end namespace Dumux
132
133#endif
Defines the index types used for grid and local indices.
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Returns the normal axis index of a unit vector (0 = x, 1 = y, 2 = z)
Definition: adapt.hh:29
static std::size_t normalAxis(const Vector &v)
Definition: normalaxis.hh:39
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:41
Definition: localintersectionindexmapper.hh:41
SmallLocalIndexType refToRealIdx(const SmallLocalIndexType localIsIdx) const
Return the intersection's local reference indexInElement given an actual local index.
Definition: localintersectionindexmapper.hh:83
void update(const GridView &gv, const Element &element)
Definition: localintersectionindexmapper.hh:56
SmallLocalIndexType realToRefIdx(const SmallLocalIndexType localIsIdx) const
Return the intersection's actual local indexInElement given a local reference index.
Definition: localintersectionindexmapper.hh:79
SmallLocalIndexType realToRefIdx(const SmallLocalIndexType localIsIdx) const
Return the intersection's actual local indexInElement given a local reference index.
Definition: localintersectionindexmapper.hh:110
SmallLocalIndexType refToRealIdx(const SmallLocalIndexType localIsIdx) const
Return the intersection's local reference indexInElement given an actual local index.
Definition: localintersectionindexmapper.hh:114
void update(const GridView &, const Element &)
Update the map for getting the corresponding local face indices in another element.
Definition: localintersectionindexmapper.hh:107