mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
NonlinearSolverEbos: put detectOscillations in compile unit
This commit is contained in:
63
opm/simulators/flow/NonlinearSolverEbos.cpp
Normal file
63
opm/simulators/flow/NonlinearSolverEbos.cpp
Normal file
@@ -0,0 +1,63 @@
|
||||
/*
|
||||
Copyright 2015 SINTEF ICT, Applied Mathematics.
|
||||
Copyright 2015 Statoil ASA.
|
||||
|
||||
This file is part of the Open Porous Media project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
#include <config.h>
|
||||
#include <opm/simulators/flow/NonlinearSolverEbos.hpp>
|
||||
|
||||
#include <cmath>
|
||||
|
||||
namespace Opm {
|
||||
namespace detail {
|
||||
|
||||
void detectOscillations(const std::vector<std::vector<double>>& residualHistory,
|
||||
const int it, const int numPhases, const double relaxRelTol,
|
||||
bool& oscillate, bool& stagnate)
|
||||
{
|
||||
// The detection of oscillation in two primary variable results in the report of the detection
|
||||
// of oscillation for the solver.
|
||||
// Only the saturations are used for oscillation detection for the black oil model.
|
||||
// Stagnate is not used for any treatment here.
|
||||
|
||||
if (it < 2) {
|
||||
oscillate = false;
|
||||
stagnate = false;
|
||||
return;
|
||||
}
|
||||
|
||||
stagnate = true;
|
||||
int oscillatePhase = 0;
|
||||
const auto& F0 = residualHistory[it];
|
||||
const auto& F1 = residualHistory[it - 1];
|
||||
const auto& F2 = residualHistory[it - 2];
|
||||
for (int p = 0; p < numPhases; ++p) {
|
||||
const double d1 = std::abs((F0[p] - F2[p]) / F0[p]);
|
||||
const double d2 = std::abs((F0[p] - F1[p]) / F0[p]);
|
||||
|
||||
oscillatePhase += (d1 < relaxRelTol) && (relaxRelTol < d2);
|
||||
|
||||
// Process is 'stagnate' unless at least one phase
|
||||
// exhibits significant residual change.
|
||||
stagnate = (stagnate && !(std::abs((F1[p] - F2[p]) / F2[p]) > 1.0e-3));
|
||||
}
|
||||
|
||||
oscillate = (oscillatePhase > 1);
|
||||
}
|
||||
|
||||
} // namespace detail
|
||||
} // namespace Opm
|
@@ -80,6 +80,15 @@ struct NewtonRelaxationType<TypeTag, TTag::FlowNonLinearSolver> {
|
||||
|
||||
namespace Opm {
|
||||
|
||||
namespace detail {
|
||||
|
||||
/// Detect oscillation or stagnation in a given residual history.
|
||||
void detectOscillations(const std::vector<std::vector<double>>& residualHistory,
|
||||
const int it, const int numPhases, const double relaxRelTol,
|
||||
bool& oscillate, bool& stagnate);
|
||||
|
||||
}
|
||||
|
||||
class WellState;
|
||||
|
||||
/// A nonlinear solver class suitable for general fully-implicit models,
|
||||
@@ -279,34 +288,8 @@ class WellState;
|
||||
void detectOscillations(const std::vector<std::vector<double>>& residualHistory,
|
||||
const int it, bool& oscillate, bool& stagnate) const
|
||||
{
|
||||
// The detection of oscillation in two primary variable results in the report of the detection
|
||||
// of oscillation for the solver.
|
||||
// Only the saturations are used for oscillation detection for the black oil model.
|
||||
// Stagnate is not used for any treatment here.
|
||||
|
||||
if ( it < 2 ) {
|
||||
oscillate = false;
|
||||
stagnate = false;
|
||||
return;
|
||||
}
|
||||
|
||||
stagnate = true;
|
||||
int oscillatePhase = 0;
|
||||
const std::vector<double>& F0 = residualHistory[it];
|
||||
const std::vector<double>& F1 = residualHistory[it - 1];
|
||||
const std::vector<double>& F2 = residualHistory[it - 2];
|
||||
for (int p= 0; p < model_->numPhases(); ++p){
|
||||
const double d1 = std::abs((F0[p] - F2[p]) / F0[p]);
|
||||
const double d2 = std::abs((F0[p] - F1[p]) / F0[p]);
|
||||
|
||||
oscillatePhase += (d1 < relaxRelTol()) && (relaxRelTol() < d2);
|
||||
|
||||
// Process is 'stagnate' unless at least one phase
|
||||
// exhibits significant residual change.
|
||||
stagnate = (stagnate && !(std::abs((F1[p] - F2[p]) / F2[p]) > 1.0e-3));
|
||||
}
|
||||
|
||||
oscillate = (oscillatePhase > 1);
|
||||
detail::detectOscillations(residualHistory, it, model_->numPhases(),
|
||||
this->relaxRelTol(), oscillate, stagnate);
|
||||
}
|
||||
|
||||
|
||||
@@ -392,6 +375,7 @@ class WellState;
|
||||
int linearIterationsLast_;
|
||||
int wellIterationsLast_;
|
||||
};
|
||||
|
||||
} // namespace Opm
|
||||
|
||||
#endif // OPM_NON_LINEAR_SOLVER_EBOS_HPP
|
||||
|
Reference in New Issue
Block a user