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