24#ifndef OPM_BLACKOILMODEL_IMPL_HEADER_INCLUDED
25#define OPM_BLACKOILMODEL_IMPL_HEADER_INCLUDED
27#ifndef OPM_BLACKOILMODEL_HEADER_INCLUDED
29#include <opm/simulators/flow/BlackoilModel.hpp>
32#include <dune/common/timer.hh>
34#include <opm/common/ErrorMacros.hpp>
35#include <opm/common/OpmLog/OpmLog.hpp>
37#include <opm/simulators/flow/countGlobalCells.hpp>
45#include <fmt/format.h>
49template <
class TypeTag>
55 : simulator_(simulator)
56 , grid_(simulator_.vanguard().grid())
60 , current_relaxation_(1.0)
61 , dx_old_(simulator_.model().numGridDof())
62 , conv_monitor_(param_.monitor_params_)
66 convergence_reports_.reserve(300);
70 OpmLog::info(
"Using Non-Linear Domain Decomposition solver (nldd).");
72 nlddSolver_ = std::make_unique<BlackoilModelNldd<TypeTag>>(*this);
75 OpmLog::info(
"Using Newton nonlinear solver.");
78 OPM_THROW(std::runtime_error,
"Unknown nonlinear solver option: " +
83template <
class TypeTag>
93 if (grid_.comm().size() > 1 && grid_.comm().max(lastStepFailed) != grid_.comm().min(lastStepFailed)) {
94 OPM_THROW(std::runtime_error,
"Misalignment of the parallel simulation run in prepareStep " +
95 "- the previous step succeeded on some ranks but failed on others.");
98 simulator_.model().updateFailed();
101 simulator_.model().advanceTimeLevel();
110 simulator_.model().newtonMethod().setIterationIndex(0);
112 simulator_.problem().beginTimeStep();
114 unsigned numDof = simulator_.model().numGridDof();
115 wasSwitched_.resize(numDof);
116 std::fill(wasSwitched_.begin(), wasSwitched_.end(),
false);
118 if (param_.update_equations_scaling_) {
119 OpmLog::error(
"Equation scaling not supported");
123 if (hasNlddSolver()) {
124 nlddSolver_->prepareStep();
127 report.pre_post_time +=
perfTimer.stop();
131 if (FluidSystem::phaseIsActive(
phaseIdx)) {
132 const unsigned sIdx = FluidSystem::solventComponentIndex(
phaseIdx);
133 return FluidSystem::canonicalToActiveCompIdx(
sIdx);
138 const auto& schedule = simulator_.vanguard().schedule();
139 auto&
rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
140 rst_conv.init(simulator_.vanguard().globalNumCells(),
142 {getIdx(FluidSystem::oilPhaseIdx),
143 getIdx(FluidSystem::gasPhaseIdx),
144 getIdx(FluidSystem::waterPhaseIdx),
152template <
class TypeTag>
166 report.total_linearizations = 1;
170 report += assembleReservoir(timer,
iteration);
171 report.assemble_time +=
perfTimer.stop();
174 report.assemble_time +=
perfTimer.stop();
175 failureReport_ += report;
187 ConvergenceReport::Severity severity =
convrep.severityOfWorstFailure();
188 convergence_reports_.back().report.push_back(std::move(
convrep));
191 if (severity == ConvergenceReport::Severity::NotANumber) {
192 failureReport_ += report;
194 }
else if (severity == ConvergenceReport::Severity::TooLarge) {
195 failureReport_ += report;
197 }
else if (severity == ConvergenceReport::Severity::ConvergenceMonitorFailure) {
198 failureReport_ += report;
201 "Total penalty count exceeded cut-off-limit of {}",
202 param_.monitor_params_.cutoff_
210template <
class TypeTag>
211template <
class NonlinearSolverType>
221 residual_norms_history_.clear();
222 conv_monitor_.reset();
223 current_relaxation_ = 1.0;
226 convergence_reports_.back().report.reserve(11);
230 if (this->param_.nonlinear_solver_ !=
"nldd")
242 auto&
rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
243 rst_conv.update(simulator_.model().linearizer().residual());
248template <
class TypeTag>
249template <
class NonlinearSolverType>
260 this->initialLinearization(report,
262 this->param_.newton_min_iter_,
263 this->param_.newton_max_iter_,
267 if (!report.converged) {
270 report.total_newton_iterations = 1;
273 unsigned nc = simulator_.model().numGridDof();
277 linear_solve_setup_time_ = 0.0;
282 wellModel().linearize(simulator().model().linearizer().jacobian(),
283 simulator().model().linearizer().residual());
286 solveJacobianSystem(x);
288 report.linear_solve_setup_time += linear_solve_setup_time_;
289 report.linear_solve_time +=
perfTimer.stop();
290 report.total_linear_iterations += linearIterationsLastSolve();
293 report.linear_solve_setup_time += linear_solve_setup_time_;
294 report.linear_solve_time +=
perfTimer.stop();
295 report.total_linear_iterations += linearIterationsLastSolve();
297 failureReport_ += report;
307 wellModel().postSolve(x);
309 if (param_.use_update_stabilization_) {
314 residual_norms_history_.size() - 1,
319 current_relaxation_ = std::max(current_relaxation_,
321 if (terminalOutputEnabled()) {
322 std::string
msg =
" Oscillating behavior detected: Relaxation set to "
323 + std::to_string(current_relaxation_);
341template <
class TypeTag>
349 simulator_.problem().endTimeStep();
350 report.pre_post_time +=
perfTimer.stop();
354template <
class TypeTag>
361 simulator_.model().newtonMethod().setIterationIndex(
iterationIdx);
362 simulator_.problem().beginIteration();
363 simulator_.model().linearizer().linearizeDomain();
364 simulator_.problem().endIteration();
365 return wellModel().lastReport();
368template <
class TypeTag>
369typename BlackoilModel<TypeTag>::Scalar
376 const auto&
elemMapper = simulator_.model().elementMapper();
377 const auto& gridView = simulator_.gridView();
387 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) &&
388 FluidSystem::numActivePhases() > 1 &&
389 priVarsNew.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
395 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
396 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
397 priVarsNew.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
399 assert(Indices::compositionSwitchIdx >= 0);
404 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
421 if (FluidSystem::numActivePhases() > 1) {
422 if (
priVarsOld.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw) {
428 if (
priVarsOld.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
430 assert(Indices::compositionSwitchIdx >= 0 );
436 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
455template <
class TypeTag>
460 auto& jacobian = simulator_.model().linearizer().jacobian().istlMatrix();
461 auto& residual = simulator_.model().linearizer().residual();
462 auto&
linSolver = simulator_.model().newtonMethod().linearSolver();
466 if (terminal_output_) {
467 OpmLog::debug(
"\nRunning speed test for comparing available linear solvers.");
476 for (
int solver = 0; solver <
numSolvers; ++solver) {
488 if (terminal_output_) {
489 OpmLog::debug(fmt::format(
"Solver time {}: {}", solver,
times[solver]));
507 linear_solve_setup_time_ =
perfTimer.stop();
517template <
class TypeTag>
523 auto& newtonMethod = simulator_.model().newtonMethod();
524 SolutionVector& solution = simulator_.model().solution(0);
526 newtonMethod.update_(solution,
536 simulator_.model().invalidateAndUpdateIntensiveQuantities(0);
540template <
class TypeTag>
541std::tuple<typename BlackoilModel<TypeTag>::Scalar,
542 typename BlackoilModel<TypeTag>::Scalar>
547 std::vector< Scalar >&
R_sum,
549 std::vector< Scalar >&
B_avg)
556 if (comm.size() > 1) {
600template <
class TypeTag>
601std::pair<typename BlackoilModel<TypeTag>::Scalar,
602 typename BlackoilModel<TypeTag>::Scalar>
606 std::vector<Scalar>&
B_avg,
612 const auto& model = simulator_.model();
613 const auto& problem = simulator_.problem();
615 const auto& residual = simulator_.model().linearizer().residual();
617 ElementContext elemCtx(simulator_);
618 const auto& gridView = simulator().gridView();
620 OPM_BEGIN_PARALLEL_TRY_CATCH();
621 for (
const auto&
elem :
elements(gridView, Dune::Partitions::interior)) {
622 elemCtx.updatePrimaryStencil(
elem);
623 elemCtx.updatePrimaryIntensiveQuantities(0);
625 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
626 const auto& intQuants = elemCtx.intensiveQuantities(0, 0);
627 const auto& fs = intQuants.fluidState();
641 OPM_END_PARALLEL_TRY_CATCH(
"BlackoilModel::localConvergenceData() failed: ", grid_.comm());
645 for (
int i = 0; i <
bSize; ++i) {
646 B_avg[i] /= Scalar(global_nc_);
652template <
class TypeTag>
653std::pair<std::vector<double>, std::vector<int>>
662 constexpr auto numPvGroups = std::vector<double>::size_type{3};
664 auto cnvPvSplit = std::pair<std::vector<double>, std::vector<int>> {
665 std::piecewise_construct,
673 std::inner_product(residual.begin(), residual.end(),
674 B_avg.begin(), Scalar{0},
675 [](
const Scalar
m,
const auto& x)
678 return std::max(m, abs(x));
679 }, std::multiplies<>{});
684 const auto& model = this->simulator().model();
685 const auto& problem = this->simulator().problem();
686 const auto& residual = model.linearizer().residual();
687 const auto& gridView = this->simulator().gridView();
691 ElementContext elemCtx(this->simulator());
693 OPM_BEGIN_PARALLEL_TRY_CATCH();
694 for (
const auto&
elem :
elements(gridView, Dune::Partitions::interior)) {
700 elemCtx.updatePrimaryStencil(
elem);
702 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
708 const auto ix = (
maxCnv > this->param_.tolerance_cnv_)
709 + (
maxCnv > this->param_.tolerance_cnv_relaxed_);
715 OPM_END_PARALLEL_TRY_CATCH(
"BlackoilModel::characteriseCnvPvSplit() failed: ",
724template <
class TypeTag>
729 this->param_.tolerance_cnv_ =
tuning.TRGCNV;
730 this->param_.tolerance_cnv_relaxed_ =
tuning.XXXCNV;
731 this->param_.tolerance_mb_ =
tuning.TRGMBE;
732 this->param_.tolerance_mb_relaxed_ =
tuning.XXXMBE;
733 this->param_.newton_max_iter_ =
tuning.NEWTMX;
734 this->param_.newton_min_iter_ =
tuning.NEWTMN;
737template <
class TypeTag>
739BlackoilModel<TypeTag>::
740getReservoirConvergence(
const double reportTime,
744 std::vector<Scalar>&
B_avg,
748 using Vector = std::vector<Scalar>;
750 ConvergenceReport report{reportTime};
763 this->convergenceReduction(this->grid_.comm(),
768 report.setCnvPoreVolSplit(this->characteriseCnvPvSplit(
B_avg,
dt),
782 (this->param_.min_strict_mb_iter_ >= 0 &&
783 iteration >= this->param_.min_strict_mb_iter_);
795 this->param_.min_strict_cnv_iter_ >= 0 &&
796 iteration >= this->param_.min_strict_cnv_iter_;
804 const auto& cnvPvSplit = report.cnvPvSplit().first;
808 return static_cast<Scalar
>(cnvPvSplit[1] + cnvPvSplit[2]) <
809 this->param_.relaxed_max_pv_fraction_ *
eligible;
817 this->terminal_output_)
819 std::string message =
820 "Number of newton iterations reached its maximum "
821 "try to continue with relaxed tolerances:";
824 message += fmt::format(
" MB: {:.1e}", param_.tolerance_mb_relaxed_);
828 message += fmt::format(
" CNV: {:.1e}", param_.tolerance_cnv_relaxed_);
831 OpmLog::debug(message);
837 const auto tol_eb =
use_relaxed_mb ? param_.tolerance_energy_balance_relaxed_ : param_.tolerance_energy_balance_;
849 using CR = ConvergenceReport;
851 const Scalar
res[2] = {
855 const CR::ReservoirFailure::Type
types[2] = {
856 CR::ReservoirFailure::Type::MassBalance,
857 CR::ReservoirFailure::Type::Cnv,
861 if (has_energy_ &&
compIdx == contiEnergyEqIdx) {
866 for (
int ii : {0, 1}) {
867 if (std::isnan(
res[
ii])) {
868 report.setReservoirFailed({
types[
ii], CR::Severity::NotANumber,
compIdx});
869 if (this->terminal_output_) {
870 OpmLog::debug(
"NaN residual for " + this->compNames_.name(
compIdx) +
" equation.");
873 else if (
res[
ii] > maxResidualAllowed()) {
874 report.setReservoirFailed({
types[
ii], CR::Severity::TooLarge,
compIdx});
875 if (this->terminal_output_) {
876 OpmLog::debug(
"Too large residual for " + this->compNames_.name(
compIdx) +
" equation.");
879 else if (
res[
ii] < 0.0) {
880 report.setReservoirFailed({
types[
ii], CR::Severity::Normal,
compIdx});
881 if (this->terminal_output_) {
882 OpmLog::debug(
"Negative residual for " + this->compNames_.name(
compIdx) +
" equation.");
886 report.setReservoirFailed({
types[
ii], CR::Severity::Normal,
compIdx});
897 if (this->terminal_output_) {
900 std::string
msg =
"Iter";
916 std::ostringstream
ss;
917 const std::streamsize
oprec =
ss.precision(3);
918 const std::ios::fmtflags
oflags =
ss.setf(std::ios::scientific);
932 OpmLog::debug(
ss.str());
938template <
class TypeTag>
947 auto&
rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
956 const auto& residual = simulator_.model().linearizer().residual();
957 const auto& gridView = this->simulator().gridView();
959 ElementContext elemCtx(this->simulator());
960 std::vector<int>
convNewt(residual.size(), 0);
961 OPM_BEGIN_PARALLEL_TRY_CATCH();
964 for (
const auto&
elem :
elements(gridView, Dune::Partitions::interior)) {
965 elemCtx.updatePrimaryStencil(
elem);
967 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
968 const auto pvValue = simulator_.problem().referencePorosity(
cell_idx, 0) *
969 simulator_.model().dofTotalVolume(
cell_idx);
980 OPM_END_PARALLEL_TRY_CATCH(
"BlackoilModel::convergencePerCell() failed: ",
985template <
class TypeTag>
995 std::vector<Scalar>
B_avg(numEq, 0.0);
1001 report += wellModel().getWellConvergence(
B_avg,
1002 report.converged());
1005 conv_monitor_.checkPenaltyCard(report,
iteration);
1010template <
class TypeTag>
1011std::vector<std::vector<typename BlackoilModel<TypeTag>::Scalar> >
1018 std::vector<std::vector<Scalar> >
regionValues(0, std::vector<Scalar>(0,0.0));
1022template <
class TypeTag>
1027 if (!hasNlddSolver()) {
1028 OPM_THROW(std::runtime_error,
"Cannot get local reports from a model without NLDD solver");
1030 return nlddSolver_->localAccumulatedReports();
1033template <
class TypeTag>
1034const std::vector<SimulatorReport>&
1039 OPM_THROW(std::runtime_error,
"Cannot get domain reports from a model without NLDD solver");
1040 return nlddSolver_->domainAccumulatedReports();
1043template <
class TypeTag>
1048 if (hasNlddSolver()) {
1049 nlddSolver_->writeNonlinearIterationsPerCell(
odir);
1053template <
class TypeTag>
1058 if (hasNlddSolver()) {
1059 nlddSolver_->writePartitions(
odir);
1063 const auto& elementMapper = this->simulator().model().elementMapper();
1064 const auto&
cartMapper = this->simulator().vanguard().cartesianIndexMapper();
1066 const auto& grid = this->simulator().vanguard().grid();
1067 const auto& comm = grid.comm();
1068 const auto nDigit = 1 +
static_cast<int>(std::floor(std::log10(comm.size())));
1070 std::ofstream
pfile {
odir / fmt::format(
"{1:0>{0}}",
nDigit, comm.rank())};
1073 pfile << comm.rank() <<
' '
1074 <<
cartMapper.cartesianIndex(elementMapper.index(cell)) <<
' '
1075 << comm.rank() <<
'\n';
1079template <
class TypeTag>
1080template<
class Flu
idState,
class Res
idual>
1082BlackoilModel<TypeTag>::
1083getMaxCoeff(
const unsigned cell_idx,
1084 const IntensiveQuantities& intQuants,
1085 const FluidState& fs,
1088 std::vector<Scalar>&
B_avg,
1089 std::vector<Scalar>&
R_sum,
1095 if (!FluidSystem::phaseIsActive(
phaseIdx)) {
1099 const unsigned sIdx = FluidSystem::solventComponentIndex(
phaseIdx);
1100 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(
sIdx);
1113 if constexpr (has_solvent_) {
1114 B_avg[contiSolventEqIdx] +=
1115 1.0 / intQuants.solventInverseFormationVolumeFactor().value();
1117 R_sum[contiSolventEqIdx] +=
R2;
1121 if constexpr (has_extbo_) {
1122 B_avg[contiZfracEqIdx] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
1124 R_sum[ contiZfracEqIdx ] +=
R2;
1128 if constexpr (has_polymer_) {
1129 B_avg[contiPolymerEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1131 R_sum[contiPolymerEqIdx] +=
R2;
1135 if constexpr (has_foam_) {
1136 B_avg[ contiFoamEqIdx ] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
1142 if constexpr (has_brine_) {
1143 B_avg[ contiBrineEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1150 if constexpr (has_polymermw_) {
1151 static_assert(has_polymer_);
1153 B_avg[contiPolymerMWEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1158 R_sum[contiPolymerMWEqIdx] +=
R2;
1163 if constexpr (has_energy_) {
1164 B_avg[contiEnergyEqIdx] += 1.0 / (4.182e1);
1166 R_sum[contiEnergyEqIdx] +=
R2;
1171 if constexpr (has_bioeffects_) {
1172 B_avg[contiMicrobialEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1174 R_sum[contiMicrobialEqIdx] +=
R1;
1177 B_avg[contiBiofilmEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1179 R_sum[contiBiofilmEqIdx] +=
R2;
1182 if constexpr (has_micp_) {
1183 B_avg[contiOxygenEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1185 R_sum[contiOxygenEqIdx] +=
R3;
1188 B_avg[contiUreaEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1193 B_avg[contiCalciteEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1195 R_sum[contiCalciteEqIdx] +=
R5;
A model implementation for three-phase black oil.
Definition BlackoilModel.hpp:61
BlackoilModel(Simulator &simulator, const ModelParameters ¶m, BlackoilWellModel< TypeTag > &well_model, const bool terminal_output)
Construct the model.
Definition BlackoilModel_impl.hpp:51
std::vector< std::vector< Scalar > > computeFluidInPlace(const T &, const std::vector< int > &fipnum) const
Wrapper required due to not following generic API.
Definition BlackoilModel.hpp:249
SimulatorReportSingle nonlinearIteration(const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Called once per nonlinear iteration.
Definition BlackoilModel_impl.hpp:214
void writeNonlinearIterationsPerCell(const std::filesystem::path &odir) const
Write the number of nonlinear iterations per cell to a file in ResInsight compatible format.
Definition BlackoilModel_impl.hpp:1046
SimulatorReportSingle assembleReservoir(const SimulatorTimerInterface &, const int iterationIdx)
Assemble the residual and Jacobian of the nonlinear system.
Definition BlackoilModel_impl.hpp:357
ConvergenceReport getConvergence(const SimulatorTimerInterface &timer, const int iteration, const int maxIter, std::vector< Scalar > &residual_norms)
Compute convergence based on total mass balance (tol_mb) and maximum residual mass balance (tol_cnv).
Definition BlackoilModel_impl.hpp:988
const std::vector< SimulatorReport > & domainAccumulatedReports() const
return the statistics of local solves accumulated for each domain on this rank
Definition BlackoilModel_impl.hpp:1036
std::pair< Scalar, Scalar > localConvergenceData(std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< Scalar > &B_avg, std::vector< int > &maxCoeffCell)
Get reservoir quantities on this process needed for convergence calculations.
Definition BlackoilModel_impl.hpp:604
const SimulatorReport & localAccumulatedReports() const
return the statistics of local solves accumulated for this rank
Definition BlackoilModel_impl.hpp:1025
SimulatorReportSingle afterStep(const SimulatorTimerInterface &)
Called once after each time step.
Definition BlackoilModel_impl.hpp:344
void updateSolution(const BVector &dx)
Apply an update to the primary variables.
Definition BlackoilModel_impl.hpp:520
long int global_nc_
The number of cells of the global grid.
Definition BlackoilModel.hpp:341
std::unique_ptr< BlackoilModelNldd< TypeTag > > nlddSolver_
Non-linear DD solver.
Definition BlackoilModel.hpp:350
std::pair< std::vector< double >, std::vector< int > > characteriseCnvPvSplit(const std::vector< Scalar > &B_avg, const double dt)
Compute pore-volume/cell count split among "converged", "relaxed converged", "unconverged" cells base...
Definition BlackoilModel_impl.hpp:655
void solveJacobianSystem(BVector &x)
Solve the Jacobian system Jx = r where J is the Jacobian and r is the residual.
Definition BlackoilModel_impl.hpp:458
SimulatorReportSingle prepareStep(const SimulatorTimerInterface &timer)
Called once before each time step.
Definition BlackoilModel_impl.hpp:86
void convergencePerCell(const std::vector< Scalar > &B_avg, const double dt, const double tol_cnv, const double tol_cnv_energy, const int iteration)
Compute the number of Newtons required by each cell in order to satisfy the solution change convergen...
Definition BlackoilModel_impl.hpp:941
Class for handling the blackoil well model.
Definition BlackoilWellModel.hpp:102
Represents the convergence status of the whole simulator, to make it possible to query and store the ...
Definition ConvergenceReport.hpp:38
Interface class for SimulatorTimer objects, to be improved.
Definition SimulatorTimerInterface.hpp:34
virtual int reportStepNum() const
Current report step number. This might differ from currentStepNum in case of sub stepping.
Definition SimulatorTimerInterface.hpp:109
virtual bool lastStepFailed() const =0
Return true if last time step failed.
virtual double currentStepLength() const =0
Current step length.
virtual double simulationTimeElapsed() const =0
Time elapsed since the start of the simulation until the beginning of the current time step [s].
virtual int currentStepNum() const =0
Current step number.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
Solver parameters for the BlackoilModel.
Definition BlackoilModelParameters.hpp:180
std::string nonlinear_solver_
Nonlinear solver type: newton or nldd.
Definition BlackoilModelParameters.hpp:336
Definition AquiferGridUtils.hpp:35
A struct for returning timing data from a simulator to its caller.
Definition SimulatorReport.hpp:34
Definition SimulatorReport.hpp:122