28#ifndef EWOMS_NCP_NEWTON_METHOD_HH
29#define EWOMS_NCP_NEWTON_METHOD_HH
31#include <opm/common/Exceptions.hpp>
42template <
class TypeTag,
class MyTypeTag>
43struct DiscNewtonMethod;
54template <
class TypeTag>
67 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
68 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
69 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
70 enum { fugacity0Idx = Indices::fugacity0Idx };
71 enum { saturation0Idx = Indices::saturation0Idx };
72 enum { pressure0Idx = Indices::pressure0Idx };
73 enum { ncp0EqIdx = Indices::ncp0EqIdx };
87 const GlobalEqVector& currentResidual)
89 const auto& constraintsMap = this->model().linearizer().constraintsMap();
90 this->lastError_ = this->error_;
95 for (
unsigned dofIdx = 0; dofIdx < currentResidual.size(); ++dofIdx) {
97 if (dofIdx >= this->model().numGridDof() || this->model().dofTotalVolume(dofIdx) <= 0.0) {
102 if (this->enableConstraints_()) {
103 if (constraintsMap.count(dofIdx) > 0) {
108 const auto& r = currentResidual[dofIdx];
109 for (
unsigned eqIdx = 0; eqIdx < r.size(); ++eqIdx) {
110 if (ncp0EqIdx <= eqIdx && eqIdx < Indices::ncp0EqIdx + numPhases) {
114 std::max(std::abs(r[eqIdx] * this->model().eqWeight(dofIdx, eqIdx)),
120 this->error_ = this->comm_.max(this->error_);
125 throw Opm::NumericalProblem(
"Newton: Error " +
std::to_string(
double(this->error_)) +
126 " is larger than maximum allowed error of " +
135 PrimaryVariables& nextValue,
136 const PrimaryVariables& currentValue,
137 const EqVector& update,
141 nextValue = currentValue;
149 Scalar sumSatDelta = 0.0;
150 Scalar maxSatDelta = 0.0;
151 for (
unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx) {
152 maxSatDelta = std::max(std::abs(update[saturation0Idx + phaseIdx]),
154 sumSatDelta += update[saturation0Idx + phaseIdx];
156 maxSatDelta = std::max(std::abs(-sumSatDelta), maxSatDelta);
158 if (maxSatDelta > 0.2) {
159 const Scalar alpha = 0.2 / maxSatDelta;
160 for (
unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx) {
161 nextValue[saturation0Idx + phaseIdx] =
162 currentValue[saturation0Idx + phaseIdx] -
163 alpha * update[saturation0Idx + phaseIdx];
168 clampValue_(nextValue[pressure0Idx],
169 currentValue[pressure0Idx] * 0.8,
170 currentValue[pressure0Idx] * 1.2);
173 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
174 Scalar& val = nextValue[fugacity0Idx + compIdx];
175 const Scalar oldVal = currentValue[fugacity0Idx + compIdx];
180 const Scalar minPhi =
181 std::max(this->problem().model().minActivityCoeff(globalDofIdx, compIdx),
182 0.001 * currentValue[pressure0Idx]);
186 const Scalar maxDelta = 0.7 * minPhi;
187 clampValue_(val, oldVal - maxDelta, oldVal + maxDelta);
190 val = std::max(val, 0.0);
195 if (this->numIterations_ < 3) {
197 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
198 Scalar& val = nextValue[fugacity0Idx + compIdx];
199 const Scalar oldVal = currentValue[fugacity0Idx + compIdx];
200 const Scalar minPhi = this->problem().model().minActivityCoeff(globalDofIdx, compIdx);
201 if (oldVal < 1.0 * minPhi && val > 1.0 * minPhi) {
204 else if (oldVal > 0.0 && val < 0.0) {
210 for (
unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx) {
211 Scalar& val = nextValue[saturation0Idx + phaseIdx];
212 const Scalar oldVal = currentValue[saturation0Idx + phaseIdx];
213 if (oldVal < 1.0 && val > 1.0) {
216 else if (oldVal > 0.0 && val < 0.0) {
224 void clampValue_(Scalar& val, Scalar minVal, Scalar maxVal)
const
225 { val = std::max(minVal, std::min(val, maxVal)); }
A Newton solver specific to the NCP model.
Definition: ncpnewtonmethod.hh:56
friend ParentType
Definition: ncpnewtonmethod.hh:83
NcpNewtonMethod(Simulator &simulator)
Definition: ncpnewtonmethod.hh:79
void updatePrimaryVariables_(unsigned globalDofIdx, PrimaryVariables &nextValue, const PrimaryVariables ¤tValue, const EqVector &update, const EqVector &)
Update a single primary variables object.
Definition: ncpnewtonmethod.hh:134
void preSolve_(const SolutionVector &, const GlobalEqVector ¤tResidual)
Definition: ncpnewtonmethod.hh:86
The multi-dimensional Newton method.
Definition: newtonmethod.hh:100
Defines the common properties required by the porous medium multi-phase models.
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:185
Definition: blackoilmodel.hh:79
Definition: blackoilboundaryratevector.hh:39
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: newtonmethodparams.hpp:31