24#ifndef OPM_BLACKOILMODEL_IMPL_HEADER_INCLUDED
25#define OPM_BLACKOILMODEL_IMPL_HEADER_INCLUDED
27#ifndef OPM_BLACKOILMODEL_HEADER_INCLUDED
32#include <dune/common/timer.hh>
34#include <opm/common/ErrorMacros.hpp>
35#include <opm/common/OpmLog/OpmLog.hpp>
55#include <fmt/format.h>
58 template <
typename TypeTag>
65 case F::STRICT:
return "Strict";
66 case F::RELAXED:
return "Relaxed";
67 case F::TUNINGDP:
return "TuningDP";
76template <
class TypeTag>
81 const bool terminal_output)
82 : simulator_(simulator)
83 , grid_(simulator_.vanguard().grid())
85 , well_model_ (well_model)
86 , terminal_output_ (terminal_output)
87 , current_relaxation_(1.0)
88 , dx_old_(simulator_.model().numGridDof())
89 , conv_monitor_(param_.monitor_params_)
96 if (terminal_output) {
97 OpmLog::info(
"Using Non-Linear Domain Decomposition solver (nldd).");
99 nlddSolver_ = std::make_unique<BlackoilModelNldd<TypeTag>>(*this);
101 if (terminal_output) {
102 OpmLog::info(
"Using Newton nonlinear solver.");
105 OPM_THROW(std::runtime_error,
"Unknown nonlinear solver option: " +
110template <
class TypeTag>
117 Dune::Timer perfTimer;
121 if (grid_.comm().size() > 1 && grid_.comm().max(lastStepFailed) != grid_.comm().min(lastStepFailed)) {
122 OPM_THROW(std::runtime_error,
"Misalignment of the parallel simulation run in prepareStep " +
123 "- the previous step succeeded on some ranks but failed on others.");
125 if (lastStepFailed) {
126 simulator_.model().updateFailed();
129 simulator_.model().advanceTimeLevel();
139 simulator_.problem().resetIterationForNewTimestep();
141 simulator_.problem().beginTimeStep();
143 unsigned numDof = simulator_.model().numGridDof();
144 wasSwitched_.resize(numDof);
145 std::fill(wasSwitched_.begin(), wasSwitched_.end(),
false);
147 if (param_.update_equations_scaling_) {
148 OpmLog::error(
"Equation scaling not supported");
152 if (hasNlddSolver()) {
153 nlddSolver_->prepareStep();
158 auto getIdx = [](
unsigned phaseIdx) ->
int
160 if (FluidSystem::phaseIsActive(phaseIdx)) {
161 const unsigned sIdx = FluidSystem::solventComponentIndex(phaseIdx);
162 return FluidSystem::canonicalToActiveCompIdx(sIdx);
167 const auto& schedule = simulator_.vanguard().schedule();
168 auto& rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
169 rst_conv.init(simulator_.vanguard().globalNumCells(),
171 {getIdx(FluidSystem::oilPhaseIdx),
172 getIdx(FluidSystem::gasPhaseIdx),
173 getIdx(FluidSystem::waterPhaseIdx),
181template <
class TypeTag>
191 Dune::Timer perfTimer;
198 report += assembleReservoir(timer);
202 simulator_.problem().markTimestepInitialized();
206 failureReport_ += report;
211 std::vector<Scalar> residual_norms;
216 auto convrep = getConvergence(timer, maxIter, residual_norms);
217 report.
converged = convrep.converged() &&
218 simulator_.problem().iterationContext().iteration() >= minIter;
220 convergence_reports_.back().report.push_back(std::move(convrep));
224 failureReport_ += report;
225 OPM_THROW_PROBLEM(NumericalProblem,
"NaN residual found!");
227 failureReport_ += report;
228 OPM_THROW_NOLOG(NumericalProblem,
"Too large residual found!");
230 failureReport_ += report;
231 OPM_THROW_PROBLEM(ConvergenceMonitorFailure,
233 "Total penalty count exceeded cut-off-limit of {}",
234 param_.monitor_params_.cutoff_
239 residual_norms_history_.push_back(residual_norms);
242template <
class TypeTag>
243template <
class NonlinearSolverType>
247 NonlinearSolverType& nonlinear_solver)
252 if (simulator_.problem().iterationContext().needsTimestepInit()) {
253 residual_norms_history_.clear();
254 conv_monitor_.reset();
255 current_relaxation_ = 1.0;
258 convergence_reports_.back().report.reserve(11);
262 if (this->param_.nonlinear_solver_ !=
"nldd")
264 result = this->nonlinearIterationNewton(timer, nonlinear_solver);
267 result = this->nlddSolver_->nonlinearIterationNldd(timer, nonlinear_solver);
270 auto& rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
271 rst_conv.update(simulator_.model().linearizer().residual());
273 simulator_.problem().advanceIteration();
277template <
class TypeTag>
278template <
class NonlinearSolverType>
282 NonlinearSolverType& nonlinear_solver)
288 Dune::Timer perfTimer;
290 this->initialLinearization(report,
291 this->param_.newton_min_iter_,
292 this->param_.newton_max_iter_,
302 unsigned nc = simulator_.model().numGridDof();
306 linear_solve_setup_time_ = 0.0;
311 wellModel().linearize(simulator().model().linearizer().jacobian(),
312 simulator().model().linearizer().residual());
315 solveJacobianSystem(x);
326 failureReport_ += report;
336 wellModel().postSolve(x);
338 if (param_.use_update_stabilization_) {
340 bool isOscillate =
false;
341 bool isStagnate =
false;
342 nonlinear_solver.detectOscillations(residual_norms_history_,
343 residual_norms_history_.size() - 1,
347 current_relaxation_ -= nonlinear_solver.relaxIncrement();
348 current_relaxation_ = std::max(current_relaxation_,
349 nonlinear_solver.relaxMax());
350 if (terminalOutputEnabled()) {
351 std::string msg =
" Oscillating behavior detected: Relaxation set to "
356 nonlinear_solver.stabilizeNonlinearUpdate(x, dx_old_, current_relaxation_);
370template <
class TypeTag>
376 simulator_.problem().beginIteration();
377 simulator_.model().linearizer().linearizeDomain();
378 simulator_.problem().endIteration();
379 return wellModel().lastReport();
382template <
class TypeTag>
390 const auto& elemMapper = simulator_.model().elementMapper();
391 const auto& gridView = simulator_.gridView();
392 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
393 unsigned globalElemIdx = elemMapper.index(elem);
394 const auto& priVarsNew = simulator_.model().solution(0)[globalElemIdx];
397 pressureNew = priVarsNew[Indices::pressureSwitchIdx];
399 Scalar saturationsNew[FluidSystem::numPhases] = { 0.0 };
400 Scalar oilSaturationNew = 1.0;
401 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) &&
402 FluidSystem::numActivePhases() > 1 &&
403 priVarsNew.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
405 saturationsNew[FluidSystem::waterPhaseIdx] = priVarsNew[Indices::waterSwitchIdx];
406 oilSaturationNew -= saturationsNew[FluidSystem::waterPhaseIdx];
409 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
410 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
411 priVarsNew.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
413 assert(Indices::compositionSwitchIdx >= 0);
414 saturationsNew[FluidSystem::gasPhaseIdx] = priVarsNew[Indices::compositionSwitchIdx];
415 oilSaturationNew -= saturationsNew[FluidSystem::gasPhaseIdx];
418 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
419 saturationsNew[FluidSystem::oilPhaseIdx] = oilSaturationNew;
422 const auto& priVarsOld = simulator_.model().solution(1)[globalElemIdx];
425 pressureOld = priVarsOld[Indices::pressureSwitchIdx];
427 Scalar saturationsOld[FluidSystem::numPhases] = { 0.0 };
428 Scalar oilSaturationOld = 1.0;
431 Scalar tmp = pressureNew - pressureOld;
432 resultDelta += tmp*tmp;
433 resultDenom += pressureNew*pressureNew;
435 if (FluidSystem::numActivePhases() > 1) {
436 if (priVarsOld.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw) {
437 saturationsOld[FluidSystem::waterPhaseIdx] =
438 priVarsOld[Indices::waterSwitchIdx];
439 oilSaturationOld -= saturationsOld[FluidSystem::waterPhaseIdx];
442 if (priVarsOld.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
444 assert(Indices::compositionSwitchIdx >= 0 );
445 saturationsOld[FluidSystem::gasPhaseIdx] =
446 priVarsOld[Indices::compositionSwitchIdx];
447 oilSaturationOld -= saturationsOld[FluidSystem::gasPhaseIdx];
450 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
451 saturationsOld[FluidSystem::oilPhaseIdx] = oilSaturationOld;
453 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++ phaseIdx) {
454 Scalar tmpSat = saturationsNew[phaseIdx] - saturationsOld[phaseIdx];
455 resultDelta += tmpSat*tmpSat;
456 resultDenom += saturationsNew[phaseIdx]*saturationsNew[phaseIdx];
457 assert(std::isfinite(resultDelta));
458 assert(std::isfinite(resultDenom));
463 resultDelta = gridView.comm().sum(resultDelta);
464 resultDenom = gridView.comm().sum(resultDenom);
466 return resultDenom > 0.0 ? resultDelta / resultDenom : 0.0;
469template <
class TypeTag>
474 auto& jacobian = simulator_.model().linearizer().jacobian().istlMatrix();
475 auto& residual = simulator_.model().linearizer().residual();
476 auto& linSolver = simulator_.model().newtonMethod().linearSolver();
478 const int numSolvers = linSolver.numAvailableSolvers();
479 if (numSolvers > 1 && (linSolver.getSolveCount() % 100 == 0)) {
480 if (terminal_output_) {
481 OpmLog::debug(
"\nRunning speed test for comparing available linear solvers.");
484 Dune::Timer perfTimer;
485 std::vector<double> times(numSolvers);
486 std::vector<double> setupTimes(numSolvers);
489 std::vector<BVector> x_trial(numSolvers, x);
490 for (
int solver = 0; solver < numSolvers; ++solver) {
492 linSolver.setActiveSolver(solver);
494 linSolver.prepare(jacobian, residual);
495 setupTimes[solver] = perfTimer.stop();
497 linSolver.setResidual(residual);
499 linSolver.solve(x_trial[solver]);
500 times[solver] = perfTimer.stop();
502 if (terminal_output_) {
503 OpmLog::debug(fmt::format(fmt::runtime(
"Solver time {}: {}"), solver, times[solver]));
507 int fastest_solver = std::ranges::min_element(times) - times.begin();
509 grid_.comm().broadcast(&fastest_solver, 1, 0);
510 linear_solve_setup_time_ = setupTimes[fastest_solver];
511 x = x_trial[fastest_solver];
512 linSolver.setActiveSolver(fastest_solver);
518 Dune::Timer perfTimer;
520 linSolver.prepare(jacobian, residual);
521 linear_solve_setup_time_ = perfTimer.stop();
522 linSolver.setResidual(residual);
531template <
class TypeTag>
536 OPM_TIMEBLOCK(updateSolution);
538 if (this->param_.tolerance_max_dp_ > 0.0 || this->param_.tolerance_max_ds_ > 0.0
539 || this->param_.tolerance_max_drs_ > 0.0 || this->param_.tolerance_max_drv_ > 0.0) {
540 prepareStoringSolutionUpdate();
543 auto& newtonMethod = simulator_.model().newtonMethod();
546 newtonMethod.update_(solution,
555 OPM_TIMEBLOCK(invalidateAndUpdateIntensiveQuantities);
556 simulator_.model().invalidateAndUpdateIntensiveQuantities(0);
560 if (this->param_.tolerance_max_dp_ > 0.0 || this->param_.tolerance_max_ds_ > 0.0
561 || this->param_.tolerance_max_drs_ > 0.0 || this->param_.tolerance_max_drv_ > 0.0) {
562 storeSolutionUpdate(dx);
566template <
class TypeTag>
572 unsigned nc = simulator_.model().numGridDof();
575 const auto& elemMapper = simulator_.model().elementMapper();
576 const auto& gridView = simulator_.gridView();
577 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
579 unsigned globalElemIdx = elemMapper.index(elem);
580 solUpd_[globalElemIdx] = simulator_.model().solution(0)[globalElemIdx];
583 std::ranges::fill(solUpd_[globalElemIdx], 0.0);
587template <
class TypeTag>
592 const auto& elemMapper = simulator_.model().elementMapper();
593 const auto& gridView = simulator_.gridView();
594 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
596 unsigned globalElemIdx = elemMapper.index(elem);
597 auto& value = solUpd_[globalElemIdx];
598 const auto& update = dx[globalElemIdx];
599 assert(value.size() == update.size());
602 std::ranges::copy(update, value.begin());
606template <
class TypeTag>
611 static constexpr bool enableSolvent = Indices::solventSaturationIdx >= 0;
612 static constexpr bool enableBrine = Indices::saltConcentrationIdx >= 0;
621 for (
const auto& ix : ixCells) {
622 const auto& value = solUpd_[ix];
623 for (
int pvIdx = 0; pvIdx < static_cast<int>(value.size()); ++pvIdx) {
624 if (pvIdx == Indices::pressureSwitchIdx) {
625 dPMax = std::max(dPMax, std::abs(value[pvIdx]));
627 else if ( (pvIdx == Indices::waterSwitchIdx
628 && value.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
629 || (pvIdx == Indices::compositionSwitchIdx
630 && value.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
631 || (enableSolvent && pvIdx == Indices::solventSaturationIdx
632 && value.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss)
633 || (enableBrine && enableSaltPrecipitation && pvIdx == Indices::saltConcentrationIdx
634 && value.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp) ) {
635 dSMax = std::max(dSMax, std::abs(value[pvIdx]));
637 else if (pvIdx == Indices::compositionSwitchIdx
638 && value.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Rs) {
639 dRsMax = std::max(dRsMax, std::abs(value[pvIdx]));
641 else if (pvIdx == Indices::compositionSwitchIdx
642 && value.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Rv) {
643 dRvMax = std::max(dRvMax, std::abs(value[pvIdx]));
649 dPMax = this->grid_.comm().max(dPMax);
650 dSMax = this->grid_.comm().max(dSMax);
651 dRsMax = this->grid_.comm().max(dRsMax);
652 dRvMax = this->grid_.comm().max(dRvMax);
654 return { dPMax, dSMax, dRsMax, dRvMax };
657template <
class TypeTag>
658std::tuple<typename BlackoilModel<TypeTag>::Scalar,
663 const Scalar numAquiferPvSumLocal,
664 std::vector< Scalar >& R_sum,
665 std::vector< Scalar >& maxCoeff,
666 std::vector< Scalar >& B_avg)
668 OPM_TIMEBLOCK(convergenceReduction);
670 Scalar pvSum = pvSumLocal;
671 Scalar numAquiferPvSum = numAquiferPvSumLocal;
673 if (comm.size() > 1) {
675 std::vector< Scalar > sumBuffer;
676 std::vector< Scalar > maxBuffer;
677 const int numComp = B_avg.size();
678 sumBuffer.reserve( 2*numComp + 2 );
679 maxBuffer.reserve( numComp );
680 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
681 sumBuffer.push_back(B_avg[compIdx]);
682 sumBuffer.push_back(R_sum[compIdx]);
683 maxBuffer.push_back(maxCoeff[ compIdx]);
687 sumBuffer.push_back( pvSum );
688 sumBuffer.push_back( numAquiferPvSum );
691 comm.sum( sumBuffer.data(), sumBuffer.size() );
694 comm.max( maxBuffer.data(), maxBuffer.size() );
697 for (
int compIdx = 0, buffIdx = 0; compIdx < numComp; ++compIdx, ++buffIdx) {
698 B_avg[compIdx] = sumBuffer[buffIdx];
701 R_sum[compIdx] = sumBuffer[buffIdx];
704 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
705 maxCoeff[compIdx] = maxBuffer[compIdx];
709 pvSum = sumBuffer[sumBuffer.size()-2];
710 numAquiferPvSum = sumBuffer.back();
714 return {pvSum, numAquiferPvSum};
717template <
class TypeTag>
718std::pair<typename BlackoilModel<TypeTag>::Scalar,
722 std::vector<Scalar>& maxCoeff,
723 std::vector<Scalar>& B_avg,
724 std::vector<int>& maxCoeffCell)
726 OPM_TIMEBLOCK(localConvergenceData);
728 Scalar numAquiferPvSumLocal = 0.0;
729 const auto& model = simulator_.model();
730 const auto& problem = simulator_.problem();
732 const auto& residual = simulator_.model().linearizer().residual();
735 const auto& gridView = simulator().gridView();
738 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
739 elemCtx.updatePrimaryStencil(elem);
740 elemCtx.updatePrimaryIntensiveQuantities(0);
742 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
743 const auto& intQuants = elemCtx.intensiveQuantities(0, 0);
744 const auto& fs = intQuants.fluidState();
746 const auto pvValue = problem.referencePorosity(cell_idx, 0) *
747 model.dofTotalVolume(cell_idx);
748 pvSumLocal += pvValue;
750 if (isNumericalAquiferCell(elem)) {
751 numAquiferPvSumLocal += pvValue;
754 this->getMaxCoeff(cell_idx, intQuants, fs, residual, pvValue,
755 B_avg, R_sum, maxCoeff, maxCoeffCell);
761 const int bSize = B_avg.size();
762 for (
int i = 0; i < bSize; ++i) {
763 B_avg[i] /=
Scalar(global_nc_);
766 return {pvSumLocal, numAquiferPvSumLocal};
769template <
class TypeTag>
774 OPM_TIMEBLOCK(computeCnvErrorPv);
779 constexpr auto numPvGroups = std::vector<double>::size_type{3};
781 auto cnvPvSplit = std::pair<std::vector<double>, std::vector<int>> {
782 std::piecewise_construct,
783 std::forward_as_tuple(numPvGroups),
784 std::forward_as_tuple(numPvGroups)
787 auto maxCNV = [&B_avg, dt](
const auto& residual,
const double pvol)
790 std::inner_product(residual.begin(), residual.end(),
792 [](
const Scalar m,
const auto& x)
795 return std::max(m, abs(x));
796 }, std::multiplies<>{});
799 auto& [splitPV, cellCntPV] = cnvPvSplit;
801 const auto& model = this->simulator().model();
802 const auto& problem = this->simulator().problem();
803 const auto& residual = model.linearizer().residual();
804 const auto& gridView = this->simulator().gridView();
810 std::vector<unsigned> ixCells;
813 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
815 if (isNumericalAquiferCell(elem)) {
819 elemCtx.updatePrimaryStencil(elem);
821 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
822 const auto pvValue = problem.referencePorosity(cell_idx, 0)
823 * model.dofTotalVolume(cell_idx);
825 const auto maxCnv = maxCNV(residual[cell_idx], pvValue);
827 const auto ix = (maxCnv > this->param_.tolerance_cnv_)
828 + (maxCnv > this->param_.tolerance_cnv_relaxed_);
830 splitPV[ix] +=
static_cast<double>(pvValue);
835 (this->param_.tolerance_max_dp_ > 0.0 || this->param_.tolerance_max_ds_ > 0.0
836 || this->param_.tolerance_max_drs_ > 0.0 || this->param_.tolerance_max_drv_ > 0.0 ) ) {
837 ixCells.push_back(cell_idx);
844 this->grid_.comm().sum(splitPV .data(), splitPV .size());
845 this->grid_.comm().sum(cellCntPV.data(), cellCntPV.size());
847 return { cnvPvSplit, ixCells };
850template <
class TypeTag>
855 this->param_.tolerance_cnv_ = tuning.TRGCNV;
856 this->param_.tolerance_cnv_relaxed_ = tuning.XXXCNV;
857 this->param_.tolerance_mb_ = tuning.TRGMBE;
858 this->param_.tolerance_mb_relaxed_ = tuning.XXXMBE;
859 this->param_.newton_max_iter_ = tuning.NEWTMX;
860 this->param_.newton_min_iter_ = tuning.NEWTMN;
863template <
class TypeTag>
869 this->param_.tolerance_max_dp_ = tuning_dp.TRGDDP;
870 this->param_.tolerance_max_ds_ = tuning_dp.TRGDDS;
871 this->param_.tolerance_max_drs_ = tuning_dp.TRGDDRS;
872 this->param_.tolerance_max_drv_ = tuning_dp.TRGDDRV;
875template <
class TypeTag>
881 std::vector<Scalar>& B_avg,
882 std::vector<Scalar>& residual_norms)
884 OPM_TIMEBLOCK(getReservoirConvergence);
885 using Vector = std::vector<Scalar>;
887 const auto& iterCtx = simulator_.problem().iterationContext();
891 const int numComp = numEq;
893 Vector R_sum(numComp,
Scalar{0});
894 Vector maxCoeff(numComp, std::numeric_limits<Scalar>::lowest());
895 std::vector<int> maxCoeffCell(numComp, -1);
897 const auto [pvSumLocal, numAquiferPvSumLocal] =
898 this->localConvergenceData(R_sum, maxCoeff, B_avg, maxCoeffCell);
901 const auto& [pvSum, numAquiferPvSum] =
902 this->convergenceReduction(this->grid_.comm(),
904 numAquiferPvSumLocal,
905 R_sum, maxCoeff, B_avg);
907 auto cnvSplitData = this->characteriseCnvPvSplit(B_avg, dt);
908 report.setCnvPoreVolSplit(cnvSplitData.cnvPvSplit,
909 pvSum - numAquiferPvSum);
918 const bool relax_final_iteration_mb =
919 this->param_.min_strict_mb_iter_ < 0 && iterCtx.iteration() == maxIter;
921 const bool relax_iter_mb = this->param_.min_strict_mb_iter_ >= 0 &&
922 iterCtx.shouldRelax(this->param_.min_strict_mb_iter_);
924 const bool use_relaxed_mb = relax_final_iteration_mb
933 const bool relax_final_iteration_cnv =
934 this->param_.min_strict_cnv_iter_ < 0 && iterCtx.iteration() == maxIter;
936 const bool relax_iter_cnv = this->param_.min_strict_cnv_iter_ >= 0 &&
937 iterCtx.shouldRelax(this->param_.min_strict_cnv_iter_);
942 const auto relax_pv_fraction_cnv =
943 [&report,
this, eligible = pvSum - numAquiferPvSum]()
945 const auto& cnvPvSplit = report.cnvPvSplit().first;
949 Scalar cnvPvSum =
static_cast<Scalar>(cnvPvSplit[1] + cnvPvSplit[2]);
950 return cnvPvSum < this->param_.relaxed_max_pv_fraction_ * eligible &&
958 const bool use_dp_tol = this->param_.tolerance_max_dp_ > 0.0;
959 const bool use_ds_tol = this->param_.tolerance_max_ds_ > 0.0;
960 const bool use_drs_tol = this->param_.tolerance_max_drs_ > 0.0;
961 const bool use_drv_tol = this->param_.tolerance_max_drv_ > 0.0;
962 const bool use_dsol_tol = use_dp_tol || use_ds_tol || use_drs_tol || use_drv_tol;
963 bool relax_dsol_cnv =
false;
964 if (!iterCtx.isFirstGlobalIteration() && use_dsol_tol) {
965 maxSolUpd = getMaxSolutionUpdate(cnvSplitData.ixCells);
967 (!use_dp_tol || (maxSolUpd.
dPMax > 0.0 && maxSolUpd.
dPMax < this->param_.tolerance_max_dp_)) &&
968 (!use_ds_tol || (maxSolUpd.
dSMax > 0.0 && maxSolUpd.
dSMax < this->param_.tolerance_max_ds_)) &&
969 (!use_drs_tol || (maxSolUpd.
dRsMax > 0.0 && maxSolUpd.
dRsMax < this->param_.tolerance_max_drs_)) &&
970 (!use_drv_tol || (maxSolUpd.
dRvMax > 0.0 && maxSolUpd.
dRvMax < this->param_.tolerance_max_drv_));
974 const bool use_relaxed_cnv = relax_final_iteration_cnv
975 || relax_pv_fraction_cnv
981 Scalar tolerance_cnv_relaxed = relax_dsol_cnv ? 1e20 : param_.tolerance_cnv_relaxed_;
983 const auto tol_cnv = use_relaxed_cnv ? tolerance_cnv_relaxed : param_.tolerance_cnv_;
984 const auto tol_mb = use_relaxed_mb ? param_.tolerance_mb_relaxed_ : param_.tolerance_mb_;
985 const auto tol_cnv_energy = use_relaxed_cnv ? param_.tolerance_cnv_energy_relaxed_ : param_.tolerance_cnv_energy_;
986 const auto tol_eb = use_relaxed_mb ? param_.tolerance_energy_balance_relaxed_ : param_.tolerance_energy_balance_;
989 std::vector<Scalar> CNV(numComp);
990 std::vector<Scalar> mass_balance_residual(numComp);
991 for (
int compIdx = 0; compIdx < numComp; ++compIdx)
993 CNV[compIdx] = B_avg[compIdx] * dt * maxCoeff[compIdx];
994 mass_balance_residual[compIdx] = std::abs(B_avg[compIdx]*R_sum[compIdx]) * dt / pvSum;
995 residual_norms.push_back(CNV[compIdx]);
999 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
1001 mass_balance_residual[compIdx], CNV[compIdx],
1004 const CR::ReservoirFailure::Type types[2] = {
1005 CR::ReservoirFailure::Type::MassBalance,
1006 CR::ReservoirFailure::Type::Cnv,
1009 Scalar tol[2] = { tol_mb, tol_cnv, };
1010 if (has_energy_ && compIdx == contiEnergyEqIdx) {
1012 tol[1] = tol_cnv_energy;
1015 for (
int ii : {0, 1}) {
1016 if (std::isnan(res[ii])) {
1018 if (this->terminal_output_) {
1019 OpmLog::debug(
"NaN residual for " + this->compNames_.name(compIdx) +
" equation.");
1022 else if (res[ii] > maxResidualAllowed()) {
1023 report.setReservoirFailed({types[ii], CR::Severity::TooLarge, compIdx});
1024 if (this->terminal_output_) {
1025 OpmLog::debug(
"Too large residual for " + this->compNames_.name(compIdx) +
" equation.");
1028 else if (res[ii] < 0.0) {
1029 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
1030 if (this->terminal_output_) {
1031 OpmLog::debug(
"Negative residual for " + this->compNames_.name(compIdx) +
" equation.");
1034 else if (res[ii] > tol[ii]) {
1035 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
1038 report.setReservoirConvergenceMetric(types[ii], compIdx, res[ii], tol[ii]);
1043 this->convergencePerCell(B_avg, dt, tol_cnv, tol_cnv_energy);
1046 if (this->terminal_output_) {
1048 if (iterCtx.isFirstGlobalIteration()) {
1049 std::string msg =
"Iter";
1050 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
1052 msg += this->compNames_.name(compIdx)[0];
1056 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
1058 msg += this->compNames_.name(compIdx)[0];
1063 msg += use_dp_tol ?
" DP " :
"";
1064 msg += use_ds_tol ?
" DS " :
"";
1065 msg += use_drs_tol ?
" DRS " :
"";
1066 msg += use_drv_tol ?
" DRV " :
"";
1075 std::ostringstream ss;
1076 const std::streamsize oprec = ss.precision(3);
1077 const std::ios::fmtflags oflags = ss.setf(std::ios::scientific);
1079 ss << std::setw(4) << iterCtx.iteration();
1080 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
1081 ss << std::setw(11) << mass_balance_residual[compIdx];
1084 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
1085 ss << std::setw(11) << CNV[compIdx];
1090 [&] (
bool use_tol,
Scalar dsol) {
1094 if (iterCtx.isFirstGlobalIteration() || dsol <= 0.0) {
1095 ss << std::string(5,
' ') <<
"-" << std::string(5,
' ');
1098 ss << std::setw(11) << dsol;
1101 print_dsol(use_dp_tol, maxSolUpd.
dPMax);
1102 print_dsol(use_ds_tol, maxSolUpd.
dSMax);
1103 print_dsol(use_drs_tol, maxSolUpd.
dRsMax);
1104 print_dsol(use_drv_tol, maxSolUpd.
dRvMax);
1107 const auto mb_flag = use_relaxed_mb
1108 ? DebugFlags::RELAXED
1109 : DebugFlags::STRICT;
1111 const auto cnv_flag = relax_dsol_cnv ?
1112 DebugFlags::TUNINGDP
1114 ? DebugFlags::RELAXED
1115 : DebugFlags::STRICT);
1117 ss << std::setw(9) << make_string<TypeTag>(mb_flag)
1118 << std::setw(9) << make_string<TypeTag>(cnv_flag);
1120 ss.precision(oprec);
1123 OpmLog::debug(ss.str());
1129template <
class TypeTag>
1134 const double tol_cnv,
1135 const double tol_cnv_energy)
1137 auto& rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
1138 if (!rst_conv.hasConv()) {
1142 if (simulator_.problem().iterationContext().isFirstGlobalIteration()) {
1143 rst_conv.prepareConv();
1146 const auto& residual = simulator_.model().linearizer().residual();
1147 const auto& gridView = this->simulator().gridView();
1150 std::vector<int> convNewt(residual.size(), 0);
1153 const int numComp = B_avg.size();
1154 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
1155 elemCtx.updatePrimaryStencil(elem);
1157 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
1158 const auto pvValue = simulator_.problem().referencePorosity(cell_idx, 0) *
1159 simulator_.model().dofTotalVolume(cell_idx);
1160 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
1161 const auto tol = (has_energy_ && compIdx == contiEnergyEqIdx) ? tol_cnv_energy : tol_cnv;
1162 const Scalar cnv = std::abs(B_avg[compIdx] * residual[cell_idx][compIdx]) * dt / pvValue;
1163 if (std::isnan(cnv) || cnv > maxResidualAllowed() || cnv < 0.0 || cnv > tol) {
1171 this->grid_.comm());
1172 rst_conv.updateNewton(convNewt);
1175template <
class TypeTag>
1180 std::vector<Scalar>& residual_norms)
1182 OPM_TIMEBLOCK(getConvergence);
1184 std::vector<Scalar> B_avg(numEq, 0.0);
1187 maxIter, B_avg, residual_norms);
1189 OPM_TIMEBLOCK(getWellConvergence);
1190 report += wellModel().getWellConvergence(B_avg,
1191 report.converged());
1194 conv_monitor_.checkPenaltyCard(report, simulator_.problem().iterationContext().iteration());
1199template <
class TypeTag>
1200std::vector<std::vector<typename BlackoilModel<TypeTag>::Scalar> >
1204 OPM_TIMEBLOCK(computeFluidInPlace);
1207 std::vector<std::vector<Scalar> > regionValues(0, std::vector<Scalar>(0,0.0));
1208 return regionValues;
1211template <
class TypeTag>
1216 if (!hasNlddSolver()) {
1217 OPM_THROW(std::runtime_error,
"Cannot get local reports from a model without NLDD solver");
1219 return nlddSolver_->localAccumulatedReports();
1222template <
class TypeTag>
1223const std::vector<SimulatorReport>&
1228 OPM_THROW(std::runtime_error,
"Cannot get domain reports from a model without NLDD solver");
1229 return nlddSolver_->domainAccumulatedReports();
1232template <
class TypeTag>
1237 if (hasNlddSolver()) {
1238 nlddSolver_->writeNonlinearIterationsPerCell(odir);
1242template <
class TypeTag>
1247 if (hasNlddSolver()) {
1248 nlddSolver_->writePartitions(odir);
1252 const auto& elementMapper = this->simulator().model().elementMapper();
1253 const auto& cartMapper = this->simulator().vanguard().cartesianIndexMapper();
1255 const auto& grid = this->simulator().vanguard().grid();
1256 const auto& comm = grid.comm();
1257 const auto nDigit = 1 +
static_cast<int>(std::floor(std::log10(comm.size())));
1259 std::ofstream pfile {odir / fmt::format(
"{1:0>{0}}", nDigit, comm.rank())};
1261 for (
const auto& cell : elements(grid.leafGridView(), Dune::Partitions::interior)) {
1262 pfile << comm.rank() <<
' '
1263 << cartMapper.cartesianIndex(elementMapper.index(cell)) <<
' '
1264 << comm.rank() <<
'\n';
1268template <
class TypeTag>
1269template<
class Flu
idState,
class Res
idual>
1274 const FluidState& fs,
1275 const Residual& modelResid,
1277 std::vector<Scalar>& B_avg,
1278 std::vector<Scalar>& R_sum,
1279 std::vector<Scalar>& maxCoeff,
1280 std::vector<int>& maxCoeffCell)
1282 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1284 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1288 const unsigned sIdx = FluidSystem::solventComponentIndex(phaseIdx);
1289 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(sIdx);
1291 B_avg[compIdx] += 1.0 / fs.invB(phaseIdx).value();
1292 const auto R2 = modelResid[cell_idx][compIdx];
1294 R_sum[compIdx] += R2;
1295 const Scalar Rval = std::abs(R2) / pvValue;
1296 if (Rval > maxCoeff[compIdx]) {
1297 maxCoeff[compIdx] = Rval;
1298 maxCoeffCell[compIdx] = cell_idx;
1302 if constexpr (has_solvent_) {
1303 B_avg[contiSolventEqIdx] +=
1304 1.0 / intQuants.solventInverseFormationVolumeFactor().value();
1305 const auto R2 = modelResid[cell_idx][contiSolventEqIdx];
1306 R_sum[contiSolventEqIdx] += R2;
1307 maxCoeff[contiSolventEqIdx] = std::max(maxCoeff[contiSolventEqIdx],
1308 std::abs(R2) / pvValue);
1310 if constexpr (has_extbo_) {
1311 B_avg[contiZfracEqIdx] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
1312 const auto R2 = modelResid[cell_idx][contiZfracEqIdx];
1313 R_sum[ contiZfracEqIdx ] += R2;
1314 maxCoeff[contiZfracEqIdx] = std::max(maxCoeff[contiZfracEqIdx],
1315 std::abs(R2) / pvValue);
1317 if constexpr (has_polymer_) {
1318 B_avg[contiPolymerEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1319 const auto R2 = modelResid[cell_idx][contiPolymerEqIdx];
1320 R_sum[contiPolymerEqIdx] += R2;
1321 maxCoeff[contiPolymerEqIdx] = std::max(maxCoeff[contiPolymerEqIdx],
1322 std::abs(R2) / pvValue);
1324 if constexpr (has_foam_) {
1325 B_avg[ contiFoamEqIdx ] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
1326 const auto R2 = modelResid[cell_idx][contiFoamEqIdx];
1327 R_sum[contiFoamEqIdx] += R2;
1328 maxCoeff[contiFoamEqIdx] = std::max(maxCoeff[contiFoamEqIdx],
1329 std::abs(R2) / pvValue);
1331 if constexpr (has_brine_) {
1332 B_avg[ contiBrineEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1333 const auto R2 = modelResid[cell_idx][contiBrineEqIdx];
1334 R_sum[contiBrineEqIdx] += R2;
1335 maxCoeff[contiBrineEqIdx] = std::max(maxCoeff[contiBrineEqIdx],
1336 std::abs(R2) / pvValue);
1339 if constexpr (has_polymermw_) {
1340 static_assert(has_polymer_);
1342 B_avg[contiPolymerMWEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1346 const auto R2 = modelResid[cell_idx][contiPolymerMWEqIdx] / 100.;
1347 R_sum[contiPolymerMWEqIdx] += R2;
1348 maxCoeff[contiPolymerMWEqIdx] = std::max(maxCoeff[contiPolymerMWEqIdx],
1349 std::abs(R2) / pvValue);
1352 if constexpr (has_energy_) {
1353 B_avg[contiEnergyEqIdx] += 1.0;
1354 const auto R2 = modelResid[cell_idx][contiEnergyEqIdx];
1355 R_sum[contiEnergyEqIdx] += R2;
1356 maxCoeff[contiEnergyEqIdx] = std::max(maxCoeff[contiEnergyEqIdx],
1357 std::abs(R2) / pvValue);
1360 if constexpr (has_bioeffects_) {
1361 B_avg[contiMicrobialEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1362 const auto R1 = modelResid[cell_idx][contiMicrobialEqIdx];
1363 R_sum[contiMicrobialEqIdx] += R1;
1364 maxCoeff[contiMicrobialEqIdx] = std::max(maxCoeff[contiMicrobialEqIdx],
1365 std::abs(R1) / pvValue);
1366 B_avg[contiBiofilmEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1367 const auto R2 = modelResid[cell_idx][contiBiofilmEqIdx];
1368 R_sum[contiBiofilmEqIdx] += R2;
1369 maxCoeff[contiBiofilmEqIdx] = std::max(maxCoeff[contiBiofilmEqIdx],
1370 std::abs(R2) / pvValue);
1371 if constexpr (has_micp_) {
1372 B_avg[contiOxygenEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1373 const auto R3 = modelResid[cell_idx][contiOxygenEqIdx];
1374 R_sum[contiOxygenEqIdx] += R3;
1375 maxCoeff[contiOxygenEqIdx] = std::max(maxCoeff[contiOxygenEqIdx],
1376 std::abs(R3) / pvValue);
1377 B_avg[contiUreaEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1378 const auto R4 = modelResid[cell_idx][contiUreaEqIdx];
1379 R_sum[contiUreaEqIdx] += R4;
1380 maxCoeff[contiUreaEqIdx] = std::max(maxCoeff[contiUreaEqIdx],
1381 std::abs(R4) / pvValue);
1382 B_avg[contiCalciteEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1383 const auto R5 = modelResid[cell_idx][contiCalciteEqIdx];
1384 R_sum[contiCalciteEqIdx] += R5;
1385 maxCoeff[contiCalciteEqIdx] = std::max(maxCoeff[contiCalciteEqIdx],
1386 std::abs(R5) / pvValue);
#define OPM_END_PARALLEL_TRY_CATCH(prefix, comm)
Catch exception and throw in a parallel try-catch clause.
Definition: DeferredLoggingErrorHelpers.hpp:192
#define OPM_BEGIN_PARALLEL_TRY_CATCH()
Macro to setup the try of a parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:158
BlackoilModel(Simulator &simulator, const ModelParameters ¶m, BlackoilWellModel< TypeTag > &well_model, const bool terminal_output)
Definition: BlackoilModel_impl.hpp:78
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:264
ModelParameters param_
Definition: BlackoilModel.hpp:354
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:1235
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: BlackoilModel.hpp:67
SimulatorReportSingle nonlinearIterationNewton(const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Definition: BlackoilModel_impl.hpp:281
const std::vector< SimulatorReport > & domainAccumulatedReports() const
return the statistics of local solves accumulated for each domain on this rank
Definition: BlackoilModel_impl.hpp:1225
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:721
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilModel.hpp:66
std::vector< StepReport > convergence_reports_
Definition: BlackoilModel.hpp:371
const SimulatorReport & localAccumulatedReports() const
return the statistics of local solves accumulated for this rank
Definition: BlackoilModel_impl.hpp:1214
std::tuple< Scalar, Scalar > convergenceReduction(Parallel::Communication comm, const Scalar pvSumLocal, const Scalar numAquiferPvSumLocal, std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< Scalar > &B_avg)
Definition: BlackoilModel_impl.hpp:661
const Grid & grid_
Definition: BlackoilModel.hpp:343
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: BlackoilModel.hpp:69
void updateSolution(const BVector &dx)
Apply an update to the primary variables.
Definition: BlackoilModel_impl.hpp:534
long int global_nc_
The number of cells of the global grid.
Definition: BlackoilModel.hpp:363
DebugFlags
Definition: BlackoilModel.hpp:127
void updateTUNING(const Tuning &tuning)
Definition: BlackoilModel_impl.hpp:853
void getMaxCoeff(const unsigned cell_idx, const IntensiveQuantities &intQuants, const FluidState &fs, const Residual &modelResid, const Scalar pvValue, std::vector< Scalar > &B_avg, std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< int > &maxCoeffCell)
Definition: BlackoilModel_impl.hpp:1272
std::unique_ptr< BlackoilModelNldd< TypeTag > > nlddSolver_
Non-linear DD solver.
Definition: BlackoilModel.hpp:374
void writePartitions(const std::filesystem::path &odir) const
Definition: BlackoilModel_impl.hpp:1245
ConvergenceReport getConvergence(const SimulatorTimerInterface &timer, const int maxIter, std::vector< Scalar > &residual_norms)
Definition: BlackoilModel_impl.hpp:1178
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilModel.hpp:64
void updateTUNINGDP(const TuningDp &tuning_dp)
Definition: BlackoilModel_impl.hpp:866
void initialLinearization(SimulatorReportSingle &report, const int minIter, const int maxIter, const SimulatorTimerInterface &timer)
Definition: BlackoilModel_impl.hpp:184
CnvPvSplitData 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:772
SimulatorReportSingle assembleReservoir(const SimulatorTimerInterface &timer)
Assemble the residual and Jacobian of the nonlinear system.
Definition: BlackoilModel_impl.hpp:373
void convergencePerCell(const std::vector< Scalar > &B_avg, const double dt, const double tol_cnv, const double tol_cnv_energy)
Compute the number of Newtons required by each cell in order to satisfy the solution change convergen...
Definition: BlackoilModel_impl.hpp:1132
SimulatorReportSingle nonlinearIteration(const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Definition: BlackoilModel_impl.hpp:246
void storeSolutionUpdate(const BVector &dx)
Definition: BlackoilModel_impl.hpp:590
void solveJacobianSystem(BVector &x)
Definition: BlackoilModel_impl.hpp:472
MaxSolutionUpdateData getMaxSolutionUpdate(const std::vector< unsigned > &ixCells)
Definition: BlackoilModel_impl.hpp:609
SimulatorReportSingle prepareStep(const SimulatorTimerInterface &timer)
Definition: BlackoilModel_impl.hpp:113
ConvergenceReport getReservoirConvergence(const double reportTime, const double dt, const int maxIter, std::vector< Scalar > &B_avg, std::vector< Scalar > &residual_norms)
Definition: BlackoilModel_impl.hpp:878
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilModel.hpp:109
void prepareStoringSolutionUpdate()
Get solution update vector as a PrimaryVariable.
Definition: BlackoilModel_impl.hpp:569
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilModel.hpp:75
Scalar relativeChange() const
Definition: BlackoilModel_impl.hpp:385
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:98
Definition: ConvergenceReport.hpp:38
Severity
Definition: ConvergenceReport.hpp:49
@ ConvergenceMonitorFailure
void setReservoirFailed(const ReservoirFailure &rf)
Definition: ConvergenceReport.hpp:266
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
virtual double simulationTimeElapsed() const =0
virtual int currentStepNum() const =0
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
std::size_t countGlobalCells(const Grid &grid)
Get the number of cells of a global grid.
Definition: countGlobalCells.hpp:80
Definition: blackoilbioeffectsmodules.hh:45
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: BlackoilModel.hpp:114
Definition: BlackoilModel.hpp:119
Scalar dRvMax
Definition: BlackoilModel.hpp:123
Scalar dPMax
Definition: BlackoilModel.hpp:120
Scalar dSMax
Definition: BlackoilModel.hpp:121
Scalar dRsMax
Definition: BlackoilModel.hpp:122
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:196
std::string nonlinear_solver_
Nonlinear solver type: newton or nldd.
Definition: BlackoilModelParameters.hpp:363
Definition: AquiferGridUtils.hpp:35
Definition: SimulatorReport.hpp:122
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34
double linear_solve_time
Definition: SimulatorReport.hpp:43
double assemble_time
Definition: SimulatorReport.hpp:39
bool converged
Definition: SimulatorReport.hpp:55
double pre_post_time
Definition: SimulatorReport.hpp:40
double linear_solve_setup_time
Definition: SimulatorReport.hpp:42
unsigned int total_newton_iterations
Definition: SimulatorReport.hpp:50
double update_time
Definition: SimulatorReport.hpp:45
unsigned int total_linearizations
Definition: SimulatorReport.hpp:49
unsigned int total_linear_iterations
Definition: SimulatorReport.hpp:51