32#include <dune/common/typetraits.hh>
33#include <dune/common/fvector.hh>
34#include <dune/common/fmatrix.hh>
35#include <dune/common/dynmatrix.hh>
36#include <dune/common/float_cmp.hh>
49template <
class Scalar>
50constexpr Scalar
arithmeticMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
52 static_assert(Dune::IsNumber<Scalar>::value,
"The arguments x, y, wx, and wy have to be numbers!");
56 return (x * wx + y * wy)/(wx + wy);
68template <
class Scalar>
69constexpr Scalar
harmonicMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
71 static_assert(Dune::IsNumber<Scalar>::value,
"The arguments x, y, wx, and wy have to be numbers!");
75 return (wx + wy) * x * y / (wy * x + wx * y);
86template <
class Scalar>
89 static_assert(Dune::IsNumber<Scalar>::value,
"The arguments x and y have to be numbers!");
94 return sqrt(x*y)*
sign(x);
108template <
class Scalar,
int m,
int n>
110 const Dune::FieldMatrix<Scalar, m, n> &Ki,
111 const Dune::FieldMatrix<Scalar, m, n> &Kj)
113 for (
int rowIdx=0; rowIdx < m; rowIdx++){
114 for (
int colIdx=0; colIdx< n; colIdx++){
115 if (Dune::FloatCmp::ne<Scalar>(Ki[rowIdx][colIdx], Kj[rowIdx][colIdx])) {
121 K[rowIdx][colIdx] = Ki[rowIdx][colIdx];
140template <
class Scalar,
class SolContainer>
168template <
class Scalar,
class SolContainer>
179 Scalar Delta = b*b - 4*a*c;
185 sol[0] = (- b + Delta)/(2*a);
186 sol[1] = (- b - Delta)/(2*a);
192 swap(sol[0], sol[1]);
198template <
class Scalar,
class SolContainer>
199void invertCubicPolynomialPostProcess_(SolContainer &sol,
201 Scalar a, Scalar b, Scalar c, Scalar d,
202 std::size_t numIterations = 1)
204 const auto eval = [&](
auto x){
return d + x*(c + x*(b + x*a)); };
205 const auto evalDeriv = [&](
auto x){
return c + x*(2*b + x*3*a); };
209 for (
int i = 0; i < numSol; ++i)
212 Scalar fCurrent = eval(x);
213 const Scalar fOld = fCurrent;
214 for (
int j = 0; j < numIterations; ++j)
216 const Scalar fPrime = evalDeriv(x);
219 x -= fCurrent/fPrime;
224 if (abs(fCurrent) < abs(fOld))
251template <
class Scalar,
class SolContainer>
253 Scalar a, Scalar b, Scalar c, Scalar d,
254 std::size_t numPostProcessIterations = 1)
267 Scalar p = c - b*b/3;
268 Scalar q = d + (2*b*b*b - 9*b*c)/27;
272 if (p == 0.0 && q == 0.0) {
274 sol[0] = sol[1] = sol[2] = 0.0 - b/3;
277 else if (p == 0.0 && q != 0.0) {
287 else if (p != 0.0 && q == 0.0) {
298 sol[0] = -sqrt(-p) - b/3;
300 sol[2] = sqrt(-p) - b/3;
335 Scalar wDisc = q*q/4 + p*p*p/27;
345 const Scalar u = [&]{
346 return q > 0 ? cbrt(-0.5*q - sqrt(wDisc)) : cbrt(-0.5*q + sqrt(wDisc));
350 sol[0] = u - p/(3*u) - b/3;
354 invertCubicPolynomialPostProcess_(sol, 1, a, b, c, d, numPostProcessIterations);
358 Scalar uCubedRe = - q/2;
360 Scalar uCubedIm = sqrt(-wDisc);
363 Scalar uAbs = cbrt(sqrt(uCubedRe*uCubedRe + uCubedIm*uCubedIm));
365 Scalar phi = atan2(uCubedIm, uCubedRe)/3;
406 for (
int i = 0; i < 3; ++i) {
408 sol[i] = cos(phi)*(uAbs - p/(3*uAbs)) - b/3;
414 invertCubicPolynomialPostProcess_(sol, 3, a, b, c, d, numPostProcessIterations);
439template <
class Scalar,
int dim>
440bool isLarger(
const Dune::FieldVector<Scalar, dim> &pos,
441 const Dune::FieldVector<Scalar, dim> &smallerVec)
443 for (
int i=0; i < dim; i++)
445 if (pos[i]<= smallerVec[i])
465template <
class Scalar,
int dim>
466bool isSmaller(
const Dune::FieldVector<Scalar, dim> &pos,
467 const Dune::FieldVector<Scalar, dim> &largerVec)
469 for (
int i=0; i < dim; i++)
471 if (pos[i]>= largerVec[i])
495template <
class Scalar,
int dim>
496bool isBetween(
const Dune::FieldVector<Scalar, dim> &pos,
497 const Dune::FieldVector<Scalar, dim> &smallerVec,
498 const Dune::FieldVector<Scalar, dim> &largerVec)
509namespace InterpolationPolicy {
struct Linear; }
525namespace InterpolationPolicy {
538 template<
class Scalar>
539 static constexpr Scalar
interpolate(Scalar ip,
const std::array<Scalar, 2>& params)
541 return params[0]*(1.0 - ip) + params[1]*ip;
558 template<
class Scalar,
class RandomAccessContainer0,
class RandomAccessContainer1>
559 static constexpr Scalar
interpolate(Scalar ip,
const RandomAccessContainer0& range,
const RandomAccessContainer1& values)
562 if (ip > range.back())
return values.back();
563 if (ip < range[0])
return values[0];
566 const auto lookUpIndex =
std::distance(range.begin(), std::lower_bound(range.begin(), range.end(), ip));
567 if (lookUpIndex == 0)
570 const auto ipLinear = (ip - range[lookUpIndex-1])/(range[lookUpIndex] - range[lookUpIndex-1]);
571 return Dumux::interpolate<Linear>(ipLinear, std::array<Scalar, 2>{{values[lookUpIndex-1], values[lookUpIndex]}});
574 template<
class Scalar,
class RandomAccessContainer>
575 static constexpr Scalar
interpolate(Scalar ip,
const std::pair<RandomAccessContainer, RandomAccessContainer>& table)
577 const auto& [range, values] = table;
593template <
class Scalar>
594std::vector<Scalar>
linspace(
const Scalar begin,
const Scalar end,
596 bool endPoint =
true)
599 samples = max(std::size_t{2}, samples);
600 const Scalar divisor = endPoint ? samples-1 : samples;
601 const Scalar delta = (end-begin)/divisor;
602 std::vector<Scalar> vec(samples);
603 for (std::size_t i = 0; i < samples; ++i)
604 vec[i] = begin + i*delta;
621template <
class Scalar>
627 const Scalar ln10 = 2.3025850929940459;
640template<
class ValueType>
641constexpr int sign(
const ValueType& value)
noexcept
643 return (ValueType(0) < value) - (value < ValueType(0));
653template <
class Scalar>
654Dune::FieldVector<Scalar, 3>
crossProduct(
const Dune::FieldVector<Scalar, 3> &vec1,
655 const Dune::FieldVector<Scalar, 3> &vec2)
657 return {vec1[1]*vec2[2]-vec1[2]*vec2[1],
658 vec1[2]*vec2[0]-vec1[0]*vec2[2],
659 vec1[0]*vec2[1]-vec1[1]*vec2[0]};
669template <
class Scalar>
671 const Dune::FieldVector<Scalar, 2> &vec2)
672{
return vec1[0]*vec2[1]-vec1[1]*vec2[0]; }
682template <
class Scalar>
684 const Dune::FieldVector<Scalar, 3> &vec2,
685 const Dune::FieldVector<Scalar, 3> &vec3)
686{
return crossProduct<Scalar>(vec1, vec2)*vec3; }
694template <
class Scalar,
int m,
int n>
695Dune::FieldMatrix<Scalar, n, m>
getTransposed(
const Dune::FieldMatrix<Scalar, m, n>& M)
697 Dune::FieldMatrix<Scalar, n, m> T;
698 for (std::size_t i = 0; i < m; ++i)
699 for (std::size_t j = 0; j < n; ++j)
711template <
class Scalar>
712Dune::DynamicMatrix<Scalar>
getTransposed(
const Dune::DynamicMatrix<Scalar>& M)
714 std::size_t rows_T = M.M();
715 std::size_t cols_T = M.N();
717 Dune::DynamicMatrix<Scalar> M_T(rows_T, cols_T, 0.0);
719 for (std::size_t i = 0; i < rows_T; ++i)
720 for (std::size_t j = 0; j < cols_T; ++j)
733template <
class Scalar>
735 const Dune::DynamicMatrix<Scalar> &M2)
737 using size_type =
typename Dune::DynamicMatrix<Scalar>::size_type;
738 const size_type rows = M1.N();
739 const size_type cols = M2.M();
741 DUNE_ASSERT_BOUNDS(M1.M() == M2.N());
743 Dune::DynamicMatrix<Scalar> result(rows, cols, 0.0);
744 for (size_type i = 0; i < rows; i++)
745 for (size_type j = 0; j < cols; j++)
746 for (size_type k = 0; k < M1.M(); k++)
747 result[i][j] += M1[i][k]*M2[k][j];
759template <
class Scalar,
int rows1,
int cols1,
int cols2>
760Dune::FieldMatrix<Scalar, rows1, cols2>
multiplyMatrices(
const Dune::FieldMatrix<Scalar, rows1, cols1> &M1,
761 const Dune::FieldMatrix<Scalar, cols1, cols2> &M2)
763 using size_type =
typename Dune::FieldMatrix<Scalar, rows1, cols2>::size_type;
765 Dune::FieldMatrix<Scalar, rows1, cols2> result(0.0);
766 for (size_type i = 0; i < rows1; i++)
767 for (size_type j = 0; j < cols2; j++)
768 for (size_type k = 0; k < cols1; k++)
769 result[i][j] += M1[i][k]*M2[k][j];
781template <
class MatrixType>
782typename Dune::DenseMatrix<MatrixType>::field_type
783trace(
const Dune::DenseMatrix<MatrixType>& M)
785 const auto rows = M.N();
786 DUNE_ASSERT_BOUNDS(rows == M.M());
788 using MatType = Dune::DenseMatrix<MatrixType>;
789 typename MatType::field_type
trace = 0.0;
791 for (
typename MatType::size_type i = 0; i < rows; ++i)
810template <
class MAT,
class V>
811typename Dune::DenseVector<V>::derived_type
812mv(
const Dune::DenseMatrix<MAT>& M,
813 const Dune::DenseVector<V>& v)
815 typename Dune::DenseVector<V>::derived_type res(v);
837template <
class FieldScalar,
class V>
838typename std::enable_if_t<Dune::IsNumber<FieldScalar>::value,
839 typename Dune::DenseVector<V>::derived_type>
840mv(
const FieldScalar m,
const Dune::DenseVector<V>& v)
842 typename Dune::DenseVector<V>::derived_type res(v);
861template <
class V1,
class MAT,
class V2>
862typename Dune::DenseMatrix<MAT>::value_type
863vtmv(
const Dune::DenseVector<V1>& v1,
864 const Dune::DenseMatrix<MAT>& M,
865 const Dune::DenseVector<V2>& v2)
889template <
class V1,
class FieldScalar,
class V2>
890typename std::enable_if_t<Dune::IsNumber<FieldScalar>::value, FieldScalar>
891vtmv(
const Dune::DenseVector<V1>& v1,
893 const Dune::DenseVector<V2>& v2)
909template <
class Scalar>
911 const std::vector<Scalar>& y)
913 if (x.size() != y.size())
914 DUNE_THROW(Dune::InvalidStateException,
"x and y array must have the same length.");
916 const Scalar averageX = std::accumulate(x.begin(), x.end(), 0.0)/x.size();
917 const Scalar averageY = std::accumulate(y.begin(), y.end(), 0.0)/y.size();
920 const Scalar numerator = std::inner_product(
921 x.begin(), x.end(), y.begin(), 0.0, std::plus<Scalar>(),
922 [&](
auto xx,
auto yy) { return (xx - averageX) * (yy - averageY); }
924 const Scalar denominator = std::inner_product(
925 x.begin(), x.end(), x.begin(), 0.0, std::plus<Scalar>(),
926 [&](
auto xx,
auto yy) { return (xx - averageX) * (yy - averageX); }
930 const Scalar slope = numerator / denominator;
931 const Scalar intercept = averageY - slope * averageX;
933 return {intercept, slope};
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:292
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:69
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:252
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:863
Dune::DynamicMatrix< Scalar > multiplyMatrices(const Dune::DynamicMatrix< Scalar > &M1, const Dune::DynamicMatrix< Scalar > &M2)
Multiply two dynamic matrices.
Definition: math.hh:734
Scalar geometricMean(Scalar x, Scalar y) noexcept
Calculate the geometric mean of two scalar values.
Definition: math.hh:87
int invertLinearPolynomial(SolContainer &sol, Scalar a, Scalar b)
Invert a linear polynomial analytically.
Definition: math.hh:141
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:496
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:109
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:654
Scalar interpolate(Scalar ip, Parameter &&... params)
a generic function to interpolate given a set of parameters and an interpolation point
Definition: math.hh:518
bool isSmaller(const Dune::FieldVector< Scalar, dim > &pos, const Dune::FieldVector< Scalar, dim > &largerVec)
Comparison of two position vectors.
Definition: math.hh:466
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:812
Dune::DenseMatrix< MatrixType >::field_type trace(const Dune::DenseMatrix< MatrixType > &M)
Trace of a dense matrix.
Definition: math.hh:783
Dune::FieldMatrix< Scalar, n, m > getTransposed(const Dune::FieldMatrix< Scalar, m, n > &M)
Transpose a FieldMatrix.
Definition: math.hh:695
bool isLarger(const Dune::FieldVector< Scalar, dim > &pos, const Dune::FieldVector< Scalar, dim > &smallerVec)
Comparison of two position vectors.
Definition: math.hh:440
int invertQuadraticPolynomial(SolContainer &sol, Scalar a, Scalar b, Scalar c)
Invert a quadratic polynomial analytically.
Definition: math.hh:169
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:910
std::vector< Scalar > linspace(const Scalar begin, const Scalar end, std::size_t samples, bool endPoint=true)
Generates linearly spaced vectors.
Definition: math.hh:594
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:683
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:641
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:622
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:50
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
interpolate linearly between two given values
Definition: math.hh:532
static constexpr Scalar interpolate(Scalar ip, const std::array< Scalar, 2 > ¶ms)
interpolate linearly between two given values
Definition: math.hh:539
interpolate linearly in a piecewise linear function (tabularized function)
Definition: math.hh:550
static constexpr Scalar interpolate(Scalar ip, const std::pair< RandomAccessContainer, RandomAccessContainer > &table)
Definition: math.hh:575
static constexpr Scalar interpolate(Scalar ip, const RandomAccessContainer0 &range, const RandomAccessContainer1 &values)
interpolate linearly in a piecewise linear function (tabularized function)
Definition: math.hh:559