version 3.10-dev
facecentered/staggered/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_FACECENTERED_ELEMENTVOLUMEVARIABLES_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_ELEMENTVOLUMEVARIABLES_HH
14
15#include <algorithm>
16#include <cassert>
17#include <vector>
18#include <utility>
19
21
22namespace Dumux {
23
28template<class GridVolumeVariables, bool cachingEnabled>
30
35template<class GVV>
36class FaceCenteredStaggeredElementVolumeVariables<GVV, /*cachingEnabled*/true>
37{
38 using ThisType = FaceCenteredStaggeredElementVolumeVariables<GVV, /*cachingEnabled*/true>;
39 using GridGeometry = std::decay_t<decltype(std::declval<GVV>().problem().gridGeometry())>;
40 using FVElementGeometry = typename GridGeometry::LocalView;
41 using SubControlVolume = typename GridGeometry::SubControlVolume;
42
43public:
46
48 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
49
51 : gridVolumeVariablesPtr_(&gridVolumeVariables)
52 , numScv_(gridVolumeVariables.problem().gridGeometry().numScv())
53 {}
54
56 const VolumeVariables& operator [](const SubControlVolume& scv) const
57 {
58 if (scv.index() < numScv_)
59 return gridVolVars().volVars(scv.index());
60 else
61 return boundaryVolumeVariables_[getLocalIdx_(scv.index())];
62 }
63
66 const VolumeVariables& operator [](const std::size_t scvIdx) const
67 {
68 if (scvIdx < numScv_)
69 return gridVolVars().volVars(scvIdx);
70 else
71 return boundaryVolumeVariables_[getLocalIdx_(scvIdx)];
72 }
73
79 template<class SolutionVector>
80 ThisType bind(const typename FVElementGeometry::Element& element,
81 const FVElementGeometry& fvGeometry,
82 const SolutionVector& sol) &&
83 {
84 this->bind(element, fvGeometry, sol);
85 return std::move(*this);
86 }
87
90 template<class SolutionVector>
91 void bind(const typename FVElementGeometry::Element& element,
92 const FVElementGeometry& fvGeometry,
93 const SolutionVector& sol) &
94 {
95 if (!fvGeometry.hasBoundaryScvf())
96 return;
97
98 clear_();
99 // upper bound of size is the number of all scvfs minus frontal scvfs
100 boundaryVolVarIndices_.reserve(fvGeometry.numScvf()-element.subEntities(1));
101 boundaryVolumeVariables_.reserve(fvGeometry.numScvf()-element.subEntities(1));
102
103 for (const auto& scvf : scvfs(fvGeometry))
104 {
105 if (!scvf.boundary() || scvf.isFrontal() || scvf.processorBoundary())
106 continue;
107
108 // check if boundary is a pure dirichlet boundary
109 const auto& problem = gridVolVars().problem();
110 const auto bcTypes = problem.boundaryTypes(element, scvf);
111
112 auto addBoundaryVolVars = [&](const auto& scvFace)
113 {
114 const auto& scvI = fvGeometry.scv(scvFace.insideScvIdx());
115 typename VolumeVariables::PrimaryVariables pv(
116 problem.dirichlet(element, scvFace)[scvI.dofAxis()]
117 );
118 const auto dirichletPriVars = elementSolution<FVElementGeometry>(pv);
119
120 VolumeVariables volVars;
121 volVars.update(dirichletPriVars, problem, element, scvI);
122
123 boundaryVolumeVariables_.emplace_back(std::move(volVars));
124 boundaryVolVarIndices_.push_back(scvFace.outsideScvIdx());
125 };
126
127 if (bcTypes.hasDirichlet())
128 {
129 addBoundaryVolVars(scvf);
130 continue;
131 }
132
133 // treat domain corners
134 if (const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf); orthogonalScvf.boundary())
135 if (problem.boundaryTypes(element, orthogonalScvf).hasDirichlet())
136 addBoundaryVolVars(scvf);
137
138 }
139
140 assert(boundaryVolumeVariables_.size() == boundaryVolVarIndices_.size());
141 }
142
148 template<class SolutionVector>
149 ThisType bindElement(const typename FVElementGeometry::Element& element,
150 const FVElementGeometry& fvGeometry,
151 const SolutionVector& sol) &&
152 {
153 this->bindElement(element, fvGeometry, sol);
154 return std::move(*this);
155 }
158 template<class SolutionVector>
159 void bindElement(const typename FVElementGeometry::Element& element,
160 const FVElementGeometry& fvGeometry,
161 const SolutionVector& sol) &
162 {}
163
164
167 { return *gridVolumeVariablesPtr_; }
168
170 bool hasVolVars(const std::size_t scvIdx) const
171 {
172 if (scvIdx < numScv_)
173 return true;
174 else
175 {
176 const auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), scvIdx);
177 return it != boundaryVolVarIndices_.end();
178 }
179 }
180
181private:
183 void clear_()
184 {
185 boundaryVolVarIndices_.clear();
186 boundaryVolumeVariables_.clear();
187 }
188
190 int getLocalIdx_(const std::size_t volVarIdx) const
191 {
192 const auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
193 assert(it != boundaryVolVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
194 return std::distance(boundaryVolVarIndices_.begin(), it);
195 }
196
197 const GridVolumeVariables* gridVolumeVariablesPtr_;
198 const std::size_t numScv_;
199 std::vector<std::size_t> boundaryVolVarIndices_;
200 std::vector<VolumeVariables> boundaryVolumeVariables_;
201};
202
207template<class GVV>
208class FaceCenteredStaggeredElementVolumeVariables<GVV, /*cachingEnabled*/false>
209{
210 using ThisType = FaceCenteredStaggeredElementVolumeVariables<GVV, /*cachingEnabled*/false>;
211 using GridGeometry = std::decay_t<decltype(std::declval<GVV>().problem().gridGeometry())>;
212 using FVElementGeometry = typename GridGeometry::LocalView;
213 using SubControlVolume = typename GridGeometry::SubControlVolume;
214
215 static constexpr auto dim = GridGeometry::GridView::dimension;
216 static constexpr auto numInsideVolVars = dim * 2;
217 static constexpr auto numOutsideVolVars = numInsideVolVars * 2 * (dim - 1);
218
219public:
222
224 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
225
227 : gridVolumeVariablesPtr_(&globalFacesVars) {}
228
230 const VolumeVariables& operator [](const SubControlVolume& scv) const
231 { return volumeVariables_[getLocalIdx_(scv.index())]; }
232
234 const VolumeVariables& operator [](const std::size_t scvIdx) const
235 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
236
238 VolumeVariables& operator [](const SubControlVolume& scv)
239 { return volumeVariables_[getLocalIdx_(scv.index())]; }
240
241 // operator for the access with an index
242 VolumeVariables& operator [](const std::size_t scvIdx)
243 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
244
250 template<class SolutionVector>
251 ThisType bind(const typename FVElementGeometry::Element& element,
252 const FVElementGeometry& fvGeometry,
253 const SolutionVector& sol) &&
254 {
255 this->bind_(element, fvGeometry, sol);
256 return std::move(*this);
257 }
258
259 template<class SolutionVector>
260 void bind(const typename FVElementGeometry::Element& element,
261 const FVElementGeometry& fvGeometry,
262 const SolutionVector& sol) &
263 {
264 this->bind_(element, fvGeometry, sol);
265 }
266
272 template<class SolutionVector>
273 ThisType bindElement(const typename FVElementGeometry::Element& element,
274 const FVElementGeometry& fvGeometry,
275 const SolutionVector& sol) &&
276 {
277 this->bindElement_(element, fvGeometry, sol);
278 return std::move(*this);
279 }
280
281 template<class SolutionVector>
282 void bindElement(const typename FVElementGeometry::Element& element,
283 const FVElementGeometry& fvGeometry,
284 const SolutionVector& sol) &
285 { this->bindElement_(element, fvGeometry, sol); }
286
289 { return *gridVolumeVariablesPtr_; }
290
292 bool hasVolVars(const std::size_t scvIdx) const
293 { return volVarsInserted_(scvIdx); }
294
295private:
298 template<class SolutionVector>
299 void bind_(const typename FVElementGeometry::Element& element,
300 const FVElementGeometry& fvGeometry,
301 const SolutionVector& sol)
302 {
303 clear_();
304
305 const auto& problem = gridVolVars().problem();
306 const auto& gridGeometry = fvGeometry.gridGeometry();
307
308 volVarIndices_.reserve(numInsideVolVars + numInsideVolVars);
309 volumeVariables_.reserve(numInsideVolVars + numInsideVolVars);
310
311 for (const auto& scv : scvs(fvGeometry))
312 {
313 for (const auto otherScvIdx : gridGeometry.connectivityMap()[scv.index()])
314 {
315 if (!volVarsInserted_(otherScvIdx))
316 {
317 const auto& otherScv = fvGeometry.scv(otherScvIdx);
318 volVarIndices_.push_back(otherScvIdx);
319 volumeVariables_.emplace_back();
320 const auto& otherElement = gridGeometry.element(otherScv.elementIndex());
321 volumeVariables_.back().update(
322 elementSolution(otherElement, sol, gridGeometry),
323 problem, otherElement, otherScv
324 );
325 }
326 }
327 }
328
329 if (fvGeometry.hasBoundaryScvf())
330 {
331 for (const auto& scvf : scvfs(fvGeometry))
332 {
333 if (!scvf.boundary() || scvf.isFrontal())
334 continue;
335
336 // check if boundary is a pure dirichlet boundary
337 const auto& problem = gridVolVars().problem();
338 const auto bcTypes = problem.boundaryTypes(element, scvf);
339
340 auto addBoundaryVolVars = [&](const auto& scvFace)
341 {
342 const auto& scvI = fvGeometry.scv(scvFace.insideScvIdx());
343 typename VolumeVariables::PrimaryVariables pv(
344 problem.dirichlet(element, scvFace)[scvI.dofAxis()]
345 );
346 const auto dirichletPriVars = elementSolution<FVElementGeometry>(pv);
347
348 VolumeVariables volVars;
349 volVars.update(dirichletPriVars,
350 problem,
351 element,
352 scvI);
353
354 volumeVariables_.emplace_back(std::move(volVars));
355 volVarIndices_.push_back(scvFace.outsideScvIdx());
356 };
357
358 if (bcTypes.hasDirichlet())
359 {
360 addBoundaryVolVars(scvf);
361 continue;
362 }
363
364 // treat domain corners
365 if (const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf); orthogonalScvf.boundary())
366 if (problem.boundaryTypes(element, orthogonalScvf).hasDirichlet())
367 addBoundaryVolVars(scvf);
368
369 }
370 }
371 }
372
375 template<class SolutionVector>
376 void bindElement_(const typename FVElementGeometry::Element& element,
377 const FVElementGeometry& fvGeometry,
378 const SolutionVector& sol)
379 {
380 clear_();
381 const auto& problem = gridVolVars().problem();
382 const auto& gridGeometry = fvGeometry.gridGeometry();
383 volVarIndices_.reserve(numInsideVolVars);
384
385 for (const auto& scv : scvs(fvGeometry))
386 {
387 volVarIndices_.push_back(scv.index());
388 volumeVariables_.emplace_back();
389 volumeVariables_.back().update(
390 elementSolution(element, sol, gridGeometry),
391 problem, element, scv
392 );
393 }
394 }
395
397 void clear_()
398 {
399 volVarIndices_.clear();
400 volumeVariables_.clear();
401 }
402
403 bool volVarsInserted_(const std::size_t scvIdx) const
404 {
405 return std::find(volVarIndices_.begin(), volVarIndices_.end(), scvIdx) != volVarIndices_.end();
406 }
407
408 int getLocalIdx_(const int scvfIdx) const
409 {
410 const auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), scvfIdx);
411 assert(it != volVarIndices_.end() && "Could not find the current face variables for scvfIdx!");
412 return std::distance(volVarIndices_.begin(), it);
413 }
414
415 const GridVolumeVariables* gridVolumeVariablesPtr_;
416 std::vector<std::size_t> volVarIndices_;
417 std::vector<VolumeVariables> volumeVariables_;
418};
419
420} // end namespace Dumux
421
422#endif
Class for the face variables vector. Specialization for the case of not storing the face variables gl...
Definition: facecentered/staggered/elementvolumevariables.hh:209
GVV GridVolumeVariables
export type of the grid volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:221
ThisType bindElement(const typename FVElementGeometry::Element &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: facecentered/staggered/elementvolumevariables.hh:273
void bind(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:260
ThisType bind(const typename FVElementGeometry::Element &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: facecentered/staggered/elementvolumevariables.hh:251
bool hasVolVars(const std::size_t scvIdx) const
Returns true if volVars exist for the given scv index.
Definition: facecentered/staggered/elementvolumevariables.hh:292
void bindElement(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:282
FaceCenteredStaggeredElementVolumeVariables(const GridVolumeVariables &globalFacesVars)
Definition: facecentered/staggered/elementvolumevariables.hh:226
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:224
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: facecentered/staggered/elementvolumevariables.hh:288
Class for the face variables vector. Specialization for the case of storing the face variables global...
Definition: facecentered/staggered/elementvolumevariables.hh:37
void bindElement(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:159
GVV GridVolumeVariables
export type of the grid volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:45
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: facecentered/staggered/elementvolumevariables.hh:166
ThisType bindElement(const typename FVElementGeometry::Element &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: facecentered/staggered/elementvolumevariables.hh:149
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:48
bool hasVolVars(const std::size_t scvIdx) const
Returns true if volVars exist for the given scv index.
Definition: facecentered/staggered/elementvolumevariables.hh:170
FaceCenteredStaggeredElementVolumeVariables(const GridVolumeVariables &gridVolumeVariables)
Definition: facecentered/staggered/elementvolumevariables.hh:50
void bind(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:91
ThisType bind(const typename FVElementGeometry::Element &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: facecentered/staggered/elementvolumevariables.hh:80
Base class for the face variables vector.
Definition: facecentered/staggered/elementvolumevariables.hh:29
Element solution classes and factory functions.
void addBoundaryVolVars(std::vector< VolumeVariables > &volVars, std::vector< IndexType > &volVarIndices, const Problem &problem, const typename FVElemGeom::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElemGeom &fvGeometry)
Adds the boundary volume variables found within the stencil to the provided containers and stores the...
Definition: cellcentered/mpfa/elementvolumevariables.hh:124
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
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