My Project
Loading...
Searching...
No Matches
flashnewtonmethod.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef OPM_FLASH_NEWTON_METHOD_HH
29#define OPM_FLASH_NEWTON_METHOD_HH
30
32
33#include <opm/common/Exceptions.hpp>
34
35#include <algorithm>
36
37namespace Opm::Properties {
38
39template <class TypeTag, class MyTypeTag>
40struct DiscNewtonMethod;
41
42} // namespace Opm::Properties
43
44namespace Opm {
51template <class TypeTag>
52class FlashNewtonMethod : public GetPropType<TypeTag, Properties::DiscNewtonMethod>
53{
55
61
62 enum { pressure0Idx = Indices::pressure0Idx };
63 enum { z0Idx = Indices::z0Idx };
64 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
65
66 static constexpr bool waterEnabled = Indices::waterEnabled;
67
68public:
72 explicit FlashNewtonMethod(Simulator& simulator) : ParentType(simulator)
73 {}
74
75protected:
76 friend ParentType;
78
82 void updatePrimaryVariables_(unsigned /* globalDofIdx */,
83 PrimaryVariables& nextValue,
84 const PrimaryVariables& currentValue,
85 const EqVector& update,
86 const EqVector& /* currentResidual */)
87 {
88 // normal Newton-Raphson update
90 nextValue -= update;
91
93 // Pressure updates
95 // limit pressure reference change to 20% of the total value per iteration
96 constexpr Scalar max_percent_change = 0.2;
97 constexpr Scalar upper_bound = 1. + max_percent_change;
98 constexpr Scalar lower_bound = 1. - max_percent_change;
99 clampValue_(nextValue[pressure0Idx],
100 currentValue[pressure0Idx] * lower_bound,
101 currentValue[pressure0Idx] * upper_bound);
102
104 // z updates
106 // restrict update to at most 0.1
107 Scalar maxDeltaZ = 0.0; // in update vector
108 Scalar sumDeltaZ = 0.0; // changes in last component (not in update vector)
109 for (unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
110 maxDeltaZ = std::max(std::abs(update[z0Idx + compIdx]), maxDeltaZ);
111 sumDeltaZ += update[z0Idx + compIdx];
112 }
113 maxDeltaZ = std::max(std::abs(sumDeltaZ), maxDeltaZ);
114
115 // if max. update is above 0.2, restrict that one to 0.2 and adjust the rest accordingly (s.t. last comp. update is sum of the changes in update vector)
116 // \Note: original code uses 0.1, while 0.1 looks like having problem make it converged. So there is some more to investigate here
117 constexpr Scalar deltaz_limit = 0.2;
118 if (maxDeltaZ > deltaz_limit) {
119 Scalar alpha = deltaz_limit / maxDeltaZ;
120 for (unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
121 nextValue[z0Idx + compIdx] = currentValue[z0Idx + compIdx] - alpha * update[z0Idx + compIdx];
122 }
123 }
124
125 // ensure that z-values are less than tol or more than 1-tol
126 Scalar tol = 1e-8;
127 for (unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
128 clampValue_(nextValue[z0Idx + compIdx], tol, 1-tol);
129 }
130
131 if constexpr (waterEnabled) {
132 // limit change in water saturation to 0.2
133 constexpr Scalar dSwMax = 0.2;
134 if (update[Indices::water0Idx] > dSwMax) {
135 nextValue[Indices::water0Idx] = currentValue[Indices::water0Idx] - dSwMax;
136 }
137 }
138 }
139private:
140 void clampValue_(Scalar& val, Scalar minVal, Scalar maxVal) const
141 {
142 val = std::clamp(val, minVal, maxVal);
143 }
144
145}; // class FlashNewtonMethod
146} // namespace Opm
147#endif
A Newton solver specific to the PTFlash model.
Definition flashnewtonmethod.hh:53
FlashNewtonMethod(Simulator &simulator)
Definition flashnewtonmethod.hh:72
void updatePrimaryVariables_(unsigned, PrimaryVariables &nextValue, const PrimaryVariables &currentValue, const EqVector &update, const EqVector &)
Update a single primary variables object.
Definition flashnewtonmethod.hh:82
The multi-dimensional Newton method.
Definition newtonmethod.hh:92
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilboundaryratevector.hh:37
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242
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
The multi-dimensional Newton method.