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::equal(this->blockVector_.begin(), this->blockVector_.end(),
379 wrapper.blockVector_.begin(), wrapper.blockVector_.end());
382 template<
class Serializer>
383 void serializeOp(Serializer& serializer)
385 serializer(blockVector_);
394 explicit FvBaseDiscretization(Simulator& simulator)
395 : simulator_(simulator)
397 , elementMapper_(gridView_, Dune::mcmgElementLayout())
398 , vertexMapper_(gridView_, Dune::mcmgVertexLayout())
399 , newtonMethod_(simulator)
400 , localLinearizer_(ThreadManager::maxThreads())
401 , linearizer_(std::make_unique<Linearizer>())
402 , enableGridAdaptation_(Parameters::Get<Parameters::EnableGridAdaptation>() )
403 , enableIntensiveQuantityCache_(Parameters::Get<Parameters::EnableIntensiveQuantityCache>())
404 , enableStorageCache_(Parameters::Get<Parameters::EnableStorageCache>())
405 , enableThermodynamicHints_(Parameters::Get<Parameters::EnableThermodynamicHints>())
406 , cachedIntensiveQuantityHistorySize_(-1)
408 const bool isEcfv = std::is_same_v<Discretization, EcfvDiscretization<TypeTag>>;
409 if (enableGridAdaptation_ && !isEcfv) {
410 throw std::invalid_argument(
"Grid adaptation currently only works for the "
411 "element-centered finite volume discretization (is: " +
412 Dune::className<Discretization>() +
")");
415 PrimaryVariables::init();
418 asImp_().registerOutputModules_();
422 FvBaseDiscretization(
const FvBaseDiscretization&) =
delete;
429 Linearizer::registerParameters();
430 LocalLinearizer::registerParameters();
431 LocalResidual::registerParameters();
432 GradientCalculator::registerParameters();
433 IntensiveQuantities::registerParameters();
434 ExtensiveQuantities::registerParameters();
436 Linearizer::registerParameters();
437 PrimaryVariables::registerParameters();
442 (
"Enable adaptive grid refinement/coarsening");
444 (
"Global switch for turning on writing VTK files");
446 (
"Enable thermodynamic hints");
448 (
"Turn on caching of intensive quantities");
450 (
"Store previous storage terms and avoid re-calculating them.");
452 (
"The directory to which result files are written");
461 const std::size_t numDof = asImp_().numGridDof();
462 dofTotalVolume_.resize(numDof);
463 std::fill(dofTotalVolume_.begin(), dofTotalVolume_.end(), 0.0);
465 ElementContext elemCtx(simulator_);
466 gridTotalVolume_ = 0.0;
469 for (
const auto& elem : elements(gridView_)) {
472 if (elem.partitionType() != Dune::InteriorEntity) {
477 elemCtx.updateStencil(elem);
478 const auto& stencil = elemCtx.stencil(0);
481 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); dofIdx++) {
483 const unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
485 const Scalar dofVolume = stencil.subControlVolume(dofIdx).volume();
486 dofTotalVolume_[globalIdx] += dofVolume;
487 gridTotalVolume_ += dofVolume;
494 isLocalDof_.resize(numDof);
495 for (
unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
496 isLocalDof_[dofIdx] = (dofTotalVolume_[dofIdx] != 0.0);
500 const auto sumHandle =
501 GridCommHandleFactory::template sumHandle<Scalar>(dofTotalVolume_,
503 gridView_.communicate(*sumHandle,
504 Dune::InteriorBorder_All_Interface,
505 Dune::ForwardCommunication);
508 gridTotalVolume_ = gridView_.comm().sum(gridTotalVolume_);
510 linearizer_->init(simulator_);
512 localLinearizer_[threadId].init(simulator_);
515 resizeAndResetIntensiveQuantitiesCache_();
517 newtonMethod_.finishInit();
524 {
return enableGridAdaptation_; }
533 SolutionVector& uCur = asImp_().solution(0);
536 ElementContext elemCtx(simulator_);
539 for (
const auto& elem : elements(gridView_)) {
542 if (elem.partitionType() != Dune::InteriorEntity) {
547 elemCtx.updateStencil(elem);
550 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); ++dofIdx) {
552 const unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
556 simulator_.problem().initial(uCur[globalIdx], elemCtx, dofIdx, 0);
557 asImp_().supplementInitialSolution_(uCur[globalIdx], elemCtx, dofIdx, 0);
558 uCur[globalIdx].checkDefined();
563 asImp_().syncOverlap();
566 for (
unsigned timeIdx = 1; timeIdx < historySize; ++timeIdx) {
573 resizeAndResetIntensiveQuantitiesCache_();
575 simulator_.problem().initialSolutionApplied();
578 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
579 const auto& sol =
solution(timeIdx);
580 for (
unsigned dofIdx = 0; dofIdx < sol.size(); ++dofIdx) {
581 sol[dofIdx].checkDefined();
600 {
return newtonMethod_; }
606 {
return newtonMethod_; }
625 if (!enableThermodynamicHints_) {
646 if (!enableIntensiveQuantityCache_ ||
647 timeIdx >= cachedIntensiveQuantityHistorySize_ ||
648 !intensiveQuantityCacheUpToDate_[timeIdx][globalIdx]) {
657 if (timeIdx > 0 && enableStorageCache_ && intensiveQuantityCache_[timeIdx].empty()) {
661 return &intensiveQuantityCache_[timeIdx][globalIdx];
674 unsigned timeIdx)
const
680 intensiveQuantityCache_[timeIdx][globalIdx] = intQuants;
681 intensiveQuantityCacheUpToDate_[timeIdx][globalIdx] =
true;
700 intensiveQuantityCacheUpToDate_[timeIdx][globalIdx] = newValue ? 1 : 0;
707 {
return cachedIntensiveQuantityHistorySize_; }
716 if (timeIdx >= cachedIntensiveQuantityHistorySize_) {
721 std::fill(intensiveQuantityCacheUpToDate_[timeIdx].begin(),
722 intensiveQuantityCacheUpToDate_[timeIdx].end(),
727 void invalidateAndUpdateIntensiveQuantities(
unsigned timeIdx)
const
737 ElementContext elemCtx(simulator_);
738 ElementIterator elemIt = threadedElemIt.beginParallel();
739 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
740 const Element& elem = *elemIt;
741 elemCtx.updatePrimaryStencil(elem);
742 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
747 template <
class Gr
idViewType>
748 void invalidateAndUpdateIntensiveQuantities(
unsigned timeIdx,
const GridViewType&
gridView)
const
751 ThreadedEntityIterator<GridViewType, 0> threadedElemIt(
gridView);
757 ElementContext elemCtx(simulator_);
758 auto elemIt = threadedElemIt.beginParallel();
759 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
760 if (elemIt->partitionType() != Dune::InteriorEntity) {
763 const Element& elem = *elemIt;
764 elemCtx.updatePrimaryStencil(elem);
766 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
767 for (
unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
768 const unsigned globalIndex = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
772 elemCtx.updatePrimaryIntensiveQuantities(0);
801 const unsigned intensiveHistorySize = cachedIntensiveQuantityHistorySize_;
802 for (
unsigned timeIdx = 0; timeIdx < intensiveHistorySize - numSlots; ++timeIdx) {
803 intensiveQuantityCache_[timeIdx + numSlots] = intensiveQuantityCache_[timeIdx];
804 intensiveQuantityCacheUpToDate_[timeIdx + numSlots] = intensiveQuantityCacheUpToDate_[timeIdx];
819 {
return enableStorageCache_; }
845 if (!enableStorageCache_ ||
846 timeIdx >= historySize ||
847 !storageCacheUpToDate_[timeIdx][globalIdx]) {
848 throw std::logic_error(
"Cached storage is not available or up to date for the requested "
849 "global index and time index. Make sure storage cache is enabled "
850 "and the entry is valid before calling this method.");
853 return storageCache_[timeIdx][globalIdx];
869 if (!enableStorageCache_ || timeIdx >= historySize) {
873 storageCache_[timeIdx][globalIdx] = value;
874 storageCacheUpToDate_[timeIdx][globalIdx] = 1;
885 if (!enableStorageCache_ || timeIdx >= historySize) {
888 return storageCacheUpToDate_[timeIdx][globalIdx] != 0;
899 if (enableStorageCache_ && timeIdx < historySize) {
900 storageCacheUpToDate_[timeIdx][globalIdx] = 0;
911 if (enableStorageCache_ && timeIdx < historySize) {
912 std::fill(storageCacheUpToDate_[timeIdx].begin(),
913 storageCacheUpToDate_[timeIdx].end(),
932 if (enableStorageCache_ && !simulator_.problem().recycleFirstIterationStorage()) {
933 for (
unsigned timeIdx = 0; timeIdx < historySize - numSlots; ++timeIdx) {
934 storageCache_[timeIdx + numSlots] = storageCache_[timeIdx];
935 storageCacheUpToDate_[timeIdx + numSlots] = storageCacheUpToDate_[timeIdx];
950 const SolutionVector& u)
const
953 const Scalar res = asImp_().globalResidual(dest);
977 ElementContext elemCtx(simulator_);
978 ElementIterator elemIt = threadedElemIt.beginParallel();
979 LocalEvalBlockVector residual, storageTerm;
981 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
982 const Element& elem = *elemIt;
983 if (elem.partitionType() != Dune::InteriorEntity) {
987 elemCtx.updateAll(elem);
988 residual.resize(elemCtx.numDof(0));
989 storageTerm.resize(elemCtx.numPrimaryDof(0));
990 asImp_().localResidual(threadId).eval(residual, elemCtx);
992 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(0);
994 for (
unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
995 const unsigned globalI = elemCtx.globalSpaceIndex(dofIdx, 0);
996 for (
unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
997 dest[globalI][eqIdx] += Toolbox::value(residual[dofIdx][eqIdx]);
1005 const auto sumHandle =
1006 GridCommHandleFactory::template sumHandle<EqVector>(dest, asImp_().
dofMapper());
1007 gridView_.communicate(*sumHandle,
1008 Dune::InteriorBorder_InteriorBorder_Interface,
1009 Dune::ForwardCommunication);
1015 return std::sqrt(asImp_().
gridView().comm().sum(dest.two_norm2()));
1037 ElementContext elemCtx(simulator_);
1038 ElementIterator elemIt = threadedElemIt.beginParallel();
1039 LocalEvalBlockVector elemStorage;
1043 elemCtx.setEnableStorageCache(
false);
1045 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
1046 const Element& elem = *elemIt;
1047 if (elem.partitionType() != Dune::InteriorEntity) {
1051 elemCtx.updateStencil(elem);
1052 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
1054 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
1055 elemStorage.resize(numPrimaryDof);
1057 localResidual(threadId).evalStorage(elemStorage, elemCtx, timeIdx);
1060 for (
unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
1061 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1062 storage[eqIdx] += Toolbox::value(elemStorage[dofIdx][eqIdx]);
1069 storage = gridView_.comm().sum(storage);
1080 [[maybe_unused]]
bool verbose =
false)
const
1083 Scalar totalBoundaryArea(0.0);
1084 Scalar totalVolume(0.0);
1085 EvalEqVector totalRate(0.0);
1089 if (tolerance <= 0) {
1091 simulator_.model().newtonMethod().tolerance() *
1092 simulator_.model().gridTotalVolume() *
1097 assert(historySize == 2);
1099 EqVector storageBeginTimeStep(0.0);
1102 EqVector storageEndTimeStep(0.0);
1106 ElementContext elemCtx(simulator_);
1107 elemCtx.setEnableStorageCache(
false);
1108 for (
const auto& elem : elements(simulator_.gridView())) {
1109 if (elem.partitionType() != Dune::InteriorEntity) {
1113 elemCtx.updateAll(elem);
1116 if (elemCtx.onBoundary()) {
1117 BoundaryContext boundaryCtx(elemCtx);
1119 for (
unsigned faceIdx = 0; faceIdx < boundaryCtx.numBoundaryFaces(0); ++faceIdx) {
1120 BoundaryRateVector values;
1121 simulator_.problem().boundary(values,
1125 Valgrind::CheckDefined(values);
1127 const unsigned dofIdx = boundaryCtx.interiorScvIndex(faceIdx, 0);
1128 const auto& insideIntQuants = elemCtx.intensiveQuantities(dofIdx, 0);
1130 const Scalar bfArea =
1131 boundaryCtx.boundarySegmentArea(faceIdx, 0) *
1132 insideIntQuants.extrusionFactor();
1134 for (
unsigned i = 0; i < values.size(); ++i) {
1135 values[i] *= bfArea;
1138 totalBoundaryArea += bfArea;
1139 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1140 totalRate[eqIdx] += values[eqIdx];
1146 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); ++dofIdx) {
1148 simulator_.problem().source(values,
1152 Valgrind::CheckDefined(values);
1154 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, 0);
1156 elemCtx.dofVolume(dofIdx, 0) *
1157 intQuants.extrusionFactor();
1158 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1159 totalRate[eqIdx] += -dofVolume*Toolbox::value(values[eqIdx]);
1161 totalVolume += dofVolume;
1166 const auto& comm = simulator_.gridView().comm();
1167 totalRate = comm.sum(totalRate);
1168 totalBoundaryArea = comm.sum(totalBoundaryArea);
1169 totalVolume = comm.sum(totalVolume);
1171 if (comm.rank() == 0) {
1172 EqVector storageRate = storageBeginTimeStep;
1173 storageRate -= storageEndTimeStep;
1174 storageRate /= simulator_.timeStepSize();
1176 std::cout <<
"storage at beginning of time step: " << storageBeginTimeStep <<
"\n";
1177 std::cout <<
"storage at end of time step: " << storageEndTimeStep <<
"\n";
1178 std::cout <<
"rate based on storage terms: " << storageRate <<
"\n";
1179 std::cout <<
"rate based on source and boundary terms: " << totalRate <<
"\n";
1180 std::cout <<
"difference in rates: ";
1181 for (
unsigned eqIdx = 0; eqIdx < EqVector::dimension; ++eqIdx) {
1182 std::cout << (storageRate[eqIdx] - Toolbox::value(totalRate[eqIdx])) <<
" ";
1186 for (
unsigned eqIdx = 0; eqIdx < EqVector::dimension; ++eqIdx) {
1188 (std::abs(storageRate[eqIdx]) + Toolbox::value(totalRate[eqIdx])) * tolerance;
1189 eps = std::max(tolerance, eps);
1190 assert(std::abs(storageRate[eqIdx] - Toolbox::value(totalRate[eqIdx])) <= eps);
1202 {
return dofTotalVolume_[globalIdx]; }
1210 {
return isLocalDof_[globalIdx]; }
1217 {
return gridTotalVolume_; }
1225 {
return solution_[timeIdx]->blockVector(); }
1231 {
return solution_[timeIdx]->blockVector(); }
1238 {
return solution_[timeIdx]->blockVector(); }
1246 {
return *linearizer_; }
1253 {
return *linearizer_; }
1264 {
return localLinearizer_[openMpThreadId]; }
1270 {
return localLinearizer_[openMpThreadId]; }
1276 {
return asImp_().localLinearizer(openMpThreadId).localResidual(); }
1282 {
return asImp_().localLinearizer(openMpThreadId).localResidual(); }
1293 const Scalar absPv = std::abs(asImp_().
solution(1)[globalDofIdx][pvIdx]);
1294 return 1.0 / std::max(absPv, 1.0);
1316 const PrimaryVariables& pv1,
1317 const PrimaryVariables& pv2)
const
1319 Scalar result = 0.0;
1320 for (
unsigned j = 0; j < numEq; ++j) {
1321 const Scalar weight = asImp_().primaryVarWeight(vertexIdx, j);
1322 const Scalar eqErr = std::abs((pv1[j] - pv2[j])*weight);
1326 result = std::max(result, eqErr);
1338 const TimerGuard prePostProcessGuard(prePostProcessTimer_);
1341 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1345 for (std::size_t i = 0; i < asImp_().solution(0).size(); ++i) {
1346 asImp_().solution(timeIdx)[i].checkDefined();
1352 prePostProcessTimer_.halt();
1353 linearizeTimer_.halt();
1355 updateTimer_.halt();
1357 prePostProcessTimer_.start();
1358 asImp_().updateBegin();
1359 prePostProcessTimer_.stop();
1361 bool converged =
false;
1364 converged = newtonMethod_.apply();
1367 prePostProcessTimer_ += newtonMethod_.prePostProcessTimer();
1368 linearizeTimer_ += newtonMethod_.linearizeTimer();
1369 solveTimer_ += newtonMethod_.solveTimer();
1370 updateTimer_ += newtonMethod_.updateTimer();
1376 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1380 for (std::size_t i = 0; i < asImp_().solution(0).size(); ++i) {
1381 asImp_().solution(timeIdx)[i].checkDefined();
1386 prePostProcessTimer_ += newtonMethod_.prePostProcessTimer();
1387 linearizeTimer_ += newtonMethod_.linearizeTimer();
1388 solveTimer_ += newtonMethod_.solveTimer();
1389 updateTimer_ += newtonMethod_.updateTimer();
1391 prePostProcessTimer_.start();
1393 asImp_().updateSuccessful();
1396 asImp_().updateFailed();
1398 prePostProcessTimer_.stop();
1401 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1405 for (std::size_t i = 0; i < asImp_().solution(0).size(); ++i) {
1406 asImp_().solution(timeIdx)[i].checkDefined();
1444 throw std::invalid_argument(
"Grid adaptation need to be implemented for "
1445 "specific settings of grid and function spaces");
1459 invalidateAndUpdateIntensiveQuantities(0);
1462 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1466 for (std::size_t i = 0; i < asImp_().solution(0).size(); ++i) {
1467 asImp_().solution(timeIdx)[i].checkDefined();
1483 if (this->enableGridAdaptation_) {
1484 asImp_().adaptGrid();
1491 asImp_().shiftStorageCache(1);
1495 asImp_().shiftIntensiveQuantityCache(1);
1505 template <
class Restarter>
1508 throw std::runtime_error(
"Not implemented: The discretization chosen for this problem "
1509 "does not support restart files. (serialize() method unimplemented)");
1519 template <
class Restarter>
1522 throw std::runtime_error(
"Not implemented: The discretization chosen for this problem "
1523 "does not support restart files. (deserialize() method unimplemented)");
1534 template <
class DofEntity>
1536 const DofEntity& dof)
1538 const unsigned dofIdx =
static_cast<unsigned>(asImp_().dofMapper().index(dof));
1541 if (!outstream.good()) {
1542 throw std::runtime_error(
"Could not serialize degree of freedom " +
1543 std::to_string(dofIdx));
1546 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1547 outstream <<
solution(0)[dofIdx][eqIdx] <<
" ";
1559 template <
class DofEntity>
1561 const DofEntity& dof)
1563 const unsigned dofIdx =
static_cast<unsigned>(asImp_().dofMapper().index(dof));
1565 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
1566 if (!instream.good()) {
1567 throw std::runtime_error(
"Could not deserialize degree of freedom " +
1568 std::to_string(dofIdx));
1570 instream >>
solution(0)[dofIdx][eqIdx];
1578 {
throw std::logic_error(
"The discretization class must implement the numGridDof() method!"); }
1585 return std::accumulate(auxEqModules_.begin(), auxEqModules_.end(),
1587 [](
const auto acc,
const auto& mod)
1588 { return acc + mod->numDofs(); });
1602 {
throw std::logic_error(
"The discretization class must implement the dofMapper() method!"); }
1608 {
return vertexMapper_; }
1614 {
return elementMapper_; }
1622 linearizer_ = std::make_unique<Linearizer>();
1623 linearizer_->init(simulator_);
1639 std::ostringstream oss;
1640 oss <<
"primary variable_" << pvIdx;
1651 std::ostringstream oss;
1652 oss <<
"equation_" << eqIdx;
1669 { outputModules_.push_back(std::move(newModule)); }
1679 template <
class VtkMultiWriter>
1681 const SolutionVector& u,
1682 const GlobalEqVector& deltaU)
const
1684 using ScalarBuffer = std::vector<double>;
1686 GlobalEqVector globalResid(u.size());
1687 asImp_().globalResidual(globalResid, u);
1690 const std::size_t numDof = asImp_().numGridDof();
1693 std::array<ScalarBuffer*, numEq> def;
1694 std::array<ScalarBuffer*, numEq> delta;
1695 std::array<ScalarBuffer*, numEq> priVars;
1696 std::array<ScalarBuffer*, numEq> priVarWeight;
1699 for (
unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
1706 Scalar minRelErr = 1e30;
1707 Scalar maxRelErr = -1e30;
1708 for (
unsigned globalIdx = 0; globalIdx < numDof; ++ globalIdx) {
1709 for (
unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
1710 (*priVars[pvIdx])[globalIdx] = u[globalIdx][pvIdx];
1711 (*priVarWeight[pvIdx])[globalIdx] = asImp_().primaryVarWeight(globalIdx, pvIdx);
1712 (*delta[pvIdx])[globalIdx] = - deltaU[globalIdx][pvIdx];
1713 (*def[pvIdx])[globalIdx] = globalResid[globalIdx][pvIdx];
1716 PrimaryVariables uOld(u[globalIdx]);
1717 PrimaryVariables uNew(uOld);
1718 uNew -= deltaU[globalIdx];
1720 const Scalar err = asImp_().relativeDofError(globalIdx, uOld, uNew);
1721 (*relError)[globalIdx] = err;
1722 (*normalizedRelError)[globalIdx] = err;
1723 minRelErr = std::min(err, minRelErr);
1724 maxRelErr = std::max(err, maxRelErr);
1728 const Scalar alpha = std::max(Scalar{1e-20},
1729 std::max(std::abs(maxRelErr),
1730 std::abs(minRelErr)));
1731 for (
unsigned globalIdx = 0; globalIdx < numDof; ++globalIdx) {
1732 (*normalizedRelError)[globalIdx] /= alpha;
1735 DiscBaseOutputModule::attachScalarDofData_(writer, *relError,
"relative error");
1736 DiscBaseOutputModule::attachScalarDofData_(writer, *normalizedRelError,
"normalized relative error");
1738 for (
unsigned i = 0; i < numEq; ++i) {
1739 std::ostringstream oss;
1740 oss.str(
""); oss <<
"priVar_" << asImp_().primaryVarName(i);
1741 DiscBaseOutputModule::attachScalarDofData_(writer,
1745 oss.str(
""); oss <<
"delta_" << asImp_().primaryVarName(i);
1746 DiscBaseOutputModule::attachScalarDofData_(writer,
1750 oss.str(
""); oss <<
"weight_" << asImp_().primaryVarName(i);
1751 DiscBaseOutputModule::attachScalarDofData_(writer,
1755 oss.str(
""); oss <<
"defect_" << asImp_().eqName(i);
1756 DiscBaseOutputModule::attachScalarDofData_(writer,
1761 asImp_().prepareOutputFields();
1762 asImp_().appendOutputFields(writer);
1771 const bool needFullContextUpdate =
1772 std::any_of(outputModules_.begin(), outputModules_.end(),
1773 [](
const auto& mod) { return mod->needExtensiveQuantities(); });
1774 std::for_each(outputModules_.begin(), outputModules_.end(),
1775 [](
auto& mod) { mod->allocBuffers(); });
1783 ElementContext elemCtx(simulator_);
1784 ElementIterator elemIt = threadedElemIt.beginParallel();
1785 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
1786 const auto& elem = *elemIt;
1787 if (elem.partitionType() != Dune::InteriorEntity) {
1792 if (needFullContextUpdate) {
1793 elemCtx.updateAll(elem);
1796 elemCtx.updatePrimaryStencil(elem);
1797 elemCtx.updatePrimaryIntensiveQuantities(0);
1800 std::for_each(outputModules_.begin(), outputModules_.end(),
1801 [&elemCtx](
auto& mod) { mod->processElement(elemCtx); });
1812 std::for_each(outputModules_.begin(), outputModules_.end(),
1813 [&writer](
auto& mod) { mod->commitBuffers(writer); });
1820 {
return gridView_; }
1836 auxEqModules_.push_back(auxMod);
1839 if (enableGridAdaptation_ && !std::is_same_v<DiscreteFunction, BlockVectorWrapper>) {
1840 throw std::invalid_argument(
"Problems which require auxiliary modules cannot be used in"
1841 " conjunction with dune-fem");
1845 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1859 auxEqModules_.clear();
1860 linearizer_->eraseMatrix();
1861 newtonMethod_.eraseMatrix();
1868 {
return auxEqModules_.size(); }
1874 {
return auxEqModules_[auxEqModIdx]; }
1880 {
return auxEqModules_[auxEqModIdx]; }
1886 {
return enableIntensiveQuantityCache_ || enableThermodynamicHints_; }
1888 const Timer& prePostProcessTimer()
const
1889 {
return prePostProcessTimer_; }
1891 const Timer& linearizeTimer()
const
1892 {
return linearizeTimer_; }
1894 const Timer& solveTimer()
const
1895 {
return solveTimer_; }
1897 const Timer& updateTimer()
const
1898 {
return updateTimer_; }
1900 template<
class Serializer>
1901 void serializeOp(Serializer& serializer)
1904 using Helper =
typename BaseDiscretization::template SerializeHelper<Serializer>;
1905 Helper::serializeOp(serializer, solution_);
1908 bool operator==(
const FvBaseDiscretization& rhs)
const
1910 return std::equal(this->solution_.begin(), this->solution_.end(),
1911 rhs.solution_.begin(), rhs.solution_.end(),
1912 [](
const auto& x,
const auto& y)
1919 void resizeAndResetIntensiveQuantitiesCache_()
1923 const std::size_t numDof = asImp_().numGridDof();
1924 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
1925 storageCache_[timeIdx].resize(numDof);
1926 storageCacheUpToDate_[timeIdx].resize(numDof, 0);
1932 const std::size_t numDof = asImp_().numGridDof();
1933 cachedIntensiveQuantityHistorySize_ = simulator_.problem().intensiveQuantityHistorySize();
1934 const unsigned intensiveHistorySize = cachedIntensiveQuantityHistorySize_;
1937 intensiveQuantityCache_.resize(intensiveHistorySize);
1938 intensiveQuantityCacheUpToDate_.resize(intensiveHistorySize);
1940 for(
unsigned timeIdx = 0; timeIdx < intensiveHistorySize; ++timeIdx) {
1941 intensiveQuantityCache_[timeIdx].resize(numDof);
1942 intensiveQuantityCacheUpToDate_[timeIdx].resize(numDof);
1948 template <
class Context>
1949 void supplementInitialSolution_(PrimaryVariables&,
1972 {
return localLinearizer_.localResidual(); }
1978 {
return gridView_.comm().rank() == 0; }
1980 Implementation& asImp_()
1981 {
return *
static_cast<Implementation*
>(
this); }
1983 const Implementation& asImp_()
const
1984 {
return *
static_cast<const Implementation*
>(
this); }
1988 Simulator& simulator_;
1994 ElementMapper elementMapper_;
1995 VertexMapper vertexMapper_;
1998 std::vector<BaseAuxiliaryModule<TypeTag>*> auxEqModules_;
2000 NewtonMethod newtonMethod_;
2002 Timer prePostProcessTimer_;
2003 Timer linearizeTimer_;
2008 std::vector<LocalLinearizer> localLinearizer_;
2011 std::unique_ptr<Linearizer> linearizer_;
2015 mutable std::vector<IntensiveQuantitiesVector> intensiveQuantityCache_;
2018 mutable std::vector<std::vector<unsigned char>> intensiveQuantityCacheUpToDate_;
2020 std::array<std::unique_ptr<DiscreteFunction>, historySize> solution_;
2022 std::list<std::unique_ptr<BaseOutputModule<TypeTag>>> outputModules_;
2024 Scalar gridTotalVolume_;
2025 std::vector<Scalar> dofTotalVolume_;
2026 std::vector<bool> isLocalDof_;
2028 mutable std::array<GlobalEqVector, historySize> storageCache_;
2031 mutable std::array<std::vector<unsigned char>, historySize> storageCacheUpToDate_;
2033 bool enableGridAdaptation_;
2034 bool enableIntensiveQuantityCache_;
2035 bool enableStorageCache_;
2036 bool enableThermodynamicHints_;
2038 mutable unsigned cachedIntensiveQuantityHistorySize_;
2046template<
class TypeTag>
2047class FvBaseDiscretizationNoAdapt :
public FvBaseDiscretization<TypeTag>
2049 using ParentType = FvBaseDiscretization<TypeTag>;
2056 template<
class Serializer>
2059 template<
class SolutionType>
2060 static void serializeOp(Serializer& serializer,
2069 explicit FvBaseDiscretizationNoAdapt(Simulator& simulator)
2070 : ParentType(simulator)
2072 if (this->enableGridAdaptation_) {
2073 throw std::invalid_argument(
"Grid adaptation need to use"
2074 " BaseDiscretization = FvBaseDiscretizationFemAdapt"
2075 " which currently requires the presence of the"
2076 " dune-fem module");
2078 const std::size_t numDof = this->asImp_().numGridDof();
2079 for (
unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
2080 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:2048
Definition fvbasediscretization.hh:350
The base class for the finite volume discretization schemes.
Definition fvbasediscretization.hh:298
LocalLinearizer & localLinearizer(unsigned openMpThreadId)
Definition fvbasediscretization.hh:1269
void prepareOutputFields() const
Prepare the quantities relevant for the current solution to be appended to the output writers.
Definition fvbasediscretization.hh:1769
void shiftIntensiveQuantityCache(unsigned numSlots=1)
Move the intensive quantities for a given time index to the back.
Definition fvbasediscretization.hh:785
void adaptGrid()
Called by the update() method when the grid should be refined.
Definition fvbasediscretization.hh:1442
void addAuxiliaryModule(BaseAuxiliaryModule< TypeTag > *auxMod)
Add a module for an auxiliary equation.
Definition fvbasediscretization.hh:1833
void setIntensiveQuantitiesCacheEntryValidity(unsigned globalIdx, unsigned timeIdx, bool newValue) const
Invalidate the cache for a given intensive quantities object.
Definition fvbasediscretization.hh:692
void finishInit()
Apply the initial conditions to the model.
Definition fvbasediscretization.hh:458
void prefetch(const Element &) const
Allows to improve the performance by prefetching all data which is associated with a given element.
Definition fvbasediscretization.hh:591
void updateSuccessful()
Called by the update() method if it was successful.
Definition fvbasediscretization.hh:1436
static std::string discretizationName()
Returns a string of discretization's human-readable name.
Definition fvbasediscretization.hh:1629
unsigned cachedIntensiveQuantityHistorySize() const
Get the cached intensive quantity history size.
Definition fvbasediscretization.hh:706
BaseAuxiliaryModule< TypeTag > * auxiliaryModule(unsigned auxEqModIdx)
Returns a given module for auxiliary equations.
Definition fvbasediscretization.hh:1873
bool isLocalDof(unsigned globalIdx) const
Returns if the overlap of the volume ofa degree of freedom is non-zero.
Definition fvbasediscretization.hh:1209
std::size_t numAuxiliaryModules() const
Returns the number of modules for auxiliary equations.
Definition fvbasediscretization.hh:1867
LocalResidual & localResidual_()
Reference to the local residal object.
Definition fvbasediscretization.hh:1971
void registerOutputModules_()
Register all output modules which make sense for the model.
Definition fvbasediscretization.hh:1962
bool enableGridAdaptation() const
Returns whether the grid ought to be adapted to the solution during the simulation.
Definition fvbasediscretization.hh:523
const NewtonMethod & newtonMethod() const
Returns the newton method object.
Definition fvbasediscretization.hh:605
SolutionVector & mutableSolution(unsigned timeIdx) const
Definition fvbasediscretization.hh:1237
void advanceTimeLevel()
Called by the problem if a time integration was successful, post processing of the solution is done a...
Definition fvbasediscretization.hh:1480
NewtonMethod & newtonMethod()
Returns the newton method object.
Definition fvbasediscretization.hh:599
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices.
Definition fvbasediscretization.hh:1607
const EqVector & cachedStorage(unsigned globalIdx, unsigned timeIdx) const
Retrieve an entry of the cache for the storage term.
Definition fvbasediscretization.hh:843
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:1662
void updateFailed()
Called by the update() method if it was unsuccessful.
Definition fvbasediscretization.hh:1453
void updateBegin()
Called by the update() method before it tries to apply the newton method.
Definition fvbasediscretization.hh:1429
Scalar globalResidual(GlobalEqVector &dest) const
Compute the global residual for the current solution vector.
Definition fvbasediscretization.hh:964
LocalResidual & localResidual(unsigned openMpThreadId)
Definition fvbasediscretization.hh:1281
void serializeEntity(std::ostream &outstream, const DofEntity &dof)
Write the current solution for a degree of freedom to a restart file.
Definition fvbasediscretization.hh:1535
void setEnableStorageCache(bool enableStorageCache)
Set the value of enable storage cache.
Definition fvbasediscretization.hh:827
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition fvbasediscretization.hh:1649
Scalar eqWeight(unsigned, unsigned) const
Returns the relative weight of an equation.
Definition fvbasediscretization.hh:1303
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition fvbasediscretization.hh:1291
void deserialize(Restarter &)
Deserializes the state of the model.
Definition fvbasediscretization.hh:1520
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:1079
Scalar gridTotalVolume() const
Returns the volume of the whole grid which represents the spatial domain.
Definition fvbasediscretization.hh:1216
const IntensiveQuantities * thermodynamicHint(unsigned globalIdx, unsigned timeIdx) const
Return the thermodynamic hint for a entity on the grid at given time.
Definition fvbasediscretization.hh:623
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:883
const Linearizer & linearizer() const
Returns the operator linearizer for the global jacobian of the problem.
Definition fvbasediscretization.hh:1245
void updateCachedStorage(unsigned globalIdx, unsigned timeIdx, const EqVector &value) const
Set an entry of the cache for the storage term.
Definition fvbasediscretization.hh:867
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:1680
Scalar relativeDofError(unsigned vertexIdx, const PrimaryVariables &pv1, const PrimaryVariables &pv2) const
Returns the relative error between two vectors of primary variables.
Definition fvbasediscretization.hh:1315
Scalar globalResidual(GlobalEqVector &dest, const SolutionVector &u) const
Compute the global residual for an arbitrary solution vector.
Definition fvbasediscretization.hh:949
SolutionVector & solution(unsigned timeIdx)
Definition fvbasediscretization.hh:1230
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:1560
const LocalLinearizer & localLinearizer(unsigned openMpThreadId) const
Returns the local jacobian which calculates the local stiffness matrix for an arbitrary element.
Definition fvbasediscretization.hh:1263
bool update()
Try to progress the model to the next timestep.
Definition fvbasediscretization.hh:1336
std::size_t numAuxiliaryDof() const
Returns the number of degrees of freedom (DOFs) of the auxiliary equations.
Definition fvbasediscretization.hh:1583
void clearAuxiliaryModules()
Causes the list of auxiliary equations to be cleared.
Definition fvbasediscretization.hh:1857
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices.
Definition fvbasediscretization.hh:1613
void syncOverlap()
Syncronize the values of the primary variables on the degrees of freedom that overlap with the neighb...
Definition fvbasediscretization.hh:1421
void appendOutputFields(BaseOutputWriter &writer) const
Append the quantities relevant for the current solution to an output writer.
Definition fvbasediscretization.hh:1810
std::size_t numTotalDof() const
Returns the total number of degrees of freedom (i.e., grid plux auxiliary DOFs).
Definition fvbasediscretization.hh:1594
void applyInitialSolution()
Applies the initial solution for all degrees of freedom to which the model applies.
Definition fvbasediscretization.hh:530
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:672
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition fvbasediscretization.hh:1637
void invalidateIntensiveQuantitiesCache(unsigned timeIdx) const
Invalidate the whole intensive quantity cache for time index.
Definition fvbasediscretization.hh:714
Linearizer & linearizer()
Returns the object which linearizes the global system of equations at the current solution.
Definition fvbasediscretization.hh:1252
const BaseAuxiliaryModule< TypeTag > * auxiliaryModule(unsigned auxEqModIdx) const
Returns a given module for auxiliary equations.
Definition fvbasediscretization.hh:1879
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:1024
Scalar dofTotalVolume(unsigned globalIdx) const
Returns the volume of a given control volume.
Definition fvbasediscretization.hh:1201
static void registerParameters()
Register all run-time parameters for the model.
Definition fvbasediscretization.hh:427
const SolutionVector & solution(unsigned timeIdx) const
Reference to the solution at a given history index as a block vector.
Definition fvbasediscretization.hh:1224
bool verbose_() const
Returns whether messages should be printed.
Definition fvbasediscretization.hh:1977
void shiftStorageCache(unsigned numSlots=1) const
Shift storage cache by a given number of time step slots.
Definition fvbasediscretization.hh:929
const LocalResidual & localResidual(unsigned openMpThreadId) const
Returns the object to calculate the local residual function.
Definition fvbasediscretization.hh:1275
std::size_t numGridDof() const
Returns the number of degrees of freedom (DOFs) for the computational grid.
Definition fvbasediscretization.hh:1577
const DofMapper & dofMapper() const
Mapper to convert the Dune entities of the discretization's degrees of freedoms are to indices.
Definition fvbasediscretization.hh:1601
bool storeIntensiveQuantities() const
Returns true if the cache for intensive quantities is enabled.
Definition fvbasediscretization.hh:1885
void invalidateStorageCache(unsigned timeIdx) const
Invalidate the whole storage cache for a given time index.
Definition fvbasediscretization.hh:909
void serialize(Restarter &)
Serializes the current state of the model.
Definition fvbasediscretization.hh:1506
void invalidateStorageCacheEntry(unsigned globalIdx, unsigned timeIdx) const
Invalidate the storage cache for a given DOF and time index.
Definition fvbasediscretization.hh:897
void resetLinearizer()
Resets the Jacobian matrix linearizer, so that the boundary types can be altered.
Definition fvbasediscretization.hh:1620
bool enableStorageCache() const
Returns true iff the storage term is cached.
Definition fvbasediscretization.hh:818
const GridView & gridView() const
Reference to the grid view of the spatial domain.
Definition fvbasediscretization.hh:1819
void addOutputModule(std::unique_ptr< BaseOutputModule< TypeTag > > newModule)
Add an module for writing visualization output after a timestep.
Definition fvbasediscretization.hh:1668
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:644
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:53
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.
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:81
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
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:2058
Definition fvbasediscretization.hh:272
The class which marks the border indices associated with the degrees of freedom on a process boundary...
Definition basicproperties.hh:125
The secondary variables of a boundary segment.
Definition fvbaseproperties.hh:143
Type of object for specifying boundary conditions.
Definition fvbaseproperties.hh:119
The secondary variables of a constraint degree of freedom.
Definition fvbaseproperties.hh:146
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:160
The discretization specific part of the intensive quantities.
Definition fvbaseproperties.hh:136
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:140
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:213
Specify whether the some degrees of fredom can be constraint.
Definition fvbaseproperties.hh:199
Specify if experimental features should be enabled or not.
Definition fvbaseproperties.hh:241
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:233
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:152
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:116
Vector containing all primary variables of the grid.
Definition fvbaseproperties.hh:126
The OpenMP threads manager.
Definition fvbaseproperties.hh:174
The history size required by the time discretization.
Definition fvbaseproperties.hh:225
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:181
Specify whether to use volumetric residuals or not.
Definition fvbaseproperties.hh:237
The mapper to find the global index of a vertex.
Definition fvbaseproperties.hh:207
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.