28#ifndef OPM_FLASH_NEWTON_METHOD_HH
29#define OPM_FLASH_NEWTON_METHOD_HH
31#include <opm/common/Exceptions.hpp>
41template <
class TypeTag,
class MyTypeTag>
42struct DiscNewtonMethod;
53template <
class TypeTag>
64 enum { pressure0Idx = Indices::pressure0Idx };
65 enum { z0Idx = Indices::z0Idx };
66 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
68 static constexpr bool waterEnabled = Indices::waterEnabled;
85 PrimaryVariables& nextValue,
86 const PrimaryVariables& currentValue,
87 const EqVector& update,
91 nextValue = currentValue;
98 constexpr Scalar max_percent_change = 0.2;
99 constexpr Scalar upper_bound = 1. + max_percent_change;
100 constexpr Scalar lower_bound = 1. - max_percent_change;
101 nextValue[pressure0Idx] = std::clamp(nextValue[pressure0Idx],
102 currentValue[pressure0Idx] * lower_bound,
103 currentValue[pressure0Idx] * upper_bound);
109 Scalar maxDeltaZ = 0.0;
110 Scalar sumDeltaZ = 0.0;
111 for (
unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
112 maxDeltaZ = std::max(std::abs(update[z0Idx + compIdx]), maxDeltaZ);
113 sumDeltaZ += update[z0Idx + compIdx];
115 maxDeltaZ = std::max(std::abs(sumDeltaZ), maxDeltaZ);
121 constexpr Scalar deltaz_limit = 0.2;
122 if (maxDeltaZ > deltaz_limit) {
123 const Scalar alpha = deltaz_limit / maxDeltaZ;
124 for (
unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
125 nextValue[z0Idx + compIdx] = currentValue[z0Idx + compIdx] - alpha * update[z0Idx + compIdx];
130 constexpr Scalar tol = 1e-8;
131 for (
unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
132 nextValue[z0Idx + compIdx] = std::clamp(nextValue[z0Idx + compIdx], tol, 1-tol);
135 if constexpr (waterEnabled) {
137 constexpr Scalar dSwMax = 0.2;
138 if (update[Indices::water0Idx] > dSwMax) {
139 nextValue[Indices::water0Idx] = currentValue[Indices::water0Idx] - dSwMax;
A Newton solver specific to the PTFlash model.
Definition: flashnewtonmethod.hh:55
FlashNewtonMethod(Simulator &simulator)
Definition: flashnewtonmethod.hh:74
void updatePrimaryVariables_(unsigned, PrimaryVariables &nextValue, const PrimaryVariables ¤tValue, const EqVector &update, const EqVector &)
Update a single primary variables object.
Definition: flashnewtonmethod.hh:84
friend ParentType
Definition: flashnewtonmethod.hh:78
The multi-dimensional Newton method.
Definition: newtonmethod.hh:100
Defines the common properties required by the porous medium multi-phase models.
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