mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Merge pull request #3592 from blattms/fix-3574
remove assert for isfinite and correctly handle exceptions
This commit is contained in:
commit
ed73f1b326
@ -335,14 +335,15 @@ namespace Opm {
|
|||||||
const int nc = UgGridHelpers::numCells(grid_);
|
const int nc = UgGridHelpers::numCells(grid_);
|
||||||
BVector x(nc);
|
BVector x(nc);
|
||||||
|
|
||||||
// apply the Schur compliment of the well model to the reservoir linearized
|
|
||||||
// equations
|
|
||||||
wellModel().linearize(ebosSimulator().model().linearizer().jacobian(),
|
|
||||||
ebosSimulator().model().linearizer().residual());
|
|
||||||
|
|
||||||
// Solve the linear system.
|
// Solve the linear system.
|
||||||
linear_solve_setup_time_ = 0.0;
|
linear_solve_setup_time_ = 0.0;
|
||||||
try {
|
try {
|
||||||
|
// apply the Schur compliment of the well model to the reservoir linearized
|
||||||
|
// equations
|
||||||
|
// Note that linearize may throw for MSwells.
|
||||||
|
wellModel().linearize(ebosSimulator().model().linearizer().jacobian(),
|
||||||
|
ebosSimulator().model().linearizer().residual());
|
||||||
|
|
||||||
solveJacobianSystem(x);
|
solveJacobianSystem(x);
|
||||||
report.linear_solve_setup_time += linear_solve_setup_time_;
|
report.linear_solve_setup_time += linear_solve_setup_time_;
|
||||||
report.linear_solve_time += perfTimer.stop();
|
report.linear_solve_time += perfTimer.stop();
|
||||||
|
@ -149,16 +149,19 @@ namespace Opm {
|
|||||||
BlackoilWellModel<TypeTag>::
|
BlackoilWellModel<TypeTag>::
|
||||||
linearize(SparseMatrixAdapter& jacobian, GlobalEqVector& res)
|
linearize(SparseMatrixAdapter& jacobian, GlobalEqVector& res)
|
||||||
{
|
{
|
||||||
if (!localWellsActive())
|
if (!param_.matrix_add_well_contributions_)
|
||||||
return;
|
{
|
||||||
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
||||||
if (!param_.matrix_add_well_contributions_) {
|
{
|
||||||
// if the well contributions are not supposed to be included explicitly in
|
// if the well contributions are not supposed to be included explicitly in
|
||||||
// the matrix, we only apply the vector part of the Schur complement here.
|
// the matrix, we only apply the vector part of the Schur complement here.
|
||||||
for (const auto& well: well_container_) {
|
for (const auto& well: well_container_) {
|
||||||
// r = r - duneC_^T * invDuneD_ * resWell_
|
// r = r - duneC_^T * invDuneD_ * resWell_
|
||||||
well->apply(res);
|
well->apply(res);
|
||||||
}
|
}
|
||||||
|
}
|
||||||
|
OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::linearize failed: ",
|
||||||
|
ebosSimulator_.gridView().comm());
|
||||||
return;
|
return;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
@ -359,7 +359,8 @@ namespace Opm
|
|||||||
DeferredLogger& deferred_logger) const;
|
DeferredLogger& deferred_logger) const;
|
||||||
|
|
||||||
void updatePrimaryVariablesNewton(const BVectorWell& dwells,
|
void updatePrimaryVariablesNewton(const BVectorWell& dwells,
|
||||||
const WellState& well_state) const;
|
const WellState& well_state,
|
||||||
|
DeferredLogger& deferred_logger) const;
|
||||||
|
|
||||||
// update extra primary vriables if there are any
|
// update extra primary vriables if there are any
|
||||||
void updateExtraPrimaryVariables(const BVectorWell& dwells) const;
|
void updateExtraPrimaryVariables(const BVectorWell& dwells) const;
|
||||||
|
@ -854,7 +854,7 @@ namespace Opm
|
|||||||
{
|
{
|
||||||
if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
|
if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
|
||||||
|
|
||||||
updatePrimaryVariablesNewton(dwells, well_state);
|
updatePrimaryVariablesNewton(dwells, well_state, deferred_logger);
|
||||||
|
|
||||||
updateWellStateFromPrimaryVariables(well_state, deferred_logger);
|
updateWellStateFromPrimaryVariables(well_state, deferred_logger);
|
||||||
Base::calculateReservoirRates(well_state.well(this->index_of_well_));
|
Base::calculateReservoirRates(well_state.well(this->index_of_well_));
|
||||||
@ -868,7 +868,8 @@ namespace Opm
|
|||||||
void
|
void
|
||||||
StandardWell<TypeTag>::
|
StandardWell<TypeTag>::
|
||||||
updatePrimaryVariablesNewton(const BVectorWell& dwells,
|
updatePrimaryVariablesNewton(const BVectorWell& dwells,
|
||||||
const WellState& /* well_state */) const
|
const WellState& /* well_state */,
|
||||||
|
DeferredLogger& deferred_logger) const
|
||||||
{
|
{
|
||||||
const double dFLimit = this->param_.dwell_fraction_max_;
|
const double dFLimit = this->param_.dwell_fraction_max_;
|
||||||
const double dBHPLimit = this->param_.dbhp_max_rel_;
|
const double dBHPLimit = this->param_.dbhp_max_rel_;
|
||||||
@ -876,11 +877,10 @@ namespace Opm
|
|||||||
|
|
||||||
updateExtraPrimaryVariables(dwells);
|
updateExtraPrimaryVariables(dwells);
|
||||||
|
|
||||||
#ifndef NDEBUG
|
|
||||||
for (double v : this->primary_variables_) {
|
for (double v : this->primary_variables_) {
|
||||||
assert(isfinite(v));
|
if(!isfinite(v))
|
||||||
|
OPM_DEFLOG_THROW(NumericalIssue, "Infinite primary variable after newton update well: " << this->name(), deferred_logger);
|
||||||
}
|
}
|
||||||
#endif
|
|
||||||
|
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -1918,11 +1918,10 @@ namespace Opm
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
#ifndef NDEBUG
|
|
||||||
for (double v : this->primary_variables_) {
|
for (double v : this->primary_variables_) {
|
||||||
assert(isfinite(v));
|
if(!isfinite(v))
|
||||||
|
OPM_DEFLOG_THROW(NumericalIssue, "Infinite primary variable after update from wellState well: " << this->name(), deferred_logger);
|
||||||
}
|
}
|
||||||
#endif
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
@ -330,8 +330,15 @@ namespace Opm
|
|||||||
const auto& summary_state = ebosSimulator.vanguard().summaryState();
|
const auto& summary_state = ebosSimulator.vanguard().summaryState();
|
||||||
const auto inj_controls = this->well_ecl_.isInjector() ? this->well_ecl_.injectionControls(summary_state) : Well::InjectionControls(0);
|
const auto inj_controls = this->well_ecl_.isInjector() ? this->well_ecl_.injectionControls(summary_state) : Well::InjectionControls(0);
|
||||||
const auto prod_controls = this->well_ecl_.isProducer() ? this->well_ecl_.productionControls(summary_state) : Well::ProductionControls(0);
|
const auto prod_controls = this->well_ecl_.isProducer() ? this->well_ecl_.productionControls(summary_state) : Well::ProductionControls(0);
|
||||||
|
bool converged = false;
|
||||||
return this->iterateWellEqWithControl(ebosSimulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
|
try {
|
||||||
|
converged = this->iterateWellEqWithControl(ebosSimulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
|
||||||
|
} catch (NumericalIssue& e ) {
|
||||||
|
const std::string msg = "Inner well iterations failed for well " + this->name() + " Treat the well as unconverged. ";
|
||||||
|
deferred_logger.warning("INNER_ITERATION_FAILED", msg);
|
||||||
|
converged = false;
|
||||||
|
}
|
||||||
|
return converged;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Loading…
Reference in New Issue
Block a user