version 3.11-dev
discretization/cvfe/hybrid/gridfluxvariablescache.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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_HYBRID_CVFE_GRID_FLUXVARSCACHE_HH
13#define DUMUX_DISCRETIZATION_HYBRID_CVFE_GRID_FLUXVARSCACHE_HH
14
15#include <memory>
16#include <ranges>
17#include <unordered_map>
18
20
21// make the local view function available whenever we use this class
25
26namespace Dumux {
27
32template<class P, class FVC>
34{
35 using Problem = P;
36 using FluxVariablesCache = FVC;
37
38 template<class GridFluxVariablesCache, bool cachingEnabled>
40};
41
47template<class Problem,
48 class FluxVariablesCache,
49 bool cachingEnabled = false,
52
59template<class P, class FVC, class Traits>
60class HybridCVFEGridFluxVariablesCache<P, FVC, true, Traits>
61{
62 using Problem = typename Traits::Problem;
64
65public:
67 using FluxVariablesCache = typename Traits::FluxVariablesCache;
68
70 static constexpr bool cachingEnabled = true;
71
73 using LocalView = typename Traits::template LocalView<ThisType, cachingEnabled>;
74
75 HybridCVFEGridFluxVariablesCache(const Problem& problem) : problemPtr_(&problem) {}
76
77 template<class GridGeometry, class GridVolumeVariables, class SolutionVector>
78 void update(const GridGeometry& gridGeometry,
79 const GridVolumeVariables& gridVolVars,
80 const SolutionVector& sol,
81 bool forceUpdate = false)
82 {
83 // Here, we do not do anything unless it is a forced update
84 if (forceUpdate)
85 {
86 fluxVarsCache_.resize(gridGeometry.gridView().size(0));
87 scvfOffset_.resize(gridGeometry.gridView().size(0));
88 elementCache_.resize(gridGeometry.gridView().size(0));
89 boundaryIntersectionCache_.resize(gridGeometry.gridView().size(0));
90 Dumux::parallelFor(gridGeometry.gridView().size(0), [&, &problem = problem()](const std::size_t eIdx)
91 {
92 // Prepare the geometries within the elements of the stencil
93 const auto element = gridGeometry.element(eIdx);
94 const auto fvGeometry = localView(gridGeometry).bind(element);
95 const auto elemVolVars = localView(gridVolVars).bind(element, fvGeometry, sol);
96
97 // Build offset vector and resize flat cache
98 scvfOffset_[eIdx].resize(fvGeometry.numScvf() + 1, 0);
99 for (const auto& scvf : scvfs(fvGeometry))
100 {
101 const auto numQps = std::ranges::size(CVFE::quadratureRule(fvGeometry, scvf));
102 scvfOffset_[eIdx][scvf.index() + 1] = numQps;
103 }
104 for (std::size_t i = 2; i < scvfOffset_[eIdx].size(); ++i)
105 scvfOffset_[eIdx][i] += scvfOffset_[eIdx][i-1];
106
107 fluxVarsCache_[eIdx].resize(scvfOffset_[eIdx].back());
108 for (const auto& scvf : scvfs(fvGeometry))
109 {
110 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scvf))
111 cache(eIdx, qpData.ipData().scvfIndex(), qpData.ipData().qpIndex()).update(problem,
112 element,
113 fvGeometry,
114 elemVolVars,
115 qpData.ipData());
116 }
117
118 // Resize element cache based on element quadrature rule
119 const auto elemQuadRule = CVFE::quadratureRule(fvGeometry, element);
120 elementCache_[eIdx].resize(std::ranges::size(elemQuadRule));
121 for (const auto& qpData : elemQuadRule)
122 elementCache_[eIdx][qpData.ipData().qpIndex()].update(problem, element, fvGeometry, elemVolVars, qpData.ipData());
123
124 // Rebuild boundary intersection cache for this element
125 if (!boundaryIntersectionCache_[eIdx])
126 boundaryIntersectionCache_[eIdx] = std::make_unique<std::unordered_map<int, std::vector<FluxVariablesCache>>>();
127 else
128 boundaryIntersectionCache_[eIdx]->clear();
129
130 for (const auto& intersection : intersections(gridGeometry.gridView(), element))
131 {
132 if (intersection.boundary())
133 {
134 const auto quadRule = CVFE::quadratureRule(fvGeometry, intersection);
135 const auto iIdx = intersection.indexInInside();
136 (*boundaryIntersectionCache_[eIdx])[iIdx].resize(std::ranges::size(quadRule));
137
138 for (const auto& qpData : quadRule)
139 (*boundaryIntersectionCache_[eIdx])[iIdx][qpData.ipData().qpIndex()].update(problem,
140 element,
141 fvGeometry,
142 elemVolVars,
143 qpData.ipData());
144 }
145 }
146 });
147 }
148 }
149
150 template<class FVElementGeometry, class ElementVolumeVariables>
151 void updateElement(const typename FVElementGeometry::Element& element,
152 const FVElementGeometry& fvGeometry,
153 const ElementVolumeVariables& elemVolVars)
154 {
155 if constexpr (FluxVariablesCache::isSolDependent)
156 {
157 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
158
159 // Build offset vector and resize flat cache
160 scvfOffset_[eIdx].resize(fvGeometry.numScvf() + 1, 0);
161 for (const auto& scvf : scvfs(fvGeometry))
162 {
163 const auto numQps = std::ranges::size(CVFE::quadratureRule(fvGeometry, scvf));
164 scvfOffset_[eIdx][scvf.index() + 1] = numQps;
165 }
166 for (std::size_t i = 2; i < scvfOffset_[eIdx].size(); ++i)
167 scvfOffset_[eIdx][i] += scvfOffset_[eIdx][i-1];
168
169 fluxVarsCache_[eIdx].resize(scvfOffset_[eIdx].back());
170 for (const auto& scvf : scvfs(fvGeometry))
171 {
172 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scvf))
173 cache(eIdx, qpData.ipData().scvfIndex(), qpData.ipData().qpIndex()).update(problem(),
174 element,
175 fvGeometry,
176 elemVolVars,
177 qpData.ipData());
178 }
179
180 // Resize element cache based on element quadrature rule
181 const auto elemQuadRule = CVFE::quadratureRule(fvGeometry, element);
182 elementCache_[eIdx].resize(std::ranges::size(elemQuadRule));
183 for (const auto& qpData : elemQuadRule)
184 elementCache_[eIdx][qpData.ipData().qpIndex()].update(problem(), element, fvGeometry, elemVolVars, qpData.ipData());
185
186 // Rebuild boundary intersection cache for this element
187 if (!boundaryIntersectionCache_[eIdx])
188 boundaryIntersectionCache_[eIdx] = std::make_unique<std::unordered_map<int, std::vector<FluxVariablesCache>>>();
189 else
190 boundaryIntersectionCache_[eIdx]->clear();
191
192 for (const auto& intersection : intersections(fvGeometry.gridGeometry().gridView(), element))
193 {
194 if (intersection.boundary())
195 {
196 const auto quadRule = CVFE::quadratureRule(fvGeometry, intersection);
197 const auto iIdx = intersection.indexInInside();
198 (*boundaryIntersectionCache_[eIdx])[iIdx].resize(std::ranges::size(quadRule));
199
200 for (const auto& qpData : quadRule)
201 (*boundaryIntersectionCache_[eIdx])[iIdx][qpData.ipData().qpIndex()].update(problem(),
202 element,
203 fvGeometry,
204 elemVolVars,
205 qpData.ipData());
206 }
207 }
208 }
209 }
210
211 const Problem& problem() const
212 { return *problemPtr_; }
213
214 // access operator
215 const FluxVariablesCache& cache(std::size_t eIdx, std::size_t scvfIdx, std::size_t qpIdx) const
216 { return fluxVarsCache_[eIdx][scvfOffset_[eIdx][scvfIdx] + qpIdx]; }
217
218 // access operator
219 FluxVariablesCache& cache(std::size_t eIdx, std::size_t scvfIdx, std::size_t qpIdx)
220 { return fluxVarsCache_[eIdx][scvfOffset_[eIdx][scvfIdx] + qpIdx]; }
221
222 // access operator
223 const FluxVariablesCache& elementCache(std::size_t eIdx, std::size_t qpIdx) const { return elementCache_[eIdx][qpIdx]; }
224 FluxVariablesCache& elementCache(std::size_t eIdx, std::size_t qpIdx) { return elementCache_[eIdx][qpIdx]; }
225
226 // access operator for boundary intersection cache
227 FluxVariablesCache& boundaryIntersectionCache(std::size_t eIdx, int iIdx, std::size_t qpIdx) { return (*boundaryIntersectionCache_[eIdx])[iIdx][qpIdx]; }
228 const FluxVariablesCache& boundaryIntersectionCache(std::size_t eIdx, int iIdx, std::size_t qpIdx) const { return (*boundaryIntersectionCache_[eIdx]).at(iIdx)[qpIdx]; }
229
230private:
231 // currently bound element
232 const Problem* problemPtr_;
233 std::vector<std::vector<FluxVariablesCache>> fluxVarsCache_;
234 std::vector<std::vector<std::size_t>> scvfOffset_;
235 std::vector<std::vector<FluxVariablesCache>> elementCache_;
236 std::vector<std::unique_ptr<std::unordered_map<int, std::vector<FluxVariablesCache>>>> boundaryIntersectionCache_;
237};
238
243template<class P, class FVC, class Traits>
244class HybridCVFEGridFluxVariablesCache<P, FVC, false, Traits>
245{
246 using Problem = typename Traits::Problem;
248
249public:
251 using FluxVariablesCache = typename Traits::FluxVariablesCache;
252
254 static constexpr bool cachingEnabled = false;
255
257 using LocalView = typename Traits::template LocalView<ThisType, cachingEnabled>;
258
259 HybridCVFEGridFluxVariablesCache(const Problem& problem) : problemPtr_(&problem) {}
260
261 template<class GridGeometry, class GridVolumeVariables, class SolutionVector>
262 void update(const GridGeometry& gridGeometry,
263 const GridVolumeVariables& gridVolVars,
264 const SolutionVector& sol,
265 bool forceUpdate = false) {}
266
267 const Problem& problem() const
268 { return *problemPtr_; }
269
270private:
271 const Problem* problemPtr_;
272};
273
274} // end namespace Dumux
275
276#endif
The flux variables caches for an element when using hybrid CVFE discretizations.
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:38
Flux variable caches on a gridview with grid caching disabled.
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:245
const Problem & problem() const
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:267
void update(const GridGeometry &gridGeometry, const GridVolumeVariables &gridVolVars, const SolutionVector &sol, bool forceUpdate=false)
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:262
typename Traits::FluxVariablesCache FluxVariablesCache
export the flux variable cache type
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:251
typename Traits::template LocalView< ThisType, cachingEnabled > LocalView
export the type of the local view
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:257
HybridCVFEGridFluxVariablesCache(const Problem &problem)
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:259
Flux variable caches on a gridview with grid caching enabled (general quadrature specialization)
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:61
const FluxVariablesCache & elementCache(std::size_t eIdx, std::size_t qpIdx) const
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:223
const Problem & problem() const
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:211
FluxVariablesCache & elementCache(std::size_t eIdx, std::size_t qpIdx)
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:224
HybridCVFEGridFluxVariablesCache(const Problem &problem)
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:75
FluxVariablesCache & boundaryIntersectionCache(std::size_t eIdx, int iIdx, std::size_t qpIdx)
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:227
const FluxVariablesCache & boundaryIntersectionCache(std::size_t eIdx, int iIdx, std::size_t qpIdx) const
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:228
FluxVariablesCache & cache(std::size_t eIdx, std::size_t scvfIdx, std::size_t qpIdx)
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:219
const FluxVariablesCache & cache(std::size_t eIdx, std::size_t scvfIdx, std::size_t qpIdx) const
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:215
typename Traits::template LocalView< ThisType, cachingEnabled > LocalView
export the type of the local view
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:73
void updateElement(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars)
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:151
typename Traits::FluxVariablesCache FluxVariablesCache
export the flux variable cache type
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:67
void update(const GridGeometry &gridGeometry, const GridVolumeVariables &gridVolVars, const SolutionVector &sol, bool forceUpdate=false)
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:78
Flux variable caches implementation on a gridview.
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:51
Global flux variable cache.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:160
Free function to get the local view of a grid cache object.
auto quadratureRule(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolume &scv, QuadratureRules::MidpointQuadrature)
Midpoint quadrature for scv.
Definition: quadraturerules.hh:148
Definition: adapt.hh:17
Parallel for loop (multithreading)
Quadrature rules over sub-control volumes and sub-control volume faces.
Flux variable caches traits.
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:34
FVC FluxVariablesCache
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:36
P Problem
Definition: discretization/cvfe/hybrid/gridfluxvariablescache.hh:35