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>
104 using Element =
typename GridView::template Codim<0>::Entity;
105 using ElementIterator =
typename GridView::template Codim<0>::Iterator;
107 using Vector = GlobalEqVector;
109 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
110 enum { historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>() };
111 enum { dimWorld = GridView::dimensionworld };
113 using MatrixBlock =
typename SparseMatrixAdapter::MatrixBlock;
114 using VectorBlock = Dune::FieldVector<Scalar, numEq>;
117 static constexpr bool linearizeNonLocalElements =
118 getPropValue<TypeTag, Properties::LinearizeNonLocalElements>();
119 static constexpr bool enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>();
120 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
121 static constexpr bool enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>();
122 static constexpr bool enableMICP = getPropValue<TypeTag, Properties::EnableMICP>();
131 simulatorPtr_ =
nullptr;
132 separateSparseSourceTerms_ = Parameters::Get<Parameters::SeparateSparseSourceTerms>();
140 Parameters::Register<Parameters::SeparateSparseSourceTerms>
141 (
"Treat well source terms all in one go, instead of on a cell by cell basis.");
155 simulatorPtr_ = &simulator;
203 catch (
const std::exception& e) {
204 std::cout <<
"rank " << simulator_().
gridView().comm().rank()
205 <<
" caught an exception while linearizing:" << e.what()
206 <<
"\n" << std::flush;
210 std::cout <<
"rank " << simulator_().
gridView().comm().rank()
211 <<
" caught an exception while linearizing"
212 <<
"\n" << std::flush;
215 succeeded = simulator_().
gridView().comm().min(succeeded);
218 throw NumericalProblem(
"A process did not succeed in linearizing the system");
232 template <
class SubDomainType>
240 initFirstIteration_();
244 if (domain.cells.size() == model_().numTotalDof()) {
256 { jacobian_->finalize(); }
264 OPM_TIMEBLOCK(linearizeAuxilaryEquations);
268 auto& model = model_();
269 const auto& comm = simulator_().
gridView().comm();
270 for (
unsigned auxModIdx = 0; auxModIdx < model.numAuxiliaryModules(); ++auxModIdx) {
271 bool succeeded =
true;
273 model.auxiliaryModule(auxModIdx)->linearize(*jacobian_, residual_);
275 catch (
const std::exception& e) {
278 std::cout <<
"rank " << simulator_().
gridView().comm().rank()
279 <<
" caught an exception while linearizing:" << e.what()
280 <<
"\n" << std::flush;
283 succeeded = comm.min(succeeded);
286 throw NumericalProblem(
"linearization of an auxiliary equation failed");
295 {
return *jacobian_; }
298 {
return *jacobian_; }
304 {
return residual_; }
307 {
return residual_; }
310 { linearizationType_ = linearizationType; }
313 {
return linearizationType_; }
321 {
return flowsInfo_; }
329 {
return floresInfo_; }
337 {
return velocityInfo_; }
341 updateStoredTransmissibilities();
346 for (
auto& bdyInfo : boundaryInfo_) {
347 const auto [type, massrateAD] = problem_().boundaryCondition(bdyInfo.cell, bdyInfo.dir);
350 VectorBlock massrate(0.0);
351 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
352 massrate[ii] = massrateAD[ii].value();
355 const auto& exFluidState = problem_().boundaryFluidState(bdyInfo.cell, bdyInfo.dir);
356 bdyInfo.bcdata.type = type;
357 bdyInfo.bcdata.massRate = massrate;
358 bdyInfo.bcdata.exFluidState = exFluidState;
371 template <
class SubDomainType>
375 initFirstIteration_();
377 for (
int globI : domain.cells) {
378 residual_[globI] = 0.0;
379 jacobian_->clearRow(globI, 0.0);
385 {
return *simulatorPtr_; }
387 const Simulator& simulator_()
const
388 {
return *simulatorPtr_; }
391 {
return simulator_().
problem(); }
393 const Problem& problem_()
const
394 {
return simulator_().
problem(); }
397 {
return simulator_().
model(); }
399 const Model& model_()
const
400 {
return simulator_().
model(); }
402 const GridView& gridView_()
const
403 {
return problem_().gridView(); }
405 void initFirstIteration_()
411 residual_.resize(model_().numTotalDof());
421 OPM_TIMEBLOCK(createMatrix);
422 if (!neighborInfo_.empty()) {
427 const auto& model = model_();
428 Stencil stencil(gridView_(), model_().dofMapper());
432 using NeighborSet = std::set<unsigned>;
433 std::vector<NeighborSet> sparsityPattern(model.numTotalDof());
434 const Scalar gravity = problem_().gravity()[dimWorld - 1];
435 unsigned numCells = model.numTotalDof();
436 neighborInfo_.reserve(numCells, 6 * numCells);
437 std::vector<NeighborInfo> loc_nbinfo;
438 for (
const auto& elem : elements(gridView_())) {
439 stencil.update(elem);
441 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
442 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
443 loc_nbinfo.resize(stencil.numDof() - 1);
445 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
446 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
447 sparsityPattern[myIdx].insert(neighborIdx);
449 const Scalar trans = problem_().transmissibility(myIdx, neighborIdx);
450 const auto scvfIdx = dofIdx - 1;
451 const auto& scvf = stencil.interiorFace(scvfIdx);
452 const Scalar area = scvf.area();
453 const Scalar Vin = problem_().model().dofTotalVolume(myIdx);
454 const Scalar Vex = problem_().model().dofTotalVolume(neighborIdx);
455 const Scalar zIn = problem_().dofCenterDepth(myIdx);
456 const Scalar zEx = problem_().dofCenterDepth(neighborIdx);
457 const Scalar dZg = (zIn - zEx)*gravity;
458 const Scalar thpres = problem_().thresholdPressure(myIdx, neighborIdx);
459 const auto dirId = scvf.dirId();
460 auto faceDir = dirId < 0 ? FaceDir::DirEnum::Unknown
461 : FaceDir::FromIntersectionIndex(dirId);
462 ResidualNBInfo nbinfo{trans, area, thpres, dZg, faceDir, Vin, Vex, {}, {}, {}, {}};
463 if constexpr (enableEnergy) {
464 nbinfo.inAlpha = problem_().thermalHalfTransmissibility(myIdx, neighborIdx);
465 nbinfo.outAlpha = problem_().thermalHalfTransmissibility(neighborIdx, myIdx);
467 if constexpr (enableDiffusion) {
468 nbinfo.diffusivity = problem_().diffusivity(myIdx, neighborIdx);
470 if constexpr (enableDispersion) {
471 nbinfo.dispersivity = problem_().dispersivity(myIdx, neighborIdx);
473 loc_nbinfo[dofIdx - 1] = NeighborInfo{neighborIdx, nbinfo,
nullptr};
476 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
477 if (problem_().nonTrivialBoundaryConditions()) {
478 for (
unsigned bfIndex = 0; bfIndex < stencil.numBoundaryFaces(); ++bfIndex) {
479 const auto& bf = stencil.boundaryFace(bfIndex);
480 const int dir_id = bf.dirId();
485 const auto [type, massrateAD] = problem_().boundaryCondition(myIdx, dir_id);
487 VectorBlock massrate(0.0);
488 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
489 massrate[ii] = massrateAD[ii].value();
491 const auto& exFluidState = problem_().boundaryFluidState(myIdx, dir_id);
492 BoundaryConditionData bcdata{type,
494 exFluidState.pvtRegionIndex(),
497 bf.integrationPos()[dimWorld - 1],
499 boundaryInfo_.push_back({myIdx, dir_id, bfIndex, bcdata});
507 const std::size_t numAuxMod = model.numAuxiliaryModules();
508 for (
unsigned auxModIdx = 0; auxModIdx < numAuxMod; ++auxModIdx) {
509 model.auxiliaryModule(auxModIdx)->addNeighbors(sparsityPattern);
513 jacobian_ = std::make_unique<SparseMatrixAdapter>(simulator_());
514 diagMatAddress_.resize(numCells);
516 jacobian_->reserve(sparsityPattern);
517 for (
unsigned globI = 0; globI < numCells; globI++) {
518 const auto& nbInfos = neighborInfo_[globI];
519 diagMatAddress_[globI] = jacobian_->blockAddress(globI, globI);
520 for (
auto& nbInfo : nbInfos) {
521 nbInfo.matBlockAddress = jacobian_->blockAddress(nbInfo.neighbor, globI);
526 fullDomain_.cells.resize(numCells);
527 std::iota(fullDomain_.cells.begin(), fullDomain_.cells.end(), 0);
541 OPM_TIMEBLOCK(createFlows);
545 const bool anyFlows = simulator_().
problem().eclWriter().outputModule().getFlows().anyFlows() ||
546 simulator_().
problem().eclWriter().outputModule().getFlows().hasBlockFlows();
547 const bool anyFlores = simulator_().
problem().eclWriter().outputModule().getFlows().anyFlores();
548 const bool dispersionActive = simulator_().
vanguard().eclState().getSimulationConfig().rock_config().dispersion();
549 if (((!anyFlows || !flowsInfo_.empty()) && (!anyFlores || !floresInfo_.empty())) && (!dispersionActive && !enableMICP)) {
552 const auto& model = model_();
553 const auto& nncOutput = simulator_().
problem().eclWriter().getOutputNnc();
554 Stencil stencil(gridView_(), model_().dofMapper());
555 const unsigned numCells = model.numTotalDof();
556 std::unordered_multimap<int, std::pair<int, int>> nncIndices;
557 std::vector<FlowInfo> loc_flinfo;
558 std::vector<VelocityInfo> loc_vlinfo;
559 unsigned int nncId = 0;
560 VectorBlock flow(0.0);
563 for (
unsigned nncIdx = 0; nncIdx < nncOutput.size(); ++nncIdx) {
564 const int ci1 = nncOutput[nncIdx].cell1;
565 const int ci2 = nncOutput[nncIdx].cell2;
566 nncIndices.emplace(ci1, std::make_pair(ci2, nncIdx));
570 flowsInfo_.reserve(numCells, 6 * numCells);
573 floresInfo_.reserve(numCells, 6 * numCells);
575 if (dispersionActive || enableMICP) {
576 velocityInfo_.reserve(numCells, 6 * numCells);
579 for (
const auto& elem : elements(gridView_())) {
580 stencil.update(elem);
581 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
582 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
583 const int numFaces = stencil.numBoundaryFaces() + stencil.numInteriorFaces();
584 loc_flinfo.resize(numFaces);
585 loc_vlinfo.resize(stencil.numDof() - 1);
587 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
588 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
590 const auto scvfIdx = dofIdx - 1;
591 const auto& scvf = stencil.interiorFace(scvfIdx);
592 int faceId = scvf.dirId();
593 const int cartMyIdx = simulator_().
vanguard().cartesianIndex(myIdx);
594 const int cartNeighborIdx = simulator_().
vanguard().cartesianIndex(neighborIdx);
595 const auto& range = nncIndices.equal_range(cartMyIdx);
596 for (
auto it = range.first; it != range.second; ++it) {
597 if (it->second.first == cartNeighborIdx){
601 nncId = it->second.second;
604 loc_flinfo[dofIdx - 1] = FlowInfo{faceId, flow, nncId};
605 loc_vlinfo[dofIdx - 1] = VelocityInfo{flow};
609 for (
unsigned bdfIdx = 0; bdfIdx < stencil.numBoundaryFaces(); ++bdfIdx) {
610 const auto& scvf = stencil.boundaryFace(bdfIdx);
611 const int faceId = scvf.dirId();
612 loc_flinfo[stencil.numInteriorFaces() + bdfIdx] = FlowInfo{faceId, flow, nncId};
616 flowsInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
619 floresInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
621 if (dispersionActive || enableMICP) {
622 velocityInfo_.appendRow(loc_vlinfo.begin(), loc_vlinfo.end());
631 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
632 res[eqIdx] = resid[eqIdx].value();
635 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
636 for (
unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
641 bMat[eqIdx][pvIdx] = resid[eqIdx].derivative(pvIdx);
648 OPM_TIMEBLOCK(updateFlows);
649 const bool enableFlows = simulator_().
problem().eclWriter().outputModule().getFlows().hasFlows() ||
650 simulator_().
problem().eclWriter().outputModule().getFlows().hasBlockFlows();
651 const bool enableFlores = simulator_().
problem().eclWriter().outputModule().getFlows().hasFlores();
652 if (!enableFlows && !enableFlores) {
655 const unsigned int numCells = model_().numTotalDof();
657#pragma omp parallel for
659 for (
unsigned globI = 0; globI < numCells; ++globI) {
660 OPM_TIMEBLOCK_LOCAL(linearizationForEachCell);
661 const auto& nbInfos = neighborInfo_[globI];
662 ADVectorBlock adres(0.0);
663 ADVectorBlock darcyFlux(0.0);
664 const IntensiveQuantities& intQuantsIn = model_().intensiveQuantities(globI, 0);
667 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell);
669 for (
const auto& nbInfo : nbInfos) {
670 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace);
671 const unsigned globJ = nbInfo.neighbor;
672 assert(globJ != globI);
675 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, 0);
676 LocalResidual::computeFlux(adres,darcyFlux, globI, globJ, intQuantsIn,
677 intQuantsEx, nbInfo.res_nbinfo, problem_().moduleParams());
678 adres *= nbInfo.res_nbinfo.faceArea;
680 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
681 flowsInfo_[globI][loc].flow[eqIdx] = adres[eqIdx].value();
685 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
686 floresInfo_[globI][loc].flow[eqIdx] = darcyFlux[eqIdx].value();
695 for (
const auto& bdyInfo : boundaryInfo_) {
700 ADVectorBlock adres(0.0);
701 const unsigned globI = bdyInfo.cell;
702 const auto& nbInfos = neighborInfo_[globI];
703 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, 0);
704 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
705 adres *= bdyInfo.bcdata.faceArea;
706 const unsigned bfIndex = bdyInfo.bfIndex;
708 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
709 flowsInfo_[globI][nbInfos.size() + bfIndex].flow[eqIdx] = adres[eqIdx].value();
717 template <
class SubDomainType>
718 void linearize_(
const SubDomainType& domain)
723 if (!problem_().recycleFirstIterationStorage()) {
724 if (!model_().storeIntensiveQuantities() && !model_().enableStorageCache()) {
725 OPM_THROW(std::runtime_error,
"Must have cached either IQs or storage when we cannot recycle.");
734 const bool dispersionActive = simulator_().
vanguard().eclState().getSimulationConfig().rock_config().dispersion();
735 const unsigned int numCells = domain.cells.size();
736 const bool on_full_domain = (numCells == model_().numTotalDof());
739#pragma omp parallel for
741 for (
unsigned ii = 0; ii < numCells; ++ii) {
742 OPM_TIMEBLOCK_LOCAL(linearizationForEachCell);
743 const unsigned globI = domain.cells[ii];
744 const auto& nbInfos = neighborInfo_[globI];
745 VectorBlock res(0.0);
746 MatrixBlock bMat(0.0);
747 ADVectorBlock adres(0.0);
748 ADVectorBlock darcyFlux(0.0);
749 const IntensiveQuantities& intQuantsIn = model_().intensiveQuantities(globI, 0);
753 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell);
755 for (
const auto& nbInfo : nbInfos) {
756 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace);
757 const unsigned globJ = nbInfo.neighbor;
758 assert(globJ != globI);
763 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, 0);
764 LocalResidual::computeFlux(adres,darcyFlux, globI, globJ, intQuantsIn, intQuantsEx,
765 nbInfo.res_nbinfo, problem_().moduleParams());
766 adres *= nbInfo.res_nbinfo.faceArea;
767 if (dispersionActive || enableMICP) {
768 for (
unsigned phaseIdx = 0; phaseIdx < numEq; ++phaseIdx) {
769 velocityInfo_[globI][loc].velocity[phaseIdx] =
770 darcyFlux[phaseIdx].value() / nbInfo.res_nbinfo.faceArea;
774 residual_[globI] += res;
776 *diagMatAddress_[globI] += bMat;
779 *nbInfo.matBlockAddress += bMat;
786 const double volume = model_().dofTotalVolume(globI);
787 const Scalar storefac = volume / dt;
790 OPM_TIMEBLOCK_LOCAL(computeStorage);
791 LocalResidual::computeStorage(adres, intQuantsIn);
795 if (model_().enableStorageCache()) {
799 model_().updateCachedStorage(globI, 0, res);
800 if (model_().newtonMethod().numIterations() == 0) {
802 if (problem_().recycleFirstIterationStorage()) {
805 if (on_full_domain) {
812 model_().updateCachedStorage(globI, 1, res);
816 Dune::FieldVector<Scalar, numEq> tmp;
817 const IntensiveQuantities intQuantOld = model_().intensiveQuantities(globI, 1);
818 LocalResidual::computeStorage(tmp, intQuantOld);
819 model_().updateCachedStorage(globI, 1, tmp);
822 res -= model_().cachedStorage(globI, 1);
825 OPM_TIMEBLOCK_LOCAL(computeStorage0);
826 Dune::FieldVector<Scalar, numEq> tmp;
827 const IntensiveQuantities intQuantOld = model_().intensiveQuantities(globI, 1);
828 LocalResidual::computeStorage(tmp, intQuantOld);
834 residual_[globI] += res;
836 *diagMatAddress_[globI] += bMat;
843 if (separateSparseSourceTerms_) {
844 LocalResidual::computeSourceDense(adres, problem_(), intQuantsIn, globI, 0);
847 LocalResidual::computeSource(adres, problem_(), intQuantsIn, globI, 0);
851 residual_[globI] += res;
853 *diagMatAddress_[globI] += bMat;
857 if (separateSparseSourceTerms_) {
858 problem_().wellModel().addReservoirSourceTerms(residual_, diagMatAddress_);
862 for (
const auto& bdyInfo : boundaryInfo_) {
867 VectorBlock res(0.0);
868 MatrixBlock bMat(0.0);
869 ADVectorBlock adres(0.0);
870 const unsigned globI = bdyInfo.cell;
871 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, 0);
872 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
873 adres *= bdyInfo.bcdata.faceArea;
875 residual_[globI] += res;
877 *diagMatAddress_[globI] += bMat;
881 void updateStoredTransmissibilities()
883 if (neighborInfo_.empty()) {
887 initFirstIteration_();
890 const unsigned numCells = model_().numTotalDof();
892#pragma omp parallel for
894 for (
unsigned globI = 0; globI < numCells; globI++) {
895 auto nbInfos = neighborInfo_[globI];
896 for (
auto& nbInfo : nbInfos) {
897 const unsigned globJ = nbInfo.neighbor;
898 nbInfo.res_nbinfo.trans = problem_().transmissibility(globI, globJ);
903 Simulator* simulatorPtr_{};
906 std::unique_ptr<SparseMatrixAdapter> jacobian_{};
909 GlobalEqVector residual_;
911 LinearizationType linearizationType_{};
913 using ResidualNBInfo =
typename LocalResidual::ResidualNBInfo;
916 unsigned int neighbor;
917 ResidualNBInfo res_nbinfo;
918 MatrixBlock* matBlockAddress;
920 SparseTable<NeighborInfo> neighborInfo_{};
921 std::vector<MatrixBlock*> diagMatAddress_{};
929 SparseTable<FlowInfo> flowsInfo_;
930 SparseTable<FlowInfo> floresInfo_;
934 VectorBlock velocity;
936 SparseTable<VelocityInfo> velocityInfo_;
938 using ScalarFluidState =
typename IntensiveQuantities::ScalarFluidState;
939 struct BoundaryConditionData
942 VectorBlock massRate;
943 unsigned pvtRegionIdx;
944 unsigned boundaryFaceIndex;
947 ScalarFluidState exFluidState;
954 unsigned int bfIndex;
955 BoundaryConditionData bcdata;
957 std::vector<BoundaryInfo> boundaryInfo_;
959 bool separateSparseSourceTerms_ =
false;
963 std::vector<int> cells;
964 std::vector<bool> interior;
966 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:227
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:328
const LinearizationType & getLinearizationType() const
Definition: tpfalinearizer.hh:312
void updateBoundaryConditionData()
Definition: tpfalinearizer.hh:344
void linearize()
Linearize the full system of non-linear equations.
Definition: tpfalinearizer.hh:180
SparseMatrixAdapter & jacobian()
Definition: tpfalinearizer.hh:297
const auto & getFlowsInfo() const
Return constant reference to the flowsInfo.
Definition: tpfalinearizer.hh:320
std::map< unsigned, Constraints > constraintsMap() const
Returns the map of constraint degrees of freedom.
Definition: tpfalinearizer.hh:368
TpfaLinearizer()
Definition: tpfalinearizer.hh:129
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:233
void finalize()
Definition: tpfalinearizer.hh:255
void init(Simulator &simulator)
Initialize the linearizer.
Definition: tpfalinearizer.hh:153
void updateFlowsInfo()
Definition: tpfalinearizer.hh:646
void setResAndJacobi(VectorBlock &res, MatrixBlock &bMat, const ADVectorBlock &resid) const
Definition: tpfalinearizer.hh:629
static void registerParameters()
Register all run-time parameters for the Jacobian linearizer.
Definition: tpfalinearizer.hh:138
void linearizeDomain()
Linearize the part of the non-linear system of equations that is associated with the spatial domain.
Definition: tpfalinearizer.hh:196
const SparseMatrixAdapter & jacobian() const
Return constant reference to global Jacobian matrix backend.
Definition: tpfalinearizer.hh:294
void setLinearizationType(LinearizationType linearizationType)
Definition: tpfalinearizer.hh:309
GlobalEqVector & residual()
Definition: tpfalinearizer.hh:306
void eraseMatrix()
Causes the Jacobian matrix to be recreated from scratch before the next iteration.
Definition: tpfalinearizer.hh:166
const auto & getVelocityInfo() const
Return constant reference to the velocityInfo.
Definition: tpfalinearizer.hh:336
void updateDiscretizationParameters()
Definition: tpfalinearizer.hh:339
void resetSystem_(const SubDomainType &domain)
Definition: tpfalinearizer.hh:372
void linearizeAuxiliaryEquations()
Linearize the part of the non-linear system of equations that is associated with the spatial domain.
Definition: tpfalinearizer.hh:262
const GlobalEqVector & residual() const
Return constant reference to global residual vector.
Definition: tpfalinearizer.hh:303
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
Definition: blackoilboundaryratevector.hh:39
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