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/simulators/timestepping/AdaptiveTimeSteppingEbos.hpp>
32 #include <opm/simulators/flow/NonlinearSolverEbos.hpp>
33 #include <opm/simulators/flow/BlackoilModelParametersEbos.hpp>
34 #include <opm/simulators/wells/BlackoilWellModel.hpp>
35 #include <opm/simulators/aquifers/BlackoilAquiferModel.hpp>
36 #include <opm/simulators/wells/WellConnectionAuxiliaryModule.hpp>
37 #include <opm/simulators/flow/countGlobalCells.hpp>
38 #include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
40 #include <opm/grid/UnstructuredGrid.h>
41 #include <opm/simulators/timestepping/SimulatorReport.hpp>
42 #include <opm/core/props/phaseUsageFromDeck.hpp>
43 #include <opm/common/ErrorMacros.hpp>
44 #include <opm/common/Exceptions.hpp>
45 #include <opm/common/OpmLog/OpmLog.hpp>
46 #include <opm/input/eclipse/Units/Units.hpp>
47 #include <opm/simulators/timestepping/SimulatorTimer.hpp>
48 #include <opm/input/eclipse/EclipseState/EclipseState.hpp>
49 #include <opm/input/eclipse/EclipseState/Tables/TableManager.hpp>
51 #include <opm/simulators/linalg/ISTLSolverEbos.hpp>
53 #include <dune/istl/owneroverlapcopy.hh>
54 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 7)
55 #include <dune/common/parallel/communication.hh>
57 #include <dune/common/parallel/collectivecommunication.hh>
59 #include <dune/common/timer.hh>
60 #include <dune/common/unused.hh>
70 namespace Opm::Properties {
78 template<
class TypeTag>
79 struct OutputDir<TypeTag, TTag::EclFlowProblem> {
80 static constexpr
auto value =
"";
82 template<
class TypeTag>
83 struct EnableDebuggingChecks<TypeTag, TTag::EclFlowProblem> {
84 static constexpr
bool value =
false;
87 template<
class TypeTag>
88 struct BlackoilConserveSurfaceVolume<TypeTag, TTag::EclFlowProblem> {
89 static constexpr
bool value =
true;
91 template<
class TypeTag>
92 struct UseVolumetricResidual<TypeTag, TTag::EclFlowProblem> {
93 static constexpr
bool value =
false;
96 template<
class TypeTag>
97 struct EclAquiferModel<TypeTag, TTag::EclFlowProblem> {
103 template<
class TypeTag>
104 struct EnablePolymer<TypeTag, TTag::EclFlowProblem> {
105 static constexpr
bool value =
false;
107 template<
class TypeTag>
108 struct EnableSolvent<TypeTag, TTag::EclFlowProblem> {
109 static constexpr
bool value =
false;
111 template<
class TypeTag>
112 struct EnableTemperature<TypeTag, TTag::EclFlowProblem> {
113 static constexpr
bool value =
true;
115 template<
class TypeTag>
116 struct EnableEnergy<TypeTag, TTag::EclFlowProblem> {
117 static constexpr
bool value =
false;
119 template<
class TypeTag>
120 struct EnableFoam<TypeTag, TTag::EclFlowProblem> {
121 static constexpr
bool value =
false;
123 template<
class TypeTag>
124 struct EnableBrine<TypeTag, TTag::EclFlowProblem> {
125 static constexpr
bool value =
false;
127 template<
class TypeTag>
128 struct EnableSaltPrecipitation<TypeTag, TTag::EclFlowProblem> {
129 static constexpr
bool value =
false;
131 template<
class TypeTag>
132 struct EnableMICP<TypeTag, TTag::EclFlowProblem> {
133 static constexpr
bool value =
false;
136 template<
class TypeTag>
140 template<
class TypeTag>
141 struct LinearSolverSplice<TypeTag, TTag::EclFlowProblem> {
154 template <
class TypeTag>
161 using Simulator = GetPropType<TypeTag, Properties::Simulator>;
162 using Grid = GetPropType<TypeTag, Properties::Grid>;
163 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
164 using SparseMatrixAdapter = GetPropType<TypeTag, Properties::SparseMatrixAdapter>;
165 using SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>;
166 using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
167 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
168 using Indices = GetPropType<TypeTag, Properties::Indices>;
169 using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
170 using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
172 typedef double Scalar;
173 static const int numEq = Indices::numEq;
174 static const int contiSolventEqIdx = Indices::contiSolventEqIdx;
175 static const int contiZfracEqIdx = Indices::contiZfracEqIdx;
176 static const int contiPolymerEqIdx = Indices::contiPolymerEqIdx;
177 static const int contiEnergyEqIdx = Indices::contiEnergyEqIdx;
178 static const int contiPolymerMWEqIdx = Indices::contiPolymerMWEqIdx;
179 static const int contiFoamEqIdx = Indices::contiFoamEqIdx;
180 static const int contiBrineEqIdx = Indices::contiBrineEqIdx;
181 static const int contiMicrobialEqIdx = Indices::contiMicrobialEqIdx;
182 static const int contiOxygenEqIdx = Indices::contiOxygenEqIdx;
183 static const int contiUreaEqIdx = Indices::contiUreaEqIdx;
184 static const int contiBiofilmEqIdx = Indices::contiBiofilmEqIdx;
185 static const int contiCalciteEqIdx = Indices::contiCalciteEqIdx;
186 static const int solventSaturationIdx = Indices::solventSaturationIdx;
187 static const int zFractionIdx = Indices::zFractionIdx;
188 static const int polymerConcentrationIdx = Indices::polymerConcentrationIdx;
189 static const int polymerMoleWeightIdx = Indices::polymerMoleWeightIdx;
190 static const int temperatureIdx = Indices::temperatureIdx;
191 static const int foamConcentrationIdx = Indices::foamConcentrationIdx;
192 static const int saltConcentrationIdx = Indices::saltConcentrationIdx;
193 static const int microbialConcentrationIdx = Indices::microbialConcentrationIdx;
194 static const int oxygenConcentrationIdx = Indices::oxygenConcentrationIdx;
195 static const int ureaConcentrationIdx = Indices::ureaConcentrationIdx;
196 static const int biofilmConcentrationIdx = Indices::biofilmConcentrationIdx;
197 static const int calciteConcentrationIdx = Indices::calciteConcentrationIdx;
199 typedef Dune::FieldVector<Scalar, numEq > VectorBlockType;
200 typedef typename SparseMatrixAdapter::MatrixBlock MatrixBlockType;
201 typedef typename SparseMatrixAdapter::IstlMatrix Mat;
202 typedef Dune::BlockVector<VectorBlockType> BVector;
222 const bool terminal_output)
223 : ebosSimulator_(ebosSimulator)
224 , grid_(ebosSimulator_.vanguard().grid())
227 , well_model_ (well_model)
229 , current_relaxation_(1.0)
230 , dx_old_(ebosSimulator_.model().numGridDof())
234 convergence_reports_.reserve(300);
237 bool isParallel()
const
238 {
return grid_.comm().size() > 1; }
240 const EclipseState& eclState()
const
241 {
return ebosSimulator_.vanguard().eclState(); }
248 Dune::Timer perfTimer;
252 ebosSimulator_.model().updateFailed();
254 ebosSimulator_.model().advanceTimeLevel();
263 ebosSimulator_.model().newtonMethod().setIterationIndex(0);
265 ebosSimulator_.problem().beginTimeStep();
267 unsigned numDof = ebosSimulator_.model().numGridDof();
268 wasSwitched_.resize(numDof);
269 std::fill(wasSwitched_.begin(), wasSwitched_.end(),
false);
272 std::cout <<
"equation scaling not supported yet" << std::endl;
275 report.pre_post_time += perfTimer.stop();
290 template <
class NonlinearSolverType>
293 NonlinearSolverType& nonlinear_solver)
297 Dune::Timer perfTimer;
300 if (iteration == 0) {
303 residual_norms_history_.clear();
304 current_relaxation_ = 1.0;
307 convergence_reports_.back().report.reserve(11);
310 report.total_linearizations = 1;
314 report.assemble_time += perfTimer.stop();
317 report.assemble_time += perfTimer.stop();
318 failureReport_ += report;
323 std::vector<double> residual_norms;
329 report.converged = convrep.converged() && iteration > nonlinear_solver.minIter();;
330 ConvergenceReport::Severity severity = convrep.severityOfWorstFailure();
331 convergence_reports_.back().report.push_back(std::move(convrep));
334 if (severity == ConvergenceReport::Severity::NotANumber) {
335 OPM_THROW(NumericalIssue,
"NaN residual found!");
336 }
else if (severity == ConvergenceReport::Severity::TooLarge) {
337 OPM_THROW_NOLOG(NumericalIssue,
"Too large residual found!");
340 report.update_time += perfTimer.stop();
341 residual_norms_history_.push_back(residual_norms);
342 if (!report.converged) {
345 report.total_newton_iterations = 1;
351 unsigned nc = ebosSimulator_.model().numGridDof();
355 linear_solve_setup_time_ = 0.0;
360 wellModel().linearize(ebosSimulator().model().linearizer().jacobian(),
361 ebosSimulator().model().linearizer().residual());
364 report.linear_solve_setup_time += linear_solve_setup_time_;
365 report.linear_solve_time += perfTimer.stop();
369 report.linear_solve_setup_time += linear_solve_setup_time_;
370 report.linear_solve_time += perfTimer.stop();
373 failureReport_ += report;
387 bool isOscillate =
false;
388 bool isStagnate =
false;
389 nonlinear_solver.detectOscillations(residual_norms_history_, iteration, isOscillate, isStagnate);
391 current_relaxation_ -= nonlinear_solver.relaxIncrement();
392 current_relaxation_ = std::max(current_relaxation_, nonlinear_solver.relaxMax());
394 std::string msg =
" Oscillating behavior detected: Relaxation set to "
395 + std::to_string(current_relaxation_);
399 nonlinear_solver.stabilizeNonlinearUpdate(x, dx_old_, current_relaxation_);
406 report.update_time += perfTimer.stop();
412 void printIf(
int c,
double x,
double y,
double eps, std::string type) {
413 if (std::abs(x-y) > eps) {
414 std::cout << type <<
" " <<c <<
": "<<x <<
" " << y << std::endl;
425 Dune::Timer perfTimer;
427 ebosSimulator_.problem().endTimeStep();
428 report.pre_post_time += perfTimer.stop();
437 const int iterationIdx)
440 ebosSimulator_.model().newtonMethod().setIterationIndex(iterationIdx);
441 ebosSimulator_.problem().beginIteration();
442 ebosSimulator_.model().linearizer().linearizeDomain();
443 ebosSimulator_.problem().endIteration();
449 double relativeChange()
const
451 Scalar resultDelta = 0.0;
452 Scalar resultDenom = 0.0;
454 const auto& elemMapper = ebosSimulator_.model().elementMapper();
455 const auto& gridView = ebosSimulator_.gridView();
456 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
457 unsigned globalElemIdx = elemMapper.index(elem);
458 const auto& priVarsNew = ebosSimulator_.model().solution(0)[globalElemIdx];
461 pressureNew = priVarsNew[Indices::pressureSwitchIdx];
463 Scalar saturationsNew[FluidSystem::numPhases] = { 0.0 };
464 Scalar oilSaturationNew = 1.0;
465 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) && FluidSystem::numActivePhases() > 1) {
466 saturationsNew[FluidSystem::waterPhaseIdx] = priVarsNew[Indices::waterSaturationIdx];
467 oilSaturationNew -= saturationsNew[FluidSystem::waterPhaseIdx];
470 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
471 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
472 priVarsNew.primaryVarsMeaning() == PrimaryVariables::Sw_po_Sg) {
473 assert(Indices::compositionSwitchIdx >= 0 );
474 saturationsNew[FluidSystem::gasPhaseIdx] = priVarsNew[Indices::compositionSwitchIdx];
475 oilSaturationNew -= saturationsNew[FluidSystem::gasPhaseIdx];
478 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
479 saturationsNew[FluidSystem::oilPhaseIdx] = oilSaturationNew;
482 const auto& priVarsOld = ebosSimulator_.model().solution(1)[globalElemIdx];
485 pressureOld = priVarsOld[Indices::pressureSwitchIdx];
487 Scalar saturationsOld[FluidSystem::numPhases] = { 0.0 };
488 Scalar oilSaturationOld = 1.0;
491 Scalar tmp = pressureNew - pressureOld;
492 resultDelta += tmp*tmp;
493 resultDenom += pressureNew*pressureNew;
495 if (FluidSystem::numActivePhases() > 1) {
496 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
497 saturationsOld[FluidSystem::waterPhaseIdx] = priVarsOld[Indices::waterSaturationIdx];
498 oilSaturationOld -= saturationsOld[FluidSystem::waterPhaseIdx];
501 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
502 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
503 priVarsOld.primaryVarsMeaning() == PrimaryVariables::Sw_po_Sg)
505 assert(Indices::compositionSwitchIdx >= 0 );
506 saturationsOld[FluidSystem::gasPhaseIdx] = priVarsOld[Indices::compositionSwitchIdx];
507 oilSaturationOld -= saturationsOld[FluidSystem::gasPhaseIdx];
510 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
511 saturationsOld[FluidSystem::oilPhaseIdx] = oilSaturationOld;
513 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++ phaseIdx) {
514 Scalar tmpSat = saturationsNew[phaseIdx] - saturationsOld[phaseIdx];
515 resultDelta += tmpSat*tmpSat;
516 resultDenom += saturationsNew[phaseIdx]*saturationsNew[phaseIdx];
517 assert(std::isfinite(resultDelta));
518 assert(std::isfinite(resultDenom));
523 resultDelta = gridView.comm().sum(resultDelta);
524 resultDenom = gridView.comm().sum(resultDenom);
526 if (resultDenom > 0.0)
527 return resultDelta/resultDenom;
535 return ebosSimulator_.model().newtonMethod().linearSolver().iterations ();
543 auto& ebosJac = ebosSimulator_.model().linearizer().jacobian();
544 auto& ebosResid = ebosSimulator_.model().linearizer().residual();
549 auto& ebosSolver = ebosSimulator_.model().newtonMethod().linearSolver();
550 Dune::Timer perfTimer;
552 ebosSolver.prepare(ebosJac, ebosResid);
553 linear_solve_setup_time_ = perfTimer.stop();
554 ebosSolver.setResidual(ebosResid);
559 ebosSolver.setMatrix(ebosJac);
568 auto& ebosNewtonMethod = ebosSimulator_.model().newtonMethod();
569 SolutionVector& solution = ebosSimulator_.model().solution(0);
571 ebosNewtonMethod.update_(solution,
579 ebosSimulator_.model().invalidateAndUpdateIntensiveQuantities(0);
588 template <
class CollectiveCommunication>
589 std::tuple<double,double> convergenceReduction(
const CollectiveCommunication& comm,
590 const double pvSumLocal,
591 const double numAquiferPvSumLocal,
592 std::vector< Scalar >& R_sum,
593 std::vector< Scalar >& maxCoeff,
594 std::vector< Scalar >& B_avg)
597 double pvSum = pvSumLocal;
598 double numAquiferPvSum = numAquiferPvSumLocal;
600 if( comm.size() > 1 )
603 std::vector< Scalar > sumBuffer;
604 std::vector< Scalar > maxBuffer;
605 const int numComp = B_avg.size();
606 sumBuffer.reserve( 2*numComp + 2 );
607 maxBuffer.reserve( numComp );
608 for(
int compIdx = 0; compIdx < numComp; ++compIdx )
610 sumBuffer.push_back( B_avg[ compIdx ] );
611 sumBuffer.push_back( R_sum[ compIdx ] );
612 maxBuffer.push_back( maxCoeff[ compIdx ] );
616 sumBuffer.push_back( pvSum );
617 sumBuffer.push_back( numAquiferPvSum );
620 comm.sum( sumBuffer.data(), sumBuffer.size() );
623 comm.max( maxBuffer.data(), maxBuffer.size() );
626 for(
int compIdx = 0, buffIdx = 0; compIdx < numComp; ++compIdx, ++buffIdx )
628 B_avg[ compIdx ] = sumBuffer[ buffIdx ];
631 R_sum[ compIdx ] = sumBuffer[ buffIdx ];
634 for(
int compIdx = 0; compIdx < numComp; ++compIdx )
636 maxCoeff[ compIdx ] = maxBuffer[ compIdx ];
640 pvSum = sumBuffer[sumBuffer.size()-2];
641 numAquiferPvSum = sumBuffer.back();
645 return {pvSum, numAquiferPvSum};
652 std::vector<Scalar>& maxCoeff,
653 std::vector<Scalar>& B_avg)
655 double pvSumLocal = 0.0;
656 double numAquiferPvSumLocal = 0.0;
657 const auto& ebosModel = ebosSimulator_.model();
658 const auto& ebosProblem = ebosSimulator_.problem();
660 const auto& ebosResid = ebosSimulator_.model().linearizer().residual();
662 ElementContext elemCtx(ebosSimulator_);
663 const auto& gridView = ebosSimulator().gridView();
664 OPM_BEGIN_PARALLEL_TRY_CATCH();
665 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
666 elemCtx.updatePrimaryStencil(elem);
667 elemCtx.updatePrimaryIntensiveQuantities(0);
668 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
669 const auto& intQuants = elemCtx.intensiveQuantities(0, 0);
670 const auto& fs = intQuants.fluidState();
672 const double pvValue = ebosProblem.referencePorosity(cell_idx, 0) * ebosModel.dofTotalVolume( cell_idx );
673 pvSumLocal += pvValue;
675 if (isNumericalAquiferCell(gridView.grid(), elem))
677 numAquiferPvSumLocal += pvValue;
680 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
682 if (!FluidSystem::phaseIsActive(phaseIdx)) {
686 const unsigned compIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
688 B_avg[ compIdx ] += 1.0 / fs.invB(phaseIdx).value();
689 const auto R2 = ebosResid[cell_idx][compIdx];
691 R_sum[ compIdx ] += R2;
692 maxCoeff[ compIdx ] = std::max( maxCoeff[ compIdx ], std::abs( R2 ) / pvValue );
695 if constexpr (has_solvent_) {
696 B_avg[ contiSolventEqIdx ] += 1.0 / intQuants.solventInverseFormationVolumeFactor().value();
697 const auto R2 = ebosResid[cell_idx][contiSolventEqIdx];
698 R_sum[ contiSolventEqIdx ] += R2;
699 maxCoeff[ contiSolventEqIdx ] = std::max( maxCoeff[ contiSolventEqIdx ], std::abs( R2 ) / pvValue );
701 if constexpr (has_extbo_) {
702 B_avg[ contiZfracEqIdx ] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
703 const auto R2 = ebosResid[cell_idx][contiZfracEqIdx];
704 R_sum[ contiZfracEqIdx ] += R2;
705 maxCoeff[ contiZfracEqIdx ] = std::max( maxCoeff[ contiZfracEqIdx ], std::abs( R2 ) / pvValue );
707 if constexpr (has_polymer_) {
708 B_avg[ contiPolymerEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
709 const auto R2 = ebosResid[cell_idx][contiPolymerEqIdx];
710 R_sum[ contiPolymerEqIdx ] += R2;
711 maxCoeff[ contiPolymerEqIdx ] = std::max( maxCoeff[ contiPolymerEqIdx ], std::abs( R2 ) / pvValue );
713 if constexpr (has_foam_) {
714 B_avg[ contiFoamEqIdx ] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
715 const auto R2 = ebosResid[cell_idx][contiFoamEqIdx];
716 R_sum[ contiFoamEqIdx ] += R2;
717 maxCoeff[ contiFoamEqIdx ] = std::max( maxCoeff[ contiFoamEqIdx ], std::abs( R2 ) / pvValue );
719 if constexpr (has_brine_) {
720 B_avg[ contiBrineEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
721 const auto R2 = ebosResid[cell_idx][contiBrineEqIdx];
722 R_sum[ contiBrineEqIdx ] += R2;
723 maxCoeff[ contiBrineEqIdx ] = std::max( maxCoeff[ contiBrineEqIdx ], std::abs( R2 ) / pvValue );
726 if constexpr (has_polymermw_) {
727 static_assert(has_polymer_);
729 B_avg[contiPolymerMWEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
733 const auto R2 = ebosResid[cell_idx][contiPolymerMWEqIdx] / 100.;
734 R_sum[contiPolymerMWEqIdx] += R2;
735 maxCoeff[contiPolymerMWEqIdx] = std::max( maxCoeff[contiPolymerMWEqIdx], std::abs( R2 ) / pvValue );
738 if constexpr (has_energy_) {
739 B_avg[ contiEnergyEqIdx ] += 1.0 / (4.182e1);
740 const auto R2 = ebosResid[cell_idx][contiEnergyEqIdx];
741 R_sum[ contiEnergyEqIdx ] += R2;
742 maxCoeff[ contiEnergyEqIdx ] = std::max( maxCoeff[ contiEnergyEqIdx ], std::abs( R2 ) / pvValue );
745 if constexpr (has_micp_) {
746 B_avg[ contiMicrobialEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
747 const auto R1 = ebosResid[cell_idx][contiMicrobialEqIdx];
748 R_sum[ contiMicrobialEqIdx ] += R1;
749 maxCoeff[ contiMicrobialEqIdx ] = std::max( maxCoeff[ contiMicrobialEqIdx ], std::abs( R1 ) / pvValue );
750 B_avg[ contiOxygenEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
751 const auto R2 = ebosResid[cell_idx][contiOxygenEqIdx];
752 R_sum[ contiOxygenEqIdx ] += R2;
753 maxCoeff[ contiOxygenEqIdx ] = std::max( maxCoeff[ contiOxygenEqIdx ], std::abs( R2 ) / pvValue );
754 B_avg[ contiUreaEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
755 const auto R3 = ebosResid[cell_idx][contiUreaEqIdx];
756 R_sum[ contiUreaEqIdx ] += R3;
757 maxCoeff[ contiUreaEqIdx ] = std::max( maxCoeff[ contiUreaEqIdx ], std::abs( R3 ) / pvValue );
758 B_avg[ contiBiofilmEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
759 const auto R4 = ebosResid[cell_idx][contiBiofilmEqIdx];
760 R_sum[ contiBiofilmEqIdx ] += R4;
761 maxCoeff[ contiBiofilmEqIdx ] = std::max( maxCoeff[ contiBiofilmEqIdx ], std::abs( R4 ) / pvValue );
762 B_avg[ contiCalciteEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
763 const auto R5 = ebosResid[cell_idx][contiCalciteEqIdx];
764 R_sum[ contiCalciteEqIdx ] += R5;
765 maxCoeff[ contiCalciteEqIdx ] = std::max( maxCoeff[ contiCalciteEqIdx ], std::abs( R5 ) / pvValue );
769 OPM_END_PARALLEL_TRY_CATCH(
"BlackoilModelEbos::localConvergenceData() failed: ", grid_.comm());
772 const int bSize = B_avg.size();
773 for (
int i = 0; i<bSize; ++i )
778 return {pvSumLocal, numAquiferPvSumLocal};
786 const auto& ebosModel = ebosSimulator_.model();
787 const auto& ebosProblem = ebosSimulator_.problem();
788 const auto& ebosResid = ebosSimulator_.model().linearizer().residual();
789 const auto& gridView = ebosSimulator().gridView();
790 ElementContext elemCtx(ebosSimulator_);
792 OPM_BEGIN_PARALLEL_TRY_CATCH();
794 for (
const auto& elem : elements(gridView, Dune::Partitions::interiorBorder))
797 if (isNumericalAquiferCell(gridView.grid(), elem))
801 elemCtx.updatePrimaryStencil(elem);
802 elemCtx.updatePrimaryIntensiveQuantities(0);
803 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
804 const double pvValue = ebosProblem.referencePorosity(cell_idx, 0) * ebosModel.dofTotalVolume( cell_idx );
805 const auto& cellResidual = ebosResid[cell_idx];
806 bool cnvViolated =
false;
808 for (
unsigned eqIdx = 0; eqIdx < cellResidual.size(); ++eqIdx)
811 Scalar CNV = cellResidual[eqIdx] * dt * B_avg[eqIdx] / pvValue;
821 OPM_END_PARALLEL_TRY_CATCH(
"BlackoilModelEbos::ComputeCnvError() failed: ", grid_.comm());
823 return grid_.comm().sum(errorPV);
828 std::vector<Scalar>& B_avg,
829 std::vector<Scalar>& residual_norms)
831 typedef std::vector< Scalar > Vector;
833 const int numComp = numEq;
834 Vector R_sum(numComp, 0.0 );
835 Vector maxCoeff(numComp, std::numeric_limits< Scalar >::lowest() );
839 const auto [ pvSum, numAquiferPvSum ] =
840 convergenceReduction(grid_.comm(), pvSumLocal,
841 numAquiferPvSumLocal,
842 R_sum, maxCoeff, B_avg);
845 cnvErrorPvFraction /= (pvSum - numAquiferPvSum);
851 const bool use_relaxed = cnvErrorPvFraction < param_.relaxed_max_pv_fraction_ && iteration >= param_.
min_strict_cnv_iter_;
855 std::vector<Scalar> CNV(numComp);
856 std::vector<Scalar> mass_balance_residual(numComp);
857 for (
int compIdx = 0; compIdx < numComp; ++compIdx )
859 CNV[compIdx] = B_avg[compIdx] * dt * maxCoeff[compIdx];
860 mass_balance_residual[compIdx] = std::abs(B_avg[compIdx]*R_sum[compIdx]) * dt / pvSum;
861 residual_norms.push_back(CNV[compIdx]);
865 static std::vector<std::string> compNames;
866 if (compNames.empty()) {
867 compNames.resize(numComp);
868 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
869 if (!FluidSystem::phaseIsActive(phaseIdx)) {
872 const unsigned canonicalCompIdx = FluidSystem::solventComponentIndex(phaseIdx);
873 const unsigned compIdx = Indices::canonicalToActiveComponentIndex(canonicalCompIdx);
874 compNames[compIdx] = FluidSystem::componentName(canonicalCompIdx);
876 if constexpr (has_solvent_) {
877 compNames[solventSaturationIdx] =
"Solvent";
879 if constexpr (has_extbo_) {
880 compNames[zFractionIdx] =
"ZFraction";
882 if constexpr (has_polymer_) {
883 compNames[polymerConcentrationIdx] =
"Polymer";
885 if constexpr (has_polymermw_) {
886 assert(has_polymer_);
887 compNames[polymerMoleWeightIdx] =
"MolecularWeightP";
889 if constexpr (has_energy_) {
890 compNames[temperatureIdx] =
"Energy";
892 if constexpr (has_foam_) {
893 compNames[foamConcentrationIdx] =
"Foam";
895 if constexpr (has_brine_) {
896 compNames[saltConcentrationIdx] =
"Brine";
898 if constexpr (has_micp_) {
899 compNames[microbialConcentrationIdx] =
"Microbes";
900 compNames[oxygenConcentrationIdx] =
"Oxygen";
901 compNames[ureaConcentrationIdx] =
"Urea";
902 compNames[biofilmConcentrationIdx] =
"Biofilm";
903 compNames[calciteConcentrationIdx] =
"Calcite";
908 ConvergenceReport report;
909 using CR = ConvergenceReport;
910 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
911 double res[2] = { mass_balance_residual[compIdx], CNV[compIdx] };
912 CR::ReservoirFailure::Type types[2] = { CR::ReservoirFailure::Type::MassBalance,
913 CR::ReservoirFailure::Type::Cnv };
914 double tol[2] = { tol_mb, tol_cnv };
915 for (
int ii : {0, 1}) {
916 if (std::isnan(res[ii])) {
917 report.setReservoirFailed({types[ii], CR::Severity::NotANumber, compIdx});
919 OpmLog::debug(
"NaN residual for " + compNames[compIdx] +
" equation.");
921 }
else if (res[ii] > maxResidualAllowed()) {
922 report.setReservoirFailed({types[ii], CR::Severity::TooLarge, compIdx});
924 OpmLog::debug(
"Too large residual for " + compNames[compIdx] +
" equation.");
926 }
else if (res[ii] < 0.0) {
927 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
929 OpmLog::debug(
"Negative residual for " + compNames[compIdx] +
" equation.");
931 }
else if (res[ii] > tol[ii]) {
932 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
941 if (iteration == 0) {
942 std::string msg =
"Iter";
943 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
945 msg += compNames[compIdx][0];
948 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
950 msg += compNames[compIdx][0];
955 std::ostringstream ss;
956 const std::streamsize oprec = ss.precision(3);
957 const std::ios::fmtflags oflags = ss.setf(std::ios::scientific);
958 ss << std::setw(4) << iteration;
959 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
960 ss << std::setw(11) << mass_balance_residual[compIdx];
962 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
963 ss << std::setw(11) << CNV[compIdx];
967 OpmLog::debug(ss.str());
980 std::vector<double>& residual_norms)
983 std::vector<Scalar> B_avg(numEq, 0.0);
984 auto report = getReservoirConvergence(timer.
currentStepLength(), iteration, B_avg, residual_norms);
985 report +=
wellModel().getWellConvergence(B_avg, report.converged());
994 return phaseUsage_.num_phases;
999 std::vector<std::vector<double> >
1006 std::vector<std::vector<double> >
1011 std::vector<std::vector<double> > regionValues(0, std::vector<double>(0,0.0));
1012 return regionValues;
1015 const Simulator& ebosSimulator()
const
1016 {
return ebosSimulator_; }
1018 Simulator& ebosSimulator()
1019 {
return ebosSimulator_; }
1023 {
return failureReport_; }
1029 std::vector<ConvergenceReport> report;
1032 const std::vector<StepReport>& stepReports()
const
1034 return convergence_reports_;
1040 Simulator& ebosSimulator_;
1043 static constexpr
bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>();
1044 static constexpr
bool has_extbo_ = getPropValue<TypeTag, Properties::EnableExtbo>();
1045 static constexpr
bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>();
1046 static constexpr
bool has_polymermw_ = getPropValue<TypeTag, Properties::EnablePolymerMW>();
1047 static constexpr
bool has_energy_ = getPropValue<TypeTag, Properties::EnableEnergy>();
1048 static constexpr
bool has_foam_ = getPropValue<TypeTag, Properties::EnableFoam>();
1049 static constexpr
bool has_brine_ = getPropValue<TypeTag, Properties::EnableBrine>();
1050 static constexpr
bool has_micp_ = getPropValue<TypeTag, Properties::EnableMICP>();
1052 ModelParameters param_;
1063 std::vector<std::vector<double>> residual_norms_history_;
1064 double current_relaxation_;
1067 std::vector<StepReport> convergence_reports_;
1074 wellModel()
const {
return well_model_; }
1076 void beginReportStep()
1078 ebosSimulator_.problem().beginEpisode();
1081 void endReportStep()
1083 ebosSimulator_.problem().endEpisode();
1088 bool isNumericalAquiferCell(
const Dune::CpGrid& grid,
const T& elem)
1090 const auto& aquiferCells = grid.sortedNumAquiferCells();
1091 if (aquiferCells.empty())
1095 auto candidate = std::lower_bound(aquiferCells.begin(), aquiferCells.end(),
1097 return candidate != aquiferCells.end() && *candidate == elem.index();
1100 template<
class G,
class T>
1101 typename std::enable_if<!std::is_same<G,Dune::CpGrid>::value,
bool>::type
1102 isNumericalAquiferCell(
const G&,
const T&)
1107 double dpMaxRel()
const {
return param_.dp_max_rel_; }
1108 double dsMax()
const {
return param_.ds_max_; }
1109 double drMaxRel()
const {
return param_.dr_max_rel_; }
1110 double maxResidualAllowed()
const {
return param_.max_residual_allowed_; }
1111 double linear_solve_setup_time_;
1113 std::vector<bool> wasSwitched_;
Class for handling the blackoil well model.
Definition: BlackoilAquiferModel.hpp:85
A model implementation for three-phase black oil.
Definition: BlackoilModelEbos.hpp:156
int linearIterationsLastSolve() const
Number of linear iterations used in last call to solveJacobianSystem().
Definition: BlackoilModelEbos.hpp:533
bool terminalOutputEnabled() const
Return true if output to cout is wanted.
Definition: BlackoilModelEbos.hpp:583
std::vector< std::vector< double > > computeFluidInPlace(const std::vector< int > &) const
Should not be called.
Definition: BlackoilModelEbos.hpp:1007
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:978
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:1000
BlackoilModelEbos(Simulator &ebosSimulator, const ModelParameters ¶m, BlackoilWellModel< TypeTag > &well_model, const bool terminal_output)
Construct the model.
Definition: BlackoilModelEbos.hpp:219
bool terminal_output_
Whether we print something to std::cout.
Definition: BlackoilModelEbos.hpp:1059
SimulatorReportSingle afterStep(const SimulatorTimerInterface &)
Called once after each time step.
Definition: BlackoilModelEbos.hpp:422
SimulatorReportSingle nonlinearIteration(const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Called once per nonlinear iteration.
Definition: BlackoilModelEbos.hpp:291
SimulatorReportSingle assembleReservoir(const SimulatorTimerInterface &, const int iterationIdx)
Assemble the residual and Jacobian of the nonlinear system.
Definition: BlackoilModelEbos.hpp:436
SimulatorReportSingle prepareStep(const SimulatorTimerInterface &timer)
Called once before each time step.
Definition: BlackoilModelEbos.hpp:245
long int global_nc_
The number of cells of the global grid.
Definition: BlackoilModelEbos.hpp:1061
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:783
const SimulatorReportSingle & failureReport() const
return the statistics if the nonlinearIteration() method failed
Definition: BlackoilModelEbos.hpp:1022
int numPhases() const
The number of active fluid phases in the model.
Definition: BlackoilModelEbos.hpp:992
BlackoilWellModel< TypeTag > & wellModel()
return the StandardWells object
Definition: BlackoilModelEbos.hpp:1071
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:651
void solveJacobianSystem(BVector &x)
Solve the Jacobian system Jx = r where J is the Jacobian and r is the residual.
Definition: BlackoilModelEbos.hpp:540
void updateSolution(const BVector &dx)
Apply an update to the primary variables.
Definition: BlackoilModelEbos.hpp:566
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:94
Represents the convergence status of the whole simulator, to make it possible to query and store the ...
Definition: ConvergenceReport.hpp:36
This class solves the fully implicit black-oil system by solving the reduced system (after eliminatin...
Definition: ISTLSolverEbos.hpp:197
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
Definition: BlackoilModelEbos.hpp:1026
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: BlackoilPhases.hpp:46
Definition: ISTLSolverEbos.hpp:66
Definition: BlackoilModelEbos.hpp:73
Definition: ISTLSolverEbos.hpp:60
Definition: BlackoilModelParametersEbos.hpp:31
Definition: NonlinearSolverEbos.hpp:41
Definition: AdaptiveTimeSteppingEbos.hpp:24
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:31