version 3.8
discretization/cellcentered/mpfa/elementfluxvariablescache.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_CCMPFA_ELEMENT_FLUXVARSCACHE_HH
13#define DUMUX_DISCRETIZATION_CCMPFA_ELEMENT_FLUXVARSCACHE_HH
14
15#include <algorithm>
16#include <cassert>
17#include <type_traits>
18#include <vector>
19#include <utility>
20
21#include <dune/common/exceptions.hh>
22
23namespace Dumux {
24
29template<class PrimaryIV, class PrimaryIVDataHandle,
30 class SecondaryIV, class SecondaryIVDataHandle>
32{
33 std::vector<PrimaryIV> primaryInteractionVolumes;
34 std::vector<SecondaryIV> secondaryInteractionVolumes;
35
36 std::vector<PrimaryIVDataHandle> primaryDataHandles;
37 std::vector<SecondaryIVDataHandle> secondaryDataHandles;
38};
39
49template<class GFVC, bool cachingEnabled>
51
56template<class GFVC>
58{
59public:
61 using PrimaryInteractionVolume = typename GFVC::PrimaryInteractionVolume;
62 using SecondaryInteractionVolume = typename GFVC::SecondaryInteractionVolume;
63
65 using PrimaryIvDataHandle = typename GFVC::PrimaryIvDataHandle;
66 using SecondaryIvDataHandle = typename GFVC::SecondaryIvDataHandle;
67
69 using FluxVariablesCache = typename GFVC::FluxVariablesCache;
70
72 static constexpr bool cachingEnabled = true;
73
76
77private:
79 using FluxVariablesCacheFiller = typename GFVC::Traits::FluxVariablesCacheFiller;
80
82 class BoundaryCacheData
83 {
84 // allow the element flux variables class access to private members
86
87 public:
89 template<class SubControlVolumeFace>
90 const FluxVariablesCache& operator [](const SubControlVolumeFace& scvf) const
91 { return fluxVarCaches_[getLocalIdx_(scvf.index())]; }
92
93 template<class SubControlVolumeFace>
94 FluxVariablesCache& operator [](const SubControlVolumeFace& scvf)
95 { return fluxVarCaches_[getLocalIdx_(scvf.index())]; }
96
98 void clear()
99 {
100 fluxVarCaches_.clear();
101 cacheScvfIndices_.clear();
102 ivDataStorage_.primaryInteractionVolumes.clear();
103 ivDataStorage_.secondaryInteractionVolumes.clear();
104 ivDataStorage_.primaryDataHandles.clear();
105 ivDataStorage_.secondaryDataHandles.clear();
106 }
107
108 private:
110 int getLocalIdx_(const int scvfIdx) const
111 {
112 auto it = std::find(cacheScvfIndices_.begin(), cacheScvfIndices_.end(), scvfIdx);
113 assert(it != cacheScvfIndices_.end() && "Could not find the local idx for the given scvf idx!");
114 return std::distance(cacheScvfIndices_.begin(), it);
115 }
116
117 std::vector<std::size_t> cacheScvfIndices_;
118 std::vector<FluxVariablesCache> fluxVarCaches_;
119
120 // stored boundary interaction volumes and handles
125 IVDataStorage ivDataStorage_;
126 };
127
128public:
131 : gridFluxVarsCachePtr_(&global)
132 {}
133
139 template<class FVElementGeometry, class ElementVolumeVariables>
140 CCMpfaElementFluxVariablesCache bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
141 const FVElementGeometry& fvGeometry,
142 const ElementVolumeVariables& elemVolVars) &&
143 {
144 this->bindElement_(element, fvGeometry, elemVolVars);
145 return std::move(*this);
146 }
147
149 template<class FVElementGeometry, class ElementVolumeVariables>
150 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
151 const FVElementGeometry& fvGeometry,
152 const ElementVolumeVariables& elemVolVars) &
153 { this->bindElement_(element, fvGeometry, elemVolVars); }
154
160 template<class FVElementGeometry, class ElementVolumeVariables>
161 CCMpfaElementFluxVariablesCache bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
162 const FVElementGeometry& fvGeometry,
163 const ElementVolumeVariables& elemVolVars) &&
164 {
165 this->bind_(element, fvGeometry, elemVolVars);
166 return std::move(*this);
167 }
168
170 template<class FVElementGeometry, class ElementVolumeVariables>
171 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
172 const FVElementGeometry& fvGeometry,
173 const ElementVolumeVariables& elemVolVars) &
174 { this->bind_(element, fvGeometry, elemVolVars); }
175
181 template<class FVElementGeometry, class ElementVolumeVariables>
182 CCMpfaElementFluxVariablesCache bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
183 const FVElementGeometry& fvGeometry,
184 const ElementVolumeVariables& elemVolVars,
185 const typename FVElementGeometry::SubControlVolumeFace& scvf) &&
186 {
187 this->bindScvf_(element, fvGeometry, elemVolVars, scvf);
188 return std::move(*this);
189 }
190
192 template<class FVElementGeometry, class ElementVolumeVariables>
193 void bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
194 const FVElementGeometry& fvGeometry,
195 const ElementVolumeVariables& elemVolVars,
196 const typename FVElementGeometry::SubControlVolumeFace& scvf) &
197 { this->bindScvf_(element, fvGeometry, elemVolVars, scvf); }
198
200 template<class FVElementGeometry, class ElementVolumeVariables>
201 void update(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
202 const FVElementGeometry& fvGeometry,
203 const ElementVolumeVariables& elemVolVars)
204 {
205 // Update only if the filler puts solution-dependent stuff into the caches
206 if (FluxVariablesCacheFiller::isSolDependent)
207 {
208 // helper class to fill flux variables caches
209 FluxVariablesCacheFiller filler(gridFluxVarsCachePtr_->problem());
210
211 // first, set all the caches to "outdated"
212 for (auto& cache : boundaryCacheData_.fluxVarCaches_)
213 cache.setUpdateStatus(false);
214
215 // go through the caches maybe update them
216 for (auto scvfIdx : boundaryCacheData_.cacheScvfIndices_)
217 {
218 const auto& scvf = fvGeometry.scvf(scvfIdx);
219 auto& scvfCache = boundaryCacheData_[scvf];
220 if (!scvfCache.isUpdated())
221 filler.fill(boundaryCacheData_, scvfCache, boundaryCacheData_.ivDataStorage_,
222 fvGeometry, elemVolVars, scvf);
223 }
224 }
225 }
226
228 template<class SubControlVolumeFace>
229 const FluxVariablesCache& operator [](const SubControlVolumeFace& scvf) const
230 { return !isEmbeddedInBoundaryIV_(scvf) ? (*gridFluxVarsCachePtr_)[scvf] : boundaryCacheData_[scvf]; }
231
233 template<class SubControlVolumeFace>
234 const PrimaryInteractionVolume& primaryInteractionVolume(const SubControlVolumeFace& scvf) const
235 {
236 return isEmbeddedInBoundaryIV_(scvf)
237 ? boundaryCacheData_.ivDataStorage_.primaryInteractionVolumes[ (*this)[scvf].ivIndexInContainer() ]
238 : gridFluxVarsCachePtr_->primaryInteractionVolume(scvf);
239 }
240
242 template<class SubControlVolumeFace>
243 const PrimaryIvDataHandle& primaryDataHandle(const SubControlVolumeFace& scvf) const
244 {
245 return isEmbeddedInBoundaryIV_(scvf)
246 ? boundaryCacheData_.ivDataStorage_.primaryDataHandles[ (*this)[scvf].ivIndexInContainer() ]
247 : gridFluxVarsCachePtr_->primaryDataHandle(scvf);
248 }
249
251 template<class SubControlVolumeFace>
252 const SecondaryInteractionVolume& secondaryInteractionVolume(const SubControlVolumeFace& scvf) const
253 {
254 return isEmbeddedInBoundaryIV_(scvf)
255 ? boundaryCacheData_.ivDataStorage_.secondaryInteractionVolumes[ (*this)[scvf].ivIndexInContainer() ]
256 : gridFluxVarsCachePtr_->secondaryInteractionVolume(scvf);
257 }
258
260 template<class SubControlVolumeFace>
261 const SecondaryIvDataHandle& secondaryDataHandle(const SubControlVolumeFace& scvf) const
262 {
263 return isEmbeddedInBoundaryIV_(scvf)
264 ? boundaryCacheData_.ivDataStorage_.secondaryDataHandles[ (*this)[scvf].ivIndexInContainer() ]
265 : gridFluxVarsCachePtr_->secondaryDataHandle(scvf);
266 }
267
270 { return *gridFluxVarsCachePtr_; }
271
272private:
273
274
276 template<class FVElementGeometry, class ElementVolumeVariables>
277 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
278 const FVElementGeometry& fvGeometry,
279 const ElementVolumeVariables& elemVolVars)
280 { DUNE_THROW(Dune::NotImplemented, "Local element binding of the flux variables cache in mpfa schemes"); }
281
283 template<class FVElementGeometry, class ElementVolumeVariables>
284 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
285 const FVElementGeometry& fvGeometry,
286 const ElementVolumeVariables& elemVolVars)
287 {
288 boundaryCacheData_.clear();
289
290 // find out how much memory needs to be reserved
291 std::size_t numPrimaryIv; numPrimaryIv = 0;
292 std::size_t numSecondaryIv; numSecondaryIv = 0;
293 std::size_t numCaches; numCaches = 0;
294
295 const auto& gridGeometry = fvGeometry.gridGeometry();
296 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
297
298 // lambda to check if a scvf was handled already
299 auto scvfHandled = [&] (auto idx)
300 {
301 return std::find(boundaryCacheData_.cacheScvfIndices_.begin(),
302 boundaryCacheData_.cacheScvfIndices_.end(),
303 idx) != boundaryCacheData_.cacheScvfIndices_.end();
304 };
305
306 // lambda to increase counters for a given scvf
307 auto handleScvf = [&] (const auto& scvf, const auto& indexSet, bool isSecondary)
308 {
309 const auto& scvfIndices = indexSet.gridScvfIndices();
310 if ( indexSet.nodalIndexSet().numBoundaryScvfs() > 0
311 && !std::any_of(scvfIndices.begin(), scvfIndices.end(), scvfHandled) )
312 {
313 boundaryCacheData_.cacheScvfIndices_.insert(boundaryCacheData_.cacheScvfIndices_.end(),
314 scvfIndices.begin(),
315 scvfIndices.end());
316 numCaches += scvfIndices.size();
317 if (isSecondary) numSecondaryIv++;
318 else numPrimaryIv++;
319 }
320 };
321
322 // search for ivs at boundary vertices
323 for (const auto& scvf : scvfs(fvGeometry))
324 gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()) ?
325 handleScvf(scvf, gridIvIndexSets.secondaryIndexSet(scvf), true) :
326 handleScvf(scvf, gridIvIndexSets.primaryIndexSet(scvf), false) ;
327
328 // skip the rest if there are no boundary caches to be created
329 if (numCaches > 0)
330 {
331 const auto& assemblyMapI = gridGeometry.connectivityMap()[gridGeometry.elementMapper().index(element)];
332
333 for (const auto& dataJ : assemblyMapI)
334 {
335 for (const auto& scvfJIdx : dataJ.scvfsJ)
336 {
337 const auto& scvfJ = fvGeometry.scvf(scvfJIdx);
338 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvfJ.vertexIndex()))
339 handleScvf(scvfJ, gridIvIndexSets.secondaryIndexSet(scvfJ), true);
340 else
341 handleScvf(scvfJ, gridIvIndexSets.primaryIndexSet(scvfJ), false);
342 }
343 }
344
345 // now prepare the caches of all scvfs that have been found to be handled
346 boundaryCacheData_.ivDataStorage_.primaryInteractionVolumes.reserve(numPrimaryIv);
347 boundaryCacheData_.ivDataStorage_.secondaryInteractionVolumes.reserve(numSecondaryIv);
348 boundaryCacheData_.ivDataStorage_.primaryDataHandles.reserve(numPrimaryIv);
349 boundaryCacheData_.ivDataStorage_.secondaryDataHandles.reserve(numSecondaryIv);
350
351 boundaryCacheData_.fluxVarCaches_.resize(numCaches);
352 for (auto& cache : boundaryCacheData_.fluxVarCaches_)
353 cache.setUpdateStatus(false);
354
355 FluxVariablesCacheFiller filler(gridFluxVarsCachePtr_->problem());
356 for (auto scvfIdx : boundaryCacheData_.cacheScvfIndices_)
357 {
358 const auto& scvf = fvGeometry.scvf(scvfIdx);
359 auto& cache = boundaryCacheData_[scvf];
360 if (!cache.isUpdated())
361 filler.fill(boundaryCacheData_, cache, boundaryCacheData_.ivDataStorage_,
362 fvGeometry, elemVolVars, scvf, /*forceUpdate*/true);
363 }
364 }
365 }
366
368 template<class FVElementGeometry, class ElementVolumeVariables>
369 void bindScvf_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
370 const FVElementGeometry& fvGeometry,
371 const ElementVolumeVariables& elemVolVars,
372 const typename FVElementGeometry::SubControlVolumeFace& scvf)
373 { DUNE_THROW(Dune::NotImplemented, "Scvf-local binding of the flux variables cache in mpfa schemes"); }
374
376 template<class SubControlVolumeFace>
377 bool isEmbeddedInBoundaryIV_(const SubControlVolumeFace& scvf) const
378 {
379 const auto& gridGeometry = gridFluxVarsCachePtr_->problem().gridGeometry();
380 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
381 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
382 return gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs() > 0;
383 else
384 return gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs() > 0;
385 }
386
387 const GridFluxVariablesCache* gridFluxVarsCachePtr_;
388
389 // we store those caches that touch the boundary locally here
390 // for the case that the boundary conditions change, which would
391 // leave the grid-wide cache outdated.
392 BoundaryCacheData boundaryCacheData_;
393};
394
399template<class GFVC>
401{
403 using FluxVariablesCacheFiller = typename GFVC::Traits::FluxVariablesCacheFiller;
404
405public:
407 using PrimaryInteractionVolume = typename GFVC::PrimaryInteractionVolume;
408 using SecondaryInteractionVolume = typename GFVC::SecondaryInteractionVolume;
409
411 using PrimaryIvDataHandle = typename GFVC::PrimaryIvDataHandle;
412 using SecondaryIvDataHandle = typename GFVC::SecondaryIvDataHandle;
413
415 using FluxVariablesCache = typename GFVC::FluxVariablesCache;
416
417
419 static constexpr bool cachingEnabled = false;
420
423
425 : gridFluxVarsCachePtr_(&global) {}
426
432 template<class FVElementGeometry, class ElementVolumeVariables>
433 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
434 const FVElementGeometry& fvGeometry,
435 const ElementVolumeVariables& elemVolVars) &
436 {
437 // For mpfa schemes we will have to prepare the caches of all scvfs that are
438 // embedded in the interaction volumes in which the element-local scvfs are embedded
439 DUNE_THROW(Dune::NotImplemented, "Local element binding of the flux variables cache in mpfa schemes");
440 }
441
447 template<class FVElementGeometry, class ElementVolumeVariables>
448 CCMpfaElementFluxVariablesCache bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
449 const FVElementGeometry& fvGeometry,
450 const ElementVolumeVariables& elemVolVars) &&
451 {
452 this->bindElement(element, fvGeometry, elemVolVars);
453 return std::move(*this);
454 }
455
461 template<class FVElementGeometry, class ElementVolumeVariables>
462 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
463 const FVElementGeometry& fvGeometry,
464 const ElementVolumeVariables& elemVolVars) &
465 {
466 clear_();
467
468 // some references for convenience
469 const auto& problem = gridFluxVarsCache().problem();
470 const auto& gridGeometry = fvGeometry.gridGeometry();
471
472 // the assembly map of the given element
473 const auto& assemblyMapI = gridGeometry.connectivityMap()[gridGeometry.elementMapper().index(element)];
474
475 // reserve memory for scvf index container
476 unsigned int numNeighborScvfs = 0;
477 for (const auto& dataJ : assemblyMapI)
478 numNeighborScvfs += dataJ.scvfsJ.size();
479 globalScvfIndices_.resize(fvGeometry.numScvf() + numNeighborScvfs);
480
481 // set the scvf indices in scvf index container
482 unsigned int i = 0;
483 for (const auto& scvf : scvfs(fvGeometry))
484 globalScvfIndices_[i++] = scvf.index();
485 for (const auto& dataJ : assemblyMapI)
486 for (auto scvfIdx : dataJ.scvfsJ)
487 globalScvfIndices_[i++] = scvfIdx;
488
489 // Reserve memory (over-) estimate for interaction volumes and corresponding data.
490 // The overestimate doesn't hurt as we are not in a memory-limited configuration.
491 // We need to avoid reallocation because in the caches we store pointers to the data handles.
492 // Default -> each facet has two neighbors (local adaption) and all scvfs belongs to different ivs.
493 // If you want to use higher local differences change the parameter below.
494 constexpr auto numIvEstimate = FVElementGeometry::maxNumElementScvfs
495 * GridFluxVariablesCache::Traits::maxLocalElementLevelDifference();
496 ivDataStorage_.primaryInteractionVolumes.reserve(numIvEstimate);
497 ivDataStorage_.secondaryInteractionVolumes.reserve(numIvEstimate);
498 ivDataStorage_.primaryDataHandles.reserve(numIvEstimate);
499 ivDataStorage_.secondaryDataHandles.reserve(numIvEstimate);
500
501 // helper class to fill flux variables caches
502 FluxVariablesCacheFiller filler(problem);
503
504 // resize the cache container
505 fluxVarsCache_.resize(globalScvfIndices_.size());
506
507 // go through the caches and fill them
508 i = 0;
509 for (const auto& scvf : scvfs(fvGeometry))
510 {
511 auto& scvfCache = fluxVarsCache_[i++];
512 if (!scvfCache.isUpdated())
513 filler.fill(*this, scvfCache, ivDataStorage_, fvGeometry, elemVolVars, scvf, true);
514 }
515
516 for (const auto& dataJ : assemblyMapI)
517 {
518 for (const auto scvfIdx : dataJ.scvfsJ)
519 {
520 auto& scvfCache = fluxVarsCache_[i++];
521 if (!scvfCache.isUpdated())
522 filler.fill(*this, scvfCache, ivDataStorage_, fvGeometry, elemVolVars, fvGeometry.scvf(scvfIdx), true);
523 }
524 }
525 }
526
532 template<class FVElementGeometry, class ElementVolumeVariables>
533 CCMpfaElementFluxVariablesCache bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
534 const FVElementGeometry& fvGeometry,
535 const ElementVolumeVariables& elemVolVars) &&
536 {
537 this->bind(element, fvGeometry, elemVolVars);
538 return std::move(*this);
539 }
540
546 template<class FVElementGeometry, class ElementVolumeVariables>
547 void bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
548 const FVElementGeometry& fvGeometry,
549 const ElementVolumeVariables& elemVolVars,
550 const typename FVElementGeometry::SubControlVolumeFace& scvf) &
551 {
552 // For mpfa schemes we will have to prepare the caches of all
553 // scvfs that are embedded in the interaction volumes this scvf is embedded
554 DUNE_THROW(Dune::NotImplemented, "Scvf-local binding of the flux variables cache in mpfa schemes");
555 }
556
562 template<class FVElementGeometry, class ElementVolumeVariables>
563 CCMpfaElementFluxVariablesCache bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
564 const FVElementGeometry& fvGeometry,
565 const ElementVolumeVariables& elemVolVars,
566 const typename FVElementGeometry::SubControlVolumeFace& scvf) &&
567 {
568 this->bindScvf(element, fvGeometry, elemVolVars, scvf);
569 return std::move(*this);
570 }
571
576 template<class FVElementGeometry, class ElementVolumeVariables>
577 void update(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
578 const FVElementGeometry& fvGeometry,
579 const ElementVolumeVariables& elemVolVars)
580 {
581 // Update only if the filler puts
582 // solution-dependent stuff into the caches
583 if (FluxVariablesCacheFiller::isSolDependent)
584 {
585 const auto& problem = gridFluxVarsCache().problem();
586 const auto& gridGeometry = fvGeometry.gridGeometry();
587 const auto& assemblyMapI = gridGeometry.connectivityMap()[gridGeometry.elementMapper().index(element)];
588
589 // helper class to fill flux variables caches
590 FluxVariablesCacheFiller filler(problem);
591
592 // set all the caches to "outdated"
593 for (auto& cache : fluxVarsCache_)
594 cache.setUpdateStatus(false);
595
596 // go through the caches maybe update them
597 unsigned int i = 0;
598 for (const auto& scvf : scvfs(fvGeometry))
599 {
600 auto& scvfCache = fluxVarsCache_[i++];
601 if (!scvfCache.isUpdated())
602 filler.fill(*this, scvfCache, ivDataStorage_, fvGeometry, elemVolVars, scvf);
603 }
604
605 for (const auto& dataJ : assemblyMapI)
606 {
607 for (const auto scvfIdx : dataJ.scvfsJ)
608 {
609 auto& scvfCache = fluxVarsCache_[i++];
610 if (!scvfCache.isUpdated())
611 filler.fill(*this, scvfCache, ivDataStorage_, fvGeometry, elemVolVars, fvGeometry.scvf(scvfIdx));
612 }
613 }
614 }
615 }
616
618 template<class SubControlVolumeFace, typename std::enable_if_t<!std::is_integral<SubControlVolumeFace>::value, int> = 0>
619 const FluxVariablesCache& operator [](const SubControlVolumeFace& scvf) const
620 { return fluxVarsCache_[getLocalScvfIdx_(scvf.index())]; }
621
623 const FluxVariablesCache& operator [](const std::size_t scvfIdx) const
624 { return fluxVarsCache_[getLocalScvfIdx_(scvfIdx)]; }
625
627 template<class SubControlVolumeFace, typename std::enable_if_t<!std::is_integral<SubControlVolumeFace>::value, int> = 0>
628 FluxVariablesCache& operator [](const SubControlVolumeFace& scvf)
629 { return fluxVarsCache_[getLocalScvfIdx_(scvf.index())]; }
630
632 FluxVariablesCache& operator [](const std::size_t scvfIdx)
633 { return fluxVarsCache_[getLocalScvfIdx_(scvfIdx)]; }
634
636 template<class SubControlVolumeFace>
637 const PrimaryInteractionVolume& primaryInteractionVolume(const SubControlVolumeFace& scvf) const
638 { return ivDataStorage_.primaryInteractionVolumes[ (*this)[scvf].ivIndexInContainer() ]; }
639
641 template<class SubControlVolumeFace>
642 const PrimaryIvDataHandle& primaryDataHandle(const SubControlVolumeFace& scvf) const
643 { return ivDataStorage_.primaryDataHandles[ (*this)[scvf].ivIndexInContainer() ]; }
644
646 template<class SubControlVolumeFace>
647 const SecondaryInteractionVolume& secondaryInteractionVolume(const SubControlVolumeFace& scvf) const
648 { return ivDataStorage_.secondaryInteractionVolumes[ (*this)[scvf].ivIndexInContainer() ]; }
649
651 template<class SubControlVolumeFace>
652 const SecondaryIvDataHandle& secondaryDataHandle(const SubControlVolumeFace& scvf) const
653 { return ivDataStorage_.secondaryDataHandles[ (*this)[scvf].ivIndexInContainer() ]; }
654
657 { return *gridFluxVarsCachePtr_; }
658
659private:
660 const GridFluxVariablesCache* gridFluxVarsCachePtr_;
661
663 void clear_()
664 {
665 fluxVarsCache_.clear();
666 globalScvfIndices_.clear();
667 ivDataStorage_.primaryInteractionVolumes.clear();
668 ivDataStorage_.secondaryInteractionVolumes.clear();
669 ivDataStorage_.primaryDataHandles.clear();
670 ivDataStorage_.secondaryDataHandles.clear();
671 }
672
674 unsigned int getLocalScvfIdx_(const int scvfIdx) const
675 {
676 auto it = std::find(globalScvfIndices_.begin(), globalScvfIndices_.end(), scvfIdx);
677 assert(it != globalScvfIndices_.end() && "Could not find the flux vars cache for scvfIdx");
678 return std::distance(globalScvfIndices_.begin(), it);
679 }
680
681 // the local flux vars caches and corresponding indices
682 std::vector<FluxVariablesCache> fluxVarsCache_;
683 std::vector<std::size_t> globalScvfIndices_;
684
685 // stored interaction volumes and handles
686 using IVDataStorage = InteractionVolumeDataStorage<PrimaryInteractionVolume,
687 PrimaryIvDataHandle,
688 SecondaryInteractionVolume,
689 SecondaryIvDataHandle>;
690 IVDataStorage ivDataStorage_;
691};
692
693} // end namespace
694
695#endif
const SecondaryInteractionVolume & secondaryInteractionVolume(const SubControlVolumeFace &scvf) const
access to the interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:647
CCMpfaElementFluxVariablesCache bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:563
const PrimaryIvDataHandle & primaryDataHandle(const SubControlVolumeFace &scvf) const
access to the data handle of an interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:642
void update(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars)
Update the transmissibilities if the volume variables have changed.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:577
const PrimaryInteractionVolume & primaryInteractionVolume(const SubControlVolumeFace &scvf) const
access to the interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:637
CCMpfaElementFluxVariablesCache bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:448
typename GFVC::SecondaryIvDataHandle SecondaryIvDataHandle
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:412
typename GFVC::PrimaryIvDataHandle PrimaryIvDataHandle
export the data handle types used
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:411
const SecondaryIvDataHandle & secondaryDataHandle(const SubControlVolumeFace &scvf) const
access to the data handle of an interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:652
typename GFVC::PrimaryInteractionVolume PrimaryInteractionVolume
export the interaction volume types
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:407
const GridFluxVariablesCache & gridFluxVarsCache() const
The global object we are a restriction of.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:656
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Prepares the transmissibilities of the scv faces in an element.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:433
CCMpfaElementFluxVariablesCache bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:533
CCMpfaElementFluxVariablesCache(const GridFluxVariablesCache &global)
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:424
typename GFVC::FluxVariablesCache FluxVariablesCache
export the flux variable cache type
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:415
GFVC GridFluxVariablesCache
export the type of the grid flux variables
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:422
typename GFVC::SecondaryInteractionVolume SecondaryInteractionVolume
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:408
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Prepares the transmissibilities of the scv faces in the stencil of an element.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:462
void bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &
Prepares the transmissibilities of a single scv face.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:547
The flux variables caches for an element with caching enabled.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:58
GFVC GridFluxVariablesCache
export the type of the grid flux variables
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:75
typename GFVC::FluxVariablesCache FluxVariablesCache
export the flux variable cache type
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:69
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Specialization for the global caching being enabled - do nothing here.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:171
CCMpfaElementFluxVariablesCache bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:161
typename GFVC::PrimaryInteractionVolume PrimaryInteractionVolume
export the interaction volume types
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:61
void bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &
Specialization for the global caching being enabled - do nothing here.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:193
const GridFluxVariablesCache & gridFluxVarsCache() const
The global object we are a restriction of.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:269
const PrimaryIvDataHandle & primaryDataHandle(const SubControlVolumeFace &scvf) const
access to the data handle of an interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:243
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Specialization for the global caching being enabled - do nothing here.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:150
const SecondaryIvDataHandle & secondaryDataHandle(const SubControlVolumeFace &scvf) const
access to the data handle of an interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:261
typename GFVC::SecondaryIvDataHandle SecondaryIvDataHandle
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:66
CCMpfaElementFluxVariablesCache bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:140
const SecondaryInteractionVolume & secondaryInteractionVolume(const SubControlVolumeFace &scvf) const
access to the interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:252
CCMpfaElementFluxVariablesCache bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:182
const PrimaryInteractionVolume & primaryInteractionVolume(const SubControlVolumeFace &scvf) const
access to the interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:234
typename GFVC::PrimaryIvDataHandle PrimaryIvDataHandle
export the data handle types used
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:65
void update(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars)
Specialization for the global caching being enabled - do nothing here.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:201
CCMpfaElementFluxVariablesCache(const GridFluxVariablesCache &global)
The constructor.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:130
typename GFVC::SecondaryInteractionVolume SecondaryInteractionVolume
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:62
The flux variables caches for an element.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:50
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
Structure to store interaction volumes and data handles.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:32
std::vector< PrimaryIVDataHandle > primaryDataHandles
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:36
std::vector< SecondaryIV > secondaryInteractionVolumes
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:34
std::vector< PrimaryIV > primaryInteractionVolumes
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:33
std::vector< SecondaryIVDataHandle > secondaryDataHandles
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:37