mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Merge pull request #182 from atgeirr/bugfix-updatestate
Bugfix for updateState()
This commit is contained in:
commit
fed1653e44
@ -37,6 +37,7 @@
|
||||
|
||||
#include "reenable_warning_pragmas.h"
|
||||
|
||||
#include <opm/core/utility/ErrorMacros.hpp>
|
||||
|
||||
namespace Opm
|
||||
{
|
||||
@ -161,8 +162,8 @@ namespace Opm
|
||||
// Solve system.
|
||||
Dune::InverseOperatorResult result;
|
||||
linsolve.apply(x, de, result);
|
||||
if (result.converged) {
|
||||
// std::cout << "solveElliptic() successful!" << std::endl;
|
||||
if (!result.converged) {
|
||||
OPM_THROW(std::runtime_error, "CPRPreconditioner failed to solve elliptic subsystem.");
|
||||
}
|
||||
return x;
|
||||
}
|
||||
|
@ -642,6 +642,15 @@ namespace {
|
||||
// Create the primary variables.
|
||||
SolutionState state = variableState(x, xw);
|
||||
|
||||
// DISKVAL(state.pressure);
|
||||
// DISKVAL(state.saturation[0]);
|
||||
// DISKVAL(state.saturation[1]);
|
||||
// DISKVAL(state.saturation[2]);
|
||||
// DISKVAL(state.rs);
|
||||
// DISKVAL(state.rv);
|
||||
// DISKVAL(state.qs);
|
||||
// DISKVAL(state.bhp);
|
||||
|
||||
// -------- Mass balance equations --------
|
||||
|
||||
// Compute b_p and the accumulation term b_p*s_p for each phase,
|
||||
@ -1345,7 +1354,7 @@ namespace {
|
||||
for (int c = 0; c < nc; ++c) {
|
||||
if (ixw[c]) {
|
||||
so[c] = so[c] / (1-sw[c]);
|
||||
sg[c] = sg[c] / (1-so[c]);
|
||||
sg[c] = sg[c] / (1-sw[c]);
|
||||
sw[c] = 0;
|
||||
}
|
||||
}
|
||||
@ -1563,8 +1572,10 @@ namespace {
|
||||
for (; quantityIt != endQuantityIt; ++quantityIt) {
|
||||
const double quantityResid = (*quantityIt).value().matrix().norm();
|
||||
if (!std::isfinite(quantityResid)) {
|
||||
const int trouble_phase = quantityIt - residual_.material_balance_eq.begin();
|
||||
OPM_THROW(Opm::NumericalProblem,
|
||||
"Encountered a non-finite residual");
|
||||
"Encountered a non-finite residual in material balance equation "
|
||||
<< trouble_phase);
|
||||
}
|
||||
globalNorm = std::max(globalNorm, quantityResid);
|
||||
}
|
||||
|
Loading…
Reference in New Issue
Block a user