24#ifndef OPM_BLACKOILMODELEBOS_HEADER_INCLUDED
25#define OPM_BLACKOILMODELEBOS_HEADER_INCLUDED
27#include <ebos/eclproblem.hh>
28#include <opm/models/utils/start.hh>
30#include <opm/common/ErrorMacros.hpp>
31#include <opm/common/Exceptions.hpp>
32#include <opm/common/OpmLog/OpmLog.hpp>
34#include <opm/core/props/phaseUsageFromDeck.hpp>
36#include <opm/grid/UnstructuredGrid.h>
38#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
39#include <opm/input/eclipse/EclipseState/Tables/TableManager.hpp>
41#include <opm/simulators/aquifers/BlackoilAquiferModel.hpp>
42#include <opm/simulators/flow/countGlobalCells.hpp>
43#include <opm/simulators/flow/NonlinearSolverEbos.hpp>
44#include <opm/simulators/flow/BlackoilModelParametersEbos.hpp>
45#include <opm/simulators/linalg/ISTLSolverEbos.hpp>
46#include <opm/simulators/timestepping/AdaptiveTimeSteppingEbos.hpp>
47#include <opm/simulators/timestepping/ConvergenceReport.hpp>
48#include <opm/simulators/timestepping/SimulatorReport.hpp>
49#include <opm/simulators/timestepping/SimulatorTimer.hpp>
51#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
52#include <opm/simulators/wells/BlackoilWellModel.hpp>
53#include <opm/simulators/wells/WellConnectionAuxiliaryModule.hpp>
55#include <dune/istl/owneroverlapcopy.hh>
56#include <dune/common/parallel/communication.hh>
57#include <dune/common/timer.hh>
58#include <dune/common/unused.hh>
69namespace Opm::Properties {
77template<
class TypeTag>
78struct OutputDir<TypeTag, TTag::EclFlowProblem> {
79 static constexpr auto value =
"";
81template<
class TypeTag>
82struct EnableDebuggingChecks<TypeTag, TTag::EclFlowProblem> {
83 static constexpr bool value =
false;
86template<
class TypeTag>
87struct BlackoilConserveSurfaceVolume<TypeTag, TTag::EclFlowProblem> {
88 static constexpr bool value =
true;
90template<
class TypeTag>
91struct UseVolumetricResidual<TypeTag, TTag::EclFlowProblem> {
92 static constexpr bool value =
false;
95template<
class TypeTag>
96struct EclAquiferModel<TypeTag, TTag::EclFlowProblem> {
102template<
class TypeTag>
103struct EnablePolymer<TypeTag, TTag::EclFlowProblem> {
104 static constexpr bool value =
false;
106template<
class TypeTag>
107struct EnableSolvent<TypeTag, TTag::EclFlowProblem> {
108 static constexpr bool value =
false;
110template<
class TypeTag>
111struct EnableTemperature<TypeTag, TTag::EclFlowProblem> {
112 static constexpr bool value =
true;
114template<
class TypeTag>
115struct EnableEnergy<TypeTag, TTag::EclFlowProblem> {
116 static constexpr bool value =
false;
118template<
class TypeTag>
119struct EnableFoam<TypeTag, TTag::EclFlowProblem> {
120 static constexpr bool value =
false;
122template<
class TypeTag>
123struct EnableBrine<TypeTag, TTag::EclFlowProblem> {
124 static constexpr bool value =
false;
126template<
class TypeTag>
127struct EnableSaltPrecipitation<TypeTag, TTag::EclFlowProblem> {
128 static constexpr bool value =
false;
130template<
class TypeTag>
131struct EnableMICP<TypeTag, TTag::EclFlowProblem> {
132 static constexpr bool value =
false;
135template<
class TypeTag>
139template<
class TypeTag>
140struct LinearSolverSplice<TypeTag, TTag::EclFlowProblem> {
153 template <
class TypeTag>
160 using Simulator = GetPropType<TypeTag, Properties::Simulator>;
161 using Grid = GetPropType<TypeTag, Properties::Grid>;
162 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
163 using SparseMatrixAdapter = GetPropType<TypeTag, Properties::SparseMatrixAdapter>;
164 using SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>;
165 using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
166 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
167 using Indices = GetPropType<TypeTag, Properties::Indices>;
168 using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
169 using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
171 typedef double Scalar;
172 static const int numEq = Indices::numEq;
173 static const int contiSolventEqIdx = Indices::contiSolventEqIdx;
174 static const int contiZfracEqIdx = Indices::contiZfracEqIdx;
175 static const int contiPolymerEqIdx = Indices::contiPolymerEqIdx;
176 static const int contiEnergyEqIdx = Indices::contiEnergyEqIdx;
177 static const int contiPolymerMWEqIdx = Indices::contiPolymerMWEqIdx;
178 static const int contiFoamEqIdx = Indices::contiFoamEqIdx;
179 static const int contiBrineEqIdx = Indices::contiBrineEqIdx;
180 static const int contiMicrobialEqIdx = Indices::contiMicrobialEqIdx;
181 static const int contiOxygenEqIdx = Indices::contiOxygenEqIdx;
182 static const int contiUreaEqIdx = Indices::contiUreaEqIdx;
183 static const int contiBiofilmEqIdx = Indices::contiBiofilmEqIdx;
184 static const int contiCalciteEqIdx = Indices::contiCalciteEqIdx;
185 static const int solventSaturationIdx = Indices::solventSaturationIdx;
186 static const int zFractionIdx = Indices::zFractionIdx;
187 static const int polymerConcentrationIdx = Indices::polymerConcentrationIdx;
188 static const int polymerMoleWeightIdx = Indices::polymerMoleWeightIdx;
189 static const int temperatureIdx = Indices::temperatureIdx;
190 static const int foamConcentrationIdx = Indices::foamConcentrationIdx;
191 static const int saltConcentrationIdx = Indices::saltConcentrationIdx;
192 static const int microbialConcentrationIdx = Indices::microbialConcentrationIdx;
193 static const int oxygenConcentrationIdx = Indices::oxygenConcentrationIdx;
194 static const int ureaConcentrationIdx = Indices::ureaConcentrationIdx;
195 static const int biofilmConcentrationIdx = Indices::biofilmConcentrationIdx;
196 static const int calciteConcentrationIdx = Indices::calciteConcentrationIdx;
198 typedef Dune::FieldVector<Scalar, numEq > VectorBlockType;
199 typedef typename SparseMatrixAdapter::MatrixBlock MatrixBlockType;
200 typedef typename SparseMatrixAdapter::IstlMatrix Mat;
201 typedef Dune::BlockVector<VectorBlockType> BVector;
211 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
212 if (!FluidSystem::phaseIsActive(phaseIdx)) {
216 const unsigned canonicalCompIdx = FluidSystem::solventComponentIndex(phaseIdx);
217 names_[Indices::canonicalToActiveComponentIndex(canonicalCompIdx)]
218 = FluidSystem::componentName(canonicalCompIdx);
221 if constexpr (has_solvent_) {
222 names_[solventSaturationIdx] =
"Solvent";
225 if constexpr (has_extbo_) {
226 names_[zFractionIdx] =
"ZFraction";
229 if constexpr (has_polymer_) {
230 names_[polymerConcentrationIdx] =
"Polymer";
233 if constexpr (has_polymermw_) {
234 assert(has_polymer_);
235 names_[polymerMoleWeightIdx] =
"MolecularWeightP";
238 if constexpr (has_energy_) {
239 names_[temperatureIdx] =
"Energy";
242 if constexpr (has_foam_) {
243 names_[foamConcentrationIdx] =
"Foam";
246 if constexpr (has_brine_) {
247 names_[saltConcentrationIdx] =
"Brine";
250 if constexpr (has_micp_) {
251 names_[microbialConcentrationIdx] =
"Microbes";
252 names_[oxygenConcentrationIdx] =
"Oxygen";
253 names_[ureaConcentrationIdx] =
"Urea";
254 names_[biofilmConcentrationIdx] =
"Biofilm";
255 names_[calciteConcentrationIdx] =
"Calcite";
259 const std::string& name(
const int compIdx)
const
261 return this->names_[compIdx];
265 std::vector<std::string> names_{};
285 const bool terminal_output)
286 : ebosSimulator_(ebosSimulator)
287 , grid_(ebosSimulator_.vanguard().grid())
290 , well_model_ (well_model)
292 , current_relaxation_(1.0)
293 , dx_old_(ebosSimulator_.model().numGridDof())
297 convergence_reports_.reserve(300);
300 bool isParallel()
const
301 {
return grid_.comm().size() > 1; }
303 const EclipseState& eclState()
const
304 {
return ebosSimulator_.vanguard().eclState(); }
311 Dune::Timer perfTimer;
315 ebosSimulator_.model().updateFailed();
317 ebosSimulator_.model().advanceTimeLevel();
326 ebosSimulator_.model().newtonMethod().setIterationIndex(0);
328 ebosSimulator_.problem().beginTimeStep();
330 unsigned numDof = ebosSimulator_.model().numGridDof();
331 wasSwitched_.resize(numDof);
332 std::fill(wasSwitched_.begin(), wasSwitched_.end(),
false);
335 std::cout <<
"equation scaling not supported yet" << std::endl;
338 report.pre_post_time += perfTimer.stop();
353 template <
class NonlinearSolverType>
356 NonlinearSolverType& nonlinear_solver)
360 Dune::Timer perfTimer;
363 if (iteration == 0) {
366 residual_norms_history_.clear();
367 current_relaxation_ = 1.0;
370 convergence_reports_.back().report.reserve(11);
373 report.total_linearizations = 1;
377 report.assemble_time += perfTimer.stop();
380 report.assemble_time += perfTimer.stop();
381 failureReport_ += report;
386 std::vector<double> residual_norms;
392 report.converged = convrep.converged() && iteration > nonlinear_solver.minIter();;
393 ConvergenceReport::Severity severity = convrep.severityOfWorstFailure();
394 convergence_reports_.back().report.push_back(std::move(convrep));
397 if (severity == ConvergenceReport::Severity::NotANumber) {
398 OPM_THROW(NumericalProblem,
"NaN residual found!");
399 }
else if (severity == ConvergenceReport::Severity::TooLarge) {
400 OPM_THROW_NOLOG(NumericalProblem,
"Too large residual found!");
403 report.update_time += perfTimer.stop();
404 residual_norms_history_.push_back(residual_norms);
405 if (!report.converged) {
408 report.total_newton_iterations = 1;
414 unsigned nc = ebosSimulator_.model().numGridDof();
418 linear_solve_setup_time_ = 0.0;
423 wellModel().linearize(ebosSimulator().model().linearizer().jacobian(),
424 ebosSimulator().model().linearizer().residual());
427 report.linear_solve_setup_time += linear_solve_setup_time_;
428 report.linear_solve_time += perfTimer.stop();
432 report.linear_solve_setup_time += linear_solve_setup_time_;
433 report.linear_solve_time += perfTimer.stop();
436 failureReport_ += report;
450 bool isOscillate =
false;
451 bool isStagnate =
false;
452 nonlinear_solver.detectOscillations(residual_norms_history_, iteration, isOscillate, isStagnate);
454 current_relaxation_ -= nonlinear_solver.relaxIncrement();
455 current_relaxation_ = std::max(current_relaxation_, nonlinear_solver.relaxMax());
457 std::string msg =
" Oscillating behavior detected: Relaxation set to "
458 + std::to_string(current_relaxation_);
462 nonlinear_solver.stabilizeNonlinearUpdate(x, dx_old_, current_relaxation_);
469 report.update_time += perfTimer.stop();
475 void printIf(
int c,
double x,
double y,
double eps, std::string type) {
476 if (std::abs(x-y) > eps) {
477 std::cout << type <<
" " <<c <<
": "<<x <<
" " << y << std::endl;
488 Dune::Timer perfTimer;
490 ebosSimulator_.problem().endTimeStep();
491 report.pre_post_time += perfTimer.stop();
500 const int iterationIdx)
503 ebosSimulator_.model().newtonMethod().setIterationIndex(iterationIdx);
504 ebosSimulator_.problem().beginIteration();
505 ebosSimulator_.model().linearizer().linearizeDomain();
506 ebosSimulator_.problem().endIteration();
512 double relativeChange()
const
514 Scalar resultDelta = 0.0;
515 Scalar resultDenom = 0.0;
517 const auto& elemMapper = ebosSimulator_.model().elementMapper();
518 const auto& gridView = ebosSimulator_.gridView();
519 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
520 unsigned globalElemIdx = elemMapper.index(elem);
521 const auto& priVarsNew = ebosSimulator_.model().solution(0)[globalElemIdx];
524 pressureNew = priVarsNew[Indices::pressureSwitchIdx];
526 Scalar saturationsNew[FluidSystem::numPhases] = { 0.0 };
527 Scalar oilSaturationNew = 1.0;
528 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) &&
529 FluidSystem::numActivePhases() > 1 &&
530 priVarsNew.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw) {
531 saturationsNew[FluidSystem::waterPhaseIdx] = priVarsNew[Indices::waterSwitchIdx];
532 oilSaturationNew -= saturationsNew[FluidSystem::waterPhaseIdx];
535 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
536 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
537 priVarsNew.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg) {
538 assert(Indices::compositionSwitchIdx >= 0 );
539 saturationsNew[FluidSystem::gasPhaseIdx] = priVarsNew[Indices::compositionSwitchIdx];
540 oilSaturationNew -= saturationsNew[FluidSystem::gasPhaseIdx];
543 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
544 saturationsNew[FluidSystem::oilPhaseIdx] = oilSaturationNew;
547 const auto& priVarsOld = ebosSimulator_.model().solution(1)[globalElemIdx];
550 pressureOld = priVarsOld[Indices::pressureSwitchIdx];
552 Scalar saturationsOld[FluidSystem::numPhases] = { 0.0 };
553 Scalar oilSaturationOld = 1.0;
556 Scalar tmp = pressureNew - pressureOld;
557 resultDelta += tmp*tmp;
558 resultDenom += pressureNew*pressureNew;
560 if (FluidSystem::numActivePhases() > 1) {
561 if (priVarsOld.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw) {
562 saturationsOld[FluidSystem::waterPhaseIdx] = priVarsOld[Indices::waterSwitchIdx];
563 oilSaturationOld -= saturationsOld[FluidSystem::waterPhaseIdx];
566 if (priVarsOld.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
568 assert(Indices::compositionSwitchIdx >= 0 );
569 saturationsOld[FluidSystem::gasPhaseIdx] = priVarsOld[Indices::compositionSwitchIdx];
570 oilSaturationOld -= saturationsOld[FluidSystem::gasPhaseIdx];
573 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
574 saturationsOld[FluidSystem::oilPhaseIdx] = oilSaturationOld;
576 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++ phaseIdx) {
577 Scalar tmpSat = saturationsNew[phaseIdx] - saturationsOld[phaseIdx];
578 resultDelta += tmpSat*tmpSat;
579 resultDenom += saturationsNew[phaseIdx]*saturationsNew[phaseIdx];
580 assert(std::isfinite(resultDelta));
581 assert(std::isfinite(resultDenom));
586 resultDelta = gridView.comm().sum(resultDelta);
587 resultDenom = gridView.comm().sum(resultDenom);
589 if (resultDenom > 0.0)
590 return resultDelta/resultDenom;
598 return ebosSimulator_.model().newtonMethod().linearSolver().iterations ();
606 auto& ebosJac = ebosSimulator_.model().linearizer().jacobian();
607 auto& ebosResid = ebosSimulator_.model().linearizer().residual();
612 auto& ebosSolver = ebosSimulator_.model().newtonMethod().linearSolver();
613 Dune::Timer perfTimer;
615 ebosSolver.prepare(ebosJac, ebosResid);
616 linear_solve_setup_time_ = perfTimer.stop();
617 ebosSolver.setResidual(ebosResid);
622 ebosSolver.setMatrix(ebosJac);
632 auto& ebosNewtonMethod = ebosSimulator_.model().newtonMethod();
633 SolutionVector& solution = ebosSimulator_.model().solution(0);
635 ebosNewtonMethod.update_(solution,
644 OPM_TIMEBLOCK(invalidateAndUpdateIntensiveQuantities);
645 ebosSimulator_.model().invalidateAndUpdateIntensiveQuantities(0);
646 ebosSimulator_.problem().eclWriter()->mutableEclOutputModule().invalidateLocalData();
656 template <
class CollectiveCommunication>
657 std::tuple<double,double> convergenceReduction(
const CollectiveCommunication& comm,
658 const double pvSumLocal,
659 const double numAquiferPvSumLocal,
660 std::vector< Scalar >& R_sum,
661 std::vector< Scalar >& maxCoeff,
662 std::vector< Scalar >& B_avg)
664 OPM_TIMEBLOCK(convergenceReduction);
666 double pvSum = pvSumLocal;
667 double numAquiferPvSum = numAquiferPvSumLocal;
669 if( comm.size() > 1 )
672 std::vector< Scalar > sumBuffer;
673 std::vector< Scalar > maxBuffer;
674 const int numComp = B_avg.size();
675 sumBuffer.reserve( 2*numComp + 2 );
676 maxBuffer.reserve( numComp );
677 for(
int compIdx = 0; compIdx < numComp; ++compIdx )
679 sumBuffer.push_back( B_avg[ compIdx ] );
680 sumBuffer.push_back( R_sum[ compIdx ] );
681 maxBuffer.push_back( maxCoeff[ compIdx ] );
685 sumBuffer.push_back( pvSum );
686 sumBuffer.push_back( numAquiferPvSum );
689 comm.sum( sumBuffer.data(), sumBuffer.size() );
692 comm.max( maxBuffer.data(), maxBuffer.size() );
695 for(
int compIdx = 0, buffIdx = 0; compIdx < numComp; ++compIdx, ++buffIdx )
697 B_avg[ compIdx ] = sumBuffer[ buffIdx ];
700 R_sum[ compIdx ] = sumBuffer[ buffIdx ];
703 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};
721 std::vector<Scalar>& maxCoeff,
722 std::vector<Scalar>& B_avg)
725 double pvSumLocal = 0.0;
726 double numAquiferPvSumLocal = 0.0;
727 const auto& ebosModel = ebosSimulator_.model();
728 const auto& ebosProblem = ebosSimulator_.problem();
730 const auto& ebosResid = ebosSimulator_.model().linearizer().residual();
732 ElementContext elemCtx(ebosSimulator_);
733 const auto& gridView = ebosSimulator().gridView();
734 OPM_BEGIN_PARALLEL_TRY_CATCH();
735 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
736 elemCtx.updatePrimaryStencil(elem);
737 elemCtx.updatePrimaryIntensiveQuantities(0);
738 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
739 const auto& intQuants = elemCtx.intensiveQuantities(0, 0);
740 const auto& fs = intQuants.fluidState();
742 const double pvValue = ebosProblem.referencePorosity(cell_idx, 0) * ebosModel.dofTotalVolume( cell_idx );
743 pvSumLocal += pvValue;
745 if (isNumericalAquiferCell(gridView.grid(), elem))
747 numAquiferPvSumLocal += pvValue;
750 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
752 if (!FluidSystem::phaseIsActive(phaseIdx)) {
756 const unsigned compIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
758 B_avg[ compIdx ] += 1.0 / fs.invB(phaseIdx).value();
759 const auto R2 = ebosResid[cell_idx][compIdx];
761 R_sum[ compIdx ] += R2;
762 maxCoeff[ compIdx ] = std::max( maxCoeff[ compIdx ], std::abs( R2 ) / pvValue );
765 if constexpr (has_solvent_) {
766 B_avg[ contiSolventEqIdx ] += 1.0 / intQuants.solventInverseFormationVolumeFactor().value();
767 const auto R2 = ebosResid[cell_idx][contiSolventEqIdx];
768 R_sum[ contiSolventEqIdx ] += R2;
769 maxCoeff[ contiSolventEqIdx ] = std::max( maxCoeff[ contiSolventEqIdx ], std::abs( R2 ) / pvValue );
771 if constexpr (has_extbo_) {
772 B_avg[ contiZfracEqIdx ] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
773 const auto R2 = ebosResid[cell_idx][contiZfracEqIdx];
774 R_sum[ contiZfracEqIdx ] += R2;
775 maxCoeff[ contiZfracEqIdx ] = std::max( maxCoeff[ contiZfracEqIdx ], std::abs( R2 ) / pvValue );
777 if constexpr (has_polymer_) {
778 B_avg[ contiPolymerEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
779 const auto R2 = ebosResid[cell_idx][contiPolymerEqIdx];
780 R_sum[ contiPolymerEqIdx ] += R2;
781 maxCoeff[ contiPolymerEqIdx ] = std::max( maxCoeff[ contiPolymerEqIdx ], std::abs( R2 ) / pvValue );
783 if constexpr (has_foam_) {
784 B_avg[ contiFoamEqIdx ] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
785 const auto R2 = ebosResid[cell_idx][contiFoamEqIdx];
786 R_sum[ contiFoamEqIdx ] += R2;
787 maxCoeff[ contiFoamEqIdx ] = std::max( maxCoeff[ contiFoamEqIdx ], std::abs( R2 ) / pvValue );
789 if constexpr (has_brine_) {
790 B_avg[ contiBrineEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
791 const auto R2 = ebosResid[cell_idx][contiBrineEqIdx];
792 R_sum[ contiBrineEqIdx ] += R2;
793 maxCoeff[ contiBrineEqIdx ] = std::max( maxCoeff[ contiBrineEqIdx ], std::abs( R2 ) / pvValue );
796 if constexpr (has_polymermw_) {
797 static_assert(has_polymer_);
799 B_avg[contiPolymerMWEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
803 const auto R2 = ebosResid[cell_idx][contiPolymerMWEqIdx] / 100.;
804 R_sum[contiPolymerMWEqIdx] += R2;
805 maxCoeff[contiPolymerMWEqIdx] = std::max( maxCoeff[contiPolymerMWEqIdx], std::abs( R2 ) / pvValue );
808 if constexpr (has_energy_) {
809 B_avg[ contiEnergyEqIdx ] += 1.0 / (4.182e1);
810 const auto R2 = ebosResid[cell_idx][contiEnergyEqIdx];
811 R_sum[ contiEnergyEqIdx ] += R2;
812 maxCoeff[ contiEnergyEqIdx ] = std::max( maxCoeff[ contiEnergyEqIdx ], std::abs( R2 ) / pvValue );
815 if constexpr (has_micp_) {
816 B_avg[ contiMicrobialEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
817 const auto R1 = ebosResid[cell_idx][contiMicrobialEqIdx];
818 R_sum[ contiMicrobialEqIdx ] += R1;
819 maxCoeff[ contiMicrobialEqIdx ] = std::max( maxCoeff[ contiMicrobialEqIdx ], std::abs( R1 ) / pvValue );
820 B_avg[ contiOxygenEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
821 const auto R2 = ebosResid[cell_idx][contiOxygenEqIdx];
822 R_sum[ contiOxygenEqIdx ] += R2;
823 maxCoeff[ contiOxygenEqIdx ] = std::max( maxCoeff[ contiOxygenEqIdx ], std::abs( R2 ) / pvValue );
824 B_avg[ contiUreaEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
825 const auto R3 = ebosResid[cell_idx][contiUreaEqIdx];
826 R_sum[ contiUreaEqIdx ] += R3;
827 maxCoeff[ contiUreaEqIdx ] = std::max( maxCoeff[ contiUreaEqIdx ], std::abs( R3 ) / pvValue );
828 B_avg[ contiBiofilmEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
829 const auto R4 = ebosResid[cell_idx][contiBiofilmEqIdx];
830 R_sum[ contiBiofilmEqIdx ] += R4;
831 maxCoeff[ contiBiofilmEqIdx ] = std::max( maxCoeff[ contiBiofilmEqIdx ], std::abs( R4 ) / pvValue );
832 B_avg[ contiCalciteEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
833 const auto R5 = ebosResid[cell_idx][contiCalciteEqIdx];
834 R_sum[ contiCalciteEqIdx ] += R5;
835 maxCoeff[ contiCalciteEqIdx ] = std::max( maxCoeff[ contiCalciteEqIdx ], std::abs( R5 ) / pvValue );
839 OPM_END_PARALLEL_TRY_CATCH(
"BlackoilModelEbos::localConvergenceData() failed: ", grid_.comm());
842 const int bSize = B_avg.size();
843 for (
int i = 0; i<bSize; ++i )
848 return {pvSumLocal, numAquiferPvSumLocal};
857 const auto& ebosModel = ebosSimulator_.model();
858 const auto& ebosProblem = ebosSimulator_.problem();
859 const auto& ebosResid = ebosSimulator_.model().linearizer().residual();
860 const auto& gridView = ebosSimulator().gridView();
861 ElementContext elemCtx(ebosSimulator_);
863 OPM_BEGIN_PARALLEL_TRY_CATCH();
865 for (
const auto& elem : elements(gridView, Dune::Partitions::interiorBorder))
868 if (isNumericalAquiferCell(gridView.grid(), elem))
872 elemCtx.updatePrimaryStencil(elem);
873 elemCtx.updatePrimaryIntensiveQuantities(0);
874 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
875 const double pvValue = ebosProblem.referencePorosity(cell_idx, 0) * ebosModel.dofTotalVolume( cell_idx );
876 const auto& cellResidual = ebosResid[cell_idx];
877 bool cnvViolated =
false;
879 for (
unsigned eqIdx = 0; eqIdx < cellResidual.size(); ++eqIdx)
882 Scalar CNV = cellResidual[eqIdx] * dt * B_avg[eqIdx] / pvValue;
892 OPM_END_PARALLEL_TRY_CATCH(
"BlackoilModelEbos::ComputeCnvError() failed: ", grid_.comm());
894 return grid_.comm().sum(errorPV);
900 std::vector<Scalar>& B_avg,
901 std::vector<Scalar>& residual_norms)
903 OPM_TIMEBLOCK(getReservoirConvergence);
904 typedef std::vector< Scalar > Vector;
906 const int numComp = numEq;
907 Vector R_sum(numComp, 0.0 );
908 Vector maxCoeff(numComp, std::numeric_limits< Scalar >::lowest() );
912 const auto [ pvSum, numAquiferPvSum ] =
913 convergenceReduction(grid_.comm(), pvSumLocal,
914 numAquiferPvSumLocal,
915 R_sum, maxCoeff, B_avg);
918 cnvErrorPvFraction /= (pvSum - numAquiferPvSum);
924 const bool use_relaxed = cnvErrorPvFraction < param_.relaxed_max_pv_fraction_ && iteration >= param_.
min_strict_cnv_iter_;
928 std::vector<Scalar> CNV(numComp);
929 std::vector<Scalar> mass_balance_residual(numComp);
930 for (
int compIdx = 0; compIdx < numComp; ++compIdx )
932 CNV[compIdx] = B_avg[compIdx] * dt * maxCoeff[compIdx];
933 mass_balance_residual[compIdx] = std::abs(B_avg[compIdx]*R_sum[compIdx]) * dt / pvSum;
934 residual_norms.push_back(CNV[compIdx]);
938 ConvergenceReport report{reportTime};
939 using CR = ConvergenceReport;
940 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
941 double res[2] = { mass_balance_residual[compIdx], CNV[compIdx] };
942 CR::ReservoirFailure::Type types[2] = { CR::ReservoirFailure::Type::MassBalance,
943 CR::ReservoirFailure::Type::Cnv };
944 double tol[2] = { tol_mb, tol_cnv };
945 for (
int ii : {0, 1}) {
946 if (std::isnan(res[ii])) {
947 report.setReservoirFailed({types[ii], CR::Severity::NotANumber, compIdx});
949 OpmLog::debug(
"NaN residual for " + this->compNames_.name(compIdx) +
" equation.");
951 }
else if (res[ii] > maxResidualAllowed()) {
952 report.setReservoirFailed({types[ii], CR::Severity::TooLarge, compIdx});
954 OpmLog::debug(
"Too large residual for " + this->compNames_.name(compIdx) +
" equation.");
956 }
else if (res[ii] < 0.0) {
957 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
959 OpmLog::debug(
"Negative residual for " + this->compNames_.name(compIdx) +
" equation.");
961 }
else if (res[ii] > tol[ii]) {
962 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
964 report.setReservoirConvergenceMetric(types[ii], compIdx, res[ii]);
972 if (iteration == 0) {
973 std::string msg =
"Iter";
974 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
976 msg += this->compNames_.name(compIdx)[0];
979 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
981 msg += this->compNames_.name(compIdx)[0];
986 std::ostringstream ss;
987 const std::streamsize oprec = ss.precision(3);
988 const std::ios::fmtflags oflags = ss.setf(std::ios::scientific);
989 ss << std::setw(4) << iteration;
990 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
991 ss << std::setw(11) << mass_balance_residual[compIdx];
993 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
994 ss << std::setw(11) << CNV[compIdx];
998 OpmLog::debug(ss.str());
1010 const int iteration,
1011 std::vector<double>& residual_norms)
1015 std::vector<Scalar> B_avg(numEq, 0.0);
1018 iteration, B_avg, residual_norms);
1020 OPM_TIMEBLOCK(getWellConvergence);
1021 report +=
wellModel().getWellConvergence(B_avg, report.converged());
1030 return phaseUsage_.num_phases;
1035 std::vector<std::vector<double> >
1042 std::vector<std::vector<double> >
1048 std::vector<std::vector<double> > regionValues(0, std::vector<double>(0,0.0));
1049 return regionValues;
1052 const Simulator& ebosSimulator()
const
1053 {
return ebosSimulator_; }
1055 Simulator& ebosSimulator()
1056 {
return ebosSimulator_; }
1060 {
return failureReport_; }
1062 const std::vector<StepReport>& stepReports()
const
1064 return convergence_reports_;
1067 std::vector<StepReport> getStepReportsDestructively()
const
1069 return std::move(this->convergence_reports_);
1075 Simulator& ebosSimulator_;
1077 const PhaseUsage phaseUsage_;
1078 static constexpr bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>();
1079 static constexpr bool has_extbo_ = getPropValue<TypeTag, Properties::EnableExtbo>();
1080 static constexpr bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>();
1081 static constexpr bool has_polymermw_ = getPropValue<TypeTag, Properties::EnablePolymerMW>();
1082 static constexpr bool has_energy_ = getPropValue<TypeTag, Properties::EnableEnergy>();
1083 static constexpr bool has_foam_ = getPropValue<TypeTag, Properties::EnableFoam>();
1084 static constexpr bool has_brine_ = getPropValue<TypeTag, Properties::EnableBrine>();
1085 static constexpr bool has_micp_ = getPropValue<TypeTag, Properties::EnableMICP>();
1087 ModelParameters param_;
1088 SimulatorReportSingle failureReport_;
1091 BlackoilWellModel<TypeTag>& well_model_;
1098 std::vector<std::vector<double>> residual_norms_history_;
1099 double current_relaxation_;
1102 std::vector<StepReport> convergence_reports_;
1111 wellModel()
const {
return well_model_; }
1113 void beginReportStep()
1115 ebosSimulator_.problem().beginEpisode();
1118 void endReportStep()
1120 ebosSimulator_.problem().endEpisode();
1125 bool isNumericalAquiferCell(
const Dune::CpGrid& grid,
const T& elem)
1127 const auto& aquiferCells = grid.sortedNumAquiferCells();
1128 if (aquiferCells.empty())
1132 auto candidate = std::lower_bound(aquiferCells.begin(), aquiferCells.end(),
1134 return candidate != aquiferCells.end() && *candidate == elem.index();
1137 template<
class G,
class T>
1138 typename std::enable_if<!std::is_same<G,Dune::CpGrid>::value,
bool>::type
1139 isNumericalAquiferCell(
const G&,
const T&)
1144 double dpMaxRel()
const {
return param_.dp_max_rel_; }
1145 double dsMax()
const {
return param_.ds_max_; }
1146 double drMaxRel()
const {
return param_.dr_max_rel_; }
1147 double maxResidualAllowed()
const {
return param_.max_residual_allowed_; }
1148 double linear_solve_setup_time_;
1151 std::vector<bool> wasSwitched_;
Class for handling the blackoil well model.
Definition: BlackoilAquiferModel.hpp:86
Definition: BlackoilModelEbos.hpp:206
A model implementation for three-phase black oil.
Definition: BlackoilModelEbos.hpp:155
BlackoilWellModel< TypeTag > & wellModel()
return the StandardWells object
Definition: BlackoilModelEbos.hpp:1108
std::vector< std::vector< double > > computeFluidInPlace(const std::vector< int > &) const
Should not be called.
Definition: BlackoilModelEbos.hpp:1043
int linearIterationsLastSolve() const
Number of linear iterations used in last call to solveJacobianSystem().
Definition: BlackoilModelEbos.hpp:596
bool terminalOutputEnabled() const
Return true if output to cout is wanted.
Definition: BlackoilModelEbos.hpp:651
const SimulatorReportSingle & failureReport() const
return the statistics if the nonlinearIteration() method failed
Definition: BlackoilModelEbos.hpp:1059
ConvergenceReport getConvergence(const SimulatorTimerInterface &timer, const int iteration, std::vector< double > &residual_norms)
Compute convergence based on total mass balance (tol_mb) and maximum residual mass balance (tol_cnv).
Definition: BlackoilModelEbos.hpp:1009
BlackoilModelEbos(Simulator &ebosSimulator, const ModelParameters ¶m, BlackoilWellModel< TypeTag > &well_model, const bool terminal_output)
Construct the model.
Definition: BlackoilModelEbos.hpp:282
bool terminal_output_
Whether we print something to std::cout.
Definition: BlackoilModelEbos.hpp:1094
SimulatorReportSingle afterStep(const SimulatorTimerInterface &)
Called once after each time step.
Definition: BlackoilModelEbos.hpp:485
std::tuple< double, double > localConvergenceData(std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< Scalar > &B_avg)
Get reservoir quantities on this process needed for convergence calculations.
Definition: BlackoilModelEbos.hpp:720
SimulatorReportSingle nonlinearIteration(const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Called once per nonlinear iteration.
Definition: BlackoilModelEbos.hpp:354
SimulatorReportSingle assembleReservoir(const SimulatorTimerInterface &, const int iterationIdx)
Assemble the residual and Jacobian of the nonlinear system.
Definition: BlackoilModelEbos.hpp:499
SimulatorReportSingle prepareStep(const SimulatorTimerInterface &timer)
Called once before each time step.
Definition: BlackoilModelEbos.hpp:308
long int global_nc_
The number of cells of the global grid.
Definition: BlackoilModelEbos.hpp:1096
double computeCnvErrorPv(const std::vector< Scalar > &B_avg, double dt)
Compute the total pore volume of cells violating CNV that are not part of a numerical aquifer.
Definition: BlackoilModelEbos.hpp:853
int numPhases() const
The number of active fluid phases in the model.
Definition: BlackoilModelEbos.hpp:1028
std::vector< std::vector< double > > computeFluidInPlace(const T &, const std::vector< int > &fipnum) const
Wrapper required due to not following generic API.
Definition: BlackoilModelEbos.hpp:1036
void solveJacobianSystem(BVector &x)
Solve the Jacobian system Jx = r where J is the Jacobian and r is the residual.
Definition: BlackoilModelEbos.hpp:603
void updateSolution(const BVector &dx)
Apply an update to the primary variables.
Definition: BlackoilModelEbos.hpp:629
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:92
Represents the convergence status of the whole simulator, to make it possible to query and store the ...
Definition: ConvergenceReport.hpp:38
This class solves the fully implicit black-oil system by solving the reduced system (after eliminatin...
Definition: ISTLSolverEbos.hpp:200
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:50
virtual bool lastStepFailed() const =0
Return true if last time step failed.
virtual double currentStepLength() const =0
Current step length.
virtual double simulationTimeElapsed() const =0
Time elapsed since the start of the simulation until the beginning of the current time step [s].
virtual int currentStepNum() const =0
Current step number.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27
PhaseUsage phaseUsageFromDeck(const EclipseState &eclipseState)
Looks at presence of WATER, OIL and GAS keywords in state object to determine active phases.
Definition: phaseUsageFromDeck.cpp:145
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParametersEbos.hpp:327
double tolerance_cnv_relaxed_
Relaxed local convergence tolerance (can be used when iter >= min_strict_cnv_iter_ && cnvViolatedPV <...
Definition: BlackoilModelParametersEbos.hpp:346
bool use_update_stabilization_
Try to detect oscillation or stagnation.
Definition: BlackoilModelParametersEbos.hpp:401
double tolerance_cnv_
Local convergence tolerance (max of local saturation errors).
Definition: BlackoilModelParametersEbos.hpp:344
bool update_equations_scaling_
Update scaling factors for mass balance equations.
Definition: BlackoilModelParametersEbos.hpp:398
double tolerance_mb_
Relative mass balance tolerance (total mass balance error).
Definition: BlackoilModelParametersEbos.hpp:342
int min_strict_cnv_iter_
Minimum number of Newton iterations before we can use relaxed CNV convergence criterion.
Definition: BlackoilModelParametersEbos.hpp:392
Definition: ISTLSolverEbos.hpp:68
Definition: BlackoilModelEbos.hpp:72
Definition: ISTLSolverEbos.hpp:62
Definition: BlackoilModelParametersEbos.hpp:31
Definition: NonlinearSolverEbos.hpp:41
Definition: AdaptiveTimeSteppingEbos.hpp:49
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34