fvbasediscretization.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_FV_BASE_DISCRETIZATION_HH
29#define EWOMS_FV_BASE_DISCRETIZATION_HH
30
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
33#include <dune/common/version.hh>
34#include <dune/istl/bvector.hh>
35
36#include <opm/material/common/MathToolbox.hpp>
37#include <opm/material/common/Valgrind.hpp>
38#include <opm/material/densead/Math.hpp>
39
55
57
60
65
68
69#include <algorithm>
70#include <array>
71#include <cstddef>
72#include <list>
73#include <memory>
74#include <stdexcept>
75#include <sstream>
76#include <string>
77#include <type_traits>
78#include <utility>
79#include <vector>
80
81namespace Opm {
82
83template<class TypeTag>
84class FvBaseDiscretizationNoAdapt;
85
86template<class TypeTag>
87class FvBaseDiscretization;
88
89} // namespace Opm
90
91namespace Opm::Properties {
92
94template<class TypeTag>
95struct Simulator<TypeTag, TTag::FvBaseDiscretization>
97
99template<class TypeTag>
100struct VertexMapper<TypeTag, TTag::FvBaseDiscretization>
101{ using type = Dune::MultipleCodimMultipleGeomTypeMapper<GetPropType<TypeTag, Properties::GridView>>; };
102
104template<class TypeTag>
106{ using type = Dune::MultipleCodimMultipleGeomTypeMapper<GetPropType<TypeTag, Properties::GridView>>; };
107
109template<class TypeTag>
111{
114
115public:
117};
118
119template<class TypeTag>
122
123template<class TypeTag>
126
127template<class TypeTag>
130
132template<class TypeTag>
135
139template<class TypeTag>
140struct EqVector<TypeTag, TTag::FvBaseDiscretization>
141{
142 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
143 getPropValue<TypeTag, Properties::NumEq>()>;
144};
145
151template<class TypeTag>
152struct RateVector<TypeTag, TTag::FvBaseDiscretization>
154
158template<class TypeTag>
161
165template<class TypeTag>
166struct Constraints<TypeTag, TTag::FvBaseDiscretization>
168
172template<class TypeTag>
174{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::EqVector>>; };
175
179template<class TypeTag>
181{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::EqVector>>; };
182
186template<class TypeTag>
189
193template<class TypeTag>
195{ using type = Dune::BlockVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; };
196
202template<class TypeTag>
205
209template<class TypeTag>
212
213template<class TypeTag>
216
217template<class TypeTag>
220
224template<class TypeTag>
227
228template<class TypeTag>
230{ static constexpr bool value = true; };
231
235template<class TypeTag>
236struct Linearizer<TypeTag, TTag::FvBaseDiscretization>
238
240template<class TypeTag>
242{ static constexpr auto value = Dune::VTK::ascii; };
243
244// disable constraints by default
245template<class TypeTag>
247{ static constexpr bool value = false; };
248
250template<class TypeTag>
252{ static constexpr int value = 2; };
253
256template<class TypeTag>
258{ static constexpr bool value = false; };
259
260// use volumetric residuals is default
261template<class TypeTag>
263{ static constexpr bool value = true; };
264
267template<class TypeTag>
269{ static constexpr bool value = true; };
270
271template <class TypeTag, class MyTypeTag>
273
274#if !HAVE_DUNE_FEM
275template<class TypeTag>
278
279template<class TypeTag>
281{
284};
285#endif
286
287} // namespace Opm::Properties
288
289namespace Opm {
290
296template<class TypeTag>
298{
299 using Implementation = GetPropType<TypeTag, Properties::Model>;
326
329
330 enum {
331 numEq = getPropValue<TypeTag, Properties::NumEq>(),
332 historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>(),
333 };
334
335 using IntensiveQuantitiesVector = std::vector<IntensiveQuantities,
336 aligned_allocator<IntensiveQuantities,
337 alignof(IntensiveQuantities)>>;
338
339 using Element = typename GridView::template Codim<0>::Entity;
340 using ElementIterator = typename GridView::template Codim<0>::Iterator;
341
342 using Toolbox = MathToolbox<Evaluation>;
343 using VectorBlock = Dune::FieldVector<Evaluation, numEq>;
344 using EvalEqVector = Dune::FieldVector<Evaluation, numEq>;
345
346 using LocalEvalBlockVector = typename LocalResidual::LocalEvalBlockVector;
347
348public:
350 {
351 protected:
352 SolutionVector blockVector_;
353 public:
354 BlockVectorWrapper(const std::string&, const std::size_t size)
355 : blockVector_(size)
356 {}
357
359
361 {
362 BlockVectorWrapper result("dummy", 3);
363 result.blockVector_[0] = 1.0;
364 result.blockVector_[1] = 2.0;
365 result.blockVector_[2] = 3.0;
366
367 return result;
368 }
369
370 SolutionVector& blockVector()
371 { return blockVector_; }
372
373 const SolutionVector& blockVector() const
374 { return blockVector_; }
375
376 bool operator==(const BlockVectorWrapper& wrapper) const
377 {
378 return std::ranges::equal(this->blockVector_, wrapper.blockVector_);
379 }
380
381 template<class Serializer>
382 void serializeOp(Serializer& serializer)
383 {
384 serializer(blockVector_);
385 }
386 };
387
388private:
391
392public:
393 explicit FvBaseDiscretization(Simulator& simulator)
394 : simulator_(simulator)
395 , gridView_(simulator.gridView())
396 , elementMapper_(gridView_, Dune::mcmgElementLayout())
397 , vertexMapper_(gridView_, Dune::mcmgVertexLayout())
398 , newtonMethod_(simulator)
399 , localLinearizer_(ThreadManager::maxThreads())
400 , linearizer_(std::make_unique<Linearizer>())
401 , enableGridAdaptation_(Parameters::Get<Parameters::EnableGridAdaptation>() )
402 , enableIntensiveQuantityCache_(Parameters::Get<Parameters::EnableIntensiveQuantityCache>())
403 , enableStorageCache_(Parameters::Get<Parameters::EnableStorageCache>())
404 , enableThermodynamicHints_(Parameters::Get<Parameters::EnableThermodynamicHints>())
405 , cachedIntensiveQuantityHistorySize_(static_cast<unsigned>(-1))
406 {
407 const bool isEcfv = std::is_same_v<Discretization, EcfvDiscretization<TypeTag>>;
408 if (enableGridAdaptation_ && !isEcfv) {
409 throw std::invalid_argument("Grid adaptation currently only works for the "
410 "element-centered finite volume discretization (is: " +
411 Dune::className<Discretization>() + ")");
412 }
413
414 PrimaryVariables::init();
415 // Setting up the intensive quantities cache and storage cache is done in finishInit()
416 // and applyInitialSolution() to ensure that the history size is correct.
417 asImp_().registerOutputModules_();
418 }
419
420 // copying a discretization object is not a good idea
422
426 static void registerParameters()
427 {
428 Linearizer::registerParameters();
429 LocalLinearizer::registerParameters();
430 LocalResidual::registerParameters();
431 GradientCalculator::registerParameters();
432 IntensiveQuantities::registerParameters();
433 ExtensiveQuantities::registerParameters();
435 Linearizer::registerParameters();
436 PrimaryVariables::registerParameters();
437 // register runtime parameters of the output modules
439
440 Parameters::Register<Parameters::EnableGridAdaptation>
441 ("Enable adaptive grid refinement/coarsening");
442 Parameters::Register<Parameters::EnableVtkOutput>
443 ("Global switch for turning on writing VTK files");
444 Parameters::Register<Parameters::EnableThermodynamicHints>
445 ("Enable thermodynamic hints");
446 Parameters::Register<Parameters::EnableIntensiveQuantityCache>
447 ("Turn on caching of intensive quantities");
448 Parameters::Register<Parameters::EnableStorageCache>
449 ("Store previous storage terms and avoid re-calculating them.");
450 Parameters::Register<Parameters::OutputDir>
451 ("The directory to which result files are written");
452 }
453
458 {
459 // initialize the volume of the finite volumes to zero
460 const std::size_t numDof = asImp_().numGridDof();
461 dofTotalVolume_.resize(numDof);
462 std::ranges::fill(dofTotalVolume_, 0.0);
463
464 ElementContext elemCtx(simulator_);
465 gridTotalVolume_ = 0.0;
466
467 // iterate through the grid and evaluate the initial condition
468 for (const auto& elem : elements(gridView_)) {
469 // ignore everything which is not in the interior if the
470 // current process' piece of the grid
471 if (elem.partitionType() != Dune::InteriorEntity) {
472 continue;
473 }
474
475 // deal with the current element
476 elemCtx.updateStencil(elem);
477 const auto& stencil = elemCtx.stencil(/*timeIdx=*/0);
478
479 // loop over all element vertices, i.e. sub control volumes
480 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); dofIdx++) {
481 // map the local degree of freedom index to the global one
482 const unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
483
484 const Scalar dofVolume = stencil.subControlVolume(dofIdx).volume();
485 dofTotalVolume_[globalIdx] += dofVolume;
486 gridTotalVolume_ += dofVolume;
487 }
488 }
489
490 // determine which DOFs should be considered to lie fully in the interior of the
491 // local process grid partition: those which do not have a non-zero volume
492 // before taking the peer processes into account...
493 isLocalDof_.resize(numDof);
494 for (unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
495 isLocalDof_[dofIdx] = (dofTotalVolume_[dofIdx] != 0.0);
496 }
497
498 // add the volumes of the DOFs on the process boundaries
499 const auto sumHandle =
500 GridCommHandleFactory::template sumHandle<Scalar>(dofTotalVolume_,
501 asImp_().dofMapper());
502 gridView_.communicate(*sumHandle,
503 Dune::InteriorBorder_All_Interface,
504 Dune::ForwardCommunication);
505
506 // sum up the volumes of the grid partitions
508
509 linearizer_->init(simulator_);
510 for (unsigned threadId = 0; threadId < ThreadManager::maxThreads(); ++threadId) {
511 localLinearizer_[threadId].init(simulator_);
512 }
513
515
516 newtonMethod_.finishInit();
517 }
518
523 { return enableGridAdaptation_; }
524
530 {
531 // first set the whole domain to zero
532 SolutionVector& uCur = asImp_().solution(/*timeIdx=*/0);
533 uCur = Scalar(0.0);
534
535 ElementContext elemCtx(simulator_);
536
537 // iterate through the grid and evaluate the initial condition
538 for (const auto& elem : elements(gridView_)) {
539 // ignore everything which is not in the interior if the
540 // current process' piece of the grid
541 if (elem.partitionType() != Dune::InteriorEntity) {
542 continue;
543 }
544
545 // deal with the current element
546 elemCtx.updateStencil(elem);
547
548 // loop over all element vertices, i.e. sub control volumes
549 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
550 // map the local degree of freedom index to the global one
551 const unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
552
553 // let the problem do the dirty work of nailing down
554 // the initial solution.
555 simulator_.problem().initial(uCur[globalIdx], elemCtx, dofIdx, /*timeIdx=*/0);
556 asImp_().supplementInitialSolution_(uCur[globalIdx], elemCtx, dofIdx, /*timeIdx=*/0);
557 uCur[globalIdx].checkDefined();
558 }
559 }
560
561 // synchronize the ghost DOFs (if necessary)
562 asImp_().syncOverlap();
563
564 // also set the solutions of the "previous" time steps to the initial solution.
565 for (unsigned timeIdx = 1; timeIdx < historySize; ++timeIdx) {
566 solution(timeIdx) = solution(/*timeIdx=*/0);
567 }
568
569 // Initialize intensive quantities cache now that all problem-specific parameters are available.
570 // This ensures intensiveQuantityHistorySize is correct based on recycleFirstIterationStorage().
571 // TODO: Where this is done should perhaps be changed once finishInit() is refactored.
573
574 simulator_.problem().initialSolutionApplied();
575
576#ifndef NDEBUG
577 for (unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
578 const auto& sol = solution(timeIdx);
579 for (unsigned dofIdx = 0; dofIdx < sol.size(); ++dofIdx) {
580 sol[dofIdx].checkDefined();
581 }
582 }
583#endif // NDEBUG
584 }
585
590 void prefetch(const Element&) const
591 {
592 // do nothing by default
593 }
594
598 NewtonMethod& newtonMethod()
599 { return newtonMethod_; }
600
604 const NewtonMethod& newtonMethod() const
605 { return newtonMethod_; }
606
622 const IntensiveQuantities* thermodynamicHint(unsigned globalIdx, unsigned timeIdx) const
623 {
625 return 0;
626 }
627
628 // the intensive quantities cache doubles as thermodynamic hint
629 return cachedIntensiveQuantities(globalIdx, timeIdx);
630 }
631
643 const IntensiveQuantities* cachedIntensiveQuantities(unsigned globalIdx, unsigned timeIdx) const
644 {
647 !intensiveQuantityCacheUpToDate_[timeIdx][globalIdx]) {
648 return nullptr;
649 }
650
651 // With the storage cache enabled, usually only the
652 // intensive quantities for the most recent time step are
653 // cached. However, this may be false for some Problem
654 // variants, so we should check if the cache exists for
655 // the timeIdx in question.
656 if (timeIdx > 0 && enableStorageCache_ && intensiveQuantityCache_[timeIdx].empty()) {
657 return nullptr;
658 }
659
660 return &intensiveQuantityCache_[timeIdx][globalIdx];
661 }
662
671 void updateCachedIntensiveQuantities(const IntensiveQuantities& intQuants,
672 unsigned globalIdx,
673 unsigned timeIdx) const
674 {
676 return;
677 }
678
679 intensiveQuantityCache_[timeIdx][globalIdx] = intQuants;
680 intensiveQuantityCacheUpToDate_[timeIdx][globalIdx] = true;
681 }
682
692 unsigned timeIdx,
693 bool newValue) const
694 {
696 return;
697 }
698
699 intensiveQuantityCacheUpToDate_[timeIdx][globalIdx] = newValue ? 1 : 0;
700 }
701
707
713 void invalidateIntensiveQuantitiesCache(unsigned timeIdx) const
714 {
716 return;
717 }
718
720 std::ranges::fill(intensiveQuantityCacheUpToDate_[timeIdx], /*value=*/0);
721 }
722 }
723
724 void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx) const
725 {
727
728 // loop over all elements...
729 ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(gridView_);
730#ifdef _OPENMP
731#pragma omp parallel
732#endif
733 {
734 ElementContext elemCtx(simulator_);
735 ElementIterator elemIt = threadedElemIt.beginParallel();
736 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
737 const Element& elem = *elemIt;
738 elemCtx.updatePrimaryStencil(elem);
739 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
740 }
741 }
742 }
743
744 template <class GridViewType>
745 void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx, const GridViewType& gridView) const
746 {
747 // loop over all elements...
748 ThreadedEntityIterator<GridViewType, /*codim=*/0> threadedElemIt(gridView);
749#ifdef _OPENMP
750#pragma omp parallel
751#endif
752 {
753
754 ElementContext elemCtx(simulator_);
755 auto elemIt = threadedElemIt.beginParallel();
756 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
757 if (elemIt->partitionType() != Dune::InteriorEntity) {
758 continue;
759 }
760 const Element& elem = *elemIt;
761 elemCtx.updatePrimaryStencil(elem);
762 // Mark cache for this element as invalid.
763 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
764 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
765 const unsigned globalIndex = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
766 setIntensiveQuantitiesCacheEntryValidity(globalIndex, timeIdx, false);
767 }
768 // Update for this element.
769 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
770 }
771 }
772 }
773
782 void shiftIntensiveQuantityCache(unsigned numSlots = 1)
783 {
784 if (!storeIntensiveQuantities() || numSlots <= 0) {
785 return;
786 }
787
788 if (enableStorageCache() && simulator_.problem().recycleFirstIterationStorage()) {
789 // If the storage term is cached, the intensive quantities of the previous
790 // time steps do not need to be accessed, and we can thus spare ourselves to
791 // copy the objects for the intensive quantities.
792 // However, if the storage term at the start of the timestep cannot be deduced
793 // from the primary variables, we must calculate it from the old intensive
794 // quantities, and need to shift them.
795 return;
796 }
797
798 const unsigned intensiveHistorySize = cachedIntensiveQuantityHistorySize_;
799 for (unsigned timeIdx = 0; timeIdx < intensiveHistorySize - numSlots; ++timeIdx) {
800 intensiveQuantityCache_[timeIdx + numSlots] = intensiveQuantityCache_[timeIdx];
802 }
803
804 // the cache for the most recent time indices do not need to be invalidated
805 // because the solution for them did not change (TODO: that assumes that there is
806 // no post-processing of the solution after a time step! fix it?)
807 }
808
816 { return enableStorageCache_; }
817
826
840 const EqVector& cachedStorage(unsigned globalIdx, unsigned timeIdx) const
841 {
842 if (!enableStorageCache_ ||
843 timeIdx >= historySize ||
844 !storageCacheUpToDate_[timeIdx][globalIdx]) {
845 throw std::logic_error("Cached storage is not available or up to date for the requested "
846 "global index and time index. Make sure storage cache is enabled "
847 "and the entry is valid before calling this method.");
848 }
849
850 return storageCache_[timeIdx][globalIdx];
851 }
852
864 void updateCachedStorage(unsigned globalIdx, unsigned timeIdx, const EqVector& value) const
865 {
866 if (!enableStorageCache_ || timeIdx >= historySize) {
867 return;
868 }
869
870 storageCache_[timeIdx][globalIdx] = value;
871 storageCacheUpToDate_[timeIdx][globalIdx] = 1;
872 }
873
880 bool storageCacheIsUpToDate(unsigned globalIdx, unsigned timeIdx) const
881 {
882 if (!enableStorageCache_ || timeIdx >= historySize) {
883 return false;
884 }
885 return storageCacheUpToDate_[timeIdx][globalIdx] != 0;
886 }
887
894 void invalidateStorageCacheEntry(unsigned globalIdx, unsigned timeIdx) const
895 {
896 if (enableStorageCache_ && timeIdx < historySize) {
897 storageCacheUpToDate_[timeIdx][globalIdx] = 0;
898 }
899 }
900
906 void invalidateStorageCache(unsigned timeIdx) const
907 {
908 if (enableStorageCache_ && timeIdx < historySize) {
909 std::ranges::fill(storageCacheUpToDate_[timeIdx], /*value=*/0);
910 }
911 }
912
924 void shiftStorageCache(unsigned numSlots = 1) const
925 {
926 // If we cannot recycle first iteration storage, it does not make sense to shift the storage cache.
927 if (enableStorageCache_ && !simulator_.problem().recycleFirstIterationStorage()) {
928 for (unsigned timeIdx = 0; timeIdx < historySize - numSlots; ++timeIdx) {
929 storageCache_[timeIdx + numSlots] = storageCache_[timeIdx];
930 storageCacheUpToDate_[timeIdx + numSlots] = storageCacheUpToDate_[timeIdx];
931 }
932
933 // should we invalidate the cache for the most recent time indices? (see shiftIntensiveQuantityCache)
934 }
935 }
936
944 Scalar globalResidual(GlobalEqVector& dest,
945 const SolutionVector& u) const
946 {
947 mutableSolution(/*timeIdx=*/0) = u;
948 const Scalar res = asImp_().globalResidual(dest);
949 mutableSolution(/*timeIdx=*/0) = asImp_().solution(/*timeIdx=*/0);
950 return res;
951 }
952
959 Scalar globalResidual(GlobalEqVector& dest) const
960 {
961 dest = 0;
962
963 std::mutex mutex;
964 ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(gridView_);
965#ifdef _OPENMP
966#pragma omp parallel
967#endif
968 {
969 // Attention: the variables below are thread specific and thus cannot be
970 // moved in front of the #pragma!
971 const unsigned threadId = ThreadManager::threadId();
972 ElementContext elemCtx(simulator_);
973 ElementIterator elemIt = threadedElemIt.beginParallel();
974 LocalEvalBlockVector residual, storageTerm;
975
976 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
977 const Element& elem = *elemIt;
978 if (elem.partitionType() != Dune::InteriorEntity) {
979 continue;
980 }
981
982 elemCtx.updateAll(elem);
983 residual.resize(elemCtx.numDof(/*timeIdx=*/0));
984 storageTerm.resize(elemCtx.numPrimaryDof(/*timeIdx=*/0));
985 asImp_().localResidual(threadId).eval(residual, elemCtx);
986
987 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
988 mutex.lock();
989 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
990 const unsigned globalI = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
991 for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
992 dest[globalI][eqIdx] += Toolbox::value(residual[dofIdx][eqIdx]);
993 }
994 }
995 mutex.unlock();
996 }
997 }
998
999 // add up the residuals on the process borders
1000 const auto sumHandle =
1001 GridCommHandleFactory::template sumHandle<EqVector>(dest, asImp_().dofMapper());
1002 gridView_.communicate(*sumHandle,
1003 Dune::InteriorBorder_InteriorBorder_Interface,
1004 Dune::ForwardCommunication);
1005
1006 // calculate the square norm of the residual. this is not
1007 // entirely correct, since the residual for the finite volumes
1008 // which are on the boundary are counted once for every
1009 // process. As often in life: shit happens (, we don't care)...
1010 return std::sqrt(asImp_().gridView().comm().sum(dest.two_norm2()));
1011 }
1012
1019 void globalStorage(EqVector& storage, unsigned timeIdx = 0) const
1020 {
1021 storage = 0;
1022
1023 std::mutex mutex;
1024 ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(gridView());
1025#ifdef _OPENMP
1026#pragma omp parallel
1027#endif
1028 {
1029 // Attention: the variables below are thread specific and thus cannot be
1030 // moved in front of the #pragma!
1031 const unsigned threadId = ThreadManager::threadId();
1032 ElementContext elemCtx(simulator_);
1033 ElementIterator elemIt = threadedElemIt.beginParallel();
1034 LocalEvalBlockVector elemStorage;
1035
1036 // in this method, we need to disable the storage cache because we want to
1037 // evaluate the storage term for other time indices than the most recent one
1038 elemCtx.setEnableStorageCache(false);
1039
1040 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
1041 const Element& elem = *elemIt;
1042 if (elem.partitionType() != Dune::InteriorEntity) {
1043 continue; // ignore ghost and overlap elements
1044 }
1045
1046 elemCtx.updateStencil(elem);
1047 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
1048
1049 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
1050 elemStorage.resize(numPrimaryDof);
1051
1052 localResidual(threadId).evalStorage(elemStorage, elemCtx, timeIdx);
1053
1054 mutex.lock();
1055 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
1056 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1057 storage[eqIdx] += Toolbox::value(elemStorage[dofIdx][eqIdx]);
1058 }
1059 }
1060 mutex.unlock();
1061 }
1062 }
1063
1064 storage = gridView_.comm().sum(storage);
1065 }
1066
1074 void checkConservativeness([[maybe_unused]] Scalar tolerance = -1,
1075 [[maybe_unused]] bool verbose = false) const
1076 {
1077#ifndef NDEBUG
1078 Scalar totalBoundaryArea(0.0);
1079 Scalar totalVolume(0.0);
1080 EvalEqVector totalRate(0.0);
1081
1082 // take the newton tolerance times the total volume of the grid if we're not
1083 // given an explicit tolerance...
1084 if (tolerance <= 0) {
1085 tolerance =
1086 simulator_.model().newtonMethod().tolerance() *
1087 simulator_.model().gridTotalVolume() *
1088 1000;
1089 }
1090
1091 // we assume the implicit Euler time discretization for now...
1092 assert(historySize == 2);
1093
1094 EqVector storageBeginTimeStep(0.0);
1095 globalStorage(storageBeginTimeStep, /*timeIdx=*/1);
1096
1097 EqVector storageEndTimeStep(0.0);
1098 globalStorage(storageEndTimeStep, /*timeIdx=*/0);
1099
1100 // calculate the rate at the boundary and the source rate
1101 ElementContext elemCtx(simulator_);
1102 elemCtx.setEnableStorageCache(false);
1103 for (const auto& elem : elements(simulator_.gridView())) {
1104 if (elem.partitionType() != Dune::InteriorEntity) {
1105 continue; // ignore ghost and overlap elements
1106 }
1107
1108 elemCtx.updateAll(elem);
1109
1110 // handle the boundary terms
1111 if (elemCtx.onBoundary()) {
1112 BoundaryContext boundaryCtx(elemCtx);
1113
1114 for (unsigned faceIdx = 0; faceIdx < boundaryCtx.numBoundaryFaces(/*timeIdx=*/0); ++faceIdx) {
1115 BoundaryRateVector values;
1116 simulator_.problem().boundary(values,
1117 boundaryCtx,
1118 faceIdx,
1119 /*timeIdx=*/0);
1120 Valgrind::CheckDefined(values);
1121
1122 const unsigned dofIdx = boundaryCtx.interiorScvIndex(faceIdx, /*timeIdx=*/0);
1123 const auto& insideIntQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
1124
1125 const Scalar bfArea =
1126 boundaryCtx.boundarySegmentArea(faceIdx, /*timeIdx=*/0) *
1127 insideIntQuants.extrusionFactor();
1128
1129 for (unsigned i = 0; i < values.size(); ++i) {
1130 values[i] *= bfArea;
1131 }
1132
1133 totalBoundaryArea += bfArea;
1134 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1135 totalRate[eqIdx] += values[eqIdx];
1136 }
1137 }
1138 }
1139
1140 // deal with the source terms
1141 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
1142 RateVector values;
1143 simulator_.problem().source(values,
1144 elemCtx,
1145 dofIdx,
1146 /*timeIdx=*/0);
1147 Valgrind::CheckDefined(values);
1148
1149 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
1150 Scalar dofVolume =
1151 elemCtx.dofVolume(dofIdx, /*timeIdx=*/0) *
1152 intQuants.extrusionFactor();
1153 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1154 totalRate[eqIdx] += -dofVolume*Toolbox::value(values[eqIdx]);
1155 }
1156 totalVolume += dofVolume;
1157 }
1158 }
1159
1160 // summarize everything over all processes
1161 const auto& comm = simulator_.gridView().comm();
1162 totalRate = comm.sum(totalRate);
1163 totalBoundaryArea = comm.sum(totalBoundaryArea);
1164 totalVolume = comm.sum(totalVolume);
1165
1166 if (comm.rank() == 0) {
1167 EqVector storageRate = storageBeginTimeStep;
1168 storageRate -= storageEndTimeStep;
1169 storageRate /= simulator_.timeStepSize();
1170 if (verbose) {
1171 std::cout << "storage at beginning of time step: " << storageBeginTimeStep << "\n";
1172 std::cout << "storage at end of time step: " << storageEndTimeStep << "\n";
1173 std::cout << "rate based on storage terms: " << storageRate << "\n";
1174 std::cout << "rate based on source and boundary terms: " << totalRate << "\n";
1175 std::cout << "difference in rates: ";
1176 for (unsigned eqIdx = 0; eqIdx < EqVector::dimension; ++eqIdx) {
1177 std::cout << (storageRate[eqIdx] - Toolbox::value(totalRate[eqIdx])) << " ";
1178 }
1179 std::cout << "\n";
1180 }
1181 for (unsigned eqIdx = 0; eqIdx < EqVector::dimension; ++eqIdx) {
1182 Scalar eps =
1183 (std::abs(storageRate[eqIdx]) + Toolbox::value(totalRate[eqIdx])) * tolerance;
1184 eps = std::max(tolerance, eps);
1185 assert(std::abs(storageRate[eqIdx] - Toolbox::value(totalRate[eqIdx])) <= eps);
1186 }
1187 }
1188#endif // NDEBUG
1189 }
1190
1196 Scalar dofTotalVolume(unsigned globalIdx) const
1197 { return dofTotalVolume_[globalIdx]; }
1198
1204 bool isLocalDof(unsigned globalIdx) const
1205 { return isLocalDof_[globalIdx]; }
1206
1211 Scalar gridTotalVolume() const
1212 { return gridTotalVolume_; }
1213
1219 const SolutionVector& solution(unsigned timeIdx) const
1220 { return solution_[timeIdx]->blockVector(); }
1221
1225 SolutionVector& solution(unsigned timeIdx)
1226 { return solution_[timeIdx]->blockVector(); }
1227
1228 protected:
1232 SolutionVector& mutableSolution(unsigned timeIdx) const
1233 { return solution_[timeIdx]->blockVector(); }
1234
1235 public:
1240 const Linearizer& linearizer() const
1241 { return *linearizer_; }
1242
1247 Linearizer& linearizer()
1248 { return *linearizer_; }
1249
1258 const LocalLinearizer& localLinearizer(unsigned openMpThreadId) const
1259 { return localLinearizer_[openMpThreadId]; }
1260
1264 LocalLinearizer& localLinearizer(unsigned openMpThreadId)
1265 { return localLinearizer_[openMpThreadId]; }
1266
1270 const LocalResidual& localResidual(unsigned openMpThreadId) const
1271 { return asImp_().localLinearizer(openMpThreadId).localResidual(); }
1272
1276 LocalResidual& localResidual(unsigned openMpThreadId)
1277 { return asImp_().localLinearizer(openMpThreadId).localResidual(); }
1278
1286 Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
1287 {
1288 const Scalar absPv = std::abs(asImp_().solution(/*timeIdx=*/1)[globalDofIdx][pvIdx]);
1289 return 1.0 / std::max(absPv, 1.0);
1290 }
1291
1298 Scalar eqWeight(unsigned, unsigned) const
1299 { return 1.0; }
1300
1310 Scalar relativeDofError(unsigned vertexIdx,
1311 const PrimaryVariables& pv1,
1312 const PrimaryVariables& pv2) const
1313 {
1314 Scalar result = 0.0;
1315 for (unsigned j = 0; j < numEq; ++j) {
1316 const Scalar weight = asImp_().primaryVarWeight(vertexIdx, j);
1317 const Scalar eqErr = std::abs((pv1[j] - pv2[j])*weight);
1318 //Scalar eqErr = std::abs(pv1[j] - pv2[j]);
1319 //eqErr *= std::max<Scalar>(1.0, std::abs(pv1[j] + pv2[j])/2);
1320
1321 result = std::max(result, eqErr);
1322 }
1323 return result;
1324 }
1325
1331 bool update()
1332 {
1333 const TimerGuard prePostProcessGuard(prePostProcessTimer_);
1334
1335#ifndef NDEBUG
1336 for (unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1337 // Make sure that the primary variables are defined. Note that because of padding
1338 // bytes, we can't just simply ask valgrind to check the whole solution vectors
1339 // for definedness...
1340 for (std::size_t i = 0; i < asImp_().solution(/*timeIdx=*/0).size(); ++i) {
1341 asImp_().solution(timeIdx)[i].checkDefined();
1342 }
1343 }
1344#endif // NDEBUG
1345
1346 // make sure all timers are prestine
1349 solveTimer_.halt();
1351
1353 asImp_().updateBegin();
1355
1356 bool converged = false;
1357
1358 try {
1359 converged = newtonMethod_.apply();
1360 }
1361 catch(...) {
1362 prePostProcessTimer_ += newtonMethod_.prePostProcessTimer();
1363 linearizeTimer_ += newtonMethod_.linearizeTimer();
1364 solveTimer_ += newtonMethod_.solveTimer();
1365 updateTimer_ += newtonMethod_.updateTimer();
1366
1367 throw;
1368 }
1369
1370#ifndef NDEBUG
1371 for (unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1372 // Make sure that the primary variables are defined. Note that because of padding
1373 // bytes, we can't just simply ask valgrind to check the whole solution vectors
1374 // for definedness...
1375 for (std::size_t i = 0; i < asImp_().solution(/*timeIdx=*/0).size(); ++i) {
1376 asImp_().solution(timeIdx)[i].checkDefined();
1377 }
1378 }
1379#endif // NDEBUG
1380
1381 prePostProcessTimer_ += newtonMethod_.prePostProcessTimer();
1382 linearizeTimer_ += newtonMethod_.linearizeTimer();
1383 solveTimer_ += newtonMethod_.solveTimer();
1384 updateTimer_ += newtonMethod_.updateTimer();
1385
1387 if (converged) {
1388 asImp_().updateSuccessful();
1389 }
1390 else {
1391 asImp_().updateFailed();
1392 }
1394
1395#ifndef NDEBUG
1396 for (unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1397 // Make sure that the primary variables are defined. Note that because of padding
1398 // bytes, we can't just simply ask valgrind to check the whole solution vectors
1399 // for definedness...
1400 for (std::size_t i = 0; i < asImp_().solution(/*timeIdx=*/0).size(); ++i) {
1401 asImp_().solution(timeIdx)[i].checkDefined();
1402 }
1403 }
1404#endif // NDEBUG
1405
1406 return converged;
1407 }
1408
1417 {}
1418
1425 {}
1426
1432 {}
1433
1438 {
1439 throw std::invalid_argument("Grid adaptation need to be implemented for "
1440 "specific settings of grid and function spaces");
1441 }
1442
1449 {
1450 // Reset the current solution to the one of the
1451 // previous time step so that we can start the next
1452 // update at a physically meaningful solution.
1453 solution(/*timeIdx=*/0) = solution(/*timeIdx=*/1);
1455
1456#ifndef NDEBUG
1457 for (unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1458 // Make sure that the primary variables are defined. Note that because of padding
1459 // bytes, we can't just simply ask valgrind to check the whole solution vectors
1460 // for definedness...
1461 for (std::size_t i = 0; i < asImp_().solution(/*timeIdx=*/0).size(); ++i) {
1462 asImp_().solution(timeIdx)[i].checkDefined();
1463 }
1464 }
1465#endif // NDEBUG
1466 }
1467
1476 {
1477 // at this point we can adapt the grid
1478 if (this->enableGridAdaptation_) {
1479 asImp_().adaptGrid();
1480 }
1481
1482 // make the current solution the previous one.
1483 solution(/*timeIdx=*/1) = solution(/*timeIdx=*/0);
1484
1485 // shift the storage cache by one position in the history
1486 asImp_().shiftStorageCache(/*numSlots=*/1);
1487
1488 // shift the intensive quantities cache by one position in the
1489 // history
1490 asImp_().shiftIntensiveQuantityCache(/*numSlots=*/1);
1491 }
1492
1500 template <class Restarter>
1501 void serialize(Restarter&)
1502 {
1503 throw std::runtime_error("Not implemented: The discretization chosen for this problem "
1504 "does not support restart files. (serialize() method unimplemented)");
1505 }
1506
1514 template <class Restarter>
1515 void deserialize(Restarter&)
1516 {
1517 throw std::runtime_error("Not implemented: The discretization chosen for this problem "
1518 "does not support restart files. (deserialize() method unimplemented)");
1519 }
1520
1529 template <class DofEntity>
1530 void serializeEntity(std::ostream& outstream,
1531 const DofEntity& dof)
1532 {
1533 const unsigned dofIdx = static_cast<unsigned>(asImp_().dofMapper().index(dof));
1534
1535 // write phase state
1536 if (!outstream.good()) {
1537 throw std::runtime_error("Could not serialize degree of freedom " +
1538 std::to_string(dofIdx));
1539 }
1540
1541 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1542 outstream << solution(/*timeIdx=*/0)[dofIdx][eqIdx] << " ";
1543 }
1544 }
1545
1554 template <class DofEntity>
1555 void deserializeEntity(std::istream& instream,
1556 const DofEntity& dof)
1557 {
1558 const unsigned dofIdx = static_cast<unsigned>(asImp_().dofMapper().index(dof));
1559
1560 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1561 if (!instream.good()) {
1562 throw std::runtime_error("Could not deserialize degree of freedom " +
1563 std::to_string(dofIdx));
1564 }
1565 instream >> solution(/*timeIdx=*/0)[dofIdx][eqIdx];
1566 }
1567 }
1568
1572 std::size_t numGridDof() const
1573 { throw std::logic_error("The discretization class must implement the numGridDof() method!"); }
1574
1578 std::size_t numAuxiliaryDof() const
1579 {
1580 return std::accumulate(auxEqModules_.begin(), auxEqModules_.end(),
1581 std::size_t{0},
1582 [](const auto acc, const auto& mod)
1583 { return acc + mod->numDofs(); });
1584 }
1585
1589 std::size_t numTotalDof() const
1590 { return asImp_().numGridDof() + numAuxiliaryDof(); }
1591
1596 const DofMapper& dofMapper() const
1597 { throw std::logic_error("The discretization class must implement the dofMapper() method!"); }
1598
1602 const VertexMapper& vertexMapper() const
1603 { return vertexMapper_; }
1604
1608 const ElementMapper& elementMapper() const
1609 { return elementMapper_; }
1610
1616 {
1617 linearizer_ = std::make_unique<Linearizer>();
1618 linearizer_->init(simulator_);
1619 }
1620
1624 static std::string discretizationName()
1625 { return ""; }
1626
1632 std::string primaryVarName(unsigned pvIdx) const
1633 {
1634 std::ostringstream oss;
1635 oss << "primary variable_" << pvIdx;
1636 return oss.str();
1637 }
1638
1644 std::string eqName(unsigned eqIdx) const
1645 {
1646 std::ostringstream oss;
1647 oss << "equation_" << eqIdx;
1648 return oss.str();
1649 }
1650
1657 void updatePVWeights(const ElementContext&) const
1658 {}
1659
1663 void addOutputModule(std::unique_ptr<BaseOutputModule<TypeTag>> newModule)
1664 { outputModules_.push_back(std::move(newModule)); }
1665
1674 template <class VtkMultiWriter>
1676 const SolutionVector& u,
1677 const GlobalEqVector& deltaU) const
1678 {
1679 using ScalarBuffer = std::vector<double>;
1680
1681 GlobalEqVector globalResid(u.size());
1682 asImp_().globalResidual(globalResid, u);
1683
1684 // create the required scalar fields
1685 const std::size_t numDof = asImp_().numGridDof();
1686
1687 // global defect of the two auxiliary equations
1688 std::array<ScalarBuffer*, numEq> def;
1689 std::array<ScalarBuffer*, numEq> delta;
1690 std::array<ScalarBuffer*, numEq> priVars;
1691 std::array<ScalarBuffer*, numEq> priVarWeight;
1692 ScalarBuffer* relError = writer.allocateManagedScalarBuffer(numDof);
1693 ScalarBuffer* normalizedRelError = writer.allocateManagedScalarBuffer(numDof);
1694 for (unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
1695 priVars[pvIdx] = writer.allocateManagedScalarBuffer(numDof);
1696 priVarWeight[pvIdx] = writer.allocateManagedScalarBuffer(numDof);
1697 delta[pvIdx] = writer.allocateManagedScalarBuffer(numDof);
1698 def[pvIdx] = writer.allocateManagedScalarBuffer(numDof);
1699 }
1700
1701 Scalar minRelErr = 1e30;
1702 Scalar maxRelErr = -1e30;
1703 for (unsigned globalIdx = 0; globalIdx < numDof; ++ globalIdx) {
1704 for (unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
1705 (*priVars[pvIdx])[globalIdx] = u[globalIdx][pvIdx];
1706 (*priVarWeight[pvIdx])[globalIdx] = asImp_().primaryVarWeight(globalIdx, pvIdx);
1707 (*delta[pvIdx])[globalIdx] = - deltaU[globalIdx][pvIdx];
1708 (*def[pvIdx])[globalIdx] = globalResid[globalIdx][pvIdx];
1709 }
1710
1711 PrimaryVariables uOld(u[globalIdx]);
1712 PrimaryVariables uNew(uOld);
1713 uNew -= deltaU[globalIdx];
1714
1715 const Scalar err = asImp_().relativeDofError(globalIdx, uOld, uNew);
1716 (*relError)[globalIdx] = err;
1717 (*normalizedRelError)[globalIdx] = err;
1718 minRelErr = std::min(err, minRelErr);
1719 maxRelErr = std::max(err, maxRelErr);
1720 }
1721
1722 // do the normalization of the relative error
1723 const Scalar alpha = std::max(Scalar{1e-20},
1724 std::max(std::abs(maxRelErr),
1725 std::abs(minRelErr)));
1726 for (unsigned globalIdx = 0; globalIdx < numDof; ++globalIdx) {
1727 (*normalizedRelError)[globalIdx] /= alpha;
1728 }
1729
1730 DiscBaseOutputModule::attachScalarDofData_(writer, *relError, "relative error");
1731 DiscBaseOutputModule::attachScalarDofData_(writer, *normalizedRelError, "normalized relative error");
1732
1733 for (unsigned i = 0; i < numEq; ++i) {
1734 std::ostringstream oss;
1735 oss.str(""); oss << "priVar_" << asImp_().primaryVarName(i);
1736 DiscBaseOutputModule::attachScalarDofData_(writer,
1737 *priVars[i],
1738 oss.str());
1739
1740 oss.str(""); oss << "delta_" << asImp_().primaryVarName(i);
1741 DiscBaseOutputModule::attachScalarDofData_(writer,
1742 *delta[i],
1743 oss.str());
1744
1745 oss.str(""); oss << "weight_" << asImp_().primaryVarName(i);
1746 DiscBaseOutputModule::attachScalarDofData_(writer,
1747 *priVarWeight[i],
1748 oss.str());
1749
1750 oss.str(""); oss << "defect_" << asImp_().eqName(i);
1751 DiscBaseOutputModule::attachScalarDofData_(writer,
1752 *def[i],
1753 oss.str());
1754 }
1755
1756 asImp_().prepareOutputFields();
1757 asImp_().appendOutputFields(writer);
1758 }
1759
1765 {
1766 const bool needFullContextUpdate =
1767 std::ranges::any_of(outputModules_,
1768 [](const auto& mod)
1769 { return mod->needExtensiveQuantities(); });
1770 std::ranges::for_each(outputModules_,
1771 [](auto& mod) { mod->allocBuffers(); });
1772
1773 // iterate over grid
1774 ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(gridView());
1775#ifdef _OPENMP
1776#pragma omp parallel
1777#endif
1778 {
1779 ElementContext elemCtx(simulator_);
1780 ElementIterator elemIt = threadedElemIt.beginParallel();
1781 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
1782 const auto& elem = *elemIt;
1783 if (elem.partitionType() != Dune::InteriorEntity) {
1784 // ignore non-interior entities
1785 continue;
1786 }
1787
1788 if (needFullContextUpdate) {
1789 elemCtx.updateAll(elem);
1790 }
1791 else {
1792 elemCtx.updatePrimaryStencil(elem);
1793 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
1794 }
1795
1796 std::ranges::for_each(outputModules_,
1797 [&elemCtx](auto& mod) { mod->processElement(elemCtx); });
1798 }
1799 }
1800 }
1801
1807 {
1808 std::ranges::for_each(outputModules_,
1809 [&writer](auto& mod) { mod->commitBuffers(writer); });
1810 }
1811
1815 const GridView& gridView() const
1816 { return gridView_; }
1817
1830 {
1831 auxMod->setDofOffset(numTotalDof());
1832 auxEqModules_.push_back(auxMod);
1833
1834 // resize the solutions
1835 if (enableGridAdaptation_ && !std::is_same_v<DiscreteFunction, BlockVectorWrapper>) {
1836 throw std::invalid_argument("Problems which require auxiliary modules cannot be used in"
1837 " conjunction with dune-fem");
1838 }
1839
1840 const std::size_t numDof = numTotalDof();
1841 for (unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1842 solution(timeIdx).resize(numDof);
1843 }
1844
1845 auxMod->applyInitial();
1846 }
1847
1854 {
1855 auxEqModules_.clear();
1856 linearizer_->eraseMatrix();
1857 newtonMethod_.eraseMatrix();
1858 }
1859
1863 std::size_t numAuxiliaryModules() const
1864 { return auxEqModules_.size(); }
1865
1870 { return auxEqModules_[auxEqModIdx]; }
1871
1875 const BaseAuxiliaryModule<TypeTag>* auxiliaryModule(unsigned auxEqModIdx) const
1876 { return auxEqModules_[auxEqModIdx]; }
1877
1883
1885 { return prePostProcessTimer_; }
1886
1887 const Timer& linearizeTimer() const
1888 { return linearizeTimer_; }
1889
1890 const Timer& solveTimer() const
1891 { return solveTimer_; }
1892
1893 const Timer& updateTimer() const
1894 { return updateTimer_; }
1895
1896 template<class Serializer>
1897 void serializeOp(Serializer& serializer)
1898 {
1900 using Helper = typename BaseDiscretization::template SerializeHelper<Serializer>;
1901 Helper::serializeOp(serializer, solution_);
1902 }
1903
1904 bool operator==(const FvBaseDiscretization& rhs) const
1905 {
1906 return std::ranges::equal(this->solution_, rhs.solution_,
1907 [](const auto& x, const auto& y)
1908 { return *x == *y; });
1909 }
1910
1911protected:
1913 {
1914 // allocate the storage cache
1915 if (enableStorageCache()) {
1916 const std::size_t numDof = asImp_().numGridDof();
1917 for (unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1918 storageCache_[timeIdx].resize(numDof);
1919 storageCacheUpToDate_[timeIdx].resize(numDof, /*value=*/0);
1920 }
1921 }
1922
1923 // allocate the intensive quantities cache
1925 const std::size_t numDof = asImp_().numGridDof();
1926 cachedIntensiveQuantityHistorySize_ = simulator_.problem().intensiveQuantityHistorySize();
1927 const unsigned intensiveHistorySize = cachedIntensiveQuantityHistorySize_;
1928
1929 // resize the vectors based on runtime history size
1930 intensiveQuantityCache_.resize(intensiveHistorySize);
1931 intensiveQuantityCacheUpToDate_.resize(intensiveHistorySize);
1932
1933 for(unsigned timeIdx = 0; timeIdx < intensiveHistorySize; ++timeIdx) {
1934 intensiveQuantityCache_[timeIdx].resize(numDof);
1935 intensiveQuantityCacheUpToDate_[timeIdx].resize(numDof);
1937 }
1938 }
1939 }
1940
1941 template <class Context>
1942 void supplementInitialSolution_(PrimaryVariables&,
1943 const Context&,
1944 unsigned,
1945 unsigned)
1946 {}
1947
1956 {
1957 // add the output modules available on all model
1958 this->outputModules_.push_back(std::make_unique<VtkPrimaryVarsModule<TypeTag>>(simulator_));
1959 }
1960
1964 LocalResidual& localResidual_()
1965 { return localLinearizer_.localResidual(); }
1966
1970 bool verbose_() const
1971 { return gridView_.comm().rank() == 0; }
1972
1973 Implementation& asImp_()
1974 { return *static_cast<Implementation*>(this); }
1975
1976 const Implementation& asImp_() const
1977 { return *static_cast<const Implementation*>(this); }
1978
1979 // the problem we want to solve. defines the constitutive
1980 // relations, matxerial laws, etc.
1981 Simulator& simulator_;
1982
1983 // the representation of the spatial domain of the problem
1984 GridView gridView_;
1985
1986 // the mappers for element and vertex entities to global indices
1987 ElementMapper elementMapper_;
1988 VertexMapper vertexMapper_;
1989
1990 // a vector with all auxiliary equations to be considered
1991 std::vector<BaseAuxiliaryModule<TypeTag>*> auxEqModules_;
1992
1993 NewtonMethod newtonMethod_;
1994
1999
2000 // calculates the local jacobian matrix for a given element
2001 std::vector<LocalLinearizer> localLinearizer_;
2002 // Linearizes the problem at the current time step using the
2003 // local jacobian
2004 std::unique_ptr<Linearizer> linearizer_;
2005
2006 // cur is the current iterative solution, prev the converged
2007 // solution of the previous time step
2008 mutable std::vector<IntensiveQuantitiesVector> intensiveQuantityCache_;
2009
2010 // while these are logically bools, concurrent writes to vector<bool> are not thread safe.
2011 mutable std::vector<std::vector<unsigned char>> intensiveQuantityCacheUpToDate_;
2012
2013 std::array<std::unique_ptr<DiscreteFunction>, historySize> solution_;
2014
2015 std::list<std::unique_ptr<BaseOutputModule<TypeTag>>> outputModules_;
2016
2018 std::vector<Scalar> dofTotalVolume_;
2019 std::vector<bool> isLocalDof_;
2020
2021 mutable std::array<GlobalEqVector, historySize> storageCache_;
2022
2023 // while these are logically bools, concurrent writes to vector<bool> are not thread safe.
2024 mutable std::array<std::vector<unsigned char>, historySize> storageCacheUpToDate_;
2025
2030
2032};
2033
2039template<class TypeTag>
2041{
2045
2046 static constexpr unsigned historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>();
2047
2048public:
2049 template<class Serializer>
2051 {
2052 template<class SolutionType>
2053 static void serializeOp(Serializer& serializer,
2054 SolutionType& solution)
2055 {
2056 for (auto& sol : solution) {
2057 serializer(*sol);
2058 }
2059 }
2060 };
2061
2062 explicit FvBaseDiscretizationNoAdapt(Simulator& simulator)
2063 : ParentType(simulator)
2064 {
2065 if (this->enableGridAdaptation_) {
2066 throw std::invalid_argument("Grid adaptation need to use"
2067 " BaseDiscretization = FvBaseDiscretizationFemAdapt"
2068 " which currently requires the presence of the"
2069 " dune-fem module");
2070 }
2071 const std::size_t numDof = this->asImp_().numGridDof();
2072 for (unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
2073 this->solution_[timeIdx] = std::make_unique<DiscreteFunction>("solution", numDof);
2074 }
2075 }
2076};
2077
2078} // namespace Opm
2079
2080#endif // EWOMS_FV_BASE_DISCRETIZATION_HH
This is a stand-alone version of boost::alignment::aligned_allocator from Boost 1....
Base class for specifying auxiliary equations.
Definition: baseauxiliarymodule.hh:56
virtual void applyInitial()=0
Set the initial condition of the auxiliary module in the solution vector.
void setDofOffset(int value)
Set the offset in the global system of equations for the first degree of freedom of this auxiliary mo...
Definition: baseauxiliarymodule.hh:78
The base class for writer modules.
Definition: baseoutputmodule.hh:68
The base class for all output writers.
Definition: baseoutputwriter.hh:46
Represents all quantities which available on boundary segments.
Definition: fvbaseboundarycontext.hh:46
Represents all quantities which available for calculating constraints.
Definition: fvbaseconstraintscontext.hh:44
Class to specify constraints for a finite volume spatial discretization.
Definition: fvbaseconstraints.hh:48
Definition: fvbasediscretization.hh:350
const SolutionVector & blockVector() const
Definition: fvbasediscretization.hh:373
SolutionVector blockVector_
Definition: fvbasediscretization.hh:352
void serializeOp(Serializer &serializer)
Definition: fvbasediscretization.hh:382
static BlockVectorWrapper serializationTestObject()
Definition: fvbasediscretization.hh:360
SolutionVector & blockVector()
Definition: fvbasediscretization.hh:370
bool operator==(const BlockVectorWrapper &wrapper) const
Definition: fvbasediscretization.hh:376
BlockVectorWrapper(const std::string &, const std::size_t size)
Definition: fvbasediscretization.hh:354
The base class for the finite volume discretization schemes without adaptation.
Definition: fvbasediscretization.hh:2041
FvBaseDiscretizationNoAdapt(Simulator &simulator)
Definition: fvbasediscretization.hh:2062
The base class for the finite volume discretization schemes.
Definition: fvbasediscretization.hh:298
Timer linearizeTimer_
Definition: fvbasediscretization.hh:1996
std::vector< IntensiveQuantitiesVector > intensiveQuantityCache_
Definition: fvbasediscretization.hh:2008
LocalLinearizer & localLinearizer(unsigned openMpThreadId)
Definition: fvbasediscretization.hh:1264
void prepareOutputFields() const
Prepare the quantities relevant for the current solution to be appended to the output writers.
Definition: fvbasediscretization.hh:1764
void shiftIntensiveQuantityCache(unsigned numSlots=1)
Move the intensive quantities for a given time index to the back.
Definition: fvbasediscretization.hh:782
void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx) const
Definition: fvbasediscretization.hh:724
void adaptGrid()
Called by the update() method when the grid should be refined.
Definition: fvbasediscretization.hh:1437
std::vector< BaseAuxiliaryModule< TypeTag > * > auxEqModules_
Definition: fvbasediscretization.hh:1991
const Implementation & asImp_() const
Definition: fvbasediscretization.hh:1976
void addAuxiliaryModule(BaseAuxiliaryModule< TypeTag > *auxMod)
Add a module for an auxiliary equation.
Definition: fvbasediscretization.hh:1829
void setIntensiveQuantitiesCacheEntryValidity(unsigned globalIdx, unsigned timeIdx, bool newValue) const
Invalidate the cache for a given intensive quantities object.
Definition: fvbasediscretization.hh:691
void finishInit()
Apply the initial conditions to the model.
Definition: fvbasediscretization.hh:457
void prefetch(const Element &) const
Allows to improve the performance by prefetching all data which is associated with a given element.
Definition: fvbasediscretization.hh:590
bool enableStorageCache_
Definition: fvbasediscretization.hh:2028
void resizeAndResetIntensiveQuantitiesCache_()
Definition: fvbasediscretization.hh:1912
std::vector< Scalar > dofTotalVolume_
Definition: fvbasediscretization.hh:2018
void updateSuccessful()
Called by the update() method if it was successful.
Definition: fvbasediscretization.hh:1431
static std::string discretizationName()
Returns a string of discretization's human-readable name.
Definition: fvbasediscretization.hh:1624
unsigned cachedIntensiveQuantityHistorySize() const
Get the cached intensive quantity history size.
Definition: fvbasediscretization.hh:705
BaseAuxiliaryModule< TypeTag > * auxiliaryModule(unsigned auxEqModIdx)
Returns a given module for auxiliary equations.
Definition: fvbasediscretization.hh:1869
bool isLocalDof(unsigned globalIdx) const
Returns if the overlap of the volume ofa degree of freedom is non-zero.
Definition: fvbasediscretization.hh:1204
std::size_t numAuxiliaryModules() const
Returns the number of modules for auxiliary equations.
Definition: fvbasediscretization.hh:1863
bool operator==(const FvBaseDiscretization &rhs) const
Definition: fvbasediscretization.hh:1904
void serializeOp(Serializer &serializer)
Definition: fvbasediscretization.hh:1897
std::vector< bool > isLocalDof_
Definition: fvbasediscretization.hh:2019
LocalResidual & localResidual_()
Reference to the local residal object.
Definition: fvbasediscretization.hh:1964
void registerOutputModules_()
Register all output modules which make sense for the model.
Definition: fvbasediscretization.hh:1955
bool enableGridAdaptation() const
Returns whether the grid ought to be adapted to the solution during the simulation.
Definition: fvbasediscretization.hh:522
const NewtonMethod & newtonMethod() const
Returns the newton method object.
Definition: fvbasediscretization.hh:604
SolutionVector & mutableSolution(unsigned timeIdx) const
Definition: fvbasediscretization.hh:1232
void advanceTimeLevel()
Called by the problem if a time integration was successful, post processing of the solution is done a...
Definition: fvbasediscretization.hh:1475
NewtonMethod & newtonMethod()
Returns the newton method object.
Definition: fvbasediscretization.hh:598
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices.
Definition: fvbasediscretization.hh:1602
const EqVector & cachedStorage(unsigned globalIdx, unsigned timeIdx) const
Retrieve an entry of the cache for the storage term.
Definition: fvbasediscretization.hh:840
void updatePVWeights(const ElementContext &) const
Update the weights of all primary variables within an element given the complete set of intensive qua...
Definition: fvbasediscretization.hh:1657
const Timer & solveTimer() const
Definition: fvbasediscretization.hh:1890
void updateFailed()
Called by the update() method if it was unsuccessful. This is primary a hook which the actual model c...
Definition: fvbasediscretization.hh:1448
void updateBegin()
Called by the update() method before it tries to apply the newton method. This is primary a hook whic...
Definition: fvbasediscretization.hh:1424
FvBaseDiscretization(const FvBaseDiscretization &)=delete
Scalar globalResidual(GlobalEqVector &dest) const
Compute the global residual for the current solution vector.
Definition: fvbasediscretization.hh:959
LocalResidual & localResidual(unsigned openMpThreadId)
Definition: fvbasediscretization.hh:1276
void serializeEntity(std::ostream &outstream, const DofEntity &dof)
Write the current solution for a degree of freedom to a restart file.
Definition: fvbasediscretization.hh:1530
void setEnableStorageCache(bool enableStorageCache)
Set the value of enable storage cache.
Definition: fvbasediscretization.hh:824
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition: fvbasediscretization.hh:1644
Scalar eqWeight(unsigned, unsigned) const
Returns the relative weight of an equation.
Definition: fvbasediscretization.hh:1298
std::array< std::vector< unsigned char >, historySize > storageCacheUpToDate_
Definition: fvbasediscretization.hh:2024
const Timer & prePostProcessTimer() const
Definition: fvbasediscretization.hh:1884
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: fvbasediscretization.hh:1286
void deserialize(Restarter &)
Deserializes the state of the model.
Definition: fvbasediscretization.hh:1515
void checkConservativeness(Scalar tolerance=-1, bool verbose=false) const
Ensure that the difference between the storage terms of the last and of the current time step is cons...
Definition: fvbasediscretization.hh:1074
Scalar gridTotalVolume() const
Returns the volume of the whole grid which represents the spatial domain.
Definition: fvbasediscretization.hh:1211
Timer updateTimer_
Definition: fvbasediscretization.hh:1998
FvBaseDiscretization(Simulator &simulator)
Definition: fvbasediscretization.hh:393
const IntensiveQuantities * thermodynamicHint(unsigned globalIdx, unsigned timeIdx) const
Return the thermodynamic hint for a entity on the grid at given time.
Definition: fvbasediscretization.hh:622
bool storageCacheIsUpToDate(unsigned globalIdx, unsigned timeIdx) const
Returns true if the storage cache entry for a given DOF and time index is up to date.
Definition: fvbasediscretization.hh:880
const Timer & updateTimer() const
Definition: fvbasediscretization.hh:1893
const Linearizer & linearizer() const
Returns the operator linearizer for the global jacobian of the problem.
Definition: fvbasediscretization.hh:1240
void updateCachedStorage(unsigned globalIdx, unsigned timeIdx, const EqVector &value) const
Set an entry of the cache for the storage term.
Definition: fvbasediscretization.hh:864
void addConvergenceVtkFields(VtkMultiWriter &writer, const SolutionVector &u, const GlobalEqVector &deltaU) const
Add the vector fields for analysing the convergence of the newton method to the a VTK writer.
Definition: fvbasediscretization.hh:1675
Scalar relativeDofError(unsigned vertexIdx, const PrimaryVariables &pv1, const PrimaryVariables &pv2) const
Returns the relative error between two vectors of primary variables.
Definition: fvbasediscretization.hh:1310
bool enableGridAdaptation_
Definition: fvbasediscretization.hh:2026
Scalar globalResidual(GlobalEqVector &dest, const SolutionVector &u) const
Compute the global residual for an arbitrary solution vector.
Definition: fvbasediscretization.hh:944
SolutionVector & solution(unsigned timeIdx)
Definition: fvbasediscretization.hh:1225
std::array< GlobalEqVector, historySize > storageCache_
Definition: fvbasediscretization.hh:2021
Timer prePostProcessTimer_
Definition: fvbasediscretization.hh:1995
void deserializeEntity(std::istream &instream, const DofEntity &dof)
Reads the current solution variables for a degree of freedom from a restart file.
Definition: fvbasediscretization.hh:1555
Timer solveTimer_
Definition: fvbasediscretization.hh:1997
void supplementInitialSolution_(PrimaryVariables &, const Context &, unsigned, unsigned)
Definition: fvbasediscretization.hh:1942
const LocalLinearizer & localLinearizer(unsigned openMpThreadId) const
Returns the local jacobian which calculates the local stiffness matrix for an arbitrary element.
Definition: fvbasediscretization.hh:1258
std::array< std::unique_ptr< DiscreteFunction >, historySize > solution_
Definition: fvbasediscretization.hh:2013
unsigned cachedIntensiveQuantityHistorySize_
Definition: fvbasediscretization.hh:2031
const Timer & linearizeTimer() const
Definition: fvbasediscretization.hh:1887
void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx, const GridViewType &gridView) const
Definition: fvbasediscretization.hh:745
bool update()
Try to progress the model to the next timestep.
Definition: fvbasediscretization.hh:1331
std::size_t numAuxiliaryDof() const
Returns the number of degrees of freedom (DOFs) of the auxiliary equations.
Definition: fvbasediscretization.hh:1578
void clearAuxiliaryModules()
Causes the list of auxiliary equations to be cleared.
Definition: fvbasediscretization.hh:1853
bool enableIntensiveQuantityCache_
Definition: fvbasediscretization.hh:2027
std::unique_ptr< Linearizer > linearizer_
Definition: fvbasediscretization.hh:2004
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices.
Definition: fvbasediscretization.hh:1608
void syncOverlap()
Syncronize the values of the primary variables on the degrees of freedom that overlap with the neighb...
Definition: fvbasediscretization.hh:1416
void appendOutputFields(BaseOutputWriter &writer) const
Append the quantities relevant for the current solution to an output writer.
Definition: fvbasediscretization.hh:1806
std::list< std::unique_ptr< BaseOutputModule< TypeTag > > > outputModules_
Definition: fvbasediscretization.hh:2015
ElementMapper elementMapper_
Definition: fvbasediscretization.hh:1987
NewtonMethod newtonMethod_
Definition: fvbasediscretization.hh:1993
std::size_t numTotalDof() const
Returns the total number of degrees of freedom (i.e., grid plux auxiliary DOFs)
Definition: fvbasediscretization.hh:1589
void applyInitialSolution()
Applies the initial solution for all degrees of freedom to which the model applies.
Definition: fvbasediscretization.hh:529
Implementation & asImp_()
Definition: fvbasediscretization.hh:1973
void updateCachedIntensiveQuantities(const IntensiveQuantities &intQuants, unsigned globalIdx, unsigned timeIdx) const
Update the intensive quantity cache for a entity on the grid at given time.
Definition: fvbasediscretization.hh:671
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition: fvbasediscretization.hh:1632
void invalidateIntensiveQuantitiesCache(unsigned timeIdx) const
Invalidate the whole intensive quantity cache for time index.
Definition: fvbasediscretization.hh:713
Linearizer & linearizer()
Returns the object which linearizes the global system of equations at the current solution.
Definition: fvbasediscretization.hh:1247
const BaseAuxiliaryModule< TypeTag > * auxiliaryModule(unsigned auxEqModIdx) const
Returns a given module for auxiliary equations.
Definition: fvbasediscretization.hh:1875
std::vector< std::vector< unsigned char > > intensiveQuantityCacheUpToDate_
Definition: fvbasediscretization.hh:2011
void globalStorage(EqVector &storage, unsigned timeIdx=0) const
Compute the integral over the domain of the storage terms of all conservation quantities.
Definition: fvbasediscretization.hh:1019
GridView gridView_
Definition: fvbasediscretization.hh:1984
Scalar dofTotalVolume(unsigned globalIdx) const
Returns the volume of a given control volume.
Definition: fvbasediscretization.hh:1196
static void registerParameters()
Register all run-time parameters for the model.
Definition: fvbasediscretization.hh:426
const SolutionVector & solution(unsigned timeIdx) const
Reference to the solution at a given history index as a block vector.
Definition: fvbasediscretization.hh:1219
bool verbose_() const
Returns whether messages should be printed.
Definition: fvbasediscretization.hh:1970
void shiftStorageCache(unsigned numSlots=1) const
Shift storage cache by a given number of time step slots.
Definition: fvbasediscretization.hh:924
Simulator & simulator_
Definition: fvbasediscretization.hh:1981
const LocalResidual & localResidual(unsigned openMpThreadId) const
Returns the object to calculate the local residual function.
Definition: fvbasediscretization.hh:1270
std::size_t numGridDof() const
Returns the number of degrees of freedom (DOFs) for the computational grid.
Definition: fvbasediscretization.hh:1572
Scalar gridTotalVolume_
Definition: fvbasediscretization.hh:2017
const DofMapper & dofMapper() const
Mapper to convert the Dune entities of the discretization's degrees of freedoms are to indices.
Definition: fvbasediscretization.hh:1596
bool storeIntensiveQuantities() const
Returns true if the cache for intensive quantities is enabled.
Definition: fvbasediscretization.hh:1881
std::vector< LocalLinearizer > localLinearizer_
Definition: fvbasediscretization.hh:2001
void invalidateStorageCache(unsigned timeIdx) const
Invalidate the whole storage cache for a given time index.
Definition: fvbasediscretization.hh:906
void serialize(Restarter &)
Serializes the current state of the model.
Definition: fvbasediscretization.hh:1501
bool enableThermodynamicHints_
Definition: fvbasediscretization.hh:2029
void invalidateStorageCacheEntry(unsigned globalIdx, unsigned timeIdx) const
Invalidate the storage cache for a given DOF and time index.
Definition: fvbasediscretization.hh:894
VertexMapper vertexMapper_
Definition: fvbasediscretization.hh:1988
void resetLinearizer()
Resets the Jacobian matrix linearizer, so that the boundary types can be altered.
Definition: fvbasediscretization.hh:1615
bool enableStorageCache() const
Returns true iff the storage term is cached.
Definition: fvbasediscretization.hh:815
const GridView & gridView() const
Reference to the grid view of the spatial domain.
Definition: fvbasediscretization.hh:1815
void addOutputModule(std::unique_ptr< BaseOutputModule< TypeTag > > newModule)
Add an module for writing visualization output after a timestep.
Definition: fvbasediscretization.hh:1663
const IntensiveQuantities * cachedIntensiveQuantities(unsigned globalIdx, unsigned timeIdx) const
Return the cached intensive quantities for a entity on the grid at given time.
Definition: fvbasediscretization.hh:643
This class stores an array of IntensiveQuantities objects, one intensive quantities object for each o...
Definition: fvbaseelementcontext.hh:55
Provide the properties at a face which make sense independently of the conserved quantities.
Definition: fvbaseextensivequantities.hh:48
This class calculates gradients of arbitrary quantities at flux integration points using the two-poin...
Definition: fvbasegradientcalculator.hh:52
Base class for the model specific class which provides access to all intensive (i....
Definition: fvbaseintensivequantities.hh:45
The common code for the linearizers of non-linear systems of equations.
Definition: fvbaselinearizer.hh:78
Element-wise caculation of the residual matrix for models based on a finite volume spatial discretiza...
Definition: fvbaselocalresidual.hh:63
Represents the primary variables used by the a model.
Definition: fvbaseprimaryvariables.hh:54
This is a grid manager which does not create any border list.
Definition: nullborderlistmanager.hh:44
static void registerParameters()
Register all run-time parameters for the Newton method.
Definition: newtonmethod.hh:136
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:84
Simplifies multi-threaded capabilities.
Definition: threadmanager.hpp:36
static unsigned maxThreads()
Return the maximum number of threads of the current process.
Definition: threadmanager.hpp:66
static unsigned threadId()
Return the index of the current OpenMP thread.
Provides an STL-iterator like interface to iterate over the enties of a GridView in OpenMP threaded a...
Definition: threadedentityiterator.hh:42
bool isFinished(const EntityIterator &it) const
Definition: threadedentityiterator.hh:67
EntityIterator increment()
Definition: threadedentityiterator.hh:80
EntityIterator beginParallel()
Definition: threadedentityiterator.hh:54
A simple class which makes sure that a timer gets stopped if an exception is thrown.
Definition: timerguard.hh:42
Provides an encapsulation to measure the system time.
Definition: timer.hpp:46
void start()
Start counting the time resources used by the simulation.
void halt()
Stop the measurement reset all timing values.
double stop()
Stop counting the time resources.
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:65
ScalarBuffer * allocateManagedScalarBuffer(std::size_t numEntities)
Allocate a managed buffer for a scalar field.
Definition: vtkmultiwriter.hh:192
VTK output module for the fluid composition.
Definition: vtkprimaryvarsmodule.hpp:48
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkprimaryvarsmodule.hpp:74
Definition: alignedallocator.hh:97
Declare the properties used by the infrastructure code of the finite volume discretizations.
Provides data handles for parallel communication which operate on DOFs.
Declares the parameters for the black oil model.
Definition: fvbaseprimaryvariables.hh:161
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:187
Definition: blackoilmodel.hh:80
Definition: blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: fvbasediscretization.hh:2051
static void serializeOp(Serializer &serializer, SolutionType &solution)
Definition: fvbasediscretization.hh:2053
Definition: fvbasediscretization.hh:272
GetPropType< TypeTag, Properties::GridView > GridView
Definition: fvbasediscretization.hh:113
GetPropType< TypeTag, Properties::DofMapper > DofMapper
Definition: fvbasediscretization.hh:112
The class which marks the border indices associated with the degrees of freedom on a process boundary...
Definition: basicproperties.hh:129
The secondary variables of a boundary segment.
Definition: fvbaseproperties.hh:147
GetPropType< TypeTag, Properties::RateVector > type
Definition: fvbasediscretization.hh:160
Type of object for specifying boundary conditions.
Definition: fvbaseproperties.hh:119
The secondary variables of a constraint degree of freedom.
Definition: fvbaseproperties.hh:150
The class which represents a constraint degree of freedom.
Definition: fvbaseproperties.hh:122
The part of the extensive quantities which is specific to the spatial discretization.
Definition: fvbaseproperties.hh:164
Definition: fvbaseproperties.hh:140
The discretization specific part of the local residual.
Definition: fvbaseproperties.hh:91
typename BaseDiscretization::BlockVectorWrapper type
Definition: fvbasediscretization.hh:283
Definition: fvbaseproperties.hh:77
The secondary variables of all degrees of freedom in an element's stencil.
Definition: fvbaseproperties.hh:144
Dune::BlockVector< GetPropType< TypeTag, Properties::EqVector > > type
Definition: fvbasediscretization.hh:174
A vector of holding a quantity for each equation for each DOF of an element.
Definition: fvbaseproperties.hh:112
Dune::MultipleCodimMultipleGeomTypeMapper< GetPropType< TypeTag, Properties::GridView > > type
Definition: fvbasediscretization.hh:106
The mapper to find the global index of an element.
Definition: fvbaseproperties.hh:217
Specify whether the some degrees of fredom can be constraint.
Definition: fvbaseproperties.hh:203
Specify if experimental features should be enabled or not.
Definition: fvbaseproperties.hh:245
Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, getPropValue< TypeTag, Properties::NumEq >()> type
Definition: fvbasediscretization.hh:143
A vector of holding a quantity for each equation (usually at a given spatial location)
Definition: fvbaseproperties.hh:109
Specify whether the storage terms use extensive quantities or not.
Definition: fvbaseproperties.hh:237
Dune::BlockVector< GetPropType< TypeTag, Properties::EqVector > > type
Definition: fvbasediscretization.hh:181
Vector containing a quantity of for equation for each DOF of the whole grid.
Definition: linalgproperties.hh:54
Calculates gradients of arbitrary quantities at flux integration points.
Definition: fvbaseproperties.hh:156
The secondary variables within a sub-control volume.
Definition: fvbaseproperties.hh:133
The class which linearizes the non-linear system of equations.
Definition: newtonmethodproperties.hh:36
A vector of primary variables within a sub-control volume.
Definition: fvbaseproperties.hh:130
GetPropType< TypeTag, Properties::EqVector > type
Definition: fvbasediscretization.hh:153
Vector containing volumetric or areal rates of quantities.
Definition: fvbaseproperties.hh:116
Manages the simulation time.
Definition: basicproperties.hh:120
Dune::BlockVector< GetPropType< TypeTag, Properties::PrimaryVariables > > type
Definition: fvbasediscretization.hh:195
Vector containing all primary variables of the grid.
Definition: fvbaseproperties.hh:126
The OpenMP threads manager.
Definition: fvbaseproperties.hh:178
The history size required by the time discretization.
Definition: fvbaseproperties.hh:229
a tag to mark properties as undefined
Definition: propertysystem.hh:38
Definition: fvbaseproperties.hh:185
Specify whether to use volumetric residuals or not.
Definition: fvbaseproperties.hh:241
Dune::MultipleCodimMultipleGeomTypeMapper< GetPropType< TypeTag, Properties::GridView > > type
Definition: fvbasediscretization.hh:101
The mapper to find the global index of a vertex.
Definition: fvbaseproperties.hh:211
Specify the format the VTK output is written to disk.
Definition: fvbaseproperties.hh:199