3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_DISCRETIZATION_FACECENTERED_ELEMENTVOLUMEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_FACECENTERED_ELEMENTVOLUMEVARIABLES_HH
26
27#include <algorithm>
28#include <cassert>
29#include <vector>
30#include <utility>
31
33
34namespace Dumux {
35
40template<class GridVolumeVariables, bool cachingEnabled>
42
47template<class GVV>
48class FaceCenteredStaggeredElementVolumeVariables<GVV, /*cachingEnabled*/true>
49{
50 using ThisType = FaceCenteredStaggeredElementVolumeVariables<GVV, /*cachingEnabled*/true>;
51 using GridGeometry = std::decay_t<decltype(std::declval<GVV>().problem().gridGeometry())>;
52 using FVElementGeometry = typename GridGeometry::LocalView;
53 using SubControlVolume = typename GridGeometry::SubControlVolume;
54
55public:
58
60 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
61
63 : gridVolumeVariablesPtr_(&gridVolumeVariables)
64 , numScv_(gridVolumeVariables.problem().gridGeometry().numScv())
65 {}
66
68 const VolumeVariables& operator [](const SubControlVolume& scv) const
69 {
70 if (scv.index() < numScv_)
71 return gridVolVars().volVars(scv.index());
72 else
73 return boundaryVolumeVariables_[getLocalIdx_(scv.index())];
74 }
75
78 const VolumeVariables& operator [](const std::size_t scvIdx) const
79 {
80 if (scvIdx < numScv_)
81 return gridVolVars().volVars(scvIdx);
82 else
83 return boundaryVolumeVariables_[getLocalIdx_(scvIdx)];
84 }
85
91 template<class SolutionVector>
92 ThisType bind(const typename FVElementGeometry::Element& element,
93 const FVElementGeometry& fvGeometry,
94 const SolutionVector& sol) &&
95 {
96 this->bind(element, fvGeometry, sol);
97 return std::move(*this);
98 }
99
102 template<class SolutionVector>
103 void bind(const typename FVElementGeometry::Element& element,
104 const FVElementGeometry& fvGeometry,
105 const SolutionVector& sol) &
106 {
107 if (!fvGeometry.hasBoundaryScvf())
108 return;
109
110 clear_();
111 boundaryVolVarIndices_.reserve(fvGeometry.gridGeometry().numBoundaryScv());
112 boundaryVolumeVariables_.reserve(fvGeometry.gridGeometry().numBoundaryScv());
113
114 for (const auto& scvf : scvfs(fvGeometry))
115 {
116 if (!scvf.boundary() || scvf.isFrontal() || scvf.processorBoundary())
117 continue;
118
119 // check if boundary is a pure dirichlet boundary
120 const auto& problem = gridVolVars().problem();
121 const auto bcTypes = problem.boundaryTypes(element, scvf);
122
123 auto addBoundaryVolVars = [&](const auto& scvFace)
124 {
125 const auto& scvI = fvGeometry.scv(scvFace.insideScvIdx());
126 typename VolumeVariables::PrimaryVariables pv(
127 problem.dirichlet(element, scvFace)[scvI.dofAxis()]
128 );
129 const auto dirichletPriVars = elementSolution<FVElementGeometry>(pv);
130
131 VolumeVariables volVars;
132 volVars.update(dirichletPriVars, problem, element, scvI);
133
134 boundaryVolumeVariables_.emplace_back(std::move(volVars));
135 boundaryVolVarIndices_.push_back(scvFace.outsideScvIdx());
136 };
137
138 if (bcTypes.hasDirichlet())
139 {
140 addBoundaryVolVars(scvf);
141 continue;
142 }
143
144 // treat domain corners
145 if (const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf); orthogonalScvf.boundary())
146 if (problem.boundaryTypes(element, orthogonalScvf).hasDirichlet())
147 addBoundaryVolVars(scvf);
148
149 }
150
151 assert(boundaryVolumeVariables_.size() == boundaryVolVarIndices_.size());
152 }
153
159 template<class SolutionVector>
160 ThisType bindElement(const typename FVElementGeometry::Element& element,
161 const FVElementGeometry& fvGeometry,
162 const SolutionVector& sol) &&
163 {
164 this->bindElement(element, fvGeometry, sol);
165 return std::move(*this);
166 }
169 template<class SolutionVector>
170 void bindElement(const typename FVElementGeometry::Element& element,
171 const FVElementGeometry& fvGeometry,
172 const SolutionVector& sol) &
173 {}
174
175
178 { return *gridVolumeVariablesPtr_; }
179
181 bool hasVolVars(const std::size_t scvIdx) const
182 {
183 if (scvIdx < numScv_)
184 return true;
185 else
186 {
187 const auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), scvIdx);
188 return it != boundaryVolVarIndices_.end();
189 }
190 }
191
192private:
194 void clear_()
195 {
196 boundaryVolVarIndices_.clear();
197 boundaryVolumeVariables_.clear();
198 }
199
201 int getLocalIdx_(const std::size_t volVarIdx) const
202 {
203 const auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
204 assert(it != boundaryVolVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
205 return std::distance(boundaryVolVarIndices_.begin(), it);
206 }
207
208 const GridVolumeVariables* gridVolumeVariablesPtr_;
209 const std::size_t numScv_;
210 std::vector<std::size_t> boundaryVolVarIndices_;
211 std::vector<VolumeVariables> boundaryVolumeVariables_;
212};
213
218template<class GVV>
219class FaceCenteredStaggeredElementVolumeVariables<GVV, /*cachingEnabled*/false>
220{
221 using ThisType = FaceCenteredStaggeredElementVolumeVariables<GVV, /*cachingEnabled*/false>;
222 using GridGeometry = std::decay_t<decltype(std::declval<GVV>().problem().gridGeometry())>;
223 using FVElementGeometry = typename GridGeometry::LocalView;
224 using SubControlVolume = typename GridGeometry::SubControlVolume;
225
226 static constexpr auto dim = GridGeometry::GridView::dimension;
227 static constexpr auto numInsideVolVars = dim * 2;
228 static constexpr auto numOutsideVolVars = numInsideVolVars * 2 * (dim - 1);
229
230public:
233
235 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
236
238 : gridVolumeVariablesPtr_(&globalFacesVars) {}
239
241 const VolumeVariables& operator [](const SubControlVolume& scv) const
242 { return volumeVariables_[getLocalIdx_(scv.index())]; }
243
245 const VolumeVariables& operator [](const std::size_t scvIdx) const
246 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
247
249 VolumeVariables& operator [](const SubControlVolume& scv)
250 { return volumeVariables_[getLocalIdx_(scv.index())]; }
251
252 // operator for the access with an index
253 VolumeVariables& operator [](const std::size_t scvIdx)
254 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
255
261 template<class SolutionVector>
262 ThisType bind(const typename FVElementGeometry::Element& element,
263 const FVElementGeometry& fvGeometry,
264 const SolutionVector& sol) &&
265 {
266 this->bind_(element, fvGeometry, sol);
267 return std::move(*this);
268 }
269
270 template<class SolutionVector>
271 void bind(const typename FVElementGeometry::Element& element,
272 const FVElementGeometry& fvGeometry,
273 const SolutionVector& sol) &
274 {
275 this->bind_(element, fvGeometry, sol);
276 }
277
283 template<class SolutionVector>
284 ThisType bindElement(const typename FVElementGeometry::Element& element,
285 const FVElementGeometry& fvGeometry,
286 const SolutionVector& sol) &&
287 {
288 this->bindElement_(element, fvGeometry, sol);
289 return std::move(*this);
290 }
291
292 template<class SolutionVector>
293 void bindElement(const typename FVElementGeometry::Element& element,
294 const FVElementGeometry& fvGeometry,
295 const SolutionVector& sol) &
296 { this->bindElement_(element, fvGeometry, sol); }
297
300 { return *gridVolumeVariablesPtr_; }
301
303 bool hasVolVars(const std::size_t scvIdx) const
304 { return volVarsInserted_(scvIdx); }
305
306private:
309 template<class SolutionVector>
310 void bind_(const typename FVElementGeometry::Element& element,
311 const FVElementGeometry& fvGeometry,
312 const SolutionVector& sol)
313 {
314 clear_();
315
316 const auto& problem = gridVolVars().problem();
317 const auto& gridGeometry = fvGeometry.gridGeometry();
318
319 volVarIndices_.reserve(numInsideVolVars + numInsideVolVars);
320 volumeVariables_.reserve(numInsideVolVars + numInsideVolVars);
321
322 for (const auto& scv : scvs(fvGeometry))
323 {
324 for (const auto otherScvIdx : gridGeometry.connectivityMap()[scv.index()])
325 {
326 if (!volVarsInserted_(otherScvIdx))
327 {
328 const auto& otherScv = fvGeometry.scv(otherScvIdx);
329 volVarIndices_.push_back(otherScvIdx);
330 volumeVariables_.emplace_back();
331 const auto& otherElement = gridGeometry.element(otherScv.elementIndex());
332 volumeVariables_.back().update(
333 elementSolution(otherElement, sol, gridGeometry),
334 problem, otherElement, otherScv
335 );
336 }
337 }
338 }
339
340 if (fvGeometry.hasBoundaryScvf())
341 {
342 for (const auto& scvf : scvfs(fvGeometry))
343 {
344 if (!scvf.boundary() || scvf.isFrontal())
345 continue;
346
347 // check if boundary is a pure dirichlet boundary
348 const auto& problem = gridVolVars().problem();
349 const auto bcTypes = problem.boundaryTypes(element, scvf);
350
351 auto addBoundaryVolVars = [&](const auto& scvFace)
352 {
353 const auto& scvI = fvGeometry.scv(scvFace.insideScvIdx());
354 typename VolumeVariables::PrimaryVariables pv(
355 problem.dirichlet(element, scvFace)[scvI.dofAxis()]
356 );
357 const auto dirichletPriVars = elementSolution<FVElementGeometry>(pv);
358
359 VolumeVariables volVars;
360 volVars.update(dirichletPriVars,
361 problem,
362 element,
363 scvI);
364
365 volumeVariables_.emplace_back(std::move(volVars));
366 volVarIndices_.push_back(scvFace.outsideScvIdx());
367 };
368
369 if (bcTypes.hasDirichlet())
370 {
371 addBoundaryVolVars(scvf);
372 continue;
373 }
374
375 // treat domain corners
376 if (const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf); orthogonalScvf.boundary())
377 if (problem.boundaryTypes(element, orthogonalScvf).hasDirichlet())
378 addBoundaryVolVars(scvf);
379
380 }
381 }
382 }
383
386 template<class SolutionVector>
387 void bindElement_(const typename FVElementGeometry::Element& element,
388 const FVElementGeometry& fvGeometry,
389 const SolutionVector& sol)
390 {
391 clear_();
392 const auto& problem = gridVolVars().problem();
393 const auto& gridGeometry = fvGeometry.gridGeometry();
394 volVarIndices_.reserve(numInsideVolVars);
395
396 for (const auto& scv : scvs(fvGeometry))
397 {
398 volVarIndices_.push_back(scv.index());
399 volumeVariables_.emplace_back();
400 volumeVariables_.back().update(
401 elementSolution(element, sol, gridGeometry),
402 problem, element, scv
403 );
404 }
405 }
406
408 void clear_()
409 {
410 volVarIndices_.clear();
411 volumeVariables_.clear();
412 }
413
414 bool volVarsInserted_(const std::size_t scvIdx) const
415 {
416 return std::find(volVarIndices_.begin(), volVarIndices_.end(), scvIdx) != volVarIndices_.end();
417 }
418
419 int getLocalIdx_(const int scvfIdx) const
420 {
421 const auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), scvfIdx);
422 assert(it != volVarIndices_.end() && "Could not find the current face variables for scvfIdx!");
423 return std::distance(volVarIndices_.begin(), it);
424 }
425
426 const GridVolumeVariables* gridVolumeVariablesPtr_;
427 std::vector<std::size_t> volVarIndices_;
428 std::vector<VolumeVariables> volumeVariables_;
429};
430
431} // end namespace Dumux
432
433#endif
Element solution classes and factory functions.
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:294
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
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:136
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Base class for the face variables vector.
Definition: facecentered/staggered/elementvolumevariables.hh:41
Class for the face variables vector. Specialization for the case of storing the face variables global...
Definition: facecentered/staggered/elementvolumevariables.hh:49
void bindElement(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:170
GVV GridVolumeVariables
export type of the grid volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:57
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: facecentered/staggered/elementvolumevariables.hh:177
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:160
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:60
bool hasVolVars(const std::size_t scvIdx) const
Returns true if volVars exist for the given scv index.
Definition: facecentered/staggered/elementvolumevariables.hh:181
FaceCenteredStaggeredElementVolumeVariables(const GridVolumeVariables &gridVolumeVariables)
Definition: facecentered/staggered/elementvolumevariables.hh:62
void bind(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:103
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:92
Class for the face variables vector. Specialization for the case of not storing the face variables gl...
Definition: facecentered/staggered/elementvolumevariables.hh:220
GVV GridVolumeVariables
export type of the grid volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:232
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:284
void bind(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:271
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:262
bool hasVolVars(const std::size_t scvIdx) const
Returns true if volVars exist for the given scv index.
Definition: facecentered/staggered/elementvolumevariables.hh:303
void bindElement(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: facecentered/staggered/elementvolumevariables.hh:293
FaceCenteredStaggeredElementVolumeVariables(const GridVolumeVariables &globalFacesVars)
Definition: facecentered/staggered/elementvolumevariables.hh:237
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: facecentered/staggered/elementvolumevariables.hh:235
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: facecentered/staggered/elementvolumevariables.hh:299