ncpnewtonmethod.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 EWOMS_NCP_NEWTON_METHOD_HH
29#define EWOMS_NCP_NEWTON_METHOD_HH
30
31#include <opm/common/Exceptions.hpp>
32
35
36#include <algorithm>
37#include <cmath>
38#include <string>
39
40namespace Opm::Properties {
41
42template <class TypeTag, class MyTypeTag>
43struct DiscNewtonMethod;
44
45} // namespace Opm::Properties
46
47namespace Opm {
48
54template <class TypeTag>
55class NcpNewtonMethod : public GetPropType<TypeTag, Properties::DiscNewtonMethod>
56{
58
66
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 };
74
75public:
79 explicit NcpNewtonMethod(Simulator& simulator) : ParentType(simulator)
80 {}
81
82protected:
83 friend ParentType;
85
86 void preSolve_(const SolutionVector&,
87 const GlobalEqVector& currentResidual)
88 {
89 const auto& constraintsMap = this->model().linearizer().constraintsMap();
90 this->lastError_ = this->error_;
91
92 // calculate the error as the maximum weighted tolerance of
93 // the solution's residual
94 this->error_ = 0;
95 for (unsigned dofIdx = 0; dofIdx < currentResidual.size(); ++dofIdx) {
96 // do not consider auxiliary DOFs for the error
97 if (dofIdx >= this->model().numGridDof() || this->model().dofTotalVolume(dofIdx) <= 0.0) {
98 continue;
99 }
100
101 // also do not consider DOFs which are constraint
102 if (this->enableConstraints_()) {
103 if (constraintsMap.count(dofIdx) > 0) {
104 continue;
105 }
106 }
107
108 const auto& r = currentResidual[dofIdx];
109 for (unsigned eqIdx = 0; eqIdx < r.size(); ++eqIdx) {
110 if (ncp0EqIdx <= eqIdx && eqIdx < Indices::ncp0EqIdx + numPhases) {
111 continue;
112 }
113 this->error_ =
114 std::max(std::abs(r[eqIdx] * this->model().eqWeight(dofIdx, eqIdx)),
115 this->error_);
116 }
117 }
118
119 // take the other processes into account
120 this->error_ = this->comm_.max(this->error_);
121
122 // make sure that the error never grows beyond the maximum
123 // allowed one
125 throw Opm::NumericalProblem("Newton: Error " + std::to_string(double(this->error_)) +
126 " is larger than maximum allowed error of " +
128 }
129 }
130
134 void updatePrimaryVariables_(unsigned globalDofIdx,
135 PrimaryVariables& nextValue,
136 const PrimaryVariables& currentValue,
137 const EqVector& update,
138 const EqVector&)
139 {
140 // normal Newton-Raphson update
141 nextValue = currentValue;
142 nextValue -= update;
143
145 // put crash barriers along the update path
147
148 // saturations: limit the change of any saturation to at most 20%
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]),
153 maxSatDelta);
154 sumSatDelta += update[saturation0Idx + phaseIdx];
155 }
156 maxSatDelta = std::max(std::abs(-sumSatDelta), maxSatDelta);
157
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];
164 }
165 }
166
167 // limit pressure reference change to 20% of the total value per iteration
168 clampValue_(nextValue[pressure0Idx],
169 currentValue[pressure0Idx] * 0.8,
170 currentValue[pressure0Idx] * 1.2);
171
172 // fugacities
173 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
174 Scalar& val = nextValue[fugacity0Idx + compIdx];
175 const Scalar oldVal = currentValue[fugacity0Idx + compIdx];
176
177 // get the minimum activity coefficient for the component (i.e., the activity
178 // coefficient of the phase for which the component has the highest affinity)
179 // Make sure that the activity coefficient does not get too small.
180 const Scalar minPhi =
181 std::max(this->problem().model().minActivityCoeff(globalDofIdx, compIdx),
182 0.001 * currentValue[pressure0Idx]);
183
184 // allow the mole fraction of the component to change at most 70% in any
185 // phase (assuming composition independent fugacity coefficients).
186 const Scalar maxDelta = 0.7 * minPhi;
187 clampValue_(val, oldVal - maxDelta, oldVal + maxDelta);
188
189 // make sure that fugacities do not become negative
190 val = std::max(val, 0.0);
191 }
192
193 // do not become grossly unphysical in a single iteration for the first few
194 // iterations of a time step
195 if (this->numIterations_ < 3) {
196 // fugacities
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) {
202 val = 1.0 * minPhi;
203 }
204 else if (oldVal > 0.0 && val < 0.0) {
205 val = 0.0;
206 }
207 }
208
209 // saturations
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) {
214 val = 1.0;
215 }
216 else if (oldVal > 0.0 && val < 0.0) {
217 val = 0.0;
218 }
219 }
220 }
221 }
222
223private:
224 void clampValue_(Scalar& val, Scalar minVal, Scalar maxVal) const
225 { val = std::max(minVal, std::min(val, maxVal)); }
226};
227} // namespace Opm
228
229#endif
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 &currentValue, const EqVector &update, const EqVector &)
Update a single primary variables object.
Definition: ncpnewtonmethod.hh:134
void preSolve_(const SolutionVector &, const GlobalEqVector &currentResidual)
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