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
66public:
70 FlashNewtonMethod(Simulator& simulator) : ParentType(simulator)
71 {}
72
73protected:
74 friend ParentType;
76
80 void updatePrimaryVariables_(unsigned /* globalDofIdx */,
81 PrimaryVariables& nextValue,
82 const PrimaryVariables& currentValue,
83 const EqVector& update,
84 const EqVector& /* currentResidual */)
85 {
86 // normal Newton-Raphson update
87 nextValue = currentValue;
88 nextValue -= update;
89
91 // Pressure updates
93 // limit pressure reference change to 20% of the total value per iteration
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);
100
102 // z updates
104 // restrict update to at most 0.1
105 Scalar maxDeltaZ = 0.0; // in update vector
106 Scalar sumDeltaZ = 0.0; // changes in last component (not in update vector)
107 for (unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
108 maxDeltaZ = std::max(std::abs(update[z0Idx + compIdx]), maxDeltaZ);
109 sumDeltaZ += update[z0Idx + compIdx];
110 }
111 maxDeltaZ = std::max(std::abs(sumDeltaZ), maxDeltaZ);
112
113 // 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)
114 // \Note: original code uses 0.1, while 0.1 looks like having problem make it converged. So there is some more to investigate here
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];
120 }
121 }
122
123 // ensure that z-values are less than tol or more than 1-tol
124 Scalar tol = 1e-8;
125 for (unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
126 clampValue_(nextValue[z0Idx + compIdx], tol, 1-tol);
127 }
128 }
129private:
130 void clampValue_(Scalar& val, Scalar minVal, Scalar maxVal) const
131 {
132 val = std::clamp(val, minVal, maxVal);
133 }
134
135}; // class FlashNewtonMethod
136} // namespace Opm
137#endif
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 &currentValue, 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