version 3.10-dev
math.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_MATH_HH
13#define DUMUX_MATH_HH
14
15#include <algorithm>
16#include <numeric>
17#include <cmath>
18#include <utility>
19
20#include <dune/common/typetraits.hh>
21#include <dune/common/fvector.hh>
22#include <dune/common/fmatrix.hh>
23#include <dune/common/dynmatrix.hh>
24#include <dune/common/float_cmp.hh>
25
26namespace Dumux {
27
37template <class Scalar>
38constexpr Scalar arithmeticMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
39{
40 static_assert(Dune::IsNumber<Scalar>::value, "The arguments x, y, wx, and wy have to be numbers!");
41
42 if (x*y <= 0)
43 return 0;
44 return (x * wx + y * wy)/(wx + wy);
45}
46
56template <class Scalar>
57constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
58{
59 static_assert(Dune::IsNumber<Scalar>::value, "The arguments x, y, wx, and wy have to be numbers!");
60
61 if (x*y <= 0)
62 return 0;
63 return (wx + wy) * x * y / (wy * x + wx * y);
64}
65
74template <class Scalar>
75Scalar geometricMean(Scalar x, Scalar y) noexcept
76{
77 static_assert(Dune::IsNumber<Scalar>::value, "The arguments x and y have to be numbers!");
78
79 if (x*y <= 0)
80 return 0;
81 using std::sqrt;
82 return sqrt(x*y)*sign(x);
83}
84
96template <class Scalar, int m, int n>
97void harmonicMeanMatrix(Dune::FieldMatrix<Scalar, m, n> &K,
98 const Dune::FieldMatrix<Scalar, m, n> &Ki,
99 const Dune::FieldMatrix<Scalar, m, n> &Kj)
100{
101 for (int rowIdx=0; rowIdx < m; rowIdx++){
102 for (int colIdx=0; colIdx< n; colIdx++){
103 if (Dune::FloatCmp::ne<Scalar>(Ki[rowIdx][colIdx], Kj[rowIdx][colIdx])) {
104 K[rowIdx][colIdx] =
105 harmonicMean(Ki[rowIdx][colIdx],
106 Kj[rowIdx][colIdx]);
107 }
108 else
109 K[rowIdx][colIdx] = Ki[rowIdx][colIdx];
110 }
111 }
112}
113
123template<class Scalar>
124Scalar smoothMin(const Scalar a, const Scalar b, const Scalar k)
125{
126 using std::max; using std::min; using std::abs;
127 const auto h = max(k-abs(a-b), 0.0 )/k;
128 return min(a, b) - h*h*h*k*(1.0/6.0);
129}
130
139template<class Scalar>
140Scalar smoothMax(const Scalar a, const Scalar b, const Scalar k)
141{ return -smoothMin(-a, -b, k); }
142
157template <class Scalar, class SolContainer>
158int invertLinearPolynomial(SolContainer &sol,
159 Scalar a,
160 Scalar b)
161{
162 if (a == 0.0)
163 return 0;
164
165 sol[0] = -b/a;
166 return 1;
167}
168
185template <class Scalar, class SolContainer>
186int invertQuadraticPolynomial(SolContainer &sol,
187 Scalar a,
188 Scalar b,
189 Scalar c)
190{
191 // check for a line
192 if (a == 0.0)
193 return invertLinearPolynomial(sol, b, c);
194
195 // discriminant
196 Scalar Delta = b*b - 4*a*c;
197 if (Delta < 0)
198 return 0; // no real roots
199
200 using std::sqrt;
201 Delta = sqrt(Delta);
202 sol[0] = (- b + Delta)/(2*a);
203 sol[1] = (- b - Delta)/(2*a);
204
205 // sort the result
206 if (sol[0] > sol[1])
207 {
208 using std::swap;
209 swap(sol[0], sol[1]);
210 }
211 return 2; // two real roots
212}
213
215template <class Scalar, class SolContainer>
216void invertCubicPolynomialPostProcess_(SolContainer &sol,
217 int numSol,
218 Scalar a, Scalar b, Scalar c, Scalar d,
219 std::size_t numIterations = 1)
220{
221 const auto eval = [&](auto x){ return d + x*(c + x*(b + x*a)); };
222 const auto evalDeriv = [&](auto x){ return c + x*(2*b + x*3*a); };
223
224 // do numIterations Newton iterations on the analytic solution
225 // and update result if the precision is increased
226 for (int i = 0; i < numSol; ++i)
227 {
228 Scalar x = sol[i];
229 Scalar fCurrent = eval(x);
230 const Scalar fOld = fCurrent;
231 for (int j = 0; j < numIterations; ++j)
232 {
233 const Scalar fPrime = evalDeriv(x);
234 if (fPrime == 0.0)
235 break;
236 x -= fCurrent/fPrime;
237 fCurrent = eval(x);
238 }
239
240 using std::abs;
241 if (abs(fCurrent) < abs(fOld))
242 sol[i] = x;
243 }
244}
246
268template <class Scalar, class SolContainer>
269int invertCubicPolynomial(SolContainer *sol,
270 Scalar a, Scalar b, Scalar c, Scalar d,
271 std::size_t numPostProcessIterations = 1)
272{
273 // reduces to a quadratic polynomial
274 if (a == 0)
275 return invertQuadraticPolynomial(sol, b, c, d);
276
277 // normalize the polynomial
278 b /= a;
279 c /= a;
280 d /= a;
281 a = 1;
282
283 // get rid of the quadratic term by substituting x = t - b/3
284 Scalar p = c - b*b/3;
285 Scalar q = d + (2*b*b*b - 9*b*c)/27;
286
287 // now we are at the form t^3 + p*t + q = 0. First we handle some
288 // special cases to avoid divisions by zero later...
289 if (p == 0.0 && q == 0.0) {
290 // t^3 = 0, i.e. triple root at t = 0
291 sol[0] = sol[1] = sol[2] = 0.0 - b/3;
292 return 3;
293 }
294 else if (p == 0.0 && q != 0.0) {
295 // t^3 + q = 0,
296 //
297 // i. e. single real root at t=curt(q)
298 using std::cbrt;
299 Scalar t = cbrt(q);
300 sol[0] = t - b/3;
301
302 return 1;
303 }
304 else if (p != 0.0 && q == 0.0) {
305 // t^3 + p*t = 0 = t*(t^2 + p),
306 //
307 // i. e. roots at t = 0, t^2 + p = 0
308 if (p > 0) {
309 sol[0] = 0.0 - b/3;
310 return 1; // only a single real root at t=0
311 }
312
313 // two additional real roots at t = sqrt(-p) and t = -sqrt(-p)
314 using std::sqrt;
315 sol[0] = -sqrt(-p) - b/3;
316 sol[1] = 0.0 - b/3;
317 sol[2] = sqrt(-p) - b/3;
318
319 return 3;
320 }
321
322 // At this point
323 //
324 // t^3 + p*t + q = 0
325 //
326 // with p != 0 and q != 0 holds. Introducing the variables u and v
327 // with the properties
328 //
329 // u + v = t and 3*u*v + p = 0
330 //
331 // leads to
332 //
333 // u^3 + v^3 + q = 0 .
334 //
335 // multiplying both sides with u^3 and taking advantage of the
336 // fact that u*v = -p/3 leads to
337 //
338 // u^6 + q*u^3 - p^3/27 = 0
339 //
340 // Now, substituting u^3 = w yields
341 //
342 // w^2 + q*w - p^3/27 = 0
343 //
344 // This is a quadratic equation with the solutions
345 //
346 // w = -q/2 + sqrt(q^2/4 + p^3/27) and
347 // w = -q/2 - sqrt(q^2/4 + p^3/27)
348 //
349 // Since w is equivalent to u^3 it is sufficient to only look at
350 // one of the two cases. Then, there are still 2 cases: positive
351 // and negative discriminant.
352 Scalar wDisc = q*q/4 + p*p*p/27;
353 if (wDisc >= 0) { // the positive discriminant case:
354 // calculate the cube root of - q/2 + sqrt(q^2/4 + p^3/27)
355 using std::cbrt;
356 using std::sqrt;
357
358 // Choose the root that is safe against the loss of precision
359 // that can cause wDisc to be equal to q*q/4 despite p != 0.
360 // We do not want u to be zero in that case. Mathematically,
361 // we are happy with either root.
362 const Scalar u = [&]{
363 return q > 0 ? cbrt(-0.5*q - sqrt(wDisc)) : cbrt(-0.5*q + sqrt(wDisc));
364 }();
365 // at this point, u != 0 since p^3 = 0 is necessary in order
366 // for u = 0 to hold, so
367 sol[0] = u - p/(3*u) - b/3;
368 // does not produce a division by zero. the remaining two
369 // roots of u are rotated by +- 2/3*pi in the complex plane
370 // and thus not considered here
371 invertCubicPolynomialPostProcess_(sol, 1, a, b, c, d, numPostProcessIterations);
372 return 1;
373 }
374 else { // the negative discriminant case:
375 Scalar uCubedRe = - q/2;
376 using std::sqrt;
377 Scalar uCubedIm = sqrt(-wDisc);
378 // calculate the cube root of - q/2 + sqrt(q^2/4 + p^3/27)
379 using std::cbrt;
380 Scalar uAbs = cbrt(sqrt(uCubedRe*uCubedRe + uCubedIm*uCubedIm));
381 using std::atan2;
382 Scalar phi = atan2(uCubedIm, uCubedRe)/3;
383
384 // calculate the length and the angle of the primitive root
385
386 // with the definitions from above it follows that
387 //
388 // x = u - p/(3*u) - b/3
389 //
390 // where x and u are complex numbers. Rewritten in polar form
391 // this is equivalent to
392 //
393 // x = |u|*e^(i*phi) - p*e^(-i*phi)/(3*|u|) - b/3 .
394 //
395 // Factoring out the e^ terms and subtracting the additional
396 // terms, yields
397 //
398 // x = (e^(i*phi) + e^(-i*phi))*(|u| - p/(3*|u|)) - y - b/3
399 //
400 // with
401 //
402 // y = - |u|*e^(-i*phi) + p*e^(i*phi)/(3*|u|) .
403 //
404 // The crucial observation is the fact that y is the conjugate
405 // of - x + b/3. This means that after taking advantage of the
406 // relation
407 //
408 // e^(i*phi) + e^(-i*phi) = 2*cos(phi)
409 //
410 // the equation
411 //
412 // x = 2*cos(phi)*(|u| - p / (3*|u|)) - conj(x) - 2*b/3
413 //
414 // holds. Since |u|, p, b and cos(phi) are real numbers, it
415 // follows that Im(x) = - Im(x) and thus Im(x) = 0. This
416 // implies
417 //
418 // Re(x) = x = cos(phi)*(|u| - p / (3*|u|)) - b/3 .
419 //
420 // Considering the fact that u is a cubic root, we have three
421 // values for phi which differ by 2/3*pi. This allows to
422 // calculate the three real roots of the polynomial:
423 for (int i = 0; i < 3; ++i) {
424 using std::cos;
425 sol[i] = cos(phi)*(uAbs - p/(3*uAbs)) - b/3;
426 phi += 2*M_PI/3;
427 }
428
429 // post process the obtained solution to increase numerical
430 // precision
431 invertCubicPolynomialPostProcess_(sol, 3, a, b, c, d, numPostProcessIterations);
432
433 // sort the result
434 using std::sort;
435 sort(sol, sol + 3);
436
437 return 3;
438 }
439
440 // NOT REACHABLE!
441 return 0;
442}
443
456template <class Scalar, int dim>
457bool isLarger(const Dune::FieldVector<Scalar, dim> &pos,
458 const Dune::FieldVector<Scalar, dim> &smallerVec)
459{
460 for (int i=0; i < dim; i++)
461 {
462 if (pos[i]<= smallerVec[i])
463 {
464 return false;
465 }
466 }
467 return true;
468}
469
482template <class Scalar, int dim>
483bool isSmaller(const Dune::FieldVector<Scalar, dim> &pos,
484 const Dune::FieldVector<Scalar, dim> &largerVec)
485{
486 for (int i=0; i < dim; i++)
487 {
488 if (pos[i]>= largerVec[i])
489 {
490 return false;
491 }
492 }
493 return true;
494}
495
512template <class Scalar, int dim>
513bool isBetween(const Dune::FieldVector<Scalar, dim> &pos,
514 const Dune::FieldVector<Scalar, dim> &smallerVec,
515 const Dune::FieldVector<Scalar, dim> &largerVec)
516{
517 if (isLarger(pos, smallerVec) && isSmaller(pos, largerVec))
518 {
519 return true;
520 }
521 else
522 return false;
523}
524
526namespace InterpolationPolicy { struct Linear; }
527
534template <class Policy = InterpolationPolicy::Linear, class Scalar, class ... Parameter>
535Scalar interpolate(Scalar ip, Parameter&& ... params)
536{ return Policy::interpolate(ip, std::forward<Parameter>(params) ...); }
537
542namespace InterpolationPolicy {
543
548struct Linear
549{
555 template<class Scalar>
556 static constexpr Scalar interpolate(Scalar ip, const std::array<Scalar, 2>& params)
557 {
558 return params[0]*(1.0 - ip) + params[1]*ip;
559 }
560};
561
567{
575 template<class Scalar, class RandomAccessContainer0, class RandomAccessContainer1>
576 static constexpr Scalar interpolate(Scalar ip, const RandomAccessContainer0& range, const RandomAccessContainer1& values)
577 {
578 // check bounds
579 if (ip > range[range.size()-1]) return values[values.size()-1];
580 if (ip < range[0]) return values[0];
581
582 // if we are within bounds find the index of the lower bound
583 const auto lookUpIndex = std::distance(range.begin(), std::lower_bound(range.begin(), range.end(), ip));
584 if (lookUpIndex == 0)
585 return values[0];
586
587 const auto ipLinear = (ip - range[lookUpIndex-1])/(range[lookUpIndex] - range[lookUpIndex-1]);
588 return Dumux::interpolate<Linear>(ipLinear, std::array<Scalar, 2>{{values[lookUpIndex-1], values[lookUpIndex]}});
589 }
590
591 template<class Scalar, class RandomAccessContainer>
592 static constexpr Scalar interpolate(Scalar ip, const std::pair<RandomAccessContainer, RandomAccessContainer>& table)
593 {
594 const auto& [range, values] = table;
595 return interpolate(ip, range, values);
596 }
597};
598
599} // end namespace InterpolationPolicy
600
610template <class Scalar>
611std::vector<Scalar> linspace(const Scalar begin, const Scalar end,
612 std::size_t samples,
613 bool endPoint = true)
614{
615 using std::max;
616 samples = max(std::size_t{2}, samples); // only makes sense for 2 or more samples
617 const Scalar divisor = endPoint ? samples-1 : samples;
618 const Scalar delta = (end-begin)/divisor;
619 std::vector<Scalar> vec(samples);
620 for (std::size_t i = 0; i < samples; ++i)
621 vec[i] = begin + i*delta;
622 return vec;
623}
624
625
638template <class Scalar>
639Scalar antoine(Scalar temperature,
640 Scalar A,
641 Scalar B,
642 Scalar C)
643{
644 const Scalar ln10 = 2.3025850929940459;
645 using std::exp;
646 return exp(ln10*(A - B/(C + temperature)));
647}
648
657template<class ValueType>
658constexpr int sign(const ValueType& value) noexcept
659{
660 return (ValueType(0) < value) - (value < ValueType(0));
661}
662
670template <class Scalar>
671Dune::FieldVector<Scalar, 3> crossProduct(const Dune::FieldVector<Scalar, 3> &vec1,
672 const Dune::FieldVector<Scalar, 3> &vec2)
673{
674 return {vec1[1]*vec2[2]-vec1[2]*vec2[1],
675 vec1[2]*vec2[0]-vec1[0]*vec2[2],
676 vec1[0]*vec2[1]-vec1[1]*vec2[0]};
677}
678
686template <class Scalar>
687Scalar crossProduct(const Dune::FieldVector<Scalar, 2> &vec1,
688 const Dune::FieldVector<Scalar, 2> &vec2)
689{ return vec1[0]*vec2[1]-vec1[1]*vec2[0]; }
690
699template <class Scalar>
700Scalar tripleProduct(const Dune::FieldVector<Scalar, 3> &vec1,
701 const Dune::FieldVector<Scalar, 3> &vec2,
702 const Dune::FieldVector<Scalar, 3> &vec3)
703{ return crossProduct<Scalar>(vec1, vec2)*vec3; }
704
711template <class Scalar, int m, int n>
712Dune::FieldMatrix<Scalar, n, m> getTransposed(const Dune::FieldMatrix<Scalar, m, n>& M)
713{
714 Dune::FieldMatrix<Scalar, n, m> T;
715 for (std::size_t i = 0; i < m; ++i)
716 for (std::size_t j = 0; j < n; ++j)
717 T[j][i] = M[i][j];
718
719 return T;
720}
721
728template <class Scalar>
729Dune::DynamicMatrix<Scalar> getTransposed(const Dune::DynamicMatrix<Scalar>& M)
730{
731 std::size_t rows_T = M.M();
732 std::size_t cols_T = M.N();
733
734 Dune::DynamicMatrix<Scalar> M_T(rows_T, cols_T, 0.0);
735
736 for (std::size_t i = 0; i < rows_T; ++i)
737 for (std::size_t j = 0; j < cols_T; ++j)
738 M_T[i][j] = M[j][i];
739
740 return M_T;
741}
742
750template <class Scalar>
751Dune::DynamicMatrix<Scalar> multiplyMatrices(const Dune::DynamicMatrix<Scalar> &M1,
752 const Dune::DynamicMatrix<Scalar> &M2)
753{
754 using size_type = typename Dune::DynamicMatrix<Scalar>::size_type;
755 const size_type rows = M1.N();
756 const size_type cols = M2.M();
757
758 DUNE_ASSERT_BOUNDS(M1.M() == M2.N());
759
760 Dune::DynamicMatrix<Scalar> result(rows, cols, 0.0);
761 for (size_type i = 0; i < rows; i++)
762 for (size_type j = 0; j < cols; j++)
763 for (size_type k = 0; k < M1.M(); k++)
764 result[i][j] += M1[i][k]*M2[k][j];
765
766 return result;
767}
768
776template <class Scalar, int rows1, int cols1, int cols2>
777Dune::FieldMatrix<Scalar, rows1, cols2> multiplyMatrices(const Dune::FieldMatrix<Scalar, rows1, cols1> &M1,
778 const Dune::FieldMatrix<Scalar, cols1, cols2> &M2)
779{
780 using size_type = typename Dune::FieldMatrix<Scalar, rows1, cols2>::size_type;
781
782 Dune::FieldMatrix<Scalar, rows1, cols2> result(0.0);
783 for (size_type i = 0; i < rows1; i++)
784 for (size_type j = 0; j < cols2; j++)
785 for (size_type k = 0; k < cols1; k++)
786 result[i][j] += M1[i][k]*M2[k][j];
787
788 return result;
789}
790
791
798template <class MatrixType>
799typename Dune::DenseMatrix<MatrixType>::field_type
800trace(const Dune::DenseMatrix<MatrixType>& M)
801{
802 const auto rows = M.N();
803 DUNE_ASSERT_BOUNDS(rows == M.M()); // rows == cols
804
805 using MatType = Dune::DenseMatrix<MatrixType>;
806 typename MatType::field_type trace = 0.0;
807
808 for (typename MatType::size_type i = 0; i < rows; ++i)
809 trace += M[i][i];
810
811 return trace;
812}
813
827template <class MAT, class V>
828typename Dune::DenseVector<V>::derived_type
829mv(const Dune::DenseMatrix<MAT>& M,
830 const Dune::DenseVector<V>& v)
831{
832 typename Dune::DenseVector<V>::derived_type res(v);
833 M.mv(v, res);
834 return res;
835}
836
854template <class FieldScalar, class V>
855typename std::enable_if_t<Dune::IsNumber<FieldScalar>::value,
856 typename Dune::DenseVector<V>::derived_type>
857mv(const FieldScalar m, const Dune::DenseVector<V>& v)
858{
859 typename Dune::DenseVector<V>::derived_type res(v);
860 res *= m;
861 return res;
862}
863
878template <class V1, class MAT, class V2>
879typename Dune::DenseMatrix<MAT>::value_type
880vtmv(const Dune::DenseVector<V1>& v1,
881 const Dune::DenseMatrix<MAT>& M,
882 const Dune::DenseVector<V2>& v2)
883{
884 return v1*mv(M, v2);
885}
886
906template <class V1, class FieldScalar, class V2>
907typename std::enable_if_t<Dune::IsNumber<FieldScalar>::value, FieldScalar>
908vtmv(const Dune::DenseVector<V1>& v1,
909 const FieldScalar m,
910 const Dune::DenseVector<V2>& v2)
911{
912 return m*(v1*v2);
913}
914
926template <class Scalar>
927std::array<Scalar, 2> linearRegression(const std::vector<Scalar>& x,
928 const std::vector<Scalar>& y)
929{
930 if (x.size() != y.size())
931 DUNE_THROW(Dune::InvalidStateException, "x and y array must have the same length.");
932
933 const Scalar averageX = std::accumulate(x.begin(), x.end(), 0.0)/x.size();
934 const Scalar averageY = std::accumulate(y.begin(), y.end(), 0.0)/y.size();
935
936 // calculate temporary variables necessary for slope computation
937 const Scalar numerator = std::inner_product(
938 x.begin(), x.end(), y.begin(), 0.0, std::plus<Scalar>(),
939 [&](auto xx, auto yy) { return (xx - averageX) * (yy - averageY); }
940 );
941 const Scalar denominator = std::inner_product(
942 x.begin(), x.end(), x.begin(), 0.0, std::plus<Scalar>(),
943 [&](auto xx, auto yy) { return (xx - averageX) * (yy - averageX); }
944 );
945
946 // compute slope and intercept of the regression line
947 const Scalar slope = numerator / denominator;
948 const Scalar intercept = averageY - slope * averageX;
949
950 return {intercept, slope};
951}
952
953} // end namespace Dumux
954
955#endif
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:57
int invertCubicPolynomial(SolContainer *sol, Scalar a, Scalar b, Scalar c, Scalar d, std::size_t numPostProcessIterations=1)
Invert a cubic polynomial analytically.
Definition: math.hh:269
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:880
Dune::DynamicMatrix< Scalar > multiplyMatrices(const Dune::DynamicMatrix< Scalar > &M1, const Dune::DynamicMatrix< Scalar > &M2)
Multiply two dynamic matrices.
Definition: math.hh:751
Scalar geometricMean(Scalar x, Scalar y) noexcept
Calculate the geometric mean of two scalar values.
Definition: math.hh:75
int invertLinearPolynomial(SolContainer &sol, Scalar a, Scalar b)
Invert a linear polynomial analytically.
Definition: math.hh:158
bool isBetween(const Dune::FieldVector< Scalar, dim > &pos, const Dune::FieldVector< Scalar, dim > &smallerVec, const Dune::FieldVector< Scalar, dim > &largerVec)
Comparison of three position vectors.
Definition: math.hh:513
void harmonicMeanMatrix(Dune::FieldMatrix< Scalar, m, n > &K, const Dune::FieldMatrix< Scalar, m, n > &Ki, const Dune::FieldMatrix< Scalar, m, n > &Kj)
Calculate the harmonic mean of a fixed-size matrix.
Definition: math.hh:97
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:671
Scalar smoothMax(const Scalar a, const Scalar b, const Scalar k)
A smoothed maximum function (using cubic interpolation)
Definition: math.hh:140
Scalar smoothMin(const Scalar a, const Scalar b, const Scalar k)
A smoothed minimum function (using cubic interpolation)
Definition: math.hh:124
Scalar interpolate(Scalar ip, Parameter &&... params)
a generic function to interpolate given a set of parameters and an interpolation point
Definition: math.hh:535
bool isSmaller(const Dune::FieldVector< Scalar, dim > &pos, const Dune::FieldVector< Scalar, dim > &largerVec)
Comparison of two position vectors.
Definition: math.hh:483
Dune::DenseVector< V >::derived_type mv(const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V > &v)
Returns the result of the projection of a vector v with a Matrix M.
Definition: math.hh:829
Dune::DenseMatrix< MatrixType >::field_type trace(const Dune::DenseMatrix< MatrixType > &M)
Trace of a dense matrix.
Definition: math.hh:800
Dune::FieldMatrix< Scalar, n, m > getTransposed(const Dune::FieldMatrix< Scalar, m, n > &M)
Transpose a FieldMatrix.
Definition: math.hh:712
bool isLarger(const Dune::FieldVector< Scalar, dim > &pos, const Dune::FieldVector< Scalar, dim > &smallerVec)
Comparison of two position vectors.
Definition: math.hh:457
int invertQuadraticPolynomial(SolContainer &sol, Scalar a, Scalar b, Scalar c)
Invert a quadratic polynomial analytically.
Definition: math.hh:186
std::array< Scalar, 2 > linearRegression(const std::vector< Scalar > &x, const std::vector< Scalar > &y)
Returns the [intercept, slope] of the regression line fitted to a set of (x, y) data points.
Definition: math.hh:927
std::vector< Scalar > linspace(const Scalar begin, const Scalar end, std::size_t samples, bool endPoint=true)
Generates linearly spaced vectors.
Definition: math.hh:611
Scalar tripleProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2, const Dune::FieldVector< Scalar, 3 > &vec3)
Triple product of three vectors in three-dimensional Euclidean space retuning scalar.
Definition: math.hh:700
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:658
Scalar antoine(Scalar temperature, Scalar A, Scalar B, Scalar C)
Evaluates the Antoine equation used to calculate the vapour pressure of various liquids.
Definition: math.hh:639
constexpr Scalar arithmeticMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) arithmetic mean of two scalar values.
Definition: math.hh:38
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:282
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:39
Definition: adapt.hh:17
interpolate linearly between two given values
Definition: math.hh:549
static constexpr Scalar interpolate(Scalar ip, const std::array< Scalar, 2 > &params)
interpolate linearly between two given values
Definition: math.hh:556
interpolate linearly in a piecewise linear function (tabularized function)
Definition: math.hh:567
static constexpr Scalar interpolate(Scalar ip, const std::pair< RandomAccessContainer, RandomAccessContainer > &table)
Definition: math.hh:592
static constexpr Scalar interpolate(Scalar ip, const RandomAccessContainer0 &range, const RandomAccessContainer1 &values)
interpolate linearly in a piecewise linear function (tabularized function)
Definition: math.hh:576