version 3.8
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// 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_CCTPFA_ELEMENT_VOLUMEVARIABLES_HH
13#define DUMUX_DISCRETIZATION_CCTPFA_ELEMENT_VOLUMEVARIABLES_HH
14
15#include <algorithm>
16#include <type_traits>
17#include <vector>
18#include <utility>
19
21
22namespace Dumux {
23
31template<class GVV, bool cachingEnabled>
33{};
34
40template<class GVV>
41class CCTpfaElementVolumeVariables<GVV, /*cachingEnabled*/true>
42{
43public:
46
48 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
49
52 : gridVolVarsPtr_(&gridVolVars)
53 , numScv_(gridVolVars.problem().gridGeometry().numScv())
54 {}
55
57 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
58 const VolumeVariables& operator [](const SubControlVolume& scv) const
59 {
60 if (scv.dofIndex() < numScv_)
61 return gridVolVars().volVars(scv.dofIndex());
62 else
63 return boundaryVolumeVariables_[getLocalIdx_(scv.dofIndex())];
64 }
65
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 CCTpfaElementVolumeVariables 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 CCTpfaElementVolumeVariables 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
129 template<class FVElementGeometry, class SolutionVector>
130 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
131 const FVElementGeometry& fvGeometry,
132 const SolutionVector& sol)
133 {
134 if (!fvGeometry.hasBoundaryScvf())
135 return;
136
137 clear_();
138 boundaryVolVarIndices_.reserve(fvGeometry.numScvf());
139 boundaryVolumeVariables_.reserve(fvGeometry.numScvf());
140
141 for (const auto& scvf : scvfs(fvGeometry))
142 {
143 if (!scvf.boundary())
144 continue;
145
146 // check if boundary is a pure dirichlet boundary
147 const auto& problem = gridVolVars().problem();
148 const auto bcTypes = problem.boundaryTypes(element, scvf);
149 if (bcTypes.hasOnlyDirichlet())
150 {
151 const auto dirichletPriVars = elementSolution<FVElementGeometry>(problem.dirichlet(element, scvf));
152 auto&& scvI = fvGeometry.scv(scvf.insideScvIdx());
153
154 VolumeVariables volVars;
155 volVars.update(dirichletPriVars,
156 problem,
157 element,
158 scvI);
159
160 boundaryVolumeVariables_.emplace_back(std::move(volVars));
161 boundaryVolVarIndices_.push_back(scvf.outsideScvIdx());
162 }
163 }
164 }
165
167 template<class FVElementGeometry, class SolutionVector>
168 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
169 const FVElementGeometry& fvGeometry,
170 const SolutionVector& sol)
171 {}
172
173 const GridVolumeVariables* gridVolVarsPtr_;
174
176 int getLocalIdx_(const int volVarIdx) const
177 {
178 auto it = std::find(boundaryVolVarIndices_.begin(), boundaryVolVarIndices_.end(), volVarIdx);
179 assert(it != boundaryVolVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
180 return std::distance(boundaryVolVarIndices_.begin(), it);
181 }
182
183 std::vector<std::size_t> boundaryVolVarIndices_;
184 std::vector<VolumeVariables> boundaryVolumeVariables_;
185 const std::size_t numScv_;
186};
187
192template<class GVV>
193class CCTpfaElementVolumeVariables<GVV, /*cachingEnabled*/false>
194{
195public:
198
200 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
201
204 : gridVolVarsPtr_(&gridVolVars) {}
205
211 template<class FVElementGeometry, class SolutionVector>
212 CCTpfaElementVolumeVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
213 const FVElementGeometry& fvGeometry,
214 const SolutionVector& sol) &&
215 {
216 this->bind_(element, fvGeometry, sol);
217 return std::move(*this);
218 }
219
220 template<class FVElementGeometry, class SolutionVector>
221 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
222 const FVElementGeometry& fvGeometry,
223 const SolutionVector& sol) &
224 { this->bind_(element, fvGeometry, sol); }
225
231 template<class FVElementGeometry, class SolutionVector>
232 CCTpfaElementVolumeVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
233 const FVElementGeometry& fvGeometry,
234 const SolutionVector& sol) &&
235 {
236 this->bindElement_(element, fvGeometry, sol);
237 return std::move(*this);
238 }
239
240 template<class FVElementGeometry, class SolutionVector>
241 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
242 const FVElementGeometry& fvGeometry,
243 const SolutionVector& sol) &
244 { this->bindElement_(element, fvGeometry, sol); }
245
247 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
248 const VolumeVariables& operator [](const SubControlVolume& scv) const
249 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
250
252 template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
253 VolumeVariables& operator [](const SubControlVolume& scv)
254 { return volumeVariables_[getLocalIdx_(scv.dofIndex())]; }
255
257 const VolumeVariables& operator [](std::size_t scvIdx) const
258 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
259
261 VolumeVariables& operator [](std::size_t scvIdx)
262 { return volumeVariables_[getLocalIdx_(scvIdx)]; }
263
266 { return *gridVolVarsPtr_; }
267
268private:
270 void clear_()
271 {
272 volVarIndices_.clear();
273 volumeVariables_.clear();
274 }
275
277 template<class FVElementGeometry, class SolutionVector>
278 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
279 const FVElementGeometry& fvGeometry,
280 const SolutionVector& sol)
281 {
282 clear_();
283
284 const auto& problem = gridVolVars().problem();
285 const auto& gridGeometry = fvGeometry.gridGeometry();
286 const auto globalI = gridGeometry.elementMapper().index(element);
287 const auto& connectivityMapI = gridGeometry.connectivityMap()[globalI];
288 const auto numDofs = connectivityMapI.size() + 1;
289
290 // resize local containers to the required size (for internal elements)
291 volumeVariables_.resize(numDofs);
292 volVarIndices_.resize(numDofs);
293 int localIdx = 0;
294
295 // update the volume variables of the element at hand
296 auto&& scvI = fvGeometry.scv(globalI);
297 volumeVariables_[localIdx].update(elementSolution(element, sol, gridGeometry),
298 problem,
299 element,
300 scvI);
301 volVarIndices_[localIdx] = scvI.dofIndex();
302 ++localIdx;
303
304 // Update the volume variables of the neighboring elements
305 for (const auto& dataJ : connectivityMapI)
306 {
307 const auto& elementJ = gridGeometry.element(dataJ.globalJ);
308 auto&& scvJ = fvGeometry.scv(dataJ.globalJ);
309 volumeVariables_[localIdx].update(elementSolution(elementJ, sol, gridGeometry),
310 problem,
311 elementJ,
312 scvJ);
313 volVarIndices_[localIdx] = scvJ.dofIndex();
314 ++localIdx;
315 }
316
317 if (fvGeometry.hasBoundaryScvf())
318 {
319 // Update boundary volume variables
320 for (auto&& scvf : scvfs(fvGeometry))
321 {
322 // if we are not on a boundary, skip to the next scvf
323 if (!scvf.boundary())
324 continue;
325
326 // check if boundary is a pure dirichlet boundary
327 const auto bcTypes = problem.boundaryTypes(element, scvf);
328 if (bcTypes.hasOnlyDirichlet())
329 {
330 const auto dirichletPriVars = elementSolution<FVElementGeometry>(problem.dirichlet(element, scvf));
331
332 volumeVariables_.resize(localIdx+1);
333 volVarIndices_.resize(localIdx+1);
334 volumeVariables_[localIdx].update(dirichletPriVars,
335 problem,
336 element,
337 scvI);
338 volVarIndices_[localIdx] = scvf.outsideScvIdx();
339 ++localIdx;
340 }
341 }
342 }
343
346 // const auto& additionalDofDependencies = problem.getAdditionalDofDependencies(globalI);
347 // if (!additionalDofDependencies.empty())
348 // {
349 // volumeVariables_.resize(volumeVariables_.size() + additionalDofDependencies.size());
350 // volVarIndices_.resize(volVarIndices_.size() + additionalDofDependencies.size());
351 // for (auto globalJ : additionalDofDependencies)
352 // {
353 // const auto& elementJ = gridGeometry.element(globalJ);
354 // auto&& scvJ = fvGeometry.scv(globalJ);
355
356 // volumeVariables_[localIdx].update(elementSolution(elementJ, sol, gridGeometry),
357 // problem,
358 // elementJ,
359 // scvJ);
360 // volVarIndices_[localIdx] = scvJ.dofIndex();
361 // ++localIdx;
362 // }
363 // }
364 }
365
366 template<class FVElementGeometry, class SolutionVector>
367 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
368 const FVElementGeometry& fvGeometry,
369 const SolutionVector& sol)
370 {
371 clear_();
372
373 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
374 volumeVariables_.resize(1);
375 volVarIndices_.resize(1);
376
377 // update the volume variables of the element
378 auto&& scv = fvGeometry.scv(eIdx);
379 volumeVariables_[0].update(elementSolution(element, sol, fvGeometry.gridGeometry()),
380 gridVolVars().problem(),
381 element,
382 scv);
383 volVarIndices_[0] = scv.dofIndex();
384 }
385
386 const GridVolumeVariables* gridVolVarsPtr_;
387
389 int getLocalIdx_(const int volVarIdx) const
390 {
391 auto it = std::find(volVarIndices_.begin(), volVarIndices_.end(), volVarIdx);
392 assert(it != volVarIndices_.end() && "Could not find the current volume variables for volVarIdx!");
393 return std::distance(volVarIndices_.begin(), it);
394 }
395
396 std::vector<std::size_t> volVarIndices_;
397 std::vector<VolumeVariables> volumeVariables_;
398};
399
400} // end namespace Dumux
401
402#endif
The local element solution class for cell-centered methods.
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:221
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:241
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:200
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:232
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:197
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:212
CCTpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/tpfa/elementvolumevariables.hh:203
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/tpfa/elementvolumevariables.hh:265
CCTpfaElementVolumeVariables(const GridVolumeVariables &gridVolVars)
Constructor.
Definition: cellcentered/tpfa/elementvolumevariables.hh:51
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:101
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:81
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:110
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: cellcentered/tpfa/elementvolumevariables.hh:90
const GridVolumeVariables & gridVolVars() const
The global volume variables object we are a restriction of.
Definition: cellcentered/tpfa/elementvolumevariables.hh:116
GVV GridVolumeVariables
export type of the grid volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:45
typename GridVolumeVariables::VolumeVariables VolumeVariables
export type of the volume variables
Definition: cellcentered/tpfa/elementvolumevariables.hh:48
The local (stencil) volume variables class for cell centered tpfa models.
Definition: cellcentered/tpfa/elementvolumevariables.hh:33
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