version 3.10-dev
staggered/freeflow/elementvolumevariables.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_DISCRETIZATION_STAGGERED_ELEMENT_VOLUMEVARIABLES_HH
13#define DUMUX_DISCRETIZATION_STAGGERED_ELEMENT_VOLUMEVARIABLES_HH
14
15#include <algorithm>
16#include <cassert>
17#include <vector>
18#include <utility>
19
20#include <dune/common/exceptions.hh>
23
24namespace Dumux {
25
30template<class GVV, bool cachingEnabled>
32{};
33
39template<class GVV>
40class StaggeredElementVolumeVariables<GVV, /*cachingEnabled*/true>
41{
42public:
45
47 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
48
51 : gridVolVarsPtr_(&gridVolVars)
52 , numScv_(gridVolVars.problem().gridGeometry().numScv())
53 {}
54
56 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
57 const VolumeVariables& operator [](const SubControlVolume& scv) const
58 {
59 if (scv.dofIndex() < numScv_)
60 return gridVolVars().volVars(scv.dofIndex());
61 else
62 return boundaryVolumeVariables_[getLocalIdx_(scv.dofIndex())];
63 }
64
67 const VolumeVariables& operator [](const std::size_t scvIdx) const
68 {
69 if (scvIdx < numScv_)
70 return gridVolVars().volVars(scvIdx);
71 else
72 return boundaryVolumeVariables_[getLocalIdx_(scvIdx)];
73 }
74
80 template<class FVElementGeometry, class SolutionVector>
81 StaggeredElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
82 const FVElementGeometry& fvGeometry,
83 const SolutionVector& sol) &&
84 {
85 this->bind_(element, fvGeometry, sol);
86 return std::move(*this);
87 }
88
89 template<class FVElementGeometry, class SolutionVector>
90 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
91 const FVElementGeometry& fvGeometry,
92 const SolutionVector& sol) &
93 { this->bind_(element, fvGeometry, sol); }
94
100 template<class FVElementGeometry, class SolutionVector>
101 StaggeredElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
102 const FVElementGeometry& fvGeometry,
103 const SolutionVector& sol) &&
104 {
105 this->bindElement_(element, fvGeometry, sol);
106 return std::move(*this);
107 }
108
109 template<class FVElementGeometry, class SolutionVector>
110 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
111 const FVElementGeometry& fvGeometry,
112 const SolutionVector& sol) &
113 { this->bindElement_(element, fvGeometry, sol); }
114
117 { return *gridVolVarsPtr_; }
118
119private:
120
122 void clear_()
123 {
124 boundaryVolVarIndices_.clear();
125 boundaryVolumeVariables_.clear();
126 }
127
130 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<isMultiTypeBlockVector<SolutionVector>::value, int> = 0>
131 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
132 const FVElementGeometry& fvGeometry,
133 const SolutionVector& sol)
134 {
135 // forward to the actual method
136 bind_(element, fvGeometry, sol[FVElementGeometry::GridGeometry::cellCenterIdx()]);
137 }
138
141 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<!isMultiTypeBlockVector<SolutionVector>::value, int> = 0>
142 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
143 const FVElementGeometry& fvGeometry,
144 const SolutionVector& sol)
145 {
146 if (!fvGeometry.hasBoundaryScvf())
147 return;
148
149 clear_();
150 boundaryVolVarIndices_.reserve(fvGeometry.numScvf());
151 boundaryVolumeVariables_.reserve(fvGeometry.numScvf());
152
153 // handle the boundary volume variables
154 for (auto&& scvf : scvfs(fvGeometry))
155 {
156 // if we are not on a boundary, skip the rest
157 if (!scvf.boundary())
158 continue;
159
160 const auto& problem = gridVolVars().problem();
161 auto boundaryPriVars = gridVolVars().getBoundaryPriVars(problem, sol, element, scvf);
162 const auto elemSol = elementSolution<FVElementGeometry>(std::move(boundaryPriVars));
163 auto&& scvI = fvGeometry.scv(scvf.insideScvIdx());
164
165 VolumeVariables volVars;
166 volVars.update(elemSol,
167 problem,
168 element,
169 scvI);
170
171 boundaryVolumeVariables_.emplace_back(std::move(volVars));
172 boundaryVolVarIndices_.push_back(scvf.outsideScvIdx());
173 }
174 }
175
177 template<class FVElementGeometry, class SolutionVector>
178 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
179 const FVElementGeometry& fvGeometry,
180 const SolutionVector& sol)
181 {}
182
183 const GridVolumeVariables* gridVolVarsPtr_;
184
186 int getLocalIdx_(const int volVarIdx) const
187 {
188 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
189 assert(it != boundaryVolVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
190 return std::distance(boundaryVolVarIndices_.begin(), it);
191 }
192
193 std::vector<std::size_t> boundaryVolVarIndices_;
194 std::vector<VolumeVariables> boundaryVolumeVariables_;
195 const std::size_t numScv_;
196};
197
198
204template<class GVV>
205class StaggeredElementVolumeVariables<GVV, /*cachingEnabled*/false>
206{
207 using PrimaryVariables = typename GVV::VolumeVariables::PrimaryVariables;
208
209public:
212
214 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
215
218 : gridVolVarsPtr_(&gridVolVars) {}
219
225 template<class FVElementGeometry, class SolutionVector>
226 StaggeredElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
227 const FVElementGeometry& fvGeometry,
228 const SolutionVector& sol) &&
229 {
230 this->bind_(element, fvGeometry, sol);
231 return std::move(*this);
232 }
233
234 template<class FVElementGeometry, class SolutionVector>
235 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
236 const FVElementGeometry& fvGeometry,
237 const SolutionVector& sol) &
238 { this->bind_(element, fvGeometry, sol); }
239
245 template<class FVElementGeometry, class SolutionVector>
246 StaggeredElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
247 const FVElementGeometry& fvGeometry,
248 const SolutionVector& sol) &&
249 {
250 this->bindElement_(element, fvGeometry, sol);
251 return std::move(*this);
252 }
253
254 template<class FVElementGeometry, class SolutionVector>
255 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
256 const FVElementGeometry& fvGeometry,
257 const SolutionVector& sol) &
258 { this->bindElement_(element, fvGeometry, sol); }
259
261 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
262 const VolumeVariables& operator [](const SubControlVolume& scv) const
263 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
264
266 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
267 VolumeVariables& operator [](const SubControlVolume& scv)
268 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
269
271 const VolumeVariables& operator [](std::size_t scvIdx) const
272 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
273
275 VolumeVariables& operator [](std::size_t scvIdx)
276 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
277
280 { return *gridVolVarsPtr_; }
281
282private:
285 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<isMultiTypeBlockVector<SolutionVector>::value, int> = 0>
286 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
287 const FVElementGeometry& fvGeometry,
288 const SolutionVector& sol)
289 {
290 // forward to the actual method
291 bind_(element, fvGeometry, sol[FVElementGeometry::GridGeometry::cellCenterIdx()]);
292 }
293
296 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<!isMultiTypeBlockVector<SolutionVector>::value, int> = 0>
297 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
298 const FVElementGeometry& fvGeometry,
299 const SolutionVector& sol)
300 {
301 clear_();
302
303 const auto& problem = gridVolVars().problem();
304 const auto& gridGeometry = fvGeometry.gridGeometry();
305 const auto globalI = gridGeometry.elementMapper().index(element);
306 const auto& map = gridGeometry.connectivityMap();
307 constexpr auto cellCenterIdx = FVElementGeometry::GridGeometry::cellCenterIdx();
308 const auto& connectivityMapI = map(cellCenterIdx, cellCenterIdx, globalI);
309 const auto numDofs = connectivityMapI.size();
310
311 auto&& scvI = fvGeometry.scv(globalI);
312
313 // resize local containers to the required size (for internal elements)
314 volumeVariables_.resize(numDofs+1);
315 volVarIndices_.resize(numDofs+1);
316 int localIdx = 0;
317
318 // Lambda to update the volume variables of the given index
319 auto doVolVarUpdate = [&](int globalJ)
320 {
321 const auto& elementJ = gridGeometry.element(globalJ);
322 auto&& scvJ = fvGeometry.scv(globalJ);
323 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PrimaryVariables>(sol[globalJ]);
324 volumeVariables_[localIdx].update(elemSol,
325 problem,
326 elementJ,
327 scvJ);
328 volVarIndices_[localIdx] = scvJ.dofIndex();
329 ++localIdx;
330 };
331
332 // Update the volume variables of the element at hand
333 doVolVarUpdate(globalI);
334
335 // Update the volume variables of the neighboring elements
336 for (const auto& globalJ : connectivityMapI)
337 doVolVarUpdate(globalJ);
338
339 if (fvGeometry.hasBoundaryScvf())
340 {
341 // Update boundary volume variables
342 for (auto&& scvf : scvfs(fvGeometry))
343 {
344 // if we are not on a boundary, skip to the next scvf
345 if (!scvf.boundary())
346 continue;
347
348 volumeVariables_.resize(localIdx+1);
349 volVarIndices_.resize(localIdx+1);
350
351 auto boundaryPriVars = gridVolVars().getBoundaryPriVars(problem, sol, element, scvf);
352 auto elemSol = elementSolution<FVElementGeometry>(std::move(boundaryPriVars));
353 volumeVariables_[localIdx].update(elemSol,
354 problem,
355 element,
356 scvI);
357 volVarIndices_[localIdx] = scvf.outsideScvIdx();
358 ++localIdx;
359 }
360 }
361 }
362
365 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<isMultiTypeBlockVector<SolutionVector>::value, int> = 0>
366 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
367 const FVElementGeometry& fvGeometry,
368 const SolutionVector& sol)
369 {
370 // forward to the actual method
371 bindElement_(element, fvGeometry, sol[FVElementGeometry::GridGeometry::cellCenterIdx()]);
372 }
373
376 template<class FVElementGeometry, class SolutionVector, typename std::enable_if_t<!isMultiTypeBlockVector<SolutionVector>::value, int> = 0>
377 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
378 const FVElementGeometry& fvGeometry,
379 const SolutionVector& sol)
380 {
381 clear_();
382
383 const auto globalI = fvGeometry.gridGeometry().elementMapper().index(element);
384 volumeVariables_.resize(1);
385 volVarIndices_.resize(1);
386
387 // update the volume variables of the element
388 auto&& scv = fvGeometry.scv(globalI);
389
390 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PrimaryVariables>(sol[globalI]);
391 volumeVariables_[0].update(elemSol,
392 gridVolVars().problem(),
393 element,
394 scv);
395 volVarIndices_[0] = scv.dofIndex();
396 }
397
399 void clear_()
400 {
401 volVarIndices_.clear();
402 volumeVariables_.clear();
403 }
404
405 const GridVolumeVariables* gridVolVarsPtr_;
406
407 int getLocalIdx_(const int volVarIdx) const
408 {
409 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
410 assert(it != volVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
411 return std::distance(volVarIndices_.begin(), it);
412 }
413
414 std::vector<std::size_t> volVarIndices_;
415 std::vector<VolumeVariables> volumeVariables_;
416};
417
418} // end namespace Dumux
419
420#endif
GVV GridVolumeVariables
export type of the grid volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:211
StaggeredElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: staggered/freeflow/elementvolumevariables.hh:226
StaggeredElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: staggered/freeflow/elementvolumevariables.hh:246
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: staggered/freeflow/elementvolumevariables.hh:235
StaggeredElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: staggered/freeflow/elementvolumevariables.hh:217
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: staggered/freeflow/elementvolumevariables.hh:255
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:214
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: staggered/freeflow/elementvolumevariables.hh:279
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:47
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: staggered/freeflow/elementvolumevariables.hh:116
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: staggered/freeflow/elementvolumevariables.hh:110
GVV GridVolumeVariables
export type of the grid volume variables
Definition: staggered/freeflow/elementvolumevariables.hh:44
StaggeredElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: staggered/freeflow/elementvolumevariables.hh:101
StaggeredElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: staggered/freeflow/elementvolumevariables.hh:50
StaggeredElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: staggered/freeflow/elementvolumevariables.hh:81
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: staggered/freeflow/elementvolumevariables.hh:90
Base class for the element volume variables vector for the staggered model.
Definition: staggered/freeflow/elementvolumevariables.hh:32
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:282
Definition: adapt.hh:17
The local element solution class for staggered methods.
Type traits to be used with vector types.