28#ifndef OPM_FLASH_NEWTON_METHOD_HH
29#define OPM_FLASH_NEWTON_METHOD_HH
33#include <opm/common/Exceptions.hpp>
39template <
class TypeTag,
class MyTypeTag>
40struct DiscNewtonMethod;
51template <
class TypeTag>
62 enum { pressure0Idx = Indices::pressure0Idx };
63 enum { z0Idx = Indices::z0Idx };
64 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
81 PrimaryVariables& nextValue,
82 const PrimaryVariables& currentValue,
83 const EqVector& update,
87 nextValue = currentValue;
94 constexpr Scalar max_percent_change = 0.2;
95 constexpr Scalar upper_bound = 1. + max_percent_change;
96 constexpr Scalar lower_bound = 1. - max_percent_change;
97 clampValue_(nextValue[pressure0Idx],
98 currentValue[pressure0Idx] * lower_bound,
99 currentValue[pressure0Idx] * upper_bound);
105 Scalar maxDeltaZ = 0.0;
106 Scalar sumDeltaZ = 0.0;
107 for (
unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
108 maxDeltaZ = std::max(std::abs(update[z0Idx + compIdx]), maxDeltaZ);
109 sumDeltaZ += update[z0Idx + compIdx];
111 maxDeltaZ = std::max(std::abs(sumDeltaZ), maxDeltaZ);
115 constexpr Scalar deltaz_limit = 0.2;
116 if (maxDeltaZ > deltaz_limit) {
117 Scalar alpha = deltaz_limit / maxDeltaZ;
118 for (
unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
119 nextValue[z0Idx + compIdx] = currentValue[z0Idx + compIdx] - alpha * update[z0Idx + compIdx];
125 for (
unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
126 clampValue_(nextValue[z0Idx + compIdx], tol, 1-tol);
130 void clampValue_(Scalar& val, Scalar minVal, Scalar maxVal)
const
132 val = std::clamp(val, minVal, maxVal);
A Newton solver specific to the PTFlash model.
Definition: flashnewtonmethod.hh:53
FlashNewtonMethod(Simulator &simulator)
Definition: flashnewtonmethod.hh:70
void updatePrimaryVariables_(unsigned, PrimaryVariables &nextValue, const PrimaryVariables ¤tValue, const EqVector &update, const EqVector &)
Update a single primary variables object.
Definition: flashnewtonmethod.hh:80
friend ParentType
Definition: flashnewtonmethod.hh:74
The multi-dimensional Newton method.
Definition: newtonmethod.hh:92
Definition: blackoilmodel.hh:72
Definition: blackoilboundaryratevector.hh:37
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:235