28#ifndef EWOMS_FV_BASE_DISCRETIZATION_HH
29#define EWOMS_FV_BASE_DISCRETIZATION_HH
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
33#include <dune/common/version.hh>
34#include <dune/istl/bvector.hh>
36#include <opm/material/common/MathToolbox.hpp>
37#include <opm/material/common/Valgrind.hpp>
38#include <opm/material/densead/Math.hpp>
83template<
class TypeTag>
86template<
class TypeTag>
91namespace Opm::Properties {
94template<
class TypeTag>
99template<
class TypeTag>
101{
using type = Dune::MultipleCodimMultipleGeomTypeMapper<GetPropType<TypeTag, Properties::GridView>>; };
104template<
class TypeTag>
106{
using type = Dune::MultipleCodimMultipleGeomTypeMapper<GetPropType<TypeTag, Properties::GridView>>; };
109template<
class TypeTag>
119template<
class TypeTag>
123template<
class TypeTag>
127template<
class TypeTag>
132template<
class TypeTag>
139template<
class TypeTag>
142 using type = Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
151template<
class TypeTag>
158template<
class TypeTag>
165template<
class TypeTag>
172template<
class TypeTag>
174{
using type = Dune::BlockVector<GetPropType<TypeTag, Properties::EqVector>>; };
179template<
class TypeTag>
181{
using type = Dune::BlockVector<GetPropType<TypeTag, Properties::EqVector>>; };
186template<
class TypeTag>
193template<
class TypeTag>
195{
using type = Dune::BlockVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; };
202template<
class TypeTag>
209template<
class TypeTag>
213template<
class TypeTag>
217template<
class TypeTag>
224template<
class TypeTag>
228template<
class TypeTag>
230{
static constexpr bool value =
true; };
235template<
class TypeTag>
240template<
class TypeTag>
242{
static constexpr auto value = Dune::VTK::ascii; };
245template<
class TypeTag>
247{
static constexpr bool value =
false; };
250template<
class TypeTag>
252{
static constexpr int value = 2; };
256template<
class TypeTag>
258{
static constexpr bool value =
false; };
261template<
class TypeTag>
263{
static constexpr bool value =
true; };
267template<
class TypeTag>
269{
static constexpr bool value =
true; };
271template <
class TypeTag,
class MyTypeTag>
275template<
class TypeTag>
279template<
class TypeTag>
296template<
class TypeTag>
297class FvBaseDiscretization
335 using IntensiveQuantitiesVector = std::vector<IntensiveQuantities,
337 alignof(IntensiveQuantities)>>;
339 using Element =
typename GridView::template Codim<0>::Entity;
340 using ElementIterator =
typename GridView::template Codim<0>::Iterator;
342 using Toolbox = MathToolbox<Evaluation>;
343 using VectorBlock = Dune::FieldVector<Evaluation, numEq>;
344 using EvalEqVector = Dune::FieldVector<Evaluation, numEq>;
346 using LocalEvalBlockVector =
typename LocalResidual::LocalEvalBlockVector;
349 class BlockVectorWrapper
352 SolutionVector blockVector_;
354 BlockVectorWrapper(
const std::string&,
const std::size_t size)
358 BlockVectorWrapper() =
default;
360 static BlockVectorWrapper serializationTestObject()
362 BlockVectorWrapper result(
"dummy", 3);
363 result.blockVector_[0] = 1.0;
364 result.blockVector_[1] = 2.0;
365 result.blockVector_[2] = 3.0;
370 SolutionVector& blockVector()
371 {
return blockVector_; }
373 const SolutionVector& blockVector()
const
374 {
return blockVector_; }
376 bool operator==(
const BlockVectorWrapper& wrapper)
const
378 return std::ranges::equal(this->blockVector_, wrapper.blockVector_);
381 template<
class Serializer>
382 void serializeOp(Serializer& serializer)
384 serializer(blockVector_);
393 explicit FvBaseDiscretization(Simulator& simulator)
394 : simulator_(simulator)
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))
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>() +
")");
414 PrimaryVariables::init();
417 asImp_().registerOutputModules_();
421 FvBaseDiscretization(
const FvBaseDiscretization&) =
delete;
428 Linearizer::registerParameters();
429 LocalLinearizer::registerParameters();
430 LocalResidual::registerParameters();
431 GradientCalculator::registerParameters();
432 IntensiveQuantities::registerParameters();
433 ExtensiveQuantities::registerParameters();
435 Linearizer::registerParameters();
436 PrimaryVariables::registerParameters();
441 (
"Enable adaptive grid refinement/coarsening");
443 (
"Global switch for turning on writing VTK files");
445 (
"Enable thermodynamic hints");
447 (
"Turn on caching of intensive quantities");
449 (
"Store previous storage terms and avoid re-calculating them.");
451 (
"The directory to which result files are written");
460 const std::size_t numDof = asImp_().numGridDof();
461 dofTotalVolume_.resize(numDof);
462 std::ranges::fill(dofTotalVolume_, 0.0);
464 ElementContext elemCtx(simulator_);
465 gridTotalVolume_ = 0.0;
468 for (
const auto& elem : elements(gridView_)) {
471 if (elem.partitionType() != Dune::InteriorEntity) {
476 elemCtx.updateStencil(elem);
477 const auto& stencil = elemCtx.stencil(0);
480 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); dofIdx++) {
482 const unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
484 const Scalar dofVolume = stencil.subControlVolume(dofIdx).volume();
485 dofTotalVolume_[globalIdx] += dofVolume;
486 gridTotalVolume_ += dofVolume;
493 isLocalDof_.resize(numDof);
494 for (
unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
495 isLocalDof_[dofIdx] = (dofTotalVolume_[dofIdx] != 0.0);
499 const auto sumHandle =
500 GridCommHandleFactory::template sumHandle<Scalar>(dofTotalVolume_,
502 gridView_.communicate(*sumHandle,
503 Dune::InteriorBorder_All_Interface,
504 Dune::ForwardCommunication);
507 gridTotalVolume_ = gridView_.comm().sum(gridTotalVolume_);
509 linearizer_->init(simulator_);
511 localLinearizer_[threadId].init(simulator_);
514 resizeAndResetIntensiveQuantitiesCache_();
516 newtonMethod_.finishInit();
523 {
return enableGridAdaptation_; }
532 SolutionVector& uCur = asImp_().solution(0);
535 ElementContext elemCtx(simulator_);
538 for (
const auto& elem : elements(gridView_)) {
541 if (elem.partitionType() != Dune::InteriorEntity) {
546 elemCtx.updateStencil(elem);
549 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); ++dofIdx) {
551 const unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
555 simulator_.problem().initial(uCur[globalIdx], elemCtx, dofIdx, 0);
556 asImp_().supplementInitialSolution_(uCur[globalIdx], elemCtx, dofIdx, 0);
557 uCur[globalIdx].checkDefined();
562 asImp_().syncOverlap();
565 for (
unsigned timeIdx = 1; timeIdx < historySize; ++timeIdx) {
572 resizeAndResetIntensiveQuantitiesCache_();
574 simulator_.problem().initialSolutionApplied();
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();
599 {
return newtonMethod_; }
605 {
return newtonMethod_; }
624 if (!enableThermodynamicHints_) {
645 if (!enableIntensiveQuantityCache_ ||
646 timeIdx >= cachedIntensiveQuantityHistorySize_ ||
647 !intensiveQuantityCacheUpToDate_[timeIdx][globalIdx]) {
656 if (timeIdx > 0 && enableStorageCache_ && intensiveQuantityCache_[timeIdx].empty()) {
660 return &intensiveQuantityCache_[timeIdx][globalIdx];
673 unsigned timeIdx)
const
679 intensiveQuantityCache_[timeIdx][globalIdx] = intQuants;
680 intensiveQuantityCacheUpToDate_[timeIdx][globalIdx] =
true;
699 intensiveQuantityCacheUpToDate_[timeIdx][globalIdx] = newValue ? 1 : 0;
706 {
return cachedIntensiveQuantityHistorySize_; }
715 if (timeIdx >= cachedIntensiveQuantityHistorySize_) {
720 std::ranges::fill(intensiveQuantityCacheUpToDate_[timeIdx], 0);
724 void invalidateAndUpdateIntensiveQuantities(
unsigned timeIdx)
const
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);
744 template <
class Gr
idViewType>
745 void invalidateAndUpdateIntensiveQuantities(
unsigned timeIdx,
const GridViewType&
gridView)
const
748 ThreadedEntityIterator<GridViewType, 0> threadedElemIt(
gridView);
754 ElementContext elemCtx(simulator_);
755 auto elemIt = threadedElemIt.beginParallel();
756 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
757 if (elemIt->partitionType() != Dune::InteriorEntity) {
760 const Element& elem = *elemIt;
761 elemCtx.updatePrimaryStencil(elem);
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);
769 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
798 const unsigned intensiveHistorySize = cachedIntensiveQuantityHistorySize_;
799 for (
unsigned timeIdx = 0; timeIdx < intensiveHistorySize - numSlots; ++timeIdx) {
800 intensiveQuantityCache_[timeIdx + numSlots] = intensiveQuantityCache_[timeIdx];
801 intensiveQuantityCacheUpToDate_[timeIdx + numSlots] = intensiveQuantityCacheUpToDate_[timeIdx];
816 {
return enableStorageCache_; }
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.");
850 return storageCache_[timeIdx][globalIdx];
866 if (!enableStorageCache_ || timeIdx >= historySize) {
870 storageCache_[timeIdx][globalIdx] = value;
871 storageCacheUpToDate_[timeIdx][globalIdx] = 1;
882 if (!enableStorageCache_ || timeIdx >= historySize) {
885 return storageCacheUpToDate_[timeIdx][globalIdx] != 0;
896 if (enableStorageCache_ && timeIdx < historySize) {
897 storageCacheUpToDate_[timeIdx][globalIdx] = 0;
908 if (enableStorageCache_ && timeIdx < historySize) {
909 std::ranges::fill(storageCacheUpToDate_[timeIdx], 0);
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];
945 const SolutionVector& u)
const
948 const Scalar res = asImp_().globalResidual(dest);
972 ElementContext elemCtx(simulator_);
973 ElementIterator elemIt = threadedElemIt.beginParallel();
974 LocalEvalBlockVector residual, storageTerm;
976 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
977 const Element& elem = *elemIt;
978 if (elem.partitionType() != Dune::InteriorEntity) {
982 elemCtx.updateAll(elem);
983 residual.resize(elemCtx.numDof(0));
984 storageTerm.resize(elemCtx.numPrimaryDof(0));
985 asImp_().localResidual(threadId).eval(residual, elemCtx);
987 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(0);
989 for (
unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
990 const unsigned globalI = elemCtx.globalSpaceIndex(dofIdx, 0);
991 for (
unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
992 dest[globalI][eqIdx] += Toolbox::value(residual[dofIdx][eqIdx]);
1000 const auto sumHandle =
1001 GridCommHandleFactory::template sumHandle<EqVector>(dest, asImp_().
dofMapper());
1002 gridView_.communicate(*sumHandle,
1003 Dune::InteriorBorder_InteriorBorder_Interface,
1004 Dune::ForwardCommunication);
1010 return std::sqrt(asImp_().
gridView().comm().sum(dest.two_norm2()));
1032 ElementContext elemCtx(simulator_);
1033 ElementIterator elemIt = threadedElemIt.beginParallel();
1034 LocalEvalBlockVector elemStorage;
1038 elemCtx.setEnableStorageCache(
false);
1040 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
1041 const Element& elem = *elemIt;
1042 if (elem.partitionType() != Dune::InteriorEntity) {
1046 elemCtx.updateStencil(elem);
1047 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
1049 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
1050 elemStorage.resize(numPrimaryDof);
1052 localResidual(threadId).evalStorage(elemStorage, elemCtx, timeIdx);
1055 for (
unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
1056 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1057 storage[eqIdx] += Toolbox::value(elemStorage[dofIdx][eqIdx]);
1064 storage = gridView_.comm().sum(storage);
1075 [[maybe_unused]]
bool verbose =
false)
const
1078 Scalar totalBoundaryArea(0.0);
1079 Scalar totalVolume(0.0);
1080 EvalEqVector totalRate(0.0);
1084 if (tolerance <= 0) {
1086 simulator_.model().newtonMethod().tolerance() *
1087 simulator_.model().gridTotalVolume() *
1092 assert(historySize == 2);
1094 EqVector storageBeginTimeStep(0.0);
1097 EqVector storageEndTimeStep(0.0);
1101 ElementContext elemCtx(simulator_);
1102 elemCtx.setEnableStorageCache(
false);
1103 for (
const auto& elem : elements(simulator_.gridView())) {
1104 if (elem.partitionType() != Dune::InteriorEntity) {
1108 elemCtx.updateAll(elem);
1111 if (elemCtx.onBoundary()) {
1112 BoundaryContext boundaryCtx(elemCtx);
1114 for (
unsigned faceIdx = 0; faceIdx < boundaryCtx.numBoundaryFaces(0); ++faceIdx) {
1115 BoundaryRateVector values;
1116 simulator_.problem().boundary(values,
1120 Valgrind::CheckDefined(values);
1122 const unsigned dofIdx = boundaryCtx.interiorScvIndex(faceIdx, 0);
1123 const auto& insideIntQuants = elemCtx.intensiveQuantities(dofIdx, 0);
1125 const Scalar bfArea =
1126 boundaryCtx.boundarySegmentArea(faceIdx, 0) *
1127 insideIntQuants.extrusionFactor();
1129 for (
unsigned i = 0; i < values.size(); ++i) {
1130 values[i] *= bfArea;
1133 totalBoundaryArea += bfArea;
1134 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1135 totalRate[eqIdx] += values[eqIdx];
1141 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); ++dofIdx) {
1143 simulator_.problem().source(values,
1147 Valgrind::CheckDefined(values);
1149 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, 0);
1151 elemCtx.dofVolume(dofIdx, 0) *
1152 intQuants.extrusionFactor();
1153 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1154 totalRate[eqIdx] += -dofVolume*Toolbox::value(values[eqIdx]);
1156 totalVolume += dofVolume;
1161 const auto& comm = simulator_.gridView().comm();
1162 totalRate = comm.sum(totalRate);
1163 totalBoundaryArea = comm.sum(totalBoundaryArea);
1164 totalVolume = comm.sum(totalVolume);
1166 if (comm.rank() == 0) {
1167 EqVector storageRate = storageBeginTimeStep;
1168 storageRate -= storageEndTimeStep;
1169 storageRate /= simulator_.timeStepSize();
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])) <<
" ";
1181 for (
unsigned eqIdx = 0; eqIdx < EqVector::dimension; ++eqIdx) {
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);
1197 {
return dofTotalVolume_[globalIdx]; }
1205 {
return isLocalDof_[globalIdx]; }
1212 {
return gridTotalVolume_; }
1220 {
return solution_[timeIdx]->blockVector(); }
1226 {
return solution_[timeIdx]->blockVector(); }
1233 {
return solution_[timeIdx]->blockVector(); }
1241 {
return *linearizer_; }
1248 {
return *linearizer_; }
1259 {
return localLinearizer_[openMpThreadId]; }
1265 {
return localLinearizer_[openMpThreadId]; }
1271 {
return asImp_().localLinearizer(openMpThreadId).localResidual(); }
1277 {
return asImp_().localLinearizer(openMpThreadId).localResidual(); }
1288 const Scalar absPv = std::abs(asImp_().
solution(1)[globalDofIdx][pvIdx]);
1289 return 1.0 / std::max(absPv, 1.0);
1311 const PrimaryVariables& pv1,
1312 const PrimaryVariables& pv2)
const
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);
1321 result = std::max(result, eqErr);
1333 const TimerGuard prePostProcessGuard(prePostProcessTimer_);
1336 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1340 for (std::size_t i = 0; i < asImp_().solution(0).size(); ++i) {
1341 asImp_().solution(timeIdx)[i].checkDefined();
1347 prePostProcessTimer_.halt();
1348 linearizeTimer_.halt();
1350 updateTimer_.halt();
1352 prePostProcessTimer_.start();
1353 asImp_().updateBegin();
1354 prePostProcessTimer_.stop();
1356 bool converged =
false;
1359 converged = newtonMethod_.apply();
1362 prePostProcessTimer_ += newtonMethod_.prePostProcessTimer();
1363 linearizeTimer_ += newtonMethod_.linearizeTimer();
1364 solveTimer_ += newtonMethod_.solveTimer();
1365 updateTimer_ += newtonMethod_.updateTimer();
1371 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1375 for (std::size_t i = 0; i < asImp_().solution(0).size(); ++i) {
1376 asImp_().solution(timeIdx)[i].checkDefined();
1381 prePostProcessTimer_ += newtonMethod_.prePostProcessTimer();
1382 linearizeTimer_ += newtonMethod_.linearizeTimer();
1383 solveTimer_ += newtonMethod_.solveTimer();
1384 updateTimer_ += newtonMethod_.updateTimer();
1386 prePostProcessTimer_.start();
1388 asImp_().updateSuccessful();
1391 asImp_().updateFailed();
1393 prePostProcessTimer_.stop();
1396 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1400 for (std::size_t i = 0; i < asImp_().solution(0).size(); ++i) {
1401 asImp_().solution(timeIdx)[i].checkDefined();
1439 throw std::invalid_argument(
"Grid adaptation need to be implemented for "
1440 "specific settings of grid and function spaces");
1454 invalidateAndUpdateIntensiveQuantities(0);
1457 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1461 for (std::size_t i = 0; i < asImp_().solution(0).size(); ++i) {
1462 asImp_().solution(timeIdx)[i].checkDefined();
1478 if (this->enableGridAdaptation_) {
1479 asImp_().adaptGrid();
1486 asImp_().shiftStorageCache(1);
1490 asImp_().shiftIntensiveQuantityCache(1);
1500 template <
class Restarter>
1503 throw std::runtime_error(
"Not implemented: The discretization chosen for this problem "
1504 "does not support restart files. (serialize() method unimplemented)");
1514 template <
class Restarter>
1517 throw std::runtime_error(
"Not implemented: The discretization chosen for this problem "
1518 "does not support restart files. (deserialize() method unimplemented)");
1529 template <
class DofEntity>
1531 const DofEntity& dof)
1533 const unsigned dofIdx =
static_cast<unsigned>(asImp_().dofMapper().index(dof));
1536 if (!outstream.good()) {
1537 throw std::runtime_error(
"Could not serialize degree of freedom " +
1538 std::to_string(dofIdx));
1541 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1542 outstream <<
solution(0)[dofIdx][eqIdx] <<
" ";
1554 template <
class DofEntity>
1556 const DofEntity& dof)
1558 const unsigned dofIdx =
static_cast<unsigned>(asImp_().dofMapper().index(dof));
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));
1565 instream >>
solution(0)[dofIdx][eqIdx];
1573 {
throw std::logic_error(
"The discretization class must implement the numGridDof() method!"); }
1580 return std::accumulate(auxEqModules_.begin(), auxEqModules_.end(),
1582 [](
const auto acc,
const auto& mod)
1583 { return acc + mod->numDofs(); });
1597 {
throw std::logic_error(
"The discretization class must implement the dofMapper() method!"); }
1603 {
return vertexMapper_; }
1609 {
return elementMapper_; }
1617 linearizer_ = std::make_unique<Linearizer>();
1618 linearizer_->init(simulator_);
1634 std::ostringstream oss;
1635 oss <<
"primary variable_" << pvIdx;
1646 std::ostringstream oss;
1647 oss <<
"equation_" << eqIdx;
1664 { outputModules_.push_back(std::move(newModule)); }
1674 template <
class VtkMultiWriter>
1676 const SolutionVector& u,
1677 const GlobalEqVector& deltaU)
const
1679 using ScalarBuffer = std::vector<double>;
1681 GlobalEqVector globalResid(u.size());
1682 asImp_().globalResidual(globalResid, u);
1685 const std::size_t numDof = asImp_().numGridDof();
1688 std::array<ScalarBuffer*, numEq> def;
1689 std::array<ScalarBuffer*, numEq> delta;
1690 std::array<ScalarBuffer*, numEq> priVars;
1691 std::array<ScalarBuffer*, numEq> priVarWeight;
1694 for (
unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
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];
1711 PrimaryVariables uOld(u[globalIdx]);
1712 PrimaryVariables uNew(uOld);
1713 uNew -= deltaU[globalIdx];
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);
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;
1730 DiscBaseOutputModule::attachScalarDofData_(writer, *relError,
"relative error");
1731 DiscBaseOutputModule::attachScalarDofData_(writer, *normalizedRelError,
"normalized relative error");
1733 for (
unsigned i = 0; i < numEq; ++i) {
1734 std::ostringstream oss;
1735 oss.str(
""); oss <<
"priVar_" << asImp_().primaryVarName(i);
1736 DiscBaseOutputModule::attachScalarDofData_(writer,
1740 oss.str(
""); oss <<
"delta_" << asImp_().primaryVarName(i);
1741 DiscBaseOutputModule::attachScalarDofData_(writer,
1745 oss.str(
""); oss <<
"weight_" << asImp_().primaryVarName(i);
1746 DiscBaseOutputModule::attachScalarDofData_(writer,
1750 oss.str(
""); oss <<
"defect_" << asImp_().eqName(i);
1751 DiscBaseOutputModule::attachScalarDofData_(writer,
1756 asImp_().prepareOutputFields();
1757 asImp_().appendOutputFields(writer);
1766 const bool needFullContextUpdate =
1767 std::ranges::any_of(outputModules_,
1769 {
return mod->needExtensiveQuantities(); });
1770 std::ranges::for_each(outputModules_,
1771 [](
auto& mod) { mod->allocBuffers(); });
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) {
1788 if (needFullContextUpdate) {
1789 elemCtx.updateAll(elem);
1792 elemCtx.updatePrimaryStencil(elem);
1793 elemCtx.updatePrimaryIntensiveQuantities(0);
1796 std::ranges::for_each(outputModules_,
1797 [&elemCtx](
auto& mod) { mod->processElement(elemCtx); });
1808 std::ranges::for_each(outputModules_,
1809 [&writer](
auto& mod) { mod->commitBuffers(writer); });
1816 {
return gridView_; }
1832 auxEqModules_.push_back(auxMod);
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");
1841 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1855 auxEqModules_.clear();
1856 linearizer_->eraseMatrix();
1857 newtonMethod_.eraseMatrix();
1864 {
return auxEqModules_.size(); }
1870 {
return auxEqModules_[auxEqModIdx]; }
1876 {
return auxEqModules_[auxEqModIdx]; }
1882 {
return enableIntensiveQuantityCache_ || enableThermodynamicHints_; }
1884 const Timer& prePostProcessTimer()
const
1885 {
return prePostProcessTimer_; }
1887 const Timer& linearizeTimer()
const
1888 {
return linearizeTimer_; }
1890 const Timer& solveTimer()
const
1891 {
return solveTimer_; }
1893 const Timer& updateTimer()
const
1894 {
return updateTimer_; }
1896 template<
class Serializer>
1897 void serializeOp(Serializer& serializer)
1900 using Helper =
typename BaseDiscretization::template SerializeHelper<Serializer>;
1901 Helper::serializeOp(serializer, solution_);
1904 bool operator==(
const FvBaseDiscretization& rhs)
const
1906 return std::ranges::equal(this->solution_, rhs.solution_,
1907 [](
const auto& x,
const auto& y)
1908 { return *x == *y; });
1912 void resizeAndResetIntensiveQuantitiesCache_()
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, 0);
1925 const std::size_t numDof = asImp_().numGridDof();
1926 cachedIntensiveQuantityHistorySize_ = simulator_.problem().intensiveQuantityHistorySize();
1927 const unsigned intensiveHistorySize = cachedIntensiveQuantityHistorySize_;
1930 intensiveQuantityCache_.resize(intensiveHistorySize);
1931 intensiveQuantityCacheUpToDate_.resize(intensiveHistorySize);
1933 for(
unsigned timeIdx = 0; timeIdx < intensiveHistorySize; ++timeIdx) {
1934 intensiveQuantityCache_[timeIdx].resize(numDof);
1935 intensiveQuantityCacheUpToDate_[timeIdx].resize(numDof);
1941 template <
class Context>
1942 void supplementInitialSolution_(PrimaryVariables&,
1965 {
return localLinearizer_.localResidual(); }
1971 {
return gridView_.comm().rank() == 0; }
1973 Implementation& asImp_()
1974 {
return *
static_cast<Implementation*
>(
this); }
1976 const Implementation& asImp_()
const
1977 {
return *
static_cast<const Implementation*
>(
this); }
1981 Simulator& simulator_;
1987 ElementMapper elementMapper_;
1988 VertexMapper vertexMapper_;
1991 std::vector<BaseAuxiliaryModule<TypeTag>*> auxEqModules_;
1993 NewtonMethod newtonMethod_;
1995 Timer prePostProcessTimer_;
1996 Timer linearizeTimer_;
2001 std::vector<LocalLinearizer> localLinearizer_;
2004 std::unique_ptr<Linearizer> linearizer_;
2008 mutable std::vector<IntensiveQuantitiesVector> intensiveQuantityCache_;
2011 mutable std::vector<std::vector<unsigned char>> intensiveQuantityCacheUpToDate_;
2013 std::array<std::unique_ptr<DiscreteFunction>, historySize> solution_;
2015 std::list<std::unique_ptr<BaseOutputModule<TypeTag>>> outputModules_;
2017 Scalar gridTotalVolume_;
2018 std::vector<Scalar> dofTotalVolume_;
2019 std::vector<bool> isLocalDof_;
2021 mutable std::array<GlobalEqVector, historySize> storageCache_;
2024 mutable std::array<std::vector<unsigned char>, historySize> storageCacheUpToDate_;
2026 bool enableGridAdaptation_;
2027 bool enableIntensiveQuantityCache_;
2028 bool enableStorageCache_;
2029 bool enableThermodynamicHints_;
2031 mutable unsigned cachedIntensiveQuantityHistorySize_;
2039template<
class TypeTag>
2040class FvBaseDiscretizationNoAdapt :
public FvBaseDiscretization<TypeTag>
2042 using ParentType = FvBaseDiscretization<TypeTag>;
2049 template<
class Serializer>
2052 template<
class SolutionType>
2053 static void serializeOp(Serializer& serializer,
2062 explicit FvBaseDiscretizationNoAdapt(Simulator& simulator)
2063 : ParentType(simulator)
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");
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);
This is a stand-alone version of boost::alignment::aligned_allocator from Boost 1....
Base class for specifying auxiliary equations.
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
The base class for the finite volume discretization schemes without adaptation.
Definition fvbasediscretization.hh:2041
Definition fvbasediscretization.hh:350
The base class for the finite volume discretization schemes.
Definition fvbasediscretization.hh:298
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 adaptGrid()
Called by the update() method when the grid should be refined.
Definition fvbasediscretization.hh:1437
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
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
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
void updateFailed()
Called by the update() method if it was unsuccessful.
Definition fvbasediscretization.hh:1448
void updateBegin()
Called by the update() method before it tries to apply the newton method.
Definition fvbasediscretization.hh:1424
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
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
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 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
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
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
const LocalLinearizer & localLinearizer(unsigned openMpThreadId) const
Returns the local jacobian which calculates the local stiffness matrix for an arbitrary element.
Definition fvbasediscretization.hh:1258
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
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::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
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
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
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
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
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
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
void invalidateStorageCacheEntry(unsigned globalIdx, unsigned timeIdx) const
Invalidate the storage cache for a given DOF and time index.
Definition fvbasediscretization.hh:894
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:135
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.
Definition threadmanager.cpp:84
Provides an STL-iterator like interface to iterate over the enties of a GridView in OpenMP threaded a...
Definition threadedentityiterator.hh:42
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
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
Calculates the local residual and its Jacobian for a single element of the grid.
Represents all quantities which available on boundary segments.
Class to specify constraints for a finite volume spatial discretization.
Represents all quantities which available for calculating constraints.
This class stores an array of IntensiveQuantities objects, one intensive quantities object for each o...
Provide the properties at a face which make sense independently of the conserved quantities.
Calculates the Jacobian of the local residual for finite volume spatial discretizations using a finit...
This class calculates gradients of arbitrary quantities at flux integration points using the two-poin...
Base class for the model specific class which provides access to all intensive (i....
The common code for the linearizers of non-linear systems of equations.
Element-wise caculation of the residual matrix for models based on a finite volume spatial discretiza...
A Newton method for models using a finite volume discretization.
Represents the primary variables used by the a model.
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.
The generic type tag for problems using the immiscible multi-phase model.
Definition blackoilmodel.hh:82
This file contains a set of helper functions used by VFPProd / VFPInj.
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
This is a grid manager which does not create any border list.
void Register(const char *usageString)
Register a run-time parameter.
Definition parametersystem.hpp:292
Manages the initializing and running of time dependent problems.
Definition fvbasediscretization.hh:2051
Definition fvbasediscretization.hh:272
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
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
Definition fvbaseproperties.hh:77
The secondary variables of all degrees of freedom in an element's stencil.
Definition fvbaseproperties.hh:144
A vector of holding a quantity for each equation for each DOF of an element.
Definition fvbaseproperties.hh:112
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
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
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
Vector containing volumetric or areal rates of quantities.
Definition fvbaseproperties.hh:116
Manages the simulation time.
Definition basicproperties.hh:120
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
use locking to prevent race conditions when linearizing the global system of equations in multi-threa...
Definition fvbaseproperties.hh:185
Specify whether to use volumetric residuals or not.
Definition fvbaseproperties.hh:241
The mapper to find the global index of a vertex.
Definition fvbaseproperties.hh:211
Simplifies multi-threaded capabilities.
Provides an encapsulation to measure the system time.
A simple class which makes sure that a timer gets stopped if an exception is thrown.
VTK output module for the fluid composition.