3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
cellcentered/tpfa/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_CCTPFA_ELEMENT_VOLUMEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_CCTPFA_ELEMENT_VOLUMEVARIABLES_HH
26
27#include <algorithm>
28#include <type_traits>
29#include <vector>
30#include <utility>
31
33
34namespace Dumux {
35
43template<class GVV, bool cachingEnabled>
45{};
46
52template<class GVV>
53class CCTpfaElementVolumeVariables<GVV, /*cachingEnabled*/true>
54{
55public:
58
60 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
61
64 : gridVolVarsPtr_(&gridVolVars)
65 , numScv_(gridVolVars.problem().gridGeometry().numScv())
66 {}
67
69 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
70 const VolumeVariables& operator [](const SubControlVolume& scv) const
71 {
72 if (scv.dofIndex() < numScv_)
73 return gridVolVars().volVars(scv.dofIndex());
74 else
75 return boundaryVolumeVariables_[getLocalIdx_(scv.dofIndex())];
76 }
77
79 const VolumeVariables& operator [](const std::size_t scvIdx) const
80 {
81 if (scvIdx < numScv_)
82 return gridVolVars().volVars(scvIdx);
83 else
84 return boundaryVolumeVariables_[getLocalIdx_(scvIdx)];
85 }
86
92 template<class FVElementGeometry, class SolutionVector>
93 CCTpfaElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
94 const FVElementGeometry& fvGeometry,
95 const SolutionVector& sol) &&
96 {
97 this->bind_(element, fvGeometry, sol);
98 return std::move(*this);
99 }
100
101 template<class FVElementGeometry, class SolutionVector>
102 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
103 const FVElementGeometry& fvGeometry,
104 const SolutionVector& sol) &
105 { this->bind_(element, fvGeometry, sol); }
106
112 template<class FVElementGeometry, class SolutionVector>
113 CCTpfaElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
114 const FVElementGeometry& fvGeometry,
115 const SolutionVector& sol) &&
116 {
117 this->bindElement_(element, fvGeometry, sol);
118 return std::move(*this);
119 }
120
121 template<class FVElementGeometry, class SolutionVector>
122 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
123 const FVElementGeometry& fvGeometry,
124 const SolutionVector& sol) &
125 { this->bindElement_(element, fvGeometry, sol); }
126
129 { return *gridVolVarsPtr_; }
130
131private:
132
134 void clear_()
135 {
136 boundaryVolVarIndices_.clear();
137 boundaryVolumeVariables_.clear();
138 }
139
141 template<class FVElementGeometry, class SolutionVector>
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 for (const auto& scvf : scvfs(fvGeometry))
154 {
155 if (!scvf.boundary())
156 continue;
157
158 // check if boundary is a pure dirichlet boundary
159 const auto& problem = gridVolVars().problem();
160 const auto bcTypes = problem.boundaryTypes(element, scvf);
161 if (bcTypes.hasOnlyDirichlet())
162 {
163 const auto dirichletPriVars = elementSolution<FVElementGeometry>(problem.dirichlet(element, scvf));
164 auto&& scvI = fvGeometry.scv(scvf.insideScvIdx());
165
166 VolumeVariables volVars;
167 volVars.update(dirichletPriVars,
168 problem,
169 element,
170 scvI);
171
172 boundaryVolumeVariables_.emplace_back(std::move(volVars));
173 boundaryVolVarIndices_.push_back(scvf.outsideScvIdx());
174 }
175 }
176 }
177
179 template<class FVElementGeometry, class SolutionVector>
180 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
181 const FVElementGeometry& fvGeometry,
182 const SolutionVector& sol)
183 {}
184
185 const GridVolumeVariables* gridVolVarsPtr_;
186
188 int getLocalIdx_(const int volVarIdx) const
189 {
190 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
191 assert(it != boundaryVolVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
192 return std::distance(boundaryVolVarIndices_.begin(), it);
193 }
194
195 std::vector<std::size_t> boundaryVolVarIndices_;
196 std::vector<VolumeVariables> boundaryVolumeVariables_;
197 const std::size_t numScv_;
198};
199
204template<class GVV>
205class CCTpfaElementVolumeVariables<GVV, /*cachingEnabled*/false>
206{
207public:
210
212 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
213
216 : gridVolVarsPtr_(&gridVolVars) {}
217
223 template<class FVElementGeometry, class SolutionVector>
224 CCTpfaElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
225 const FVElementGeometry& fvGeometry,
226 const SolutionVector& sol) &&
227 {
228 this->bind_(element, fvGeometry, sol);
229 return std::move(*this);
230 }
231
232 template<class FVElementGeometry, class SolutionVector>
233 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
234 const FVElementGeometry& fvGeometry,
235 const SolutionVector& sol) &
236 { this->bind_(element, fvGeometry, sol); }
237
243 template<class FVElementGeometry, class SolutionVector>
244 CCTpfaElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
245 const FVElementGeometry& fvGeometry,
246 const SolutionVector& sol) &&
247 {
248 this->bindElement_(element, fvGeometry, sol);
249 return std::move(*this);
250 }
251
252 template<class FVElementGeometry, class SolutionVector>
253 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
254 const FVElementGeometry& fvGeometry,
255 const SolutionVector& sol) &
256 { this->bindElement_(element, fvGeometry, sol); }
257
259 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
260 const VolumeVariables& operator [](const SubControlVolume& scv) const
261 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
262
264 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
265 VolumeVariables& operator [](const SubControlVolume& scv)
266 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
267
269 const VolumeVariables& operator [](std::size_t scvIdx) const
270 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
271
273 VolumeVariables& operator [](std::size_t scvIdx)
274 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
275
278 { return *gridVolVarsPtr_; }
279
280private:
282 void clear_()
283 {
284 volVarIndices_.clear();
285 volumeVariables_.clear();
286 }
287
289 template<class FVElementGeometry, class SolutionVector>
290 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
291 const FVElementGeometry& fvGeometry,
292 const SolutionVector& sol)
293 {
294 clear_();
295
296 const auto& problem = gridVolVars().problem();
297 const auto& gridGeometry = fvGeometry.gridGeometry();
298 const auto globalI = gridGeometry.elementMapper().index(element);
299 const auto& connectivityMapI = gridGeometry.connectivityMap()[globalI];
300 const auto numDofs = connectivityMapI.size() + 1;
301
302 // resize local containers to the required size (for internal elements)
303 volumeVariables_.resize(numDofs);
304 volVarIndices_.resize(numDofs);
305 int localIdx = 0;
306
307 // update the volume variables of the element at hand
308 auto&& scvI = fvGeometry.scv(globalI);
309 volumeVariables_[localIdx].update(elementSolution(element, sol, gridGeometry),
310 problem,
311 element,
312 scvI);
313 volVarIndices_[localIdx] = scvI.dofIndex();
314 ++localIdx;
315
316 // Update the volume variables of the neighboring elements
317 for (const auto& dataJ : connectivityMapI)
318 {
319 const auto& elementJ = gridGeometry.element(dataJ.globalJ);
320 auto&& scvJ = fvGeometry.scv(dataJ.globalJ);
321 volumeVariables_[localIdx].update(elementSolution(elementJ, sol, gridGeometry),
322 problem,
323 elementJ,
324 scvJ);
325 volVarIndices_[localIdx] = scvJ.dofIndex();
326 ++localIdx;
327 }
328
329 if (fvGeometry.hasBoundaryScvf())
330 {
331 // Update boundary volume variables
332 for (auto&& scvf : scvfs(fvGeometry))
333 {
334 // if we are not on a boundary, skip to the next scvf
335 if (!scvf.boundary())
336 continue;
337
338 // check if boundary is a pure dirichlet boundary
339 const auto bcTypes = problem.boundaryTypes(element, scvf);
340 if (bcTypes.hasOnlyDirichlet())
341 {
342 const auto dirichletPriVars = elementSolution<FVElementGeometry>(problem.dirichlet(element, scvf));
343
344 volumeVariables_.resize(localIdx+1);
345 volVarIndices_.resize(localIdx+1);
346 volumeVariables_[localIdx].update(dirichletPriVars,
347 problem,
348 element,
349 scvI);
350 volVarIndices_[localIdx] = scvf.outsideScvIdx();
351 ++localIdx;
352 }
353 }
354 }
355
358 // const auto& additionalDofDependencies = problem.getAdditionalDofDependencies(globalI);
359 // if (!additionalDofDependencies.empty())
360 // {
361 // volumeVariables_.resize(volumeVariables_.size() + additionalDofDependencies.size());
362 // volVarIndices_.resize(volVarIndices_.size() + additionalDofDependencies.size());
363 // for (auto globalJ : additionalDofDependencies)
364 // {
365 // const auto& elementJ = gridGeometry.element(globalJ);
366 // auto&& scvJ = fvGeometry.scv(globalJ);
367
368 // volumeVariables_[localIdx].update(elementSolution(elementJ, sol, gridGeometry),
369 // problem,
370 // elementJ,
371 // scvJ);
372 // volVarIndices_[localIdx] = scvJ.dofIndex();
373 // ++localIdx;
374 // }
375 // }
376 }
377
378 template<class FVElementGeometry, class SolutionVector>
379 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
380 const FVElementGeometry& fvGeometry,
381 const SolutionVector& sol)
382 {
383 clear_();
384
385 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
386 volumeVariables_.resize(1);
387 volVarIndices_.resize(1);
388
389 // update the volume variables of the element
390 auto&& scv = fvGeometry.scv(eIdx);
391 volumeVariables_[0].update(elementSolution(element, sol, fvGeometry.gridGeometry()),
392 gridVolVars().problem(),
393 element,
394 scv);
395 volVarIndices_[0] = scv.dofIndex();
396 }
397
398 const GridVolumeVariables* gridVolVarsPtr_;
399
401 int getLocalIdx_(const int volVarIdx) const
402 {
403 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
404 assert(it != volVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
405 return std::distance(volVarIndices_.begin(), it);
406 }
407
408 std::vector<std::size_t> volVarIndices_;
409 std::vector<VolumeVariables> volumeVariables_;
410};
411
412} // end namespace Dumux
413
414#endif
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:292
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
Definition: adapt.hh:29
The local (stencil) volume variables class for cell centered tpfa models.
Definition: cellcentered/tpfa/elementvolumevariables.hh:45
CCTpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/tpfa/elementvolumevariables.hh:63
CCTpfaElementVolumeVariables 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: cellcentered/tpfa/elementvolumevariables.hh:113
CCTpfaElementVolumeVariables 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: cellcentered/tpfa/elementvolumevariables.hh:93
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:122
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:102
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/tpfa/elementvolumevariables.hh:128
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:57
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:60
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:233
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:253
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:212
CCTpfaElementVolumeVariables 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: cellcentered/tpfa/elementvolumevariables.hh:244
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:209
CCTpfaElementVolumeVariables 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: cellcentered/tpfa/elementvolumevariables.hh:224
CCTpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/tpfa/elementvolumevariables.hh:215
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/tpfa/elementvolumevariables.hh:277
The local element solution class for cell-centered methods.