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>
45#include <fmt/format.h>
49template <
class TypeTag>
54 const bool terminal_output)
55 : simulator_(simulator)
56 , grid_(simulator_.vanguard().grid())
58 , well_model_ (well_model)
59 , terminal_output_ (terminal_output)
60 , current_relaxation_(1.0)
61 , dx_old_(simulator_.model().numGridDof())
62 , conv_monitor_(param_.monitor_params_)
69 if (terminal_output) {
70 OpmLog::info(
"Using Non-Linear Domain Decomposition solver (nldd).");
72 nlddSolver_ = std::make_unique<BlackoilModelNldd<TypeTag>>(*this);
74 if (terminal_output) {
75 OpmLog::info(
"Using Newton nonlinear solver.");
78 OPM_THROW(std::runtime_error,
"Unknown nonlinear solver option: " +
83template <
class TypeTag>
89 Dune::Timer perfTimer;
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();
129 auto getIdx = [](
unsigned phaseIdx) ->
int
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>
163 Dune::Timer perfTimer;
170 report += assembleReservoir(timer, iteration);
175 failureReport_ += report;
180 std::vector<Scalar> residual_norms;
185 auto convrep = getConvergence(timer, iteration, maxIter, residual_norms);
186 report.
converged = convrep.converged() && iteration >= minIter;
188 convergence_reports_.back().report.push_back(std::move(convrep));
192 failureReport_ += report;
193 OPM_THROW_PROBLEM(NumericalProblem,
"NaN residual found!");
195 failureReport_ += report;
196 OPM_THROW_NOLOG(NumericalProblem,
"Too large residual found!");
198 failureReport_ += report;
199 OPM_THROW_PROBLEM(ConvergenceMonitorFailure,
201 "Total penalty count exceeded cut-off-limit of {}",
202 param_.monitor_params_.cutoff_
207 residual_norms_history_.push_back(residual_norms);
210template <
class TypeTag>
211template <
class NonlinearSolverType>
216 NonlinearSolverType& nonlinear_solver)
218 if (iteration == 0) {
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")
232 result = this->nonlinearIterationNewton(iteration,
237 result = this->nlddSolver_->nonlinearIterationNldd(iteration,
242 auto& rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
243 rst_conv.update(simulator_.model().linearizer().residual());
248template <
class TypeTag>
249template <
class NonlinearSolverType>
254 NonlinearSolverType& nonlinear_solver)
258 Dune::Timer perfTimer;
260 this->initialLinearization(report,
262 this->param_.newton_min_iter_,
263 this->param_.newton_max_iter_,
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);
297 failureReport_ += report;
307 wellModel().postSolve(x);
309 if (param_.use_update_stabilization_) {
311 bool isOscillate =
false;
312 bool isStagnate =
false;
313 nonlinear_solver.detectOscillations(residual_norms_history_,
314 residual_norms_history_.size() - 1,
318 current_relaxation_ -= nonlinear_solver.relaxIncrement();
319 current_relaxation_ = std::max(current_relaxation_,
320 nonlinear_solver.relaxMax());
321 if (terminalOutputEnabled()) {
322 std::string msg =
" Oscillating behavior detected: Relaxation set to "
327 nonlinear_solver.stabilizeNonlinearUpdate(x, dx_old_, current_relaxation_);
341template <
class TypeTag>
347 Dune::Timer perfTimer;
349 simulator_.problem().endTimeStep();
354template <
class TypeTag>
358 const int iterationIdx)
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>
376 const auto& elemMapper = simulator_.model().elementMapper();
377 const auto& gridView = simulator_.gridView();
378 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
379 unsigned globalElemIdx = elemMapper.index(elem);
380 const auto& priVarsNew = simulator_.model().solution(0)[globalElemIdx];
383 pressureNew = priVarsNew[Indices::pressureSwitchIdx];
385 Scalar saturationsNew[FluidSystem::numPhases] = { 0.0 };
386 Scalar oilSaturationNew = 1.0;
387 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) &&
388 FluidSystem::numActivePhases() > 1 &&
389 priVarsNew.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
391 saturationsNew[FluidSystem::waterPhaseIdx] = priVarsNew[Indices::waterSwitchIdx];
392 oilSaturationNew -= saturationsNew[FluidSystem::waterPhaseIdx];
395 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
396 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
397 priVarsNew.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
399 assert(Indices::compositionSwitchIdx >= 0);
400 saturationsNew[FluidSystem::gasPhaseIdx] = priVarsNew[Indices::compositionSwitchIdx];
401 oilSaturationNew -= saturationsNew[FluidSystem::gasPhaseIdx];
404 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
405 saturationsNew[FluidSystem::oilPhaseIdx] = oilSaturationNew;
408 const auto& priVarsOld = simulator_.model().solution(1)[globalElemIdx];
411 pressureOld = priVarsOld[Indices::pressureSwitchIdx];
413 Scalar saturationsOld[FluidSystem::numPhases] = { 0.0 };
414 Scalar oilSaturationOld = 1.0;
417 Scalar tmp = pressureNew - pressureOld;
418 resultDelta += tmp*tmp;
419 resultDenom += pressureNew*pressureNew;
421 if (FluidSystem::numActivePhases() > 1) {
422 if (priVarsOld.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw) {
423 saturationsOld[FluidSystem::waterPhaseIdx] =
424 priVarsOld[Indices::waterSwitchIdx];
425 oilSaturationOld -= saturationsOld[FluidSystem::waterPhaseIdx];
428 if (priVarsOld.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
430 assert(Indices::compositionSwitchIdx >= 0 );
431 saturationsOld[FluidSystem::gasPhaseIdx] =
432 priVarsOld[Indices::compositionSwitchIdx];
433 oilSaturationOld -= saturationsOld[FluidSystem::gasPhaseIdx];
436 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
437 saturationsOld[FluidSystem::oilPhaseIdx] = oilSaturationOld;
439 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++ phaseIdx) {
440 Scalar tmpSat = saturationsNew[phaseIdx] - saturationsOld[phaseIdx];
441 resultDelta += tmpSat*tmpSat;
442 resultDenom += saturationsNew[phaseIdx]*saturationsNew[phaseIdx];
443 assert(std::isfinite(resultDelta));
444 assert(std::isfinite(resultDenom));
449 resultDelta = gridView.comm().sum(resultDelta);
450 resultDenom = gridView.comm().sum(resultDenom);
452 return resultDenom > 0.0 ? resultDelta / resultDenom : 0.0;
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();
464 const int numSolvers = linSolver.numAvailableSolvers();
465 if (numSolvers > 1 && (linSolver.getSolveCount() % 100 == 0)) {
466 if (terminal_output_) {
467 OpmLog::debug(
"\nRunning speed test for comparing available linear solvers.");
470 Dune::Timer perfTimer;
471 std::vector<double> times(numSolvers);
472 std::vector<double> setupTimes(numSolvers);
475 std::vector<BVector> x_trial(numSolvers, x);
476 for (
int solver = 0; solver < numSolvers; ++solver) {
478 linSolver.setActiveSolver(solver);
480 linSolver.prepare(jacobian, residual);
481 setupTimes[solver] = perfTimer.stop();
483 linSolver.setResidual(residual);
485 linSolver.solve(x_trial[solver]);
486 times[solver] = perfTimer.stop();
488 if (terminal_output_) {
489 OpmLog::debug(fmt::format(
"Solver time {}: {}", solver, times[solver]));
493 int fastest_solver = std::min_element(times.begin(), times.end()) - times.begin();
495 grid_.comm().broadcast(&fastest_solver, 1, 0);
496 linear_solve_setup_time_ = setupTimes[fastest_solver];
497 x = x_trial[fastest_solver];
498 linSolver.setActiveSolver(fastest_solver);
504 Dune::Timer perfTimer;
506 linSolver.prepare(jacobian, residual);
507 linear_solve_setup_time_ = perfTimer.stop();
508 linSolver.setResidual(residual);
517template <
class TypeTag>
522 OPM_TIMEBLOCK(updateSolution);
523 auto& newtonMethod = simulator_.model().newtonMethod();
526 newtonMethod.update_(solution,
535 OPM_TIMEBLOCK(invalidateAndUpdateIntensiveQuantities);
536 simulator_.model().invalidateAndUpdateIntensiveQuantities(0);
540template <
class TypeTag>
541std::tuple<typename BlackoilModel<TypeTag>::Scalar,
546 const Scalar numAquiferPvSumLocal,
547 std::vector< Scalar >& R_sum,
548 std::vector< Scalar >& maxCoeff,
549 std::vector< Scalar >& B_avg)
551 OPM_TIMEBLOCK(convergenceReduction);
553 Scalar pvSum = pvSumLocal;
554 Scalar numAquiferPvSum = numAquiferPvSumLocal;
556 if (comm.size() > 1) {
558 std::vector< Scalar > sumBuffer;
559 std::vector< Scalar > maxBuffer;
560 const int numComp = B_avg.size();
561 sumBuffer.reserve( 2*numComp + 2 );
562 maxBuffer.reserve( numComp );
563 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
564 sumBuffer.push_back(B_avg[compIdx]);
565 sumBuffer.push_back(R_sum[compIdx]);
566 maxBuffer.push_back(maxCoeff[ compIdx]);
570 sumBuffer.push_back( pvSum );
571 sumBuffer.push_back( numAquiferPvSum );
574 comm.sum( sumBuffer.data(), sumBuffer.size() );
577 comm.max( maxBuffer.data(), maxBuffer.size() );
580 for (
int compIdx = 0, buffIdx = 0; compIdx < numComp; ++compIdx, ++buffIdx) {
581 B_avg[compIdx] = sumBuffer[buffIdx];
584 R_sum[compIdx] = sumBuffer[buffIdx];
587 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
588 maxCoeff[compIdx] = maxBuffer[compIdx];
592 pvSum = sumBuffer[sumBuffer.size()-2];
593 numAquiferPvSum = sumBuffer.back();
597 return {pvSum, numAquiferPvSum};
600template <
class TypeTag>
601std::pair<typename BlackoilModel<TypeTag>::Scalar,
605 std::vector<Scalar>& maxCoeff,
606 std::vector<Scalar>& B_avg,
607 std::vector<int>& maxCoeffCell)
609 OPM_TIMEBLOCK(localConvergenceData);
611 Scalar numAquiferPvSumLocal = 0.0;
612 const auto& model = simulator_.model();
613 const auto& problem = simulator_.problem();
615 const auto& residual = simulator_.model().linearizer().residual();
618 const auto& gridView = simulator().gridView();
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();
629 const auto pvValue = problem.referencePorosity(cell_idx, 0) *
630 model.dofTotalVolume(cell_idx);
631 pvSumLocal += pvValue;
633 if (isNumericalAquiferCell(elem)) {
634 numAquiferPvSumLocal += pvValue;
637 this->getMaxCoeff(cell_idx, intQuants, fs, residual, pvValue,
638 B_avg, R_sum, maxCoeff, maxCoeffCell);
644 const int bSize = B_avg.size();
645 for (
int i = 0; i < bSize; ++i) {
646 B_avg[i] /=
Scalar(global_nc_);
649 return {pvSumLocal, numAquiferPvSumLocal};
652template <
class TypeTag>
653std::pair<std::vector<double>, std::vector<int>>
657 OPM_TIMEBLOCK(computeCnvErrorPv);
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,
666 std::forward_as_tuple(numPvGroups),
667 std::forward_as_tuple(numPvGroups)
670 auto maxCNV = [&B_avg, dt](
const auto& residual,
const double pvol)
673 std::inner_product(residual.begin(), residual.end(),
675 [](
const Scalar m,
const auto& x)
678 return std::max(m, abs(x));
679 }, std::multiplies<>{});
682 auto& [splitPV, cellCntPV] = cnvPvSplit;
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();
694 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
696 if (isNumericalAquiferCell(elem)) {
700 elemCtx.updatePrimaryStencil(elem);
702 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
703 const auto pvValue = problem.referencePorosity(cell_idx, 0)
704 * model.dofTotalVolume(cell_idx);
706 const auto maxCnv = maxCNV(residual[cell_idx], pvValue);
708 const auto ix = (maxCnv > this->param_.tolerance_cnv_)
709 + (maxCnv > this->param_.tolerance_cnv_relaxed_);
711 splitPV[ix] +=
static_cast<double>(pvValue);
718 this->grid_.comm().sum(splitPV .data(), splitPV .size());
719 this->grid_.comm().sum(cellCntPV.data(), cellCntPV.size());
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>
744 std::vector<Scalar>& B_avg,
745 std::vector<Scalar>& residual_norms)
747 OPM_TIMEBLOCK(getReservoirConvergence);
748 using Vector = std::vector<Scalar>;
752 const int numComp = numEq;
754 Vector R_sum(numComp,
Scalar{0});
755 Vector maxCoeff(numComp, std::numeric_limits<Scalar>::lowest());
756 std::vector<int> maxCoeffCell(numComp, -1);
758 const auto [pvSumLocal, numAquiferPvSumLocal] =
759 this->localConvergenceData(R_sum, maxCoeff, B_avg, maxCoeffCell);
762 const auto& [pvSum, numAquiferPvSum] =
763 this->convergenceReduction(this->grid_.comm(),
765 numAquiferPvSumLocal,
766 R_sum, maxCoeff, B_avg);
768 report.setCnvPoreVolSplit(this->characteriseCnvPvSplit(B_avg, dt),
769 pvSum - numAquiferPvSum);
778 const bool relax_final_iteration_mb =
779 this->param_.min_strict_mb_iter_ < 0 && iteration == maxIter;
781 const bool use_relaxed_mb = relax_final_iteration_mb ||
782 (this->param_.min_strict_mb_iter_ >= 0 &&
783 iteration >= this->param_.min_strict_mb_iter_);
791 const bool relax_final_iteration_cnv =
792 this->param_.min_strict_cnv_iter_ < 0 && iteration == maxIter;
794 const bool relax_iter_cnv =
795 this->param_.min_strict_cnv_iter_ >= 0 &&
796 iteration >= this->param_.min_strict_cnv_iter_;
801 const auto relax_pv_fraction_cnv =
802 [&report,
this, eligible = pvSum - numAquiferPvSum]()
804 const auto& cnvPvSplit = report.cnvPvSplit().first;
808 return static_cast<Scalar>(cnvPvSplit[1] + cnvPvSplit[2]) <
809 this->param_.relaxed_max_pv_fraction_ * eligible;
812 const bool use_relaxed_cnv = relax_final_iteration_cnv
813 || relax_pv_fraction_cnv
816 if ((relax_final_iteration_mb || relax_final_iteration_cnv) &&
817 this->terminal_output_)
819 std::string message =
820 "Number of newton iterations reached its maximum "
821 "try to continue with relaxed tolerances:";
823 if (relax_final_iteration_mb) {
824 message += fmt::format(
" MB: {:.1e}", param_.tolerance_mb_relaxed_);
827 if (relax_final_iteration_cnv) {
828 message += fmt::format(
" CNV: {:.1e}", param_.tolerance_cnv_relaxed_);
831 OpmLog::debug(message);
834 const auto tol_cnv = use_relaxed_cnv ? param_.tolerance_cnv_relaxed_ : param_.tolerance_cnv_;
835 const auto tol_mb = use_relaxed_mb ? param_.tolerance_mb_relaxed_ : param_.tolerance_mb_;
836 const auto tol_cnv_energy = use_relaxed_cnv ? param_.tolerance_cnv_energy_relaxed_ : param_.tolerance_cnv_energy_;
837 const auto tol_eb = use_relaxed_mb ? param_.tolerance_energy_balance_relaxed_ : param_.tolerance_energy_balance_;
840 std::vector<Scalar> CNV(numComp);
841 std::vector<Scalar> mass_balance_residual(numComp);
842 for (
int compIdx = 0; compIdx < numComp; ++compIdx)
844 CNV[compIdx] = B_avg[compIdx] * dt * maxCoeff[compIdx];
845 mass_balance_residual[compIdx] = std::abs(B_avg[compIdx]*R_sum[compIdx]) * dt / pvSum;
846 residual_norms.push_back(CNV[compIdx]);
850 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
852 mass_balance_residual[compIdx], CNV[compIdx],
855 const CR::ReservoirFailure::Type types[2] = {
856 CR::ReservoirFailure::Type::MassBalance,
857 CR::ReservoirFailure::Type::Cnv,
860 Scalar tol[2] = { tol_mb, tol_cnv, };
861 if (has_energy_ && compIdx == contiEnergyEqIdx) {
863 tol[1] = tol_cnv_energy;
866 for (
int ii : {0, 1}) {
867 if (std::isnan(res[ii])) {
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.");
885 else if (res[ii] > tol[ii]) {
886 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
889 report.setReservoirConvergenceMetric(types[ii], compIdx, res[ii], tol[ii]);
894 this->convergencePerCell(B_avg, dt, tol_cnv, tol_cnv_energy, iteration);
897 if (this->terminal_output_) {
899 if (iteration == 0) {
900 std::string msg =
"Iter";
901 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
903 msg += this->compNames_.name(compIdx)[0];
907 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
909 msg += this->compNames_.name(compIdx)[0];
916 std::ostringstream ss;
917 const std::streamsize oprec = ss.precision(3);
918 const std::ios::fmtflags oflags = ss.setf(std::ios::scientific);
920 ss << std::setw(4) << iteration;
921 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
922 ss << std::setw(11) << mass_balance_residual[compIdx];
925 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
926 ss << std::setw(11) << CNV[compIdx];
932 OpmLog::debug(ss.str());
938template <
class TypeTag>
943 const double tol_cnv,
944 const double tol_cnv_energy,
947 auto& rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
948 if (!rst_conv.hasConv()) {
952 if (iteration == 0) {
953 rst_conv.prepareConv();
956 const auto& residual = simulator_.model().linearizer().residual();
957 const auto& gridView = this->simulator().gridView();
960 std::vector<int> convNewt(residual.size(), 0);
963 const int numComp = B_avg.size();
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);
970 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
971 const auto tol = (has_energy_ && compIdx == contiEnergyEqIdx) ? tol_cnv_energy : tol_cnv;
972 const Scalar cnv = std::abs(B_avg[compIdx] * residual[cell_idx][compIdx]) * dt / pvValue;
973 if (std::isnan(cnv) || cnv > maxResidualAllowed() || cnv < 0.0 || cnv > tol) {
982 rst_conv.updateNewton(convNewt);
985template <
class TypeTag>
991 std::vector<Scalar>& residual_norms)
993 OPM_TIMEBLOCK(getConvergence);
995 std::vector<Scalar> B_avg(numEq, 0.0);
998 iteration, maxIter, B_avg, residual_norms);
1000 OPM_TIMEBLOCK(getWellConvergence);
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> >
1015 OPM_TIMEBLOCK(computeFluidInPlace);
1018 std::vector<std::vector<Scalar> > regionValues(0, std::vector<Scalar>(0,0.0));
1019 return regionValues;
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())};
1072 for (
const auto& cell : elements(grid.leafGridView(), Dune::Partitions::interior)) {
1073 pfile << comm.rank() <<
' '
1074 << cartMapper.cartesianIndex(elementMapper.index(cell)) <<
' '
1075 << comm.rank() <<
'\n';
1079template <
class TypeTag>
1080template<
class Flu
idState,
class Res
idual>
1085 const FluidState& fs,
1086 const Residual& modelResid,
1088 std::vector<Scalar>& B_avg,
1089 std::vector<Scalar>& R_sum,
1090 std::vector<Scalar>& maxCoeff,
1091 std::vector<int>& maxCoeffCell)
1093 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1095 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1099 const unsigned sIdx = FluidSystem::solventComponentIndex(phaseIdx);
1100 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(sIdx);
1102 B_avg[compIdx] += 1.0 / fs.invB(phaseIdx).value();
1103 const auto R2 = modelResid[cell_idx][compIdx];
1105 R_sum[compIdx] += R2;
1106 const Scalar Rval = std::abs(R2) / pvValue;
1107 if (Rval > maxCoeff[compIdx]) {
1108 maxCoeff[compIdx] = Rval;
1109 maxCoeffCell[compIdx] = cell_idx;
1113 if constexpr (has_solvent_) {
1114 B_avg[contiSolventEqIdx] +=
1115 1.0 / intQuants.solventInverseFormationVolumeFactor().value();
1116 const auto R2 = modelResid[cell_idx][contiSolventEqIdx];
1117 R_sum[contiSolventEqIdx] += R2;
1118 maxCoeff[contiSolventEqIdx] = std::max(maxCoeff[contiSolventEqIdx],
1119 std::abs(R2) / pvValue);
1121 if constexpr (has_extbo_) {
1122 B_avg[contiZfracEqIdx] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
1123 const auto R2 = modelResid[cell_idx][contiZfracEqIdx];
1124 R_sum[ contiZfracEqIdx ] += R2;
1125 maxCoeff[contiZfracEqIdx] = std::max(maxCoeff[contiZfracEqIdx],
1126 std::abs(R2) / pvValue);
1128 if constexpr (has_polymer_) {
1129 B_avg[contiPolymerEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1130 const auto R2 = modelResid[cell_idx][contiPolymerEqIdx];
1131 R_sum[contiPolymerEqIdx] += R2;
1132 maxCoeff[contiPolymerEqIdx] = std::max(maxCoeff[contiPolymerEqIdx],
1133 std::abs(R2) / pvValue);
1135 if constexpr (has_foam_) {
1136 B_avg[ contiFoamEqIdx ] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
1137 const auto R2 = modelResid[cell_idx][contiFoamEqIdx];
1138 R_sum[contiFoamEqIdx] += R2;
1139 maxCoeff[contiFoamEqIdx] = std::max(maxCoeff[contiFoamEqIdx],
1140 std::abs(R2) / pvValue);
1142 if constexpr (has_brine_) {
1143 B_avg[ contiBrineEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1144 const auto R2 = modelResid[cell_idx][contiBrineEqIdx];
1145 R_sum[contiBrineEqIdx] += R2;
1146 maxCoeff[contiBrineEqIdx] = std::max(maxCoeff[contiBrineEqIdx],
1147 std::abs(R2) / pvValue);
1150 if constexpr (has_polymermw_) {
1151 static_assert(has_polymer_);
1153 B_avg[contiPolymerMWEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1157 const auto R2 = modelResid[cell_idx][contiPolymerMWEqIdx] / 100.;
1158 R_sum[contiPolymerMWEqIdx] += R2;
1159 maxCoeff[contiPolymerMWEqIdx] = std::max(maxCoeff[contiPolymerMWEqIdx],
1160 std::abs(R2) / pvValue);
1163 if constexpr (has_energy_) {
1164 B_avg[contiEnergyEqIdx] += 1.0 / (4.182e1);
1165 const auto R2 = modelResid[cell_idx][contiEnergyEqIdx];
1166 R_sum[contiEnergyEqIdx] += R2;
1167 maxCoeff[contiEnergyEqIdx] = std::max(maxCoeff[contiEnergyEqIdx],
1168 std::abs(R2) / pvValue);
1171 if constexpr (has_micp_) {
1172 B_avg[contiMicrobialEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1173 const auto R1 = modelResid[cell_idx][contiMicrobialEqIdx];
1174 R_sum[contiMicrobialEqIdx] += R1;
1175 maxCoeff[contiMicrobialEqIdx] = std::max(maxCoeff[contiMicrobialEqIdx],
1176 std::abs(R1) / pvValue);
1177 B_avg[contiOxygenEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1178 const auto R2 = modelResid[cell_idx][contiOxygenEqIdx];
1179 R_sum[contiOxygenEqIdx] += R2;
1180 maxCoeff[contiOxygenEqIdx] = std::max(maxCoeff[contiOxygenEqIdx],
1181 std::abs(R2) / pvValue);
1182 B_avg[contiUreaEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1183 const auto R3 = modelResid[cell_idx][contiUreaEqIdx];
1184 R_sum[contiUreaEqIdx] += R3;
1185 maxCoeff[contiUreaEqIdx] = std::max(maxCoeff[contiUreaEqIdx],
1186 std::abs(R3) / pvValue);
1187 B_avg[contiBiofilmEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1188 const auto R4 = modelResid[cell_idx][contiBiofilmEqIdx];
1189 R_sum[contiBiofilmEqIdx] += R4;
1190 maxCoeff[contiBiofilmEqIdx] = std::max(maxCoeff[contiBiofilmEqIdx],
1191 std::abs(R4) / pvValue);
1192 B_avg[contiCalciteEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1193 const auto R5 = modelResid[cell_idx][contiCalciteEqIdx];
1194 R_sum[contiCalciteEqIdx] += R5;
1195 maxCoeff[contiCalciteEqIdx] = std::max(maxCoeff[contiCalciteEqIdx],
1196 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: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
ModelParameters param_
Definition: BlackoilModel.hpp:331
SimulatorReportSingle nonlinearIteration(const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
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)
Definition: BlackoilModel_impl.hpp:988
ConvergenceReport getReservoirConvergence(const double reportTime, const double dt, const int iteration, const int maxIter, std::vector< Scalar > &B_avg, std::vector< Scalar > &residual_norms)
Definition: BlackoilModel_impl.hpp:740
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: BlackoilModel.hpp:67
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
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilModel.hpp:66
std::vector< StepReport > convergence_reports_
Definition: BlackoilModel.hpp:346
const SimulatorReport & localAccumulatedReports() const
return the statistics of local solves accumulated for this rank
Definition: BlackoilModel_impl.hpp:1025
SimulatorReportSingle afterStep(const SimulatorTimerInterface &)
Definition: BlackoilModel_impl.hpp:344
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:544
const Grid & grid_
Definition: BlackoilModel.hpp:321
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: BlackoilModel.hpp:69
void initialLinearization(SimulatorReportSingle &report, const int iteration, const int minIter, const int maxIter, const SimulatorTimerInterface &timer)
Definition: BlackoilModel_impl.hpp:155
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:340
void updateTUNING(const Tuning &tuning)
Definition: BlackoilModel_impl.hpp:727
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:1083
std::unique_ptr< BlackoilModelNldd< TypeTag > > nlddSolver_
Non-linear DD solver.
Definition: BlackoilModel.hpp:349
SimulatorReportSingle nonlinearIterationNewton(const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Definition: BlackoilModel_impl.hpp:252
void writePartitions(const std::filesystem::path &odir) const
Definition: BlackoilModel_impl.hpp:1056
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
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilModel.hpp:64
void solveJacobianSystem(BVector &x)
Definition: BlackoilModel_impl.hpp:458
SimulatorReportSingle prepareStep(const SimulatorTimerInterface &timer)
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
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilModel.hpp:107
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilModel.hpp:75
Scalar relativeChange() const
Definition: BlackoilModel_impl.hpp:371
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:102
Definition: ConvergenceReport.hpp:38
Severity
Definition: ConvergenceReport.hpp:49
@ ConvergenceMonitorFailure
void setReservoirFailed(const ReservoirFailure &rf)
Definition: ConvergenceReport.hpp:264
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: blackoilboundaryratevector.hh:39
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:179
std::string nonlinear_solver_
Nonlinear solver type: newton or nldd.
Definition: BlackoilModelParameters.hpp:332
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