28#ifndef TPFA_LINEARIZER_HH
29#define TPFA_LINEARIZER_HH
31#include <dune/common/version.hh>
32#include <dune/common/fvector.hh>
33#include <dune/common/fmatrix.hh>
35#include <opm/common/Exceptions.hpp>
36#include <opm/common/TimingMacros.hpp>
38#include <opm/grid/utility/SparseTable.hpp>
40#include <opm/material/common/ConditionalStorage.hpp>
42#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
43#include <opm/input/eclipse/Schedule/BCProp.hpp>
60#include <unordered_map>
72template<
class TypeTag>
84template<
class TypeTag>
105 using Element =
typename GridView::template Codim<0>::Entity;
106 using ElementIterator =
typename GridView::template Codim<0>::Iterator;
108 using Vector = GlobalEqVector;
110 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
111 enum { historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>() };
112 enum { dimWorld = GridView::dimensionworld };
114 using MatrixBlock =
typename SparseMatrixAdapter::MatrixBlock;
115 using VectorBlock = Dune::FieldVector<Scalar, numEq>;
118 static constexpr bool linearizeNonLocalElements =
119 getPropValue<TypeTag, Properties::LinearizeNonLocalElements>();
120 static constexpr bool enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>();
121 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
122 static constexpr bool enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>();
123 static const bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
132 simulatorPtr_ =
nullptr;
133 separateSparseSourceTerms_ = Parameters::Get<Parameters::SeparateSparseSourceTerms>();
141 Parameters::Register<Parameters::SeparateSparseSourceTerms>
142 (
"Treat well source terms all in one go, instead of on a cell by cell basis.");
156 simulatorPtr_ = &simulator;
204 catch (
const std::exception& e) {
205 std::cout <<
"rank " << simulator_().
gridView().comm().rank()
206 <<
" caught an exception while linearizing:" << e.what()
207 <<
"\n" << std::flush;
211 std::cout <<
"rank " << simulator_().
gridView().comm().rank()
212 <<
" caught an exception while linearizing"
213 <<
"\n" << std::flush;
216 succeeded = simulator_().
gridView().comm().min(succeeded);
219 throw NumericalProblem(
"A process did not succeed in linearizing the system");
233 template <
class SubDomainType>
241 initFirstIteration_();
245 if (domain.cells.size() == model_().numTotalDof()) {
257 { jacobian_->finalize(); }
265 OPM_TIMEBLOCK(linearizeAuxilaryEquations);
269 auto& model = model_();
270 const auto& comm = simulator_().
gridView().comm();
271 for (
unsigned auxModIdx = 0; auxModIdx < model.numAuxiliaryModules(); ++auxModIdx) {
272 bool succeeded =
true;
274 model.auxiliaryModule(auxModIdx)->linearize(*jacobian_, residual_);
276 catch (
const std::exception& e) {
279 std::cout <<
"rank " << simulator_().
gridView().comm().rank()
280 <<
" caught an exception while linearizing:" << e.what()
281 <<
"\n" << std::flush;
284 succeeded = comm.min(succeeded);
287 throw NumericalProblem(
"linearization of an auxiliary equation failed");
296 {
return *jacobian_; }
299 {
return *jacobian_; }
305 {
return residual_; }
308 {
return residual_; }
311 { linearizationType_ = linearizationType; }
314 {
return linearizationType_; }
322 {
return flowsInfo_; }
330 {
return floresInfo_; }
338 {
return velocityInfo_; }
342 updateStoredTransmissibilities();
347 for (
auto& bdyInfo : boundaryInfo_) {
348 const auto [type, massrateAD] = problem_().boundaryCondition(bdyInfo.cell, bdyInfo.dir);
351 VectorBlock massrate(0.0);
352 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
353 massrate[ii] = massrateAD[ii].value();
356 const auto& exFluidState = problem_().boundaryFluidState(bdyInfo.cell, bdyInfo.dir);
357 bdyInfo.bcdata.type = type;
358 bdyInfo.bcdata.massRate = massrate;
359 bdyInfo.bcdata.exFluidState = exFluidState;
372 template <
class SubDomainType>
376 initFirstIteration_();
378 for (
int globI : domain.cells) {
379 residual_[globI] = 0.0;
380 jacobian_->clearRow(globI, 0.0);
386 {
return *simulatorPtr_; }
388 const Simulator& simulator_()
const
389 {
return *simulatorPtr_; }
392 {
return simulator_().
problem(); }
394 const Problem& problem_()
const
395 {
return simulator_().
problem(); }
398 {
return simulator_().
model(); }
400 const Model& model_()
const
401 {
return simulator_().
model(); }
403 const GridView& gridView_()
const
404 {
return problem_().gridView(); }
406 void initFirstIteration_()
412 residual_.resize(model_().numTotalDof());
423 if (!neighborInfo_.empty()) {
428 const auto& model = model_();
429 Stencil stencil(gridView_(), model_().dofMapper());
433 using NeighborSet = std::set<unsigned>;
434 std::vector<NeighborSet> sparsityPattern(model.numTotalDof());
435 const Scalar gravity = problem_().gravity()[dimWorld - 1];
436 unsigned numCells = model.numTotalDof();
437 neighborInfo_.reserve(numCells, 6 * numCells);
438 std::vector<NeighborInfo> loc_nbinfo;
439 for (
const auto& elem : elements(gridView_())) {
440 stencil.update(elem);
442 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
443 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
444 loc_nbinfo.resize(stencil.numDof() - 1);
446 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
447 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
448 sparsityPattern[myIdx].insert(neighborIdx);
450 const Scalar trans = problem_().transmissibility(myIdx, neighborIdx);
451 const auto scvfIdx = dofIdx - 1;
452 const auto& scvf = stencil.interiorFace(scvfIdx);
453 const Scalar area = scvf.area();
454 const Scalar Vin = problem_().model().dofTotalVolume(myIdx);
455 const Scalar Vex = problem_().model().dofTotalVolume(neighborIdx);
456 const Scalar zIn = problem_().dofCenterDepth(myIdx);
457 const Scalar zEx = problem_().dofCenterDepth(neighborIdx);
458 const Scalar dZg = (zIn - zEx)*gravity;
459 const Scalar thpres = problem_().thresholdPressure(myIdx, neighborIdx);
460 const auto dirId = scvf.dirId();
461 auto faceDir = dirId < 0 ? FaceDir::DirEnum::Unknown
462 : FaceDir::FromIntersectionIndex(dirId);
463 ResidualNBInfo nbinfo{trans, area, thpres, dZg, faceDir, Vin, Vex, {}, {}, {}, {}};
464 if constexpr (enableEnergy) {
465 nbinfo.inAlpha = problem_().thermalHalfTransmissibility(myIdx, neighborIdx);
466 nbinfo.outAlpha = problem_().thermalHalfTransmissibility(neighborIdx, myIdx);
468 if constexpr (enableDiffusion) {
469 nbinfo.diffusivity = problem_().diffusivity(myIdx, neighborIdx);
471 if constexpr (enableDispersion) {
472 nbinfo.dispersivity = problem_().dispersivity(myIdx, neighborIdx);
474 loc_nbinfo[dofIdx - 1] = NeighborInfo{neighborIdx, nbinfo,
nullptr};
477 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
478 if (problem_().nonTrivialBoundaryConditions()) {
479 for (
unsigned bfIndex = 0; bfIndex < stencil.numBoundaryFaces(); ++bfIndex) {
480 const auto& bf = stencil.boundaryFace(bfIndex);
481 const int dir_id = bf.dirId();
486 const auto [type, massrateAD] = problem_().boundaryCondition(myIdx, dir_id);
488 VectorBlock massrate(0.0);
489 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
490 massrate[ii] = massrateAD[ii].value();
492 const auto& exFluidState = problem_().boundaryFluidState(myIdx, dir_id);
493 BoundaryConditionData bcdata{type,
495 exFluidState.pvtRegionIndex(),
498 bf.integrationPos()[dimWorld - 1],
500 boundaryInfo_.push_back({myIdx, dir_id, bfIndex, bcdata});
508 const std::size_t numAuxMod = model.numAuxiliaryModules();
509 for (
unsigned auxModIdx = 0; auxModIdx < numAuxMod; ++auxModIdx) {
510 model.auxiliaryModule(auxModIdx)->addNeighbors(sparsityPattern);
514 jacobian_ = std::make_unique<SparseMatrixAdapter>(simulator_());
515 diagMatAddress_.resize(numCells);
517 jacobian_->reserve(sparsityPattern);
518 for (
unsigned globI = 0; globI < numCells; globI++) {
519 const auto& nbInfos = neighborInfo_[globI];
520 diagMatAddress_[globI] = jacobian_->blockAddress(globI, globI);
521 for (
auto& nbInfo : nbInfos) {
522 nbInfo.matBlockAddress = jacobian_->blockAddress(nbInfo.neighbor, globI);
527 fullDomain_.cells.resize(numCells);
528 std::iota(fullDomain_.cells.begin(), fullDomain_.cells.end(), 0);
542 OPM_TIMEBLOCK(createFlows);
546 const bool anyFlows = simulator_().
problem().eclWriter().outputModule().getFlows().anyFlows() ||
547 simulator_().
problem().eclWriter().outputModule().getFlows().hasBlockFlows();
548 const bool anyFlores = simulator_().
problem().eclWriter().outputModule().getFlows().anyFlores();
549 const bool dispersionActive = simulator_().
vanguard().eclState().getSimulationConfig().rock_config().dispersion();
550 if (((!anyFlows || !flowsInfo_.empty()) && (!anyFlores || !floresInfo_.empty())) && (!dispersionActive && !enableBioeffects)) {
553 const auto& model = model_();
554 const auto& nncOutput = simulator_().
problem().eclWriter().getOutputNnc();
555 Stencil stencil(gridView_(), model_().dofMapper());
556 const unsigned numCells = model.numTotalDof();
557 std::unordered_multimap<int, std::pair<int, int>> nncIndices;
558 std::vector<FlowInfo> loc_flinfo;
559 std::vector<VelocityInfo> loc_vlinfo;
560 unsigned int nncId = 0;
561 VectorBlock flow(0.0);
564 for (
unsigned nncIdx = 0; nncIdx < nncOutput.size(); ++nncIdx) {
565 const int ci1 = nncOutput[nncIdx].cell1;
566 const int ci2 = nncOutput[nncIdx].cell2;
567 nncIndices.emplace(ci1, std::make_pair(ci2, nncIdx));
571 flowsInfo_.reserve(numCells, 6 * numCells);
574 floresInfo_.reserve(numCells, 6 * numCells);
576 if (dispersionActive || enableBioeffects) {
577 velocityInfo_.reserve(numCells, 6 * numCells);
580 for (
const auto& elem : elements(gridView_())) {
581 stencil.update(elem);
582 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
583 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
584 const int numFaces = stencil.numBoundaryFaces() + stencil.numInteriorFaces();
585 loc_flinfo.resize(numFaces);
586 loc_vlinfo.resize(stencil.numDof() - 1);
588 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
589 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
591 const auto scvfIdx = dofIdx - 1;
592 const auto& scvf = stencil.interiorFace(scvfIdx);
593 int faceId = scvf.dirId();
594 const int cartMyIdx = simulator_().
vanguard().cartesianIndex(myIdx);
595 const int cartNeighborIdx = simulator_().
vanguard().cartesianIndex(neighborIdx);
596 const auto& range = nncIndices.equal_range(cartMyIdx);
597 for (
auto it = range.first; it != range.second; ++it) {
598 if (it->second.first == cartNeighborIdx){
602 nncId = it->second.second;
605 loc_flinfo[dofIdx - 1] = FlowInfo{faceId, flow, nncId};
606 loc_vlinfo[dofIdx - 1] = VelocityInfo{flow};
610 for (
unsigned bdfIdx = 0; bdfIdx < stencil.numBoundaryFaces(); ++bdfIdx) {
611 const auto& scvf = stencil.boundaryFace(bdfIdx);
612 const int faceId = scvf.dirId();
613 loc_flinfo[stencil.numInteriorFaces() + bdfIdx] = FlowInfo{faceId, flow, nncId};
617 flowsInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
620 floresInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
622 if (dispersionActive || enableBioeffects) {
623 velocityInfo_.appendRow(loc_vlinfo.begin(), loc_vlinfo.end());
632 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
633 res[eqIdx] = resid[eqIdx].value();
636 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
637 for (
unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
642 bMat[eqIdx][pvIdx] = resid[eqIdx].derivative(pvIdx);
649 OPM_TIMEBLOCK(updateFlows);
650 const bool enableFlows = simulator_().
problem().eclWriter().outputModule().getFlows().hasFlows() ||
651 simulator_().
problem().eclWriter().outputModule().getFlows().hasBlockFlows();
652 const bool enableFlores = simulator_().
problem().eclWriter().outputModule().getFlows().hasFlores();
653 if (!enableFlows && !enableFlores) {
656 const unsigned int numCells = model_().numTotalDof();
658#pragma omp parallel for
660 for (
unsigned globI = 0; globI < numCells; ++globI) {
661 OPM_TIMEBLOCK_LOCAL(linearizationForEachCell, Subsystem::Assembly);
662 const auto& nbInfos = neighborInfo_[globI];
663 ADVectorBlock adres(0.0);
664 ADVectorBlock darcyFlux(0.0);
665 const IntensiveQuantities& intQuantsIn = model_().intensiveQuantities(globI, 0);
668 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell, Subsystem::Assembly);
670 for (
const auto& nbInfo : nbInfos) {
671 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace, Subsystem::Assembly);
672 const unsigned globJ = nbInfo.neighbor;
673 assert(globJ != globI);
676 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, 0);
677 LocalResidual::computeFlux(adres, darcyFlux, globI, globJ, intQuantsIn,
678 intQuantsEx, nbInfo.res_nbinfo, problem_().moduleParams());
679 adres *= nbInfo.res_nbinfo.faceArea;
681 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
682 flowsInfo_[globI][loc].flow[eqIdx] = adres[eqIdx].value();
686 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
687 floresInfo_[globI][loc].flow[eqIdx] = darcyFlux[eqIdx].value();
696 for (
const auto& bdyInfo : boundaryInfo_) {
701 ADVectorBlock adres(0.0);
702 const unsigned globI = bdyInfo.cell;
703 const auto& nbInfos = neighborInfo_[globI];
704 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, 0);
705 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
706 adres *= bdyInfo.bcdata.faceArea;
707 const unsigned bfIndex = bdyInfo.bfIndex;
709 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
710 flowsInfo_[globI][nbInfos.size() + bfIndex].flow[eqIdx] = adres[eqIdx].value();
718 template <
class SubDomainType>
719 void linearize_(
const SubDomainType& domain)
724 if (!problem_().recycleFirstIterationStorage()) {
725 if (!model_().storeIntensiveQuantities() && !model_().enableStorageCache()) {
726 OPM_THROW(std::runtime_error,
"Must have cached either IQs or storage when we cannot recycle.");
735 const bool dispersionActive = simulator_().
vanguard().eclState().getSimulationConfig().rock_config().dispersion();
736 const unsigned int numCells = domain.cells.size();
737 const bool on_full_domain = (numCells == model_().numTotalDof());
743#pragma omp parallel for
745 for (
unsigned ii = 0; ii < numCells; ++ii) {
746 OPM_TIMEBLOCK_LOCAL(linearizationForEachCell, Subsystem::Assembly);
747 const unsigned globI = domain.cells[ii];
748 const auto& nbInfos = neighborInfo_[globI];
749 VectorBlock res(0.0);
750 MatrixBlock bMat(0.0);
751 ADVectorBlock adres(0.0);
752 ADVectorBlock darcyFlux(0.0);
753 const IntensiveQuantities& intQuantsIn = model_().intensiveQuantities(globI, 0);
757 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell, Subsystem::Assembly);
759 for (
const auto& nbInfo : nbInfos) {
760 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace, Subsystem::Assembly);
761 const unsigned globJ = nbInfo.neighbor;
762 assert(globJ != globI);
767 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, 0);
768 LocalResidual::computeFlux(adres, darcyFlux, globI, globJ, intQuantsIn, intQuantsEx,
769 nbInfo.res_nbinfo, problem_().moduleParams());
770 adres *= nbInfo.res_nbinfo.faceArea;
771 if (dispersionActive || enableBioeffects) {
772 for (
unsigned phaseIdx = 0; phaseIdx < numEq; ++phaseIdx) {
773 velocityInfo_[globI][loc].velocity[phaseIdx] =
774 darcyFlux[phaseIdx].value() / nbInfo.res_nbinfo.faceArea;
778 residual_[globI] += res;
780 *diagMatAddress_[globI] += bMat;
783 *nbInfo.matBlockAddress += bMat;
789 const double volume = model_().dofTotalVolume(globI);
790 const Scalar storefac = volume / dt;
793 OPM_TIMEBLOCK_LOCAL(computeStorage, Subsystem::Assembly);
794 LocalResidual::computeStorage(adres, intQuantsIn);
798 if (model_().enableStorageCache()) {
802 model_().updateCachedStorage(globI, 0, res);
803 if (model_().newtonMethod().numIterations() == 0) {
805 if (problem_().recycleFirstIterationStorage()) {
808 if (on_full_domain) {
815 model_().updateCachedStorage(globI, 1, res);
819 Dune::FieldVector<Scalar, numEq> tmp;
820 const IntensiveQuantities intQuantOld = model_().intensiveQuantities(globI, 1);
821 LocalResidual::computeStorage(tmp, intQuantOld);
822 model_().updateCachedStorage(globI, 1, tmp);
825 res -= model_().cachedStorage(globI, 1);
828 OPM_TIMEBLOCK_LOCAL(computeStorage0, Subsystem::Assembly);
829 Dune::FieldVector<Scalar, numEq> tmp;
830 const IntensiveQuantities intQuantOld = model_().intensiveQuantities(globI, 1);
831 LocalResidual::computeStorage(tmp, intQuantOld);
837 residual_[globI] += res;
839 *diagMatAddress_[globI] += bMat;
846 if (separateSparseSourceTerms_) {
847 LocalResidual::computeSourceDense(adres, problem_(), intQuantsIn, globI, 0);
850 LocalResidual::computeSource(adres, problem_(), intQuantsIn, globI, 0);
854 residual_[globI] += res;
856 *diagMatAddress_[globI] += bMat;
860 if (separateSparseSourceTerms_) {
861 problem_().wellModel().addReservoirSourceTerms(residual_, diagMatAddress_);
865 for (
const auto& bdyInfo : boundaryInfo_) {
870 VectorBlock res(0.0);
871 MatrixBlock bMat(0.0);
872 ADVectorBlock adres(0.0);
873 const unsigned globI = bdyInfo.cell;
874 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, 0);
875 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
876 adres *= bdyInfo.bcdata.faceArea;
878 residual_[globI] += res;
880 *diagMatAddress_[globI] += bMat;
884 void updateStoredTransmissibilities()
886 if (neighborInfo_.empty()) {
890 initFirstIteration_();
893 const unsigned numCells = model_().numTotalDof();
895#pragma omp parallel for
897 for (
unsigned globI = 0; globI < numCells; globI++) {
898 auto nbInfos = neighborInfo_[globI];
899 for (
auto& nbInfo : nbInfos) {
900 const unsigned globJ = nbInfo.neighbor;
901 nbInfo.res_nbinfo.trans = problem_().transmissibility(globI, globJ);
906 Simulator* simulatorPtr_{};
909 std::unique_ptr<SparseMatrixAdapter> jacobian_{};
912 GlobalEqVector residual_;
914 LinearizationType linearizationType_{};
916 using ResidualNBInfo =
typename LocalResidual::ResidualNBInfo;
919 unsigned int neighbor;
920 ResidualNBInfo res_nbinfo;
921 MatrixBlock* matBlockAddress;
923 SparseTable<NeighborInfo> neighborInfo_{};
924 std::vector<MatrixBlock*> diagMatAddress_{};
932 SparseTable<FlowInfo> flowsInfo_;
933 SparseTable<FlowInfo> floresInfo_;
937 VectorBlock velocity;
939 SparseTable<VelocityInfo> velocityInfo_;
941 using ScalarFluidState =
typename IntensiveQuantities::ScalarFluidState;
942 struct BoundaryConditionData
945 VectorBlock massRate;
946 unsigned pvtRegionIdx;
947 unsigned boundaryFaceIndex;
950 ScalarFluidState exFluidState;
957 unsigned int bfIndex;
958 BoundaryConditionData bcdata;
960 std::vector<BoundaryInfo> boundaryInfo_;
962 bool separateSparseSourceTerms_ =
false;
966 std::vector<int> cells;
967 std::vector<bool> interior;
969 FullDomain fullDomain_;
Declares the properties required by the black oil model.
The base class for the element-centered finite-volume discretization scheme.
Definition: ecfvdiscretization.hh:160
Definition: matrixblock.hh:229
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:84
Scalar timeStepSize() const
Returns the time step length so that we don't miss the beginning of the next episode or cross the en...
Definition: simulator.hh:413
Vanguard & vanguard()
Return a reference to the grid manager of simulation.
Definition: simulator.hh:234
Problem & problem()
Return the object which specifies the pysical setup of the simulation.
Definition: simulator.hh:265
const GridView & gridView() const
Return the grid view for which the simulation is done.
Definition: simulator.hh:246
Model & model()
Return the physical model used in the simulation.
Definition: simulator.hh:252
The common code for the linearizers of non-linear systems of equations.
Definition: tpfalinearizer.hh:86
const auto & getFloresInfo() const
Return constant reference to the floresInfo.
Definition: tpfalinearizer.hh:329
const LinearizationType & getLinearizationType() const
Definition: tpfalinearizer.hh:313
void updateBoundaryConditionData()
Definition: tpfalinearizer.hh:345
void linearize()
Linearize the full system of non-linear equations.
Definition: tpfalinearizer.hh:181
SparseMatrixAdapter & jacobian()
Definition: tpfalinearizer.hh:298
const auto & getFlowsInfo() const
Return constant reference to the flowsInfo.
Definition: tpfalinearizer.hh:321
std::map< unsigned, Constraints > constraintsMap() const
Returns the map of constraint degrees of freedom.
Definition: tpfalinearizer.hh:369
TpfaLinearizer()
Definition: tpfalinearizer.hh:130
void linearizeDomain(const SubDomainType &domain)
Linearize the part of the non-linear system of equations that is associated with a part of the spatia...
Definition: tpfalinearizer.hh:234
void finalize()
Definition: tpfalinearizer.hh:256
void init(Simulator &simulator)
Initialize the linearizer.
Definition: tpfalinearizer.hh:154
void updateFlowsInfo()
Definition: tpfalinearizer.hh:647
void setResAndJacobi(VectorBlock &res, MatrixBlock &bMat, const ADVectorBlock &resid) const
Definition: tpfalinearizer.hh:630
static void registerParameters()
Register all run-time parameters for the Jacobian linearizer.
Definition: tpfalinearizer.hh:139
void linearizeDomain()
Linearize the part of the non-linear system of equations that is associated with the spatial domain.
Definition: tpfalinearizer.hh:197
const SparseMatrixAdapter & jacobian() const
Return constant reference to global Jacobian matrix backend.
Definition: tpfalinearizer.hh:295
void setLinearizationType(LinearizationType linearizationType)
Definition: tpfalinearizer.hh:310
GlobalEqVector & residual()
Definition: tpfalinearizer.hh:307
void eraseMatrix()
Causes the Jacobian matrix to be recreated from scratch before the next iteration.
Definition: tpfalinearizer.hh:167
const auto & getVelocityInfo() const
Return constant reference to the velocityInfo.
Definition: tpfalinearizer.hh:337
void updateDiscretizationParameters()
Definition: tpfalinearizer.hh:340
void resetSystem_(const SubDomainType &domain)
Definition: tpfalinearizer.hh:373
void linearizeAuxiliaryEquations()
Linearize the part of the non-linear system of equations that is associated with the spatial domain.
Definition: tpfalinearizer.hh:263
const GlobalEqVector & residual() const
Return constant reference to global residual vector.
Definition: tpfalinearizer.hh:304
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
@ NONE
Definition: DeferredLogger.hpp:46
Definition: blackoilnewtonmethodparams.hpp:31
HYPRE_IJMatrix createMatrix(HYPRE_Int N, HYPRE_Int dof_offset, const CommType &comm)
Create Hypre matrix.
Definition: HypreSetup.hpp:170
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
Definition: linearizationtype.hh:34
Definition: tpfalinearizer.hh:65
static constexpr bool value
Definition: tpfalinearizer.hh:65