3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
staggered/freeflow/connectivitymap.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_STAGGERED_FREEFLOW_CONNECTIVITY_MAP_HH
25#define DUMUX_STAGGERED_FREEFLOW_CONNECTIVITY_MAP_HH
26
27#include <vector>
29
30namespace Dumux {
31
37template<class GridGeometry>
39{
40 using GridView = typename GridGeometry::GridView;
41 using FVElementGeometry = typename GridGeometry::LocalView;
42 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
43
44 using Element = typename GridView::template Codim<0>::Entity;
45 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
46
47 using CellCenterIdxType = typename GridGeometry::DofTypeIndices::CellCenterIdx;
48 using FaceIdxType = typename GridGeometry::DofTypeIndices::FaceIdx;
49
50 using SmallLocalIndex = typename IndexTraits<GridView>::SmallLocalIndex;
51
52 using Stencil = std::vector<GridIndexType>;
53 using Map = std::vector<Stencil>;
54
55 static constexpr SmallLocalIndex upwindSchemeOrder = GridGeometry::upwindSchemeOrder;
56 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
57
58public:
59
61 void update(const GridGeometry& gridGeometry)
62 {
63 const auto numDofsCC = gridGeometry.gridView().size(0);
64 const auto numDofsFace = gridGeometry.gridView().size(1);
65 const auto numBoundaryFacets = gridGeometry.numBoundaryScvf();
66 cellCenterToCellCenterMap_.resize(numDofsCC);
67 cellCenterToFaceMap_.resize(numDofsCC);
68 faceToCellCenterMap_.resize(2*numDofsFace - numBoundaryFacets);
69 faceToFaceMap_.resize(2*numDofsFace - numBoundaryFacets);
70
71 for(auto&& element: elements(gridGeometry.gridView()))
72 {
73 // restrict the FvGeometry locally and bind to the element
74 auto fvGeometry = localView(gridGeometry);
75 fvGeometry.bindElement(element);
76
77 // loop over sub control faces
78 for (auto&& scvf : scvfs(fvGeometry))
79 {
80 // handle the cell center dof stencils first
81 const auto dofIdxCellCenter = gridGeometry.elementMapper().index(element);
82
83 // the stencil for cell center dofs w.r.t. to other cell center dofs,
84 // includes all neighboring element indices
85 if (!scvf.boundary())
86 cellCenterToCellCenterMap_[dofIdxCellCenter].push_back(scvf.outsideScvIdx());
87
88 // the stencil for cell center dofs w.r.t. face dofs, includes the face dof indices of the current element
89 cellCenterToFaceMap_[dofIdxCellCenter].push_back(scvf.dofIndex());
90
91 // handle the face dof stencils
92 const auto scvfIdx = scvf.index();
93 computeFaceToCellCenterStencil_(faceToCellCenterMap_[scvfIdx], fvGeometry, scvf);
94 computeFaceToFaceStencil_(faceToFaceMap_[scvfIdx], fvGeometry, scvf);
95 }
96 }
97 }
98
100 const Stencil& operator() (CellCenterIdxType, CellCenterIdxType, const GridIndexType globalI) const
101 {
102 return cellCenterToCellCenterMap_[globalI];
103 }
104
106 const Stencil& operator() (CellCenterIdxType, FaceIdxType, const GridIndexType globalI) const
107 {
108 return cellCenterToFaceMap_[globalI];
109 }
110
112 const Stencil& operator() (FaceIdxType, CellCenterIdxType, const GridIndexType globalI) const
113 {
114 return faceToCellCenterMap_[globalI];
115 }
116
118 const Stencil& operator() (FaceIdxType, FaceIdxType, const GridIndexType globalI) const
119 {
120 return faceToFaceMap_[globalI];
121 }
122
123private:
124
125 /*
126 * \brief Computes the stencil for face dofs w.r.t to cell center dofs.
127 * Basically, these are the dof indices of the elements adjacent to the face and those of
128 * the elements adjacent to the faces parallel to the own face.
129 */
130 void computeFaceToCellCenterStencil_(Stencil& stencil,
131 const FVElementGeometry& fvGeometry,
132 const SubControlVolumeFace& scvf)
133 {
134 const auto eIdx = scvf.insideScvIdx();
135 stencil.push_back(eIdx);
136
137 for (const auto& data : scvf.pairData())
138 {
139 auto& lateralFace = fvGeometry.scvf(eIdx, data.localLateralFaceIdx);
140 if (!lateralFace.boundary())
141 {
142 const auto firstParallelElementDofIdx = lateralFace.outsideScvIdx();
143 stencil.push_back(firstParallelElementDofIdx);
144 }
145 }
146 }
147
148 /*
149 * \brief Computes the stencil for face dofs w.r.t to face dofs.
150 * For a full description of the stencil, please see the document under dumux/doc/docextra/staggered
151 */
152 void computeFaceToFaceStencil_(Stencil& stencil,
153 const FVElementGeometry& fvGeometry,
154 const SubControlVolumeFace& scvf)
155 {
156 stencil.push_back(scvf.axisData().oppositeDof);
157 addHigherOrderInAxisDofs_(scvf, stencil, std::integral_constant<bool, useHigherOrder>{});
158
159 for (const auto& data : scvf.pairData())
160 {
161 // add normal dofs
162 stencil.push_back(data.lateralPair.first);
163 if (!scvf.boundary())
164 stencil.push_back(data.lateralPair.second);
165
166 // add parallel dofs
167 for (SmallLocalIndex i = 0; i < upwindSchemeOrder; i++)
168 {
169 if (data.hasParallelNeighbor[i])
170 stencil.push_back(data.parallelDofs[i]);
171 }
172 }
173 }
174
175 void addHigherOrderInAxisDofs_(const SubControlVolumeFace& scvf, Stencil& stencil, std::false_type) {}
176
177 void addHigherOrderInAxisDofs_(const SubControlVolumeFace& scvf, Stencil& stencil, std::true_type)
178 {
179 for (SmallLocalIndex i = 0; i < upwindSchemeOrder - 1; i++)
180 {
181 if (scvf.hasBackwardNeighbor(i))
182 stencil.push_back(scvf.axisData().inAxisBackwardDofs[i]);
183
184 if (scvf.hasForwardNeighbor(i))
185 stencil.push_back(scvf.axisData().inAxisForwardDofs[i]);
186 }
187 }
188
189 Map cellCenterToCellCenterMap_;
190 Map cellCenterToFaceMap_;
191 Map faceToCellCenterMap_;
192 Map faceToFaceMap_;
193};
194
195} // end namespace Dumux
196
197#endif // DUMUX_STAGGERED_FREEFLOW_CONNECTIVITY_MAP_HH
Defines the index types used for grid and local indices.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:41
Stores the dof indices corresponding to the neighboring cell centers and faces that contribute to the...
Definition: staggered/freeflow/connectivitymap.hh:39
void update(const GridGeometry &gridGeometry)
Update the map and prepare the stencils.
Definition: staggered/freeflow/connectivitymap.hh:61
const Stencil & operator()(CellCenterIdxType, CellCenterIdxType, const GridIndexType globalI) const
Returns the stencil of a cell center dof w.r.t. other cell center dofs.
Definition: staggered/freeflow/connectivitymap.hh:100