24#ifndef DUMUX_COMMON_SCALAR_ROOT_FINDING_HH
25#define DUMUX_COMMON_SCALAR_ROOT_FINDING_HH
46template<
class Scalar,
class ResFunc,
class DerivFunc,
47 typename std::enable_if_t<std::is_invocable_r_v<Scalar, ResFunc, Scalar> &&
48 std::is_invocable_r_v<Scalar, DerivFunc, Scalar>>...>
50 const Scalar tol = 1e-13,
const int maxIter = 200)
53 Scalar r = residual(xNew);
56 Scalar relativeShift = std::numeric_limits<Scalar>::max();
57 while (relativeShift > tol && n > 0)
59 xNew = xOld - r/derivative(xOld);
62 using std::abs;
using std::max;
63 relativeShift = abs(xOld-xNew)/max(abs(xOld), abs(xNew));
70 DUNE_THROW(
NumericalProblem,
"Residual is not finite: " << r <<
" after " << maxIter - n <<
" iterations!");
72 if (relativeShift > tol)
73 DUNE_THROW(
NumericalProblem,
"Scalar newton solver did not converge after " << maxIter <<
" iterations!");
87template<
class Scalar,
class ResFunc,
88 typename std::enable_if_t<std::is_invocable_r_v<Scalar, ResFunc, Scalar>>...>
90 const Scalar tol = 1e-13,
const int maxIter = 200)
93 auto derivative = [&](
const auto x){
return (residual(x + eps)-residual(x))/eps; };
109template<
class Scalar,
class ResFunc,
110 typename std::enable_if_t<std::is_invocable_r_v<Scalar, ResFunc, Scalar>>...>
112 const Scalar tol = 1e-13,
const int maxIter = 200)
115 Scalar fa = residual(a);
116 Scalar fb = residual(b);
122 DUNE_THROW(
NumericalProblem,
"Brent's algorithm failed: [a,b] does not contain any, or no uniquely defined root!");
125 using std::abs;
using std::swap;
126 if (abs(fa) < abs(fb))
138 for (
int i = 0; i < maxIter; ++i)
142 if (abs(b-a) < tol*max(abs(a), abs(b)))
146 if (fa != fc && fb != fc)
148 const auto fab = fa-fb;
149 const auto fac = fa-fc;
150 const auto fbc = fb-fc;
151 s = a*fb*fc/(fab*fac) - b*fa*fc/(fab*fbc) + c*fa*fb/(fac*fbc);
157 s = b - fb*(b-a)/(fb-fa);
161 if ( (s < (3*a + b)*0.25 || s > b)
162 || (flag && abs(s-b) >= abs(b-c)*0.5)
163 || (!flag && abs(s-b) >= abs(c-d)*0.5)
164 || (flag && abs(b-c) < tol*max(abs(b), abs(c)))
165 || (!flag && abs(c-d) < tol*max(abs(c), abs(d))) )
192 if (abs(fa) < abs(fb))
199 DUNE_THROW(
NumericalProblem,
"Brent's algorithm didn't converge after " << maxIter <<
" iterations!");
Some exceptions thrown in DuMux
A class for numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Scalar findScalarRootNewton(Scalar xOld, const ResFunc &residual, const DerivFunc &derivative, const Scalar tol=1e-13, const int maxIter=200)
Newton's root finding algorithm for scalar functions (secant method)
Definition: findscalarroot.hh:49
Scalar findScalarRootBrent(Scalar a, Scalar b, const ResFunc &residual, const Scalar tol=1e-13, const int maxIter=200)
Brent's root finding algorithm for scalar functions.
Definition: findscalarroot.hh:111
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
static Scalar epsilon(const Scalar value, const Scalar baseEps=1e-10)
Computes the epsilon used for numeric differentiation.
Definition: numericdifferentiation.hh:47