3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_CCMPFA_ELEMENT_FLUXVARSCACHE_HH
25#define DUMUX_DISCRETIZATION_CCMPFA_ELEMENT_FLUXVARSCACHE_HH
26
27#include <algorithm>
28#include <cassert>
29#include <type_traits>
30#include <vector>
31#include <utility>
32
33#include <dune/common/exceptions.hh>
34
35namespace Dumux {
36
41template<class PrimaryIV, class PrimaryIVDataHandle,
42 class SecondaryIV, class SecondaryIVDataHandle>
44{
45 std::vector<PrimaryIV> primaryInteractionVolumes;
46 std::vector<SecondaryIV> secondaryInteractionVolumes;
47
48 std::vector<PrimaryIVDataHandle> primaryDataHandles;
49 std::vector<SecondaryIVDataHandle> secondaryDataHandles;
50};
51
61template<class GFVC, bool cachingEnabled>
63
68template<class GFVC>
70{
71public:
73 using PrimaryInteractionVolume = typename GFVC::PrimaryInteractionVolume;
74 using SecondaryInteractionVolume = typename GFVC::SecondaryInteractionVolume;
75
77 using PrimaryIvDataHandle = typename GFVC::PrimaryIvDataHandle;
78 using SecondaryIvDataHandle = typename GFVC::SecondaryIvDataHandle;
79
81 using FluxVariablesCache = typename GFVC::FluxVariablesCache;
82
84 static constexpr bool cachingEnabled = true;
85
88
89private:
91 using FluxVariablesCacheFiller = typename GFVC::Traits::FluxVariablesCacheFiller;
92
94 class BoundaryCacheData
95 {
96 // allow the element flux variables class access to private members
98
99 public:
101 template<class SubControlVolumeFace>
102 const FluxVariablesCache& operator [](const SubControlVolumeFace& scvf) const
103 { return fluxVarCaches_[getLocalIdx_(scvf.index())]; }
104
105 template<class SubControlVolumeFace>
106 FluxVariablesCache& operator [](const SubControlVolumeFace& scvf)
107 { return fluxVarCaches_[getLocalIdx_(scvf.index())]; }
108
110 void clear()
111 {
112 fluxVarCaches_.clear();
113 cacheScvfIndices_.clear();
114 ivDataStorage_.primaryInteractionVolumes.clear();
115 ivDataStorage_.secondaryInteractionVolumes.clear();
116 ivDataStorage_.primaryDataHandles.clear();
117 ivDataStorage_.secondaryDataHandles.clear();
118 }
119
120 private:
122 int getLocalIdx_(const int scvfIdx) const
123 {
124 auto it = std::find(cacheScvfIndices_.begin(), cacheScvfIndices_.end(), scvfIdx);
125 assert(it != cacheScvfIndices_.end() && "Could not find the local idx for the given scvf idx!");
126 return std::distance(cacheScvfIndices_.begin(), it);
127 }
128
129 std::vector<std::size_t> cacheScvfIndices_;
130 std::vector<FluxVariablesCache> fluxVarCaches_;
131
132 // stored boundary interaction volumes and handles
137 IVDataStorage ivDataStorage_;
138 };
139
140public:
143 : gridFluxVarsCachePtr_(&global)
144 {}
145
151 template<class FVElementGeometry, class ElementVolumeVariables>
152 CCMpfaElementFluxVariablesCache bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
153 const FVElementGeometry& fvGeometry,
154 const ElementVolumeVariables& elemVolVars) &&
155 {
156 this->bindElement_(element, fvGeometry, elemVolVars);
157 return std::move(*this);
158 }
159
161 template<class FVElementGeometry, class ElementVolumeVariables>
162 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
163 const FVElementGeometry& fvGeometry,
164 const ElementVolumeVariables& elemVolVars) &
165 { this->bindElement_(element, fvGeometry, elemVolVars); }
166
172 template<class FVElementGeometry, class ElementVolumeVariables>
173 CCMpfaElementFluxVariablesCache bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
174 const FVElementGeometry& fvGeometry,
175 const ElementVolumeVariables& elemVolVars) &&
176 {
177 this->bind_(element, fvGeometry, elemVolVars);
178 return std::move(*this);
179 }
180
182 template<class FVElementGeometry, class ElementVolumeVariables>
183 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
184 const FVElementGeometry& fvGeometry,
185 const ElementVolumeVariables& elemVolVars) &
186 { this->bind_(element, fvGeometry, elemVolVars); }
187
193 template<class FVElementGeometry, class ElementVolumeVariables>
194 CCMpfaElementFluxVariablesCache bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
195 const FVElementGeometry& fvGeometry,
196 const ElementVolumeVariables& elemVolVars,
197 const typename FVElementGeometry::SubControlVolumeFace& scvf) &&
198 {
199 this->bindScvf_(element, fvGeometry, elemVolVars, scvf);
200 return std::move(*this);
201 }
202
204 template<class FVElementGeometry, class ElementVolumeVariables>
205 void bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
206 const FVElementGeometry& fvGeometry,
207 const ElementVolumeVariables& elemVolVars,
208 const typename FVElementGeometry::SubControlVolumeFace& scvf) &
209 { this->bindScvf_(element, fvGeometry, elemVolVars, scvf); }
210
212 template<class FVElementGeometry, class ElementVolumeVariables>
213 void update(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
214 const FVElementGeometry& fvGeometry,
215 const ElementVolumeVariables& elemVolVars)
216 {
217 // Update only if the filler puts solution-dependent stuff into the caches
218 if (FluxVariablesCacheFiller::isSolDependent)
219 {
220 // helper class to fill flux variables caches
221 FluxVariablesCacheFiller filler(gridFluxVarsCachePtr_->problem());
222
223 // first, set all the caches to "outdated"
224 for (auto& cache : boundaryCacheData_.fluxVarCaches_)
225 cache.setUpdateStatus(false);
226
227 // go through the caches maybe update them
228 for (auto scvfIdx : boundaryCacheData_.cacheScvfIndices_)
229 {
230 const auto& scvf = fvGeometry.scvf(scvfIdx);
231 auto& scvfCache = boundaryCacheData_[scvf];
232 if (!scvfCache.isUpdated())
233 filler.fill(boundaryCacheData_, scvfCache, boundaryCacheData_.ivDataStorage_,
234 fvGeometry, elemVolVars, scvf);
235 }
236 }
237 }
238
240 template<class SubControlVolumeFace>
241 const FluxVariablesCache& operator [](const SubControlVolumeFace& scvf) const
242 { return !isEmbeddedInBoundaryIV_(scvf) ? (*gridFluxVarsCachePtr_)[scvf] : boundaryCacheData_[scvf]; }
243
245 template<class SubControlVolumeFace>
246 const PrimaryInteractionVolume& primaryInteractionVolume(const SubControlVolumeFace& scvf) const
247 {
248 return isEmbeddedInBoundaryIV_(scvf)
249 ? boundaryCacheData_.ivDataStorage_.primaryInteractionVolumes[ (*this)[scvf].ivIndexInContainer() ]
250 : gridFluxVarsCachePtr_->primaryInteractionVolume(scvf);
251 }
252
254 template<class SubControlVolumeFace>
255 const PrimaryIvDataHandle& primaryDataHandle(const SubControlVolumeFace& scvf) const
256 {
257 return isEmbeddedInBoundaryIV_(scvf)
258 ? boundaryCacheData_.ivDataStorage_.primaryDataHandles[ (*this)[scvf].ivIndexInContainer() ]
259 : gridFluxVarsCachePtr_->primaryDataHandle(scvf);
260 }
261
263 template<class SubControlVolumeFace>
264 const SecondaryInteractionVolume& secondaryInteractionVolume(const SubControlVolumeFace& scvf) const
265 {
266 return isEmbeddedInBoundaryIV_(scvf)
267 ? boundaryCacheData_.ivDataStorage_.secondaryInteractionVolumes[ (*this)[scvf].ivIndexInContainer() ]
268 : gridFluxVarsCachePtr_->secondaryInteractionVolume(scvf);
269 }
270
272 template<class SubControlVolumeFace>
273 const SecondaryIvDataHandle& secondaryDataHandle(const SubControlVolumeFace& scvf) const
274 {
275 return isEmbeddedInBoundaryIV_(scvf)
276 ? boundaryCacheData_.ivDataStorage_.secondaryDataHandles[ (*this)[scvf].ivIndexInContainer() ]
277 : gridFluxVarsCachePtr_->secondaryDataHandle(scvf);
278 }
279
282 { return *gridFluxVarsCachePtr_; }
283
284private:
285
286
288 template<class FVElementGeometry, class ElementVolumeVariables>
289 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
290 const FVElementGeometry& fvGeometry,
291 const ElementVolumeVariables& elemVolVars)
292 { DUNE_THROW(Dune::NotImplemented, "Local element binding of the flux variables cache in mpfa schemes"); }
293
295 template<class FVElementGeometry, class ElementVolumeVariables>
296 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
297 const FVElementGeometry& fvGeometry,
298 const ElementVolumeVariables& elemVolVars)
299 {
300 boundaryCacheData_.clear();
301
302 // find out how much memory needs to be reserved
303 std::size_t numPrimaryIv; numPrimaryIv = 0;
304 std::size_t numSecondaryIv; numSecondaryIv = 0;
305 std::size_t numCaches; numCaches = 0;
306
307 const auto& gridGeometry = fvGeometry.gridGeometry();
308 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
309
310 // lambda to check if a scvf was handled already
311 auto scvfHandled = [&] (auto idx)
312 {
313 return std::find(boundaryCacheData_.cacheScvfIndices_.begin(),
314 boundaryCacheData_.cacheScvfIndices_.end(),
315 idx) != boundaryCacheData_.cacheScvfIndices_.end();
316 };
317
318 // lambda to increase counters for a given scvf
319 auto handleScvf = [&] (const auto& scvf, const auto& indexSet, bool isSecondary)
320 {
321 const auto& scvfIndices = indexSet.gridScvfIndices();
322 if ( indexSet.nodalIndexSet().numBoundaryScvfs() > 0
323 && !std::any_of(scvfIndices.begin(), scvfIndices.end(), scvfHandled) )
324 {
325 boundaryCacheData_.cacheScvfIndices_.insert(boundaryCacheData_.cacheScvfIndices_.end(),
326 scvfIndices.begin(),
327 scvfIndices.end());
328 numCaches += scvfIndices.size();
329 if (isSecondary) numSecondaryIv++;
330 else numPrimaryIv++;
331 }
332 };
333
334 // search for ivs at boundary vertices
335 for (const auto& scvf : scvfs(fvGeometry))
336 gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()) ?
337 handleScvf(scvf, gridIvIndexSets.secondaryIndexSet(scvf), true) :
338 handleScvf(scvf, gridIvIndexSets.primaryIndexSet(scvf), false) ;
339
340 // skip the rest if there are no boundary caches to be created
341 if (numCaches > 0)
342 {
343 const auto& assemblyMapI = gridGeometry.connectivityMap()[gridGeometry.elementMapper().index(element)];
344
345 for (const auto& dataJ : assemblyMapI)
346 {
347 for (const auto& scvfJIdx : dataJ.scvfsJ)
348 {
349 const auto& scvfJ = fvGeometry.scvf(scvfJIdx);
350 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvfJ.vertexIndex()))
351 handleScvf(scvfJ, gridIvIndexSets.secondaryIndexSet(scvfJ), true);
352 else
353 handleScvf(scvfJ, gridIvIndexSets.primaryIndexSet(scvfJ), false);
354 }
355 }
356
357 // now prepare the caches of all scvfs that have been found to be handled
358 boundaryCacheData_.ivDataStorage_.primaryInteractionVolumes.reserve(numPrimaryIv);
359 boundaryCacheData_.ivDataStorage_.secondaryInteractionVolumes.reserve(numSecondaryIv);
360 boundaryCacheData_.ivDataStorage_.primaryDataHandles.reserve(numPrimaryIv);
361 boundaryCacheData_.ivDataStorage_.secondaryDataHandles.reserve(numSecondaryIv);
362
363 boundaryCacheData_.fluxVarCaches_.resize(numCaches);
364 for (auto& cache : boundaryCacheData_.fluxVarCaches_)
365 cache.setUpdateStatus(false);
366
367 FluxVariablesCacheFiller filler(gridFluxVarsCachePtr_->problem());
368 for (auto scvfIdx : boundaryCacheData_.cacheScvfIndices_)
369 {
370 const auto& scvf = fvGeometry.scvf(scvfIdx);
371 auto& cache = boundaryCacheData_[scvf];
372 if (!cache.isUpdated())
373 filler.fill(boundaryCacheData_, cache, boundaryCacheData_.ivDataStorage_,
374 fvGeometry, elemVolVars, scvf, /*forceUpdate*/true);
375 }
376 }
377 }
378
380 template<class FVElementGeometry, class ElementVolumeVariables>
381 void bindScvf_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
382 const FVElementGeometry& fvGeometry,
383 const ElementVolumeVariables& elemVolVars,
384 const typename FVElementGeometry::SubControlVolumeFace& scvf)
385 { DUNE_THROW(Dune::NotImplemented, "Scvf-local binding of the flux variables cache in mpfa schemes"); }
386
388 template<class SubControlVolumeFace>
389 bool isEmbeddedInBoundaryIV_(const SubControlVolumeFace& scvf) const
390 {
391 const auto& gridGeometry = gridFluxVarsCachePtr_->problem().gridGeometry();
392 const auto& gridIvIndexSets = gridGeometry.gridInteractionVolumeIndexSets();
393 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
394 return gridIvIndexSets.secondaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs() > 0;
395 else
396 return gridIvIndexSets.primaryIndexSet(scvf).nodalIndexSet().numBoundaryScvfs() > 0;
397 }
398
399 const GridFluxVariablesCache* gridFluxVarsCachePtr_;
400
401 // we store those caches that touch the boundary locally here
402 // for the case that the boundary conditions change, which would
403 // leave the grid-wide cache outdated.
404 BoundaryCacheData boundaryCacheData_;
405};
406
411template<class GFVC>
413{
415 using FluxVariablesCacheFiller = typename GFVC::Traits::FluxVariablesCacheFiller;
416
417public:
419 using PrimaryInteractionVolume = typename GFVC::PrimaryInteractionVolume;
420 using SecondaryInteractionVolume = typename GFVC::SecondaryInteractionVolume;
421
423 using PrimaryIvDataHandle = typename GFVC::PrimaryIvDataHandle;
424 using SecondaryIvDataHandle = typename GFVC::SecondaryIvDataHandle;
425
427 using FluxVariablesCache = typename GFVC::FluxVariablesCache;
428
429
431 static constexpr bool cachingEnabled = false;
432
435
437 : gridFluxVarsCachePtr_(&global) {}
438
444 template<class FVElementGeometry, class ElementVolumeVariables>
445 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
446 const FVElementGeometry& fvGeometry,
447 const ElementVolumeVariables& elemVolVars) &
448 {
449 // For mpfa schemes we will have to prepare the caches of all scvfs that are
450 // embedded in the interaction volumes in which the element-local scvfs are embedded
451 DUNE_THROW(Dune::NotImplemented, "Local element binding of the flux variables cache in mpfa schemes");
452 }
453
459 template<class FVElementGeometry, class ElementVolumeVariables>
460 CCMpfaElementFluxVariablesCache bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
461 const FVElementGeometry& fvGeometry,
462 const ElementVolumeVariables& elemVolVars) &&
463 {
464 this->bindElement(element, fvGeometry, elemVolVars);
465 return std::move(*this);
466 }
467
473 template<class FVElementGeometry, class ElementVolumeVariables>
474 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
475 const FVElementGeometry& fvGeometry,
476 const ElementVolumeVariables& elemVolVars) &
477 {
478 clear_();
479
480 // some references for convenience
481 const auto& problem = gridFluxVarsCache().problem();
482 const auto& gridGeometry = fvGeometry.gridGeometry();
483
484 // the assembly map of the given element
485 const auto& assemblyMapI = gridGeometry.connectivityMap()[gridGeometry.elementMapper().index(element)];
486
487 // reserve memory for scvf index container
488 unsigned int numNeighborScvfs = 0;
489 for (const auto& dataJ : assemblyMapI)
490 numNeighborScvfs += dataJ.scvfsJ.size();
491 globalScvfIndices_.resize(fvGeometry.numScvf() + numNeighborScvfs);
492
493 // set the scvf indices in scvf index container
494 unsigned int i = 0;
495 for (const auto& scvf : scvfs(fvGeometry))
496 globalScvfIndices_[i++] = scvf.index();
497 for (const auto& dataJ : assemblyMapI)
498 for (auto scvfIdx : dataJ.scvfsJ)
499 globalScvfIndices_[i++] = scvfIdx;
500
501 // Reserve memory (over-) estimate for interaction volumes and corresponding data.
502 // The overestimate doesn't hurt as we are not in a memory-limited configuration.
503 // We need to avoid reallocation because in the caches we store pointers to the data handles.
504 // Default -> each facet has two neighbors (local adaption) and all scvfs belongs to different ivs.
505 // If you want to use higher local differences change the parameter below.
506 constexpr auto numIvEstimate = FVElementGeometry::maxNumElementScvfs
507 * GridFluxVariablesCache::Traits::maxLocalElementLevelDifference();
508 ivDataStorage_.primaryInteractionVolumes.reserve(numIvEstimate);
509 ivDataStorage_.secondaryInteractionVolumes.reserve(numIvEstimate);
510 ivDataStorage_.primaryDataHandles.reserve(numIvEstimate);
511 ivDataStorage_.secondaryDataHandles.reserve(numIvEstimate);
512
513 // helper class to fill flux variables caches
514 FluxVariablesCacheFiller filler(problem);
515
516 // resize the cache container
517 fluxVarsCache_.resize(globalScvfIndices_.size());
518
519 // go through the caches and fill them
520 i = 0;
521 for (const auto& scvf : scvfs(fvGeometry))
522 {
523 auto& scvfCache = fluxVarsCache_[i++];
524 if (!scvfCache.isUpdated())
525 filler.fill(*this, scvfCache, ivDataStorage_, fvGeometry, elemVolVars, scvf, true);
526 }
527
528 for (const auto& dataJ : assemblyMapI)
529 {
530 for (const auto scvfIdx : dataJ.scvfsJ)
531 {
532 auto& scvfCache = fluxVarsCache_[i++];
533 if (!scvfCache.isUpdated())
534 filler.fill(*this, scvfCache, ivDataStorage_, fvGeometry, elemVolVars, fvGeometry.scvf(scvfIdx), true);
535 }
536 }
537 }
538
544 template<class FVElementGeometry, class ElementVolumeVariables>
545 CCMpfaElementFluxVariablesCache bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
546 const FVElementGeometry& fvGeometry,
547 const ElementVolumeVariables& elemVolVars) &&
548 {
549 this->bind(element, fvGeometry, elemVolVars);
550 return std::move(*this);
551 }
552
558 template<class FVElementGeometry, class ElementVolumeVariables>
559 void bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
560 const FVElementGeometry& fvGeometry,
561 const ElementVolumeVariables& elemVolVars,
562 const typename FVElementGeometry::SubControlVolumeFace& scvf) &
563 {
564 // For mpfa schemes we will have to prepare the caches of all
565 // scvfs that are embedded in the interaction volumes this scvf is embedded
566 DUNE_THROW(Dune::NotImplemented, "Scvf-local binding of the flux variables cache in mpfa schemes");
567 }
568
574 template<class FVElementGeometry, class ElementVolumeVariables>
575 CCMpfaElementFluxVariablesCache bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
576 const FVElementGeometry& fvGeometry,
577 const ElementVolumeVariables& elemVolVars,
578 const typename FVElementGeometry::SubControlVolumeFace& scvf) &&
579 {
580 this->bindScvf(element, fvGeometry, elemVolVars, scvf);
581 return std::move(*this);
582 }
583
588 template<class FVElementGeometry, class ElementVolumeVariables>
589 void update(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
590 const FVElementGeometry& fvGeometry,
591 const ElementVolumeVariables& elemVolVars)
592 {
593 // Update only if the filler puts
594 // solution-dependent stuff into the caches
595 if (FluxVariablesCacheFiller::isSolDependent)
596 {
597 const auto& problem = gridFluxVarsCache().problem();
598 const auto& gridGeometry = fvGeometry.gridGeometry();
599 const auto& assemblyMapI = gridGeometry.connectivityMap()[gridGeometry.elementMapper().index(element)];
600
601 // helper class to fill flux variables caches
602 FluxVariablesCacheFiller filler(problem);
603
604 // set all the caches to "outdated"
605 for (auto& cache : fluxVarsCache_)
606 cache.setUpdateStatus(false);
607
608 // go through the caches maybe update them
609 unsigned int i = 0;
610 for (const auto& scvf : scvfs(fvGeometry))
611 {
612 auto& scvfCache = fluxVarsCache_[i++];
613 if (!scvfCache.isUpdated())
614 filler.fill(*this, scvfCache, ivDataStorage_, fvGeometry, elemVolVars, scvf);
615 }
616
617 for (const auto& dataJ : assemblyMapI)
618 {
619 for (const auto scvfIdx : dataJ.scvfsJ)
620 {
621 auto& scvfCache = fluxVarsCache_[i++];
622 if (!scvfCache.isUpdated())
623 filler.fill(*this, scvfCache, ivDataStorage_, fvGeometry, elemVolVars, fvGeometry.scvf(scvfIdx));
624 }
625 }
626 }
627 }
628
630 template<class SubControlVolumeFace, typename std::enable_if_t<!std::is_integral<SubControlVolumeFace>::value, int> = 0>
631 const FluxVariablesCache& operator [](const SubControlVolumeFace& scvf) const
632 { return fluxVarsCache_[getLocalScvfIdx_(scvf.index())]; }
633
635 const FluxVariablesCache& operator [](const std::size_t scvfIdx) const
636 { return fluxVarsCache_[getLocalScvfIdx_(scvfIdx)]; }
637
639 template<class SubControlVolumeFace, typename std::enable_if_t<!std::is_integral<SubControlVolumeFace>::value, int> = 0>
640 FluxVariablesCache& operator [](const SubControlVolumeFace& scvf)
641 { return fluxVarsCache_[getLocalScvfIdx_(scvf.index())]; }
642
644 FluxVariablesCache& operator [](const std::size_t scvfIdx)
645 { return fluxVarsCache_[getLocalScvfIdx_(scvfIdx)]; }
646
648 template<class SubControlVolumeFace>
649 const PrimaryInteractionVolume& primaryInteractionVolume(const SubControlVolumeFace& scvf) const
650 { return ivDataStorage_.primaryInteractionVolumes[ (*this)[scvf].ivIndexInContainer() ]; }
651
653 template<class SubControlVolumeFace>
654 const PrimaryIvDataHandle& primaryDataHandle(const SubControlVolumeFace& scvf) const
655 { return ivDataStorage_.primaryDataHandles[ (*this)[scvf].ivIndexInContainer() ]; }
656
658 template<class SubControlVolumeFace>
659 const SecondaryInteractionVolume& secondaryInteractionVolume(const SubControlVolumeFace& scvf) const
660 { return ivDataStorage_.secondaryInteractionVolumes[ (*this)[scvf].ivIndexInContainer() ]; }
661
663 template<class SubControlVolumeFace>
664 const SecondaryIvDataHandle& secondaryDataHandle(const SubControlVolumeFace& scvf) const
665 { return ivDataStorage_.secondaryDataHandles[ (*this)[scvf].ivIndexInContainer() ]; }
666
669 { return *gridFluxVarsCachePtr_; }
670
671private:
672 const GridFluxVariablesCache* gridFluxVarsCachePtr_;
673
675 void clear_()
676 {
677 fluxVarsCache_.clear();
678 globalScvfIndices_.clear();
679 ivDataStorage_.primaryInteractionVolumes.clear();
680 ivDataStorage_.secondaryInteractionVolumes.clear();
681 ivDataStorage_.primaryDataHandles.clear();
682 ivDataStorage_.secondaryDataHandles.clear();
683 }
684
686 unsigned int getLocalScvfIdx_(const int scvfIdx) const
687 {
688 auto it = std::find(globalScvfIndices_.begin(), globalScvfIndices_.end(), scvfIdx);
689 assert(it != globalScvfIndices_.end() && "Could not find the flux vars cache for scvfIdx");
690 return std::distance(globalScvfIndices_.begin(), it);
691 }
692
693 // the local flux vars caches and corresponding indices
694 std::vector<FluxVariablesCache> fluxVarsCache_;
695 std::vector<std::size_t> globalScvfIndices_;
696
697 // stored interaction volumes and handles
698 using IVDataStorage = InteractionVolumeDataStorage<PrimaryInteractionVolume,
699 PrimaryIvDataHandle,
700 SecondaryInteractionVolume,
701 SecondaryIvDataHandle>;
702 IVDataStorage ivDataStorage_;
703};
704
705} // end namespace
706
707#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
Definition: adapt.hh:29
Structure to store interaction volumes and data handles.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:44
std::vector< PrimaryIVDataHandle > primaryDataHandles
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:48
std::vector< SecondaryIV > secondaryInteractionVolumes
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:46
std::vector< PrimaryIV > primaryInteractionVolumes
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:45
std::vector< SecondaryIVDataHandle > secondaryDataHandles
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:49
The flux variables caches for an element.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:62
The flux variables caches for an element with caching enabled.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:70
GFVC GridFluxVariablesCache
export the type of the grid flux variables
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:87
typename GFVC::FluxVariablesCache FluxVariablesCache
export the flux variable cache type
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:81
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:183
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:173
typename GFVC::PrimaryInteractionVolume PrimaryInteractionVolume
export the interaction volume types
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:73
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:205
const GridFluxVariablesCache & gridFluxVarsCache() const
The global object we are a restriction of.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:281
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:255
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:162
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:273
typename GFVC::SecondaryIvDataHandle SecondaryIvDataHandle
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:78
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:152
const SecondaryInteractionVolume & secondaryInteractionVolume(const SubControlVolumeFace &scvf) const
access to the interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:264
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:194
const PrimaryInteractionVolume & primaryInteractionVolume(const SubControlVolumeFace &scvf) const
access to the interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:246
typename GFVC::PrimaryIvDataHandle PrimaryIvDataHandle
export the data handle types used
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:77
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:213
CCMpfaElementFluxVariablesCache(const GridFluxVariablesCache &global)
The constructor.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:142
typename GFVC::SecondaryInteractionVolume SecondaryInteractionVolume
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:74
const SecondaryInteractionVolume & secondaryInteractionVolume(const SubControlVolumeFace &scvf) const
access to the interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:659
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:575
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:654
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:589
const PrimaryInteractionVolume & primaryInteractionVolume(const SubControlVolumeFace &scvf) const
access to the interaction volume an scvf is embedded in
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:649
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:460
typename GFVC::SecondaryIvDataHandle SecondaryIvDataHandle
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:424
typename GFVC::PrimaryIvDataHandle PrimaryIvDataHandle
export the data handle types used
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:423
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:664
typename GFVC::PrimaryInteractionVolume PrimaryInteractionVolume
export the interaction volume types
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:419
const GridFluxVariablesCache & gridFluxVarsCache() const
The global object we are a restriction of.
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:668
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:445
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:545
CCMpfaElementFluxVariablesCache(const GridFluxVariablesCache &global)
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:436
typename GFVC::FluxVariablesCache FluxVariablesCache
export the flux variable cache type
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:427
GFVC GridFluxVariablesCache
export the type of the grid flux variables
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:434
typename GFVC::SecondaryInteractionVolume SecondaryInteractionVolume
Definition: discretization/cellcentered/mpfa/elementfluxvariablescache.hh:420
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:474
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:559