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>
47#include <fmt/format.h>
51template <
class TypeTag>
56 const bool terminal_output)
57 : simulator_(simulator)
58 , grid_(simulator_.vanguard().grid())
61 , well_model_ (well_model)
62 , terminal_output_ (terminal_output)
63 , current_relaxation_(1.0)
64 , dx_old_(simulator_.model().numGridDof())
65 , conv_monitor_(param_.monitor_params_)
72 if (terminal_output) {
73 OpmLog::info(
"Using Non-Linear Domain Decomposition solver (nldd).");
75 nlddSolver_ = std::make_unique<BlackoilModelNldd<TypeTag>>(*this);
77 if (terminal_output) {
78 OpmLog::info(
"Using Newton nonlinear solver.");
81 OPM_THROW(std::runtime_error,
"Unknown nonlinear solver option: " +
86template <
class TypeTag>
92 Dune::Timer perfTimer;
96 if (grid_.comm().size() > 1 && grid_.comm().max(lastStepFailed) != grid_.comm().min(lastStepFailed)) {
97 OPM_THROW(std::runtime_error,
"Misalignment of the parallel simulation run in prepareStep " +
98 "- the previous step succeeded on some ranks but failed on others.");
100 if (lastStepFailed) {
101 simulator_.model().updateFailed();
104 simulator_.model().advanceTimeLevel();
113 simulator_.model().newtonMethod().setIterationIndex(0);
115 simulator_.problem().beginTimeStep();
117 unsigned numDof = simulator_.model().numGridDof();
118 wasSwitched_.resize(numDof);
119 std::fill(wasSwitched_.begin(), wasSwitched_.end(),
false);
121 if (param_.update_equations_scaling_) {
122 OpmLog::error(
"Equation scaling not supported");
126 if (hasNlddSolver()) {
127 nlddSolver_->prepareStep();
132 auto getIdx = [](
unsigned phaseIdx) ->
int
134 if (FluidSystem::phaseIsActive(phaseIdx)) {
135 const unsigned sIdx = FluidSystem::solventComponentIndex(phaseIdx);
136 return Indices::canonicalToActiveComponentIndex(sIdx);
141 const auto& schedule = simulator_.vanguard().schedule();
142 auto& rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
143 rst_conv.init(simulator_.vanguard().globalNumCells(),
145 {getIdx(FluidSystem::oilPhaseIdx),
146 getIdx(FluidSystem::gasPhaseIdx),
147 getIdx(FluidSystem::waterPhaseIdx),
155template <
class TypeTag>
166 Dune::Timer perfTimer;
173 report += assembleReservoir(timer, iteration);
178 failureReport_ += report;
183 std::vector<Scalar> residual_norms;
188 auto convrep = getConvergence(timer, iteration, maxIter, residual_norms);
189 report.
converged = convrep.converged() && iteration >= minIter;
191 convergence_reports_.back().report.push_back(std::move(convrep));
195 failureReport_ += report;
196 OPM_THROW_PROBLEM(NumericalProblem,
"NaN residual found!");
198 failureReport_ += report;
199 OPM_THROW_NOLOG(NumericalProblem,
"Too large residual found!");
201 failureReport_ += report;
202 OPM_THROW_PROBLEM(ConvergenceMonitorFailure,
204 "Total penalty count exceeded cut-off-limit of {}",
205 param_.monitor_params_.cutoff_
210 residual_norms_history_.push_back(residual_norms);
213template <
class TypeTag>
214template <
class NonlinearSolverType>
219 NonlinearSolverType& nonlinear_solver)
221 if (iteration == 0) {
224 residual_norms_history_.clear();
225 conv_monitor_.reset();
226 current_relaxation_ = 1.0;
229 convergence_reports_.back().report.reserve(11);
233 if (this->param_.nonlinear_solver_ !=
"nldd")
235 result = this->nonlinearIterationNewton(iteration,
240 result = this->nlddSolver_->nonlinearIterationNldd(iteration,
245 auto& rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
246 rst_conv.update(simulator_.model().linearizer().residual());
251template <
class TypeTag>
252template <
class NonlinearSolverType>
257 NonlinearSolverType& nonlinear_solver)
261 Dune::Timer perfTimer;
263 this->initialLinearization(report,
265 nonlinear_solver.minIter(),
266 nonlinear_solver.maxIter(),
276 unsigned nc = simulator_.model().numGridDof();
280 linear_solve_setup_time_ = 0.0;
285 wellModel().linearize(simulator().model().linearizer().jacobian(),
286 simulator().model().linearizer().residual());
289 solveJacobianSystem(x);
300 failureReport_ += report;
310 wellModel().postSolve(x);
312 if (param_.use_update_stabilization_) {
314 bool isOscillate =
false;
315 bool isStagnate =
false;
316 nonlinear_solver.detectOscillations(residual_norms_history_,
317 residual_norms_history_.size() - 1,
321 current_relaxation_ -= nonlinear_solver.relaxIncrement();
322 current_relaxation_ = std::max(current_relaxation_,
323 nonlinear_solver.relaxMax());
324 if (terminalOutputEnabled()) {
325 std::string msg =
" Oscillating behavior detected: Relaxation set to "
330 nonlinear_solver.stabilizeNonlinearUpdate(x, dx_old_, current_relaxation_);
344template <
class TypeTag>
350 Dune::Timer perfTimer;
352 simulator_.problem().endTimeStep();
357template <
class TypeTag>
361 const int iterationIdx)
364 simulator_.model().newtonMethod().setIterationIndex(iterationIdx);
365 simulator_.problem().beginIteration();
366 simulator_.model().linearizer().linearizeDomain();
367 simulator_.problem().endIteration();
368 return wellModel().lastReport();
371template <
class TypeTag>
379 const auto& elemMapper = simulator_.model().elementMapper();
380 const auto& gridView = simulator_.gridView();
381 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
382 unsigned globalElemIdx = elemMapper.index(elem);
383 const auto& priVarsNew = simulator_.model().solution(0)[globalElemIdx];
386 pressureNew = priVarsNew[Indices::pressureSwitchIdx];
388 Scalar saturationsNew[FluidSystem::numPhases] = { 0.0 };
389 Scalar oilSaturationNew = 1.0;
390 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) &&
391 FluidSystem::numActivePhases() > 1 &&
392 priVarsNew.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
394 saturationsNew[FluidSystem::waterPhaseIdx] = priVarsNew[Indices::waterSwitchIdx];
395 oilSaturationNew -= saturationsNew[FluidSystem::waterPhaseIdx];
398 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
399 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
400 priVarsNew.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
402 assert(Indices::compositionSwitchIdx >= 0);
403 saturationsNew[FluidSystem::gasPhaseIdx] = priVarsNew[Indices::compositionSwitchIdx];
404 oilSaturationNew -= saturationsNew[FluidSystem::gasPhaseIdx];
407 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
408 saturationsNew[FluidSystem::oilPhaseIdx] = oilSaturationNew;
411 const auto& priVarsOld = simulator_.model().solution(1)[globalElemIdx];
414 pressureOld = priVarsOld[Indices::pressureSwitchIdx];
416 Scalar saturationsOld[FluidSystem::numPhases] = { 0.0 };
417 Scalar oilSaturationOld = 1.0;
420 Scalar tmp = pressureNew - pressureOld;
421 resultDelta += tmp*tmp;
422 resultDenom += pressureNew*pressureNew;
424 if (FluidSystem::numActivePhases() > 1) {
425 if (priVarsOld.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw) {
426 saturationsOld[FluidSystem::waterPhaseIdx] =
427 priVarsOld[Indices::waterSwitchIdx];
428 oilSaturationOld -= saturationsOld[FluidSystem::waterPhaseIdx];
431 if (priVarsOld.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
433 assert(Indices::compositionSwitchIdx >= 0 );
434 saturationsOld[FluidSystem::gasPhaseIdx] =
435 priVarsOld[Indices::compositionSwitchIdx];
436 oilSaturationOld -= saturationsOld[FluidSystem::gasPhaseIdx];
439 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
440 saturationsOld[FluidSystem::oilPhaseIdx] = oilSaturationOld;
442 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++ phaseIdx) {
443 Scalar tmpSat = saturationsNew[phaseIdx] - saturationsOld[phaseIdx];
444 resultDelta += tmpSat*tmpSat;
445 resultDenom += saturationsNew[phaseIdx]*saturationsNew[phaseIdx];
446 assert(std::isfinite(resultDelta));
447 assert(std::isfinite(resultDenom));
452 resultDelta = gridView.comm().sum(resultDelta);
453 resultDenom = gridView.comm().sum(resultDenom);
455 return resultDenom > 0.0 ? resultDelta / resultDenom : 0.0;
458template <
class TypeTag>
463 auto& jacobian = simulator_.model().linearizer().jacobian().istlMatrix();
464 auto& residual = simulator_.model().linearizer().residual();
465 auto& linSolver = simulator_.model().newtonMethod().linearSolver();
467 const int numSolvers = linSolver.numAvailableSolvers();
468 if (numSolvers > 1 && (linSolver.getSolveCount() % 100 == 0)) {
469 if (terminal_output_) {
470 OpmLog::debug(
"\nRunning speed test for comparing available linear solvers.");
473 Dune::Timer perfTimer;
474 std::vector<double> times(numSolvers);
475 std::vector<double> setupTimes(numSolvers);
478 std::vector<BVector> x_trial(numSolvers, x);
479 for (
int solver = 0; solver < numSolvers; ++solver) {
481 linSolver.setActiveSolver(solver);
483 linSolver.prepare(jacobian, residual);
484 setupTimes[solver] = perfTimer.stop();
486 linSolver.setResidual(residual);
488 linSolver.solve(x_trial[solver]);
489 times[solver] = perfTimer.stop();
491 if (terminal_output_) {
492 OpmLog::debug(fmt::format(
"Solver time {}: {}", solver, times[solver]));
496 int fastest_solver = std::min_element(times.begin(), times.end()) - times.begin();
498 grid_.comm().broadcast(&fastest_solver, 1, 0);
499 linear_solve_setup_time_ = setupTimes[fastest_solver];
500 x = x_trial[fastest_solver];
501 linSolver.setActiveSolver(fastest_solver);
507 Dune::Timer perfTimer;
509 linSolver.prepare(jacobian, residual);
510 linear_solve_setup_time_ = perfTimer.stop();
511 linSolver.setResidual(residual);
520template <
class TypeTag>
525 OPM_TIMEBLOCK(updateSolution);
526 auto& newtonMethod = simulator_.model().newtonMethod();
529 newtonMethod.update_(solution,
538 OPM_TIMEBLOCK(invalidateAndUpdateIntensiveQuantities);
539 simulator_.model().invalidateAndUpdateIntensiveQuantities(0);
543template <
class TypeTag>
544std::tuple<typename BlackoilModel<TypeTag>::Scalar,
549 const Scalar numAquiferPvSumLocal,
550 std::vector< Scalar >& R_sum,
551 std::vector< Scalar >& maxCoeff,
552 std::vector< Scalar >& B_avg)
554 OPM_TIMEBLOCK(convergenceReduction);
556 Scalar pvSum = pvSumLocal;
557 Scalar numAquiferPvSum = numAquiferPvSumLocal;
559 if (comm.size() > 1) {
561 std::vector< Scalar > sumBuffer;
562 std::vector< Scalar > maxBuffer;
563 const int numComp = B_avg.size();
564 sumBuffer.reserve( 2*numComp + 2 );
565 maxBuffer.reserve( numComp );
566 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
567 sumBuffer.push_back(B_avg[compIdx]);
568 sumBuffer.push_back(R_sum[compIdx]);
569 maxBuffer.push_back(maxCoeff[ compIdx]);
573 sumBuffer.push_back( pvSum );
574 sumBuffer.push_back( numAquiferPvSum );
577 comm.sum( sumBuffer.data(), sumBuffer.size() );
580 comm.max( maxBuffer.data(), maxBuffer.size() );
583 for (
int compIdx = 0, buffIdx = 0; compIdx < numComp; ++compIdx, ++buffIdx) {
584 B_avg[compIdx] = sumBuffer[buffIdx];
587 R_sum[compIdx] = sumBuffer[buffIdx];
590 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
591 maxCoeff[compIdx] = maxBuffer[compIdx];
595 pvSum = sumBuffer[sumBuffer.size()-2];
596 numAquiferPvSum = sumBuffer.back();
600 return {pvSum, numAquiferPvSum};
603template <
class TypeTag>
604std::pair<typename BlackoilModel<TypeTag>::Scalar,
608 std::vector<Scalar>& maxCoeff,
609 std::vector<Scalar>& B_avg,
610 std::vector<int>& maxCoeffCell)
612 OPM_TIMEBLOCK(localConvergenceData);
614 Scalar numAquiferPvSumLocal = 0.0;
615 const auto& model = simulator_.model();
616 const auto& problem = simulator_.problem();
618 const auto& residual = simulator_.model().linearizer().residual();
621 const auto& gridView = simulator().gridView();
624 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
625 elemCtx.updatePrimaryStencil(elem);
626 elemCtx.updatePrimaryIntensiveQuantities(0);
628 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
629 const auto& intQuants = elemCtx.intensiveQuantities(0, 0);
630 const auto& fs = intQuants.fluidState();
632 const auto pvValue = problem.referencePorosity(cell_idx, 0) *
633 model.dofTotalVolume(cell_idx);
634 pvSumLocal += pvValue;
636 if (isNumericalAquiferCell(elem)) {
637 numAquiferPvSumLocal += pvValue;
640 this->getMaxCoeff(cell_idx, intQuants, fs, residual, pvValue,
641 B_avg, R_sum, maxCoeff, maxCoeffCell);
647 const int bSize = B_avg.size();
648 for (
int i = 0; i < bSize; ++i) {
649 B_avg[i] /=
Scalar(global_nc_);
652 return {pvSumLocal, numAquiferPvSumLocal};
655template <
class TypeTag>
656std::pair<std::vector<double>, std::vector<int>>
660 OPM_TIMEBLOCK(computeCnvErrorPv);
665 constexpr auto numPvGroups = std::vector<double>::size_type{3};
667 auto cnvPvSplit = std::pair<std::vector<double>, std::vector<int>> {
668 std::piecewise_construct,
669 std::forward_as_tuple(numPvGroups),
670 std::forward_as_tuple(numPvGroups)
673 auto maxCNV = [&B_avg, dt](
const auto& residual,
const double pvol)
676 std::inner_product(residual.begin(), residual.end(),
678 [](
const Scalar m,
const auto& x)
681 return std::max(m, abs(x));
682 }, std::multiplies<>{});
685 auto& [splitPV, cellCntPV] = cnvPvSplit;
687 const auto& model = this->simulator().model();
688 const auto& problem = this->simulator().problem();
689 const auto& residual = model.linearizer().residual();
690 const auto& gridView = this->simulator().gridView();
697 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
699 if (isNumericalAquiferCell(elem)) {
703 elemCtx.updatePrimaryStencil(elem);
705 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
706 const auto pvValue = problem.referencePorosity(cell_idx, 0)
707 * model.dofTotalVolume(cell_idx);
709 const auto maxCnv = maxCNV(residual[cell_idx], pvValue);
711 const auto ix = (maxCnv > this->param_.tolerance_cnv_)
712 + (maxCnv > this->param_.tolerance_cnv_relaxed_);
714 splitPV[ix] +=
static_cast<double>(pvValue);
721 this->grid_.comm().sum(splitPV .data(), splitPV .size());
722 this->grid_.comm().sum(cellCntPV.data(), cellCntPV.size());
727template <
class TypeTag>
732 this->param_.tolerance_mb_ = tuning.XXXMBE;
734 if (terminal_output_) {
735 OpmLog::debug(fmt::format(
"Setting BlackoilModel mass "
736 "balance limit (XXXMBE) to {:.2e}",
741template <
class TypeTag>
748 std::vector<Scalar>& B_avg,
749 std::vector<Scalar>& residual_norms)
751 OPM_TIMEBLOCK(getReservoirConvergence);
752 using Vector = std::vector<Scalar>;
756 const int numComp = numEq;
758 Vector R_sum(numComp,
Scalar{0});
759 Vector maxCoeff(numComp, std::numeric_limits<Scalar>::lowest());
760 std::vector<int> maxCoeffCell(numComp, -1);
762 const auto [pvSumLocal, numAquiferPvSumLocal] =
763 this->localConvergenceData(R_sum, maxCoeff, B_avg, maxCoeffCell);
766 const auto& [pvSum, numAquiferPvSum] =
767 this->convergenceReduction(this->grid_.comm(),
769 numAquiferPvSumLocal,
770 R_sum, maxCoeff, B_avg);
772 report.setCnvPoreVolSplit(this->characteriseCnvPvSplit(B_avg, dt),
773 pvSum - numAquiferPvSum);
782 const bool relax_final_iteration_mb =
783 this->param_.min_strict_mb_iter_ < 0 && iteration == maxIter;
785 const bool use_relaxed_mb = relax_final_iteration_mb ||
786 (this->param_.min_strict_mb_iter_ >= 0 &&
787 iteration >= this->param_.min_strict_mb_iter_);
795 const bool relax_final_iteration_cnv =
796 this->param_.min_strict_cnv_iter_ < 0 && iteration == maxIter;
798 const bool relax_iter_cnv =
799 this->param_.min_strict_cnv_iter_ >= 0 &&
800 iteration >= this->param_.min_strict_cnv_iter_;
805 const auto relax_pv_fraction_cnv =
806 [&report,
this, eligible = pvSum - numAquiferPvSum]()
808 const auto& cnvPvSplit = report.cnvPvSplit().first;
812 return static_cast<Scalar>(cnvPvSplit[1] + cnvPvSplit[2]) <
813 this->param_.relaxed_max_pv_fraction_ * eligible;
816 const bool use_relaxed_cnv = relax_final_iteration_cnv
817 || relax_pv_fraction_cnv
820 if ((relax_final_iteration_mb || relax_final_iteration_cnv) &&
821 this->terminal_output_)
823 std::string message =
824 "Number of newton iterations reached its maximum "
825 "try to continue with relaxed tolerances:";
827 if (relax_final_iteration_mb) {
828 message += fmt::format(
" MB: {:.1e}", param_.tolerance_mb_relaxed_);
831 if (relax_final_iteration_cnv) {
832 message += fmt::format(
" CNV: {:.1e}", param_.tolerance_cnv_relaxed_);
835 OpmLog::debug(message);
838 const auto tol_cnv = use_relaxed_cnv ? param_.tolerance_cnv_relaxed_ : param_.tolerance_cnv_;
839 const auto tol_mb = use_relaxed_mb ? param_.tolerance_mb_relaxed_ : param_.tolerance_mb_;
840 const auto tol_cnv_energy = use_relaxed_cnv ? param_.tolerance_cnv_energy_relaxed_ : param_.tolerance_cnv_energy_;
841 const auto tol_eb = use_relaxed_mb ? param_.tolerance_energy_balance_relaxed_ : param_.tolerance_energy_balance_;
844 std::vector<Scalar> CNV(numComp);
845 std::vector<Scalar> mass_balance_residual(numComp);
846 for (
int compIdx = 0; compIdx < numComp; ++compIdx)
848 CNV[compIdx] = B_avg[compIdx] * dt * maxCoeff[compIdx];
849 mass_balance_residual[compIdx] = std::abs(B_avg[compIdx]*R_sum[compIdx]) * dt / pvSum;
850 residual_norms.push_back(CNV[compIdx]);
854 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
856 mass_balance_residual[compIdx], CNV[compIdx],
859 const CR::ReservoirFailure::Type types[2] = {
860 CR::ReservoirFailure::Type::MassBalance,
861 CR::ReservoirFailure::Type::Cnv,
864 Scalar tol[2] = { tol_mb, tol_cnv, };
865 if (has_energy_ && compIdx == contiEnergyEqIdx) {
867 tol[1] = tol_cnv_energy;
870 for (
int ii : {0, 1}) {
871 if (std::isnan(res[ii])) {
873 if (this->terminal_output_) {
874 OpmLog::debug(
"NaN residual for " + this->compNames_.name(compIdx) +
" equation.");
877 else if (res[ii] > maxResidualAllowed()) {
878 report.setReservoirFailed({types[ii], CR::Severity::TooLarge, compIdx});
879 if (this->terminal_output_) {
880 OpmLog::debug(
"Too large residual for " + this->compNames_.name(compIdx) +
" equation.");
883 else if (res[ii] < 0.0) {
884 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
885 if (this->terminal_output_) {
886 OpmLog::debug(
"Negative residual for " + this->compNames_.name(compIdx) +
" equation.");
889 else if (res[ii] > tol[ii]) {
890 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
893 report.setReservoirConvergenceMetric(types[ii], compIdx, res[ii], tol[ii]);
898 if (this->terminal_output_) {
900 if (iteration == 0) {
901 std::string msg =
"Iter";
902 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
904 msg += this->compNames_.name(compIdx)[0];
908 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
910 msg += this->compNames_.name(compIdx)[0];
917 std::ostringstream ss;
918 const std::streamsize oprec = ss.precision(3);
919 const std::ios::fmtflags oflags = ss.setf(std::ios::scientific);
921 ss << std::setw(4) << iteration;
922 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
923 ss << std::setw(11) << mass_balance_residual[compIdx];
926 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
927 ss << std::setw(11) << CNV[compIdx];
933 OpmLog::debug(ss.str());
939template <
class TypeTag>
945 std::vector<Scalar>& residual_norms)
947 OPM_TIMEBLOCK(getConvergence);
949 std::vector<Scalar> B_avg(numEq, 0.0);
952 iteration, maxIter, B_avg, residual_norms);
954 OPM_TIMEBLOCK(getWellConvergence);
955 report += wellModel().getWellConvergence(B_avg,
959 conv_monitor_.checkPenaltyCard(report, iteration);
964template <
class TypeTag>
965std::vector<std::vector<typename BlackoilModel<TypeTag>::Scalar> >
969 OPM_TIMEBLOCK(computeFluidInPlace);
972 std::vector<std::vector<Scalar> > regionValues(0, std::vector<Scalar>(0,0.0));
976template <
class TypeTag>
981 if (!hasNlddSolver()) {
982 OPM_THROW(std::runtime_error,
"Cannot get local reports from a model without NLDD solver");
984 return nlddSolver_->localAccumulatedReports();
987template <
class TypeTag>
988const std::vector<SimulatorReport>&
993 OPM_THROW(std::runtime_error,
"Cannot get domain reports from a model without NLDD solver");
994 return nlddSolver_->domainAccumulatedReports();
997template <
class TypeTag>
1002 if (hasNlddSolver()) {
1003 nlddSolver_->writeNonlinearIterationsPerCell(odir);
1007template <
class TypeTag>
1012 if (hasNlddSolver()) {
1013 nlddSolver_->writePartitions(odir);
1017 const auto& elementMapper = this->simulator().model().elementMapper();
1018 const auto& cartMapper = this->simulator().vanguard().cartesianIndexMapper();
1020 const auto& grid = this->simulator().vanguard().grid();
1021 const auto& comm = grid.comm();
1022 const auto nDigit = 1 +
static_cast<int>(std::floor(std::log10(comm.size())));
1024 std::ofstream pfile {odir / fmt::format(
"{1:0>{0}}", nDigit, comm.rank())};
1026 for (
const auto& cell : elements(grid.leafGridView(), Dune::Partitions::interior)) {
1027 pfile << comm.rank() <<
' '
1028 << cartMapper.cartesianIndex(elementMapper.index(cell)) <<
' '
1029 << comm.rank() <<
'\n';
1033template <
class TypeTag>
1034template<
class Flu
idState,
class Res
idual>
1039 const FluidState& fs,
1040 const Residual& modelResid,
1042 std::vector<Scalar>& B_avg,
1043 std::vector<Scalar>& R_sum,
1044 std::vector<Scalar>& maxCoeff,
1045 std::vector<int>& maxCoeffCell)
1047 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1049 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1053 const unsigned sIdx = FluidSystem::solventComponentIndex(phaseIdx);
1054 const unsigned compIdx = Indices::canonicalToActiveComponentIndex(sIdx);
1056 B_avg[compIdx] += 1.0 / fs.invB(phaseIdx).value();
1057 const auto R2 = modelResid[cell_idx][compIdx];
1059 R_sum[compIdx] += R2;
1060 const Scalar Rval = std::abs(R2) / pvValue;
1061 if (Rval > maxCoeff[compIdx]) {
1062 maxCoeff[compIdx] = Rval;
1063 maxCoeffCell[compIdx] = cell_idx;
1067 if constexpr (has_solvent_) {
1068 B_avg[contiSolventEqIdx] +=
1069 1.0 / intQuants.solventInverseFormationVolumeFactor().value();
1070 const auto R2 = modelResid[cell_idx][contiSolventEqIdx];
1071 R_sum[contiSolventEqIdx] += R2;
1072 maxCoeff[contiSolventEqIdx] = std::max(maxCoeff[contiSolventEqIdx],
1073 std::abs(R2) / pvValue);
1075 if constexpr (has_extbo_) {
1076 B_avg[contiZfracEqIdx] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
1077 const auto R2 = modelResid[cell_idx][contiZfracEqIdx];
1078 R_sum[ contiZfracEqIdx ] += R2;
1079 maxCoeff[contiZfracEqIdx] = std::max(maxCoeff[contiZfracEqIdx],
1080 std::abs(R2) / pvValue);
1082 if constexpr (has_polymer_) {
1083 B_avg[contiPolymerEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1084 const auto R2 = modelResid[cell_idx][contiPolymerEqIdx];
1085 R_sum[contiPolymerEqIdx] += R2;
1086 maxCoeff[contiPolymerEqIdx] = std::max(maxCoeff[contiPolymerEqIdx],
1087 std::abs(R2) / pvValue);
1089 if constexpr (has_foam_) {
1090 B_avg[ contiFoamEqIdx ] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
1091 const auto R2 = modelResid[cell_idx][contiFoamEqIdx];
1092 R_sum[contiFoamEqIdx] += R2;
1093 maxCoeff[contiFoamEqIdx] = std::max(maxCoeff[contiFoamEqIdx],
1094 std::abs(R2) / pvValue);
1096 if constexpr (has_brine_) {
1097 B_avg[ contiBrineEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1098 const auto R2 = modelResid[cell_idx][contiBrineEqIdx];
1099 R_sum[contiBrineEqIdx] += R2;
1100 maxCoeff[contiBrineEqIdx] = std::max(maxCoeff[contiBrineEqIdx],
1101 std::abs(R2) / pvValue);
1104 if constexpr (has_polymermw_) {
1105 static_assert(has_polymer_);
1107 B_avg[contiPolymerMWEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1111 const auto R2 = modelResid[cell_idx][contiPolymerMWEqIdx] / 100.;
1112 R_sum[contiPolymerMWEqIdx] += R2;
1113 maxCoeff[contiPolymerMWEqIdx] = std::max(maxCoeff[contiPolymerMWEqIdx],
1114 std::abs(R2) / pvValue);
1117 if constexpr (has_energy_) {
1118 B_avg[contiEnergyEqIdx] += 1.0 / (4.182e1);
1119 const auto R2 = modelResid[cell_idx][contiEnergyEqIdx];
1120 R_sum[contiEnergyEqIdx] += R2;
1121 maxCoeff[contiEnergyEqIdx] = std::max(maxCoeff[contiEnergyEqIdx],
1122 std::abs(R2) / pvValue);
1125 if constexpr (has_micp_) {
1126 B_avg[contiMicrobialEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1127 const auto R1 = modelResid[cell_idx][contiMicrobialEqIdx];
1128 R_sum[contiMicrobialEqIdx] += R1;
1129 maxCoeff[contiMicrobialEqIdx] = std::max(maxCoeff[contiMicrobialEqIdx],
1130 std::abs(R1) / pvValue);
1131 B_avg[contiOxygenEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1132 const auto R2 = modelResid[cell_idx][contiOxygenEqIdx];
1133 R_sum[contiOxygenEqIdx] += R2;
1134 maxCoeff[contiOxygenEqIdx] = std::max(maxCoeff[contiOxygenEqIdx],
1135 std::abs(R2) / pvValue);
1136 B_avg[contiUreaEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1137 const auto R3 = modelResid[cell_idx][contiUreaEqIdx];
1138 R_sum[contiUreaEqIdx] += R3;
1139 maxCoeff[contiUreaEqIdx] = std::max(maxCoeff[contiUreaEqIdx],
1140 std::abs(R3) / pvValue);
1141 B_avg[contiBiofilmEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1142 const auto R4 = modelResid[cell_idx][contiBiofilmEqIdx];
1143 R_sum[contiBiofilmEqIdx] += R4;
1144 maxCoeff[contiBiofilmEqIdx] = std::max(maxCoeff[contiBiofilmEqIdx],
1145 std::abs(R4) / pvValue);
1146 B_avg[contiCalciteEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1147 const auto R5 = modelResid[cell_idx][contiCalciteEqIdx];
1148 R_sum[contiCalciteEqIdx] += R5;
1149 maxCoeff[contiCalciteEqIdx] = std::max(maxCoeff[contiCalciteEqIdx],
1150 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:53
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:245
ModelParameters param_
Definition: BlackoilModel.hpp:328
SimulatorReportSingle nonlinearIteration(const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Definition: BlackoilModel_impl.hpp:217
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:1000
SimulatorReportSingle assembleReservoir(const SimulatorTimerInterface &, const int iterationIdx)
Assemble the residual and Jacobian of the nonlinear system.
Definition: BlackoilModel_impl.hpp:360
ConvergenceReport getConvergence(const SimulatorTimerInterface &timer, const int iteration, const int maxIter, std::vector< Scalar > &residual_norms)
Definition: BlackoilModel_impl.hpp:942
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:744
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: BlackoilModel.hpp:68
const std::vector< SimulatorReport > & domainAccumulatedReports() const
return the statistics of local solves accumulated for each domain on this rank
Definition: BlackoilModel_impl.hpp:990
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:607
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilModel.hpp:67
std::vector< StepReport > convergence_reports_
Definition: BlackoilModel.hpp:343
const SimulatorReport & localAccumulatedReports() const
return the statistics of local solves accumulated for this rank
Definition: BlackoilModel_impl.hpp:979
SimulatorReportSingle afterStep(const SimulatorTimerInterface &)
Definition: BlackoilModel_impl.hpp:347
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:547
const Grid & grid_
Definition: BlackoilModel.hpp:317
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: BlackoilModel.hpp:70
void initialLinearization(SimulatorReportSingle &report, const int iteration, const int minIter, const int maxIter, const SimulatorTimerInterface &timer)
Definition: BlackoilModel_impl.hpp:158
void updateSolution(const BVector &dx)
Apply an update to the primary variables.
Definition: BlackoilModel_impl.hpp:523
long int global_nc_
The number of cells of the global grid.
Definition: BlackoilModel.hpp:337
void updateTUNING(const Tuning &tuning)
Definition: BlackoilModel_impl.hpp:730
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:1037
std::unique_ptr< BlackoilModelNldd< TypeTag > > nlddSolver_
Non-linear DD solver.
Definition: BlackoilModel.hpp:346
SimulatorReportSingle nonlinearIterationNewton(const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Definition: BlackoilModel_impl.hpp:255
void writePartitions(const std::filesystem::path &odir) const
Definition: BlackoilModel_impl.hpp:1010
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:658
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilModel.hpp:65
void solveJacobianSystem(BVector &x)
Definition: BlackoilModel_impl.hpp:461
SimulatorReportSingle prepareStep(const SimulatorTimerInterface &timer)
Definition: BlackoilModel_impl.hpp:89
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilModel.hpp:108
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilModel.hpp:76
Scalar relativeChange() const
Definition: BlackoilModel_impl.hpp:374
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:82
Definition: blackoilboundaryratevector.hh:39
PhaseUsage phaseUsageFromDeck(const EclipseState &eclipseState)
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:177
std::string nonlinear_solver_
Nonlinear solver type: newton or nldd.
Definition: BlackoilModelParameters.hpp:327
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