31#include <dune/common/typetraits.hh>
32#include <dune/common/fvector.hh>
33#include <dune/common/fmatrix.hh>
34#include <dune/common/dynmatrix.hh>
35#include <dune/common/float_cmp.hh>
48template <
class Scalar>
49constexpr Scalar
arithmeticMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
51 static_assert(Dune::IsNumber<Scalar>::value,
"The arguments x, y, wx, and wy have to be numbers!");
55 return (x * wx + y * wy)/(wx + wy);
67template <
class Scalar>
68constexpr Scalar
harmonicMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
70 static_assert(Dune::IsNumber<Scalar>::value,
"The arguments x, y, wx, and wy have to be numbers!");
74 return (wx + wy) * x * y / (wy * x + wx * y);
85template <
class Scalar>
88 static_assert(Dune::IsNumber<Scalar>::value,
"The arguments x and y have to be numbers!");
93 return sqrt(x*y)*
sign(x);
107template <
class Scalar,
int m,
int n>
109 const Dune::FieldMatrix<Scalar, m, n> &Ki,
110 const Dune::FieldMatrix<Scalar, m, n> &Kj)
112 for (
int rowIdx=0; rowIdx < m; rowIdx++){
113 for (
int colIdx=0; colIdx< n; colIdx++){
114 if (Dune::FloatCmp::ne<Scalar>(Ki[rowIdx][colIdx], Kj[rowIdx][colIdx])) {
120 K[rowIdx][colIdx] = Ki[rowIdx][colIdx];
139template <
class Scalar,
class SolContainer>
167template <
class Scalar,
class SolContainer>
178 Scalar Delta = b*b - 4*a*c;
184 sol[0] = (- b + Delta)/(2*a);
185 sol[1] = (- b - Delta)/(2*a);
191 swap(sol[0], sol[1]);
197template <
class Scalar,
class SolContainer>
198void invertCubicPolynomialPostProcess_(SolContainer &sol,
207 for (
int i = 0; i < numSol; ++i) {
209 Scalar fOld = d + x*(c + x*(b + x*a));
211 Scalar fPrime = c + x*(2*b + x*3*a);
216 Scalar fNew = d + x*(c + x*(b + x*a));
218 if (abs(fNew) < abs(fOld))
241template <
class Scalar,
class SolContainer>
259 Scalar p = c - b*b/3;
260 Scalar q = d + (2*b*b*b - 9*b*c)/27;
264 if (p == 0.0 && q == 0.0) {
266 sol[0] = sol[1] = sol[2] = 0.0 - b/3;
269 else if (p == 0.0 && q != 0.0) {
279 else if (p != 0.0 && q == 0.0) {
290 sol[0] = -sqrt(-p) - b/3;
292 sol[2] = sqrt(-p) - b/3;
327 Scalar wDisc = q*q/4 + p*p*p/27;
332 Scalar u = cbrt(-q/2 + sqrt(wDisc));
336 sol[0] = u - p/(3*u) - b/3;
340 invertCubicPolynomialPostProcess_(sol, 1, a, b, c, d);
344 Scalar uCubedRe = - q/2;
346 Scalar uCubedIm = sqrt(-wDisc);
349 Scalar uAbs = cbrt(sqrt(uCubedRe*uCubedRe + uCubedIm*uCubedIm));
351 Scalar phi = atan2(uCubedIm, uCubedRe)/3;
392 for (
int i = 0; i < 3; ++i) {
394 sol[i] = cos(phi)*(uAbs - p/(3*uAbs)) - b/3;
400 invertCubicPolynomialPostProcess_(sol, 3, a, b, c, d);
425template <
class Scalar,
int dim>
426bool isLarger(
const Dune::FieldVector<Scalar, dim> &pos,
427 const Dune::FieldVector<Scalar, dim> &smallerVec)
429 for (
int i=0; i < dim; i++)
431 if (pos[i]<= smallerVec[i])
451template <
class Scalar,
int dim>
452bool isSmaller(
const Dune::FieldVector<Scalar, dim> &pos,
453 const Dune::FieldVector<Scalar, dim> &largerVec)
455 for (
int i=0; i < dim; i++)
457 if (pos[i]>= largerVec[i])
481template <
class Scalar,
int dim>
482bool isBetween(
const Dune::FieldVector<Scalar, dim> &pos,
483 const Dune::FieldVector<Scalar, dim> &smallerVec,
484 const Dune::FieldVector<Scalar, dim> &largerVec)
495namespace InterpolationPolicy {
struct Linear; }
511namespace InterpolationPolicy {
524 template<
class Scalar>
525 static constexpr Scalar
interpolate(Scalar ip,
const std::array<Scalar, 2>& params)
527 return params[0]*(1.0 - ip) + params[1]*ip;
544 template<
class Scalar,
class RandomAccessContainer0,
class RandomAccessContainer1>
545 static constexpr Scalar
interpolate(Scalar ip,
const RandomAccessContainer0& range,
const RandomAccessContainer1& values)
548 if (ip > range.back())
return values.back();
549 if (ip < range[0])
return values[0];
552 const auto lookUpIndex =
std::distance(range.begin(), std::lower_bound(range.begin(), range.end(), ip));
553 if (lookUpIndex == 0)
556 const auto ipLinear = (ip - range[lookUpIndex-1])/(range[lookUpIndex] - range[lookUpIndex-1]);
557 return Dumux::interpolate<Linear>(ipLinear, std::array<Scalar, 2>{{values[lookUpIndex-1], values[lookUpIndex]}});
560 template<
class Scalar,
class RandomAccessContainer>
561 static constexpr Scalar
interpolate(Scalar ip,
const std::pair<RandomAccessContainer, RandomAccessContainer>& table)
563 const auto& [range, values] = table;
579template <
class Scalar>
580std::vector<Scalar>
linspace(
const Scalar begin,
const Scalar end,
582 bool endPoint =
true)
585 samples = max(std::size_t{2}, samples);
586 const Scalar divisor = endPoint ? samples-1 : samples;
587 const Scalar delta = (end-begin)/divisor;
588 std::vector<Scalar> vec(samples);
589 for (std::size_t i = 0; i < samples; ++i)
590 vec[i] = begin + i*delta;
607template <
class Scalar>
613 const Scalar ln10 = 2.3025850929940459;
626template<
class ValueType>
627constexpr int sign(
const ValueType& value)
noexcept
629 return (ValueType(0) < value) - (value < ValueType(0));
639template <
class Scalar>
640Dune::FieldVector<Scalar, 3>
crossProduct(
const Dune::FieldVector<Scalar, 3> &vec1,
641 const Dune::FieldVector<Scalar, 3> &vec2)
643 return {vec1[1]*vec2[2]-vec1[2]*vec2[1],
644 vec1[2]*vec2[0]-vec1[0]*vec2[2],
645 vec1[0]*vec2[1]-vec1[1]*vec2[0]};
655template <
class Scalar>
657 const Dune::FieldVector<Scalar, 2> &vec2)
658{
return vec1[0]*vec2[1]-vec1[1]*vec2[0]; }
668template <
class Scalar>
670 const Dune::FieldVector<Scalar, 3> &vec2,
671 const Dune::FieldVector<Scalar, 3> &vec3)
672{
return crossProduct<Scalar>(vec1, vec2)*vec3; }
680template <
class Scalar,
int m,
int n>
681Dune::FieldMatrix<Scalar, n, m>
getTransposed(
const Dune::FieldMatrix<Scalar, m, n>& M)
683 Dune::FieldMatrix<Scalar, n, m> T;
684 for (std::size_t i = 0; i < m; ++i)
685 for (std::size_t j = 0; j < n; ++j)
697template <
class Scalar>
698Dune::DynamicMatrix<Scalar>
getTransposed(
const Dune::DynamicMatrix<Scalar>& M)
700 std::size_t rows_T = M.M();
701 std::size_t cols_T = M.N();
703 Dune::DynamicMatrix<Scalar> M_T(rows_T, cols_T, 0.0);
705 for (std::size_t i = 0; i < rows_T; ++i)
706 for (std::size_t j = 0; j < cols_T; ++j)
719template <
class Scalar>
721 const Dune::DynamicMatrix<Scalar> &M2)
723 using size_type =
typename Dune::DynamicMatrix<Scalar>::size_type;
724 const size_type rows = M1.N();
725 const size_type cols = M2.M();
727 DUNE_ASSERT_BOUNDS(M1.M() == M2.N());
729 Dune::DynamicMatrix<Scalar> result(rows, cols, 0.0);
730 for (size_type i = 0; i < rows; i++)
731 for (size_type j = 0; j < cols; j++)
732 for (size_type k = 0; k < M1.M(); k++)
733 result[i][j] += M1[i][k]*M2[k][j];
745template <
class Scalar,
int rows1,
int cols1,
int cols2>
746Dune::FieldMatrix<Scalar, rows1, cols2>
multiplyMatrices(
const Dune::FieldMatrix<Scalar, rows1, cols1> &M1,
747 const Dune::FieldMatrix<Scalar, cols1, cols2> &M2)
749 using size_type =
typename Dune::FieldMatrix<Scalar, rows1, cols2>::size_type;
751 Dune::FieldMatrix<Scalar, rows1, cols2> result(0.0);
752 for (size_type i = 0; i < rows1; i++)
753 for (size_type j = 0; j < cols2; j++)
754 for (size_type k = 0; k < cols1; k++)
755 result[i][j] += M1[i][k]*M2[k][j];
767template <
class MatrixType>
768typename Dune::DenseMatrix<MatrixType>::field_type
769trace(
const Dune::DenseMatrix<MatrixType>& M)
771 const auto rows = M.N();
772 DUNE_ASSERT_BOUNDS(rows == M.M());
774 using MatType = Dune::DenseMatrix<MatrixType>;
775 typename MatType::field_type
trace = 0.0;
777 for (
typename MatType::size_type i = 0; i < rows; ++i)
796template <
class MAT,
class V>
797typename Dune::DenseVector<V>::derived_type
798mv(
const Dune::DenseMatrix<MAT>& M,
799 const Dune::DenseVector<V>& v)
801 typename Dune::DenseVector<V>::derived_type res(v);
823template <
class FieldScalar,
class V>
824typename std::enable_if_t<Dune::IsNumber<FieldScalar>::value,
825 typename Dune::DenseVector<V>::derived_type>
826mv(
const FieldScalar m,
const Dune::DenseVector<V>& v)
828 typename Dune::DenseVector<V>::derived_type res(v);
847template <
class V1,
class MAT,
class V2>
848typename Dune::DenseMatrix<MAT>::value_type
849vtmv(
const Dune::DenseVector<V1>& v1,
850 const Dune::DenseMatrix<MAT>& M,
851 const Dune::DenseVector<V2>& v2)
875template <
class V1,
class FieldScalar,
class V2>
876typename std::enable_if_t<Dune::IsNumber<FieldScalar>::value, FieldScalar>
877vtmv(
const Dune::DenseVector<V1>& v1,
879 const Dune::DenseVector<V2>& v2)
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: geometry/distance.hh:138
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:68
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:849
Dune::DynamicMatrix< Scalar > multiplyMatrices(const Dune::DynamicMatrix< Scalar > &M1, const Dune::DynamicMatrix< Scalar > &M2)
Multiply two dynamic matrices.
Definition: math.hh:720
Scalar geometricMean(Scalar x, Scalar y) noexcept
Calculate the geometric mean of two scalar values.
Definition: math.hh:86
int invertLinearPolynomial(SolContainer &sol, Scalar a, Scalar b)
Invert a linear polynomial analytically.
Definition: math.hh:140
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:482
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:108
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:640
int invertCubicPolynomial(SolContainer *sol, Scalar a, Scalar b, Scalar c, Scalar d)
Invert a cubic polynomial analytically.
Definition: math.hh:242
Scalar interpolate(Scalar ip, Parameter &&... params)
a generic function to interpolate given a set of parameters and an interpolation point
Definition: math.hh:504
bool isSmaller(const Dune::FieldVector< Scalar, dim > &pos, const Dune::FieldVector< Scalar, dim > &largerVec)
Comparison of two position vectors.
Definition: math.hh:452
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:798
Dune::DenseMatrix< MatrixType >::field_type trace(const Dune::DenseMatrix< MatrixType > &M)
Trace of a dense matrix.
Definition: math.hh:769
Dune::FieldMatrix< Scalar, n, m > getTransposed(const Dune::FieldMatrix< Scalar, m, n > &M)
Transpose a FieldMatrix.
Definition: math.hh:681
bool isLarger(const Dune::FieldVector< Scalar, dim > &pos, const Dune::FieldVector< Scalar, dim > &smallerVec)
Comparison of two position vectors.
Definition: math.hh:426
int invertQuadraticPolynomial(SolContainer &sol, Scalar a, Scalar b, Scalar c)
Invert a quadratic polynomial analytically.
Definition: math.hh:168
std::vector< Scalar > linspace(const Scalar begin, const Scalar end, std::size_t samples, bool endPoint=true)
Generates linearly spaced vectors.
Definition: math.hh:580
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:669
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:627
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:608
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:49
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:518
static constexpr Scalar interpolate(Scalar ip, const std::array< Scalar, 2 > ¶ms)
interpolate linearly between two given values
Definition: math.hh:525
interpolate linearly in a piecewise linear function (tabularized function)
Definition: math.hh:536
static constexpr Scalar interpolate(Scalar ip, const std::pair< RandomAccessContainer, RandomAccessContainer > &table)
Definition: math.hh:561
static constexpr Scalar interpolate(Scalar ip, const RandomAccessContainer0 &range, const RandomAccessContainer1 &values)
interpolate linearly in a piecewise linear function (tabularized function)
Definition: math.hh:545