GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/nonlinear/findscalarroot.hh
Date: 2024-05-04 19:09:25
Exec Total Coverage
Lines: 42 60 70.0%
Functions: 5 9 55.6%
Branches: 84 213 39.4%

Line Branch Exec Source
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 // SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5 // SPDX-License-Identifier: GPL-3.0-or-later
6 //
7 /*!
8 * \file
9 * \ingroup Nonlinear
10 * \brief Root finding algorithms for scalar functions
11 */
12 #ifndef DUMUX_COMMON_SCALAR_ROOT_FINDING_HH
13 #define DUMUX_COMMON_SCALAR_ROOT_FINDING_HH
14
15 #include <cmath>
16 #include <limits>
17 #include <type_traits>
18
19 #include <dumux/common/exceptions.hh>
20 #include <dumux/common/parameters.hh>
21 #include <dumux/common/numericdifferentiation.hh>
22
23 namespace Dumux {
24
25 /*!
26 * \ingroup Nonlinear
27 * \brief Newton's root finding algorithm for scalar functions (secant method)
28 * \param xOld initial guess
29 * \param residual Residual function
30 * \param derivative Derivative of the residual
31 * \param tol Relative shift tolerance
32 * \param maxIter Maximum number of iterations
33 */
34 template<class Scalar, class ResFunc, class DerivFunc,
35 typename std::enable_if_t<std::is_invocable_r_v<Scalar, ResFunc, Scalar> &&
36 std::is_invocable_r_v<Scalar, DerivFunc, Scalar>>...>
37 Scalar findScalarRootNewton(Scalar xOld, const ResFunc& residual, const DerivFunc& derivative,
38 const Scalar tol = 1e-13, const int maxIter = 200)
39 {
40 Scalar xNew = xOld;
41 Scalar r = residual(xNew);
42
43 int n = maxIter;
44 Scalar relativeShift = std::numeric_limits<Scalar>::max();
45 while (relativeShift > tol && n > 0)
46 {
47 xNew = xOld - r/derivative(xOld);
48 r = residual(xNew);
49
50 using std::abs; using std::max;
51 relativeShift = abs(xOld-xNew)/max(abs(xOld), abs(xNew));
52 xOld = xNew;
53 n--;
54 }
55
56 using std::isfinite;
57 if (!isfinite(r))
58 DUNE_THROW(NumericalProblem, "Residual is not finite: " << r << " after " << maxIter - n << " iterations!");
59
60 if (relativeShift > tol)
61 DUNE_THROW(NumericalProblem, "Scalar newton solver did not converge after " << maxIter << " iterations!");
62
63 return xNew;
64 }
65
66 /*!
67 * \ingroup Nonlinear
68 * \brief Newton's root finding algorithm for scalar functions (secant method)
69 * \note The derivative is numerically computed. If the derivative is know use signature with derivative function.
70 * \param xOld initial guess
71 * \param residual Residual function
72 * \param tol Relative shift tolerance
73 * \param maxIter Maximum number of iterations
74 */
75 template<class Scalar, class ResFunc,
76 typename std::enable_if_t<std::is_invocable_r_v<Scalar, ResFunc, Scalar>>...>
77 Scalar findScalarRootNewton(Scalar xOld, const ResFunc& residual,
78 const Scalar tol = 1e-13, const int maxIter = 200)
79 {
80 1 const auto eps = NumericDifferentiation::epsilon(xOld);
81 auto derivative = [&](const auto x){ return (residual(x + eps)-residual(x))/eps; };
82 1 return findScalarRootNewton(xOld, residual, derivative, tol, maxIter);
83 }
84
85 /*!
86 * \ingroup Nonlinear
87 * \brief Brent's root finding algorithm for scalar functions
88 * \note Modified from pseudo-code on wikipedia: https://en.wikipedia.org/wiki/Brent%27s_method
89 * \note See also R.P. Brent "An algorithm with guaranteed convergence for finding a zero of a function", The Computer Journal (1971).
90 * \note This is usually more robust than Newton's method
91 * \param a Lower bound
92 * \param b Upper bound
93 * \param residual Residual function
94 * \param tol Relative shift tolerance
95 * \param maxIter Maximum number of iterations
96 */
97 template<class Scalar, class ResFunc,
98 typename std::enable_if_t<std::is_invocable_r_v<Scalar, ResFunc, Scalar>>...>
99 10809086 Scalar findScalarRootBrent(Scalar a, Scalar b, const ResFunc& residual,
100 const Scalar tol = 1e-13, const int maxIter = 200)
101 {
102 // precompute the residuals (minimize function evaluations)
103
0/2
✗ Branch 0 not taken.
✗ Branch 1 not taken.
10809086 Scalar fa = residual(a);
104
0/2
✗ Branch 0 not taken.
✗ Branch 1 not taken.
10809086 Scalar fb = residual(b);
105 10809086 Scalar fs = 0.0;
106
107 // check if the root is inside the given interval
108 using std::signbit;
109
4/4
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 10809085 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 10809085 times.
21618172 if (!signbit(fa*fb))
110
7/28
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✓ Branch 11 taken 1 times.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 1 times.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✓ Branch 18 taken 1 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 1 times.
✗ Branch 22 not taken.
✓ Branch 23 taken 1 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 1 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
✗ Branch 33 not taken.
10 DUNE_THROW(NumericalProblem, "Brent's algorithm failed: [a,b] does not contain any, or no uniquely defined root!");
111
112 // sort bounds
113 using std::abs; using std::swap;
114
6/6
✓ Branch 0 taken 9779696 times.
✓ Branch 1 taken 1029389 times.
✓ Branch 2 taken 9779696 times.
✓ Branch 3 taken 1029389 times.
✓ Branch 4 taken 9779696 times.
✓ Branch 5 taken 1029389 times.
32427255 if (abs(fa) < abs(fb))
115 {
116 9779696 swap(a, b);
117 9779696 swap(fa, fb);
118 }
119
120 10809085 Scalar c = a;
121 10809085 Scalar fc = fa;
122 10809085 Scalar d = 0.0;
123 10809085 Scalar s = 0.0;
124 10809085 bool flag = true;
125
126
1/2
✓ Branch 0 taken 261200795 times.
✗ Branch 1 not taken.
261200795 for (int i = 0; i < maxIter; ++i)
127 {
128 // stopping criterion
129 using std::max;
130
10/10
✓ Branch 0 taken 107460204 times.
✓ Branch 1 taken 153740591 times.
✓ Branch 2 taken 107460204 times.
✓ Branch 3 taken 153740591 times.
✓ Branch 4 taken 107460204 times.
✓ Branch 5 taken 153740591 times.
✓ Branch 6 taken 107460204 times.
✓ Branch 7 taken 153740591 times.
✓ Branch 8 taken 10809085 times.
✓ Branch 9 taken 250391710 times.
1152263384 if (abs(b-a) < tol*max(abs(a), abs(b)))
131 10809085 return b;
132
133 // inverse quadratic interpolation
134
4/4
✓ Branch 0 taken 206536472 times.
✓ Branch 1 taken 43855238 times.
✓ Branch 2 taken 34625364 times.
✓ Branch 3 taken 171911108 times.
250391710 if (fa != fc && fb != fc)
135 {
136 34625364 const auto fab = fa-fb;
137 34625364 const auto fac = fa-fc;
138 34625364 const auto fbc = fb-fc;
139 34625364 s = a*fb*fc/(fab*fac) - b*fa*fc/(fab*fbc) + c*fa*fb/(fac*fbc);
140 }
141
142 // secant method
143 else
144 {
145 215766346 s = b - fb*(b-a)/(fb-fa);
146 }
147
148 // bisection method
149
2/2
✓ Branch 0 taken 34918683 times.
✓ Branch 1 taken 12147761 times.
47066444 if ( (s < (3*a + b)*0.25 || s > b)
150
8/8
✓ Branch 0 taken 23530110 times.
✓ Branch 1 taken 11388573 times.
✓ Branch 2 taken 18645349 times.
✓ Branch 3 taken 4884761 times.
✓ Branch 4 taken 18645349 times.
✓ Branch 5 taken 4884761 times.
✓ Branch 6 taken 18645349 times.
✓ Branch 7 taken 4884761 times.
34918683 || (flag && abs(s-b) >= abs(b-c)*0.5)
151
8/8
✓ Branch 0 taken 11388573 times.
✓ Branch 1 taken 18645349 times.
✓ Branch 2 taken 11201409 times.
✓ Branch 3 taken 187164 times.
✓ Branch 4 taken 11201409 times.
✓ Branch 5 taken 187164 times.
✓ Branch 6 taken 11201409 times.
✓ Branch 7 taken 187164 times.
30033922 || (!flag && abs(s-b) >= abs(c-d)*0.5)
152
12/12
✓ Branch 0 taken 18645349 times.
✓ Branch 1 taken 11201409 times.
✓ Branch 2 taken 5854323 times.
✓ Branch 3 taken 12791026 times.
✓ Branch 4 taken 5854323 times.
✓ Branch 5 taken 12791026 times.
✓ Branch 6 taken 5854323 times.
✓ Branch 7 taken 12791026 times.
✓ Branch 8 taken 5854323 times.
✓ Branch 9 taken 12791026 times.
✓ Branch 10 taken 18478852 times.
✓ Branch 11 taken 166497 times.
35701081 || (flag && abs(b-c) < tol*max(abs(b), abs(c)))
153
14/14
✓ Branch 0 taken 47066444 times.
✓ Branch 1 taken 203325266 times.
✓ Branch 2 taken 11201409 times.
✓ Branch 3 taken 18478852 times.
✓ Branch 4 taken 2084935 times.
✓ Branch 5 taken 9116474 times.
✓ Branch 6 taken 2084935 times.
✓ Branch 7 taken 9116474 times.
✓ Branch 8 taken 2084935 times.
✓ Branch 9 taken 9116474 times.
✓ Branch 10 taken 2084935 times.
✓ Branch 11 taken 9116474 times.
✓ Branch 12 taken 75043 times.
✓ Branch 13 taken 11126366 times.
282156906 || (!flag && abs(c-d) < tol*max(abs(c), abs(d))) )
154 {
155 220786492 s = (a+b)*0.5;
156 220786492 flag = true;
157 }
158 else
159 flag = false;
160
161 // compute residual at new guess
162
0/2
✗ Branch 0 not taken.
✗ Branch 1 not taken.
250391710 fs = residual(s);
163 250391710 d = c;
164 250391710 c = b;
165 250391710 fc = fb;
166
167 // check on which side of the root s falls
168
2/2
✓ Branch 0 taken 214834303 times.
✓ Branch 1 taken 35557407 times.
250391710 if (fa*fs < 0.0)
169 {
170 b = s;
171 fb = fs;
172 }
173 else
174 {
175 214834303 a = s;
176 214834303 fa = fs;
177 }
178
179 // sort if necessary
180
6/6
✓ Branch 0 taken 34890588 times.
✓ Branch 1 taken 215501122 times.
✓ Branch 2 taken 34890588 times.
✓ Branch 3 taken 215501122 times.
✓ Branch 4 taken 34890588 times.
✓ Branch 5 taken 215501122 times.
751175130 if (abs(fa) < abs(fb))
181 {
182 34890588 swap(a, b);
183 34890588 swap(fa, fb);
184 }
185 }
186
187 DUNE_THROW(NumericalProblem, "Brent's algorithm didn't converge after " << maxIter << " iterations!");
188 }
189
190 } // end namespace Dumux
191
192 #endif
193