3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_MATH_HH
25#define DUMUX_MATH_HH
26
27#include <algorithm>
28#include <numeric>
29#include <cmath>
30#include <utility>
31
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>
37
38namespace Dumux {
39
49template <class Scalar>
50constexpr Scalar arithmeticMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
51{
52 static_assert(Dune::IsNumber<Scalar>::value, "The arguments x, y, wx, and wy have to be numbers!");
53
54 if (x*y <= 0)
55 return 0;
56 return (x * wx + y * wy)/(wx + wy);
57}
58
68template <class Scalar>
69constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
70{
71 static_assert(Dune::IsNumber<Scalar>::value, "The arguments x, y, wx, and wy have to be numbers!");
72
73 if (x*y <= 0)
74 return 0;
75 return (wx + wy) * x * y / (wy * x + wx * y);
76}
77
86template <class Scalar>
87Scalar geometricMean(Scalar x, Scalar y) noexcept
88{
89 static_assert(Dune::IsNumber<Scalar>::value, "The arguments x and y have to be numbers!");
90
91 if (x*y <= 0)
92 return 0;
93 using std::sqrt;
94 return sqrt(x*y)*sign(x);
95}
96
108template <class Scalar, int m, int n>
109void harmonicMeanMatrix(Dune::FieldMatrix<Scalar, m, n> &K,
110 const Dune::FieldMatrix<Scalar, m, n> &Ki,
111 const Dune::FieldMatrix<Scalar, m, n> &Kj)
112{
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])) {
116 K[rowIdx][colIdx] =
117 harmonicMean(Ki[rowIdx][colIdx],
118 Kj[rowIdx][colIdx]);
119 }
120 else
121 K[rowIdx][colIdx] = Ki[rowIdx][colIdx];
122 }
123 }
124}
125
140template <class Scalar, class SolContainer>
141int invertLinearPolynomial(SolContainer &sol,
142 Scalar a,
143 Scalar b)
144{
145 if (a == 0.0)
146 return 0;
147
148 sol[0] = -b/a;
149 return 1;
150}
151
168template <class Scalar, class SolContainer>
169int invertQuadraticPolynomial(SolContainer &sol,
170 Scalar a,
171 Scalar b,
172 Scalar c)
173{
174 // check for a line
175 if (a == 0.0)
176 return invertLinearPolynomial(sol, b, c);
177
178 // discriminant
179 Scalar Delta = b*b - 4*a*c;
180 if (Delta < 0)
181 return 0; // no real roots
182
183 using std::sqrt;
184 Delta = sqrt(Delta);
185 sol[0] = (- b + Delta)/(2*a);
186 sol[1] = (- b - Delta)/(2*a);
187
188 // sort the result
189 if (sol[0] > sol[1])
190 {
191 using std::swap;
192 swap(sol[0], sol[1]);
193 }
194 return 2; // two real roots
195}
196
198template <class Scalar, class SolContainer>
199void invertCubicPolynomialPostProcess_(SolContainer &sol,
200 int numSol,
201 Scalar a, Scalar b, Scalar c, Scalar d,
202 std::size_t numIterations = 1)
203{
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); };
206
207 // do numIterations Newton iterations on the analytic solution
208 // and update result if the precision is increased
209 for (int i = 0; i < numSol; ++i)
210 {
211 Scalar x = sol[i];
212 Scalar fCurrent = eval(x);
213 const Scalar fOld = fCurrent;
214 for (int j = 0; j < numIterations; ++j)
215 {
216 const Scalar fPrime = evalDeriv(x);
217 if (fPrime == 0.0)
218 break;
219 x -= fCurrent/fPrime;
220 fCurrent = eval(x);
221 }
222
223 using std::abs;
224 if (abs(fCurrent) < abs(fOld))
225 sol[i] = x;
226 }
227}
229
251template <class Scalar, class SolContainer>
252int invertCubicPolynomial(SolContainer *sol,
253 Scalar a, Scalar b, Scalar c, Scalar d,
254 std::size_t numPostProcessIterations = 1)
255{
256 // reduces to a quadratic polynomial
257 if (a == 0)
258 return invertQuadraticPolynomial(sol, b, c, d);
259
260 // normalize the polynomial
261 b /= a;
262 c /= a;
263 d /= a;
264 a = 1;
265
266 // get rid of the quadratic term by substituting x = t - b/3
267 Scalar p = c - b*b/3;
268 Scalar q = d + (2*b*b*b - 9*b*c)/27;
269
270 // now we are at the form t^3 + p*t + q = 0. First we handle some
271 // special cases to avoid divisions by zero later...
272 if (p == 0.0 && q == 0.0) {
273 // t^3 = 0, i.e. triple root at t = 0
274 sol[0] = sol[1] = sol[2] = 0.0 - b/3;
275 return 3;
276 }
277 else if (p == 0.0 && q != 0.0) {
278 // t^3 + q = 0,
279 //
280 // i. e. single real root at t=curt(q)
281 using std::cbrt;
282 Scalar t = cbrt(q);
283 sol[0] = t - b/3;
284
285 return 1;
286 }
287 else if (p != 0.0 && q == 0.0) {
288 // t^3 + p*t = 0 = t*(t^2 + p),
289 //
290 // i. e. roots at t = 0, t^2 + p = 0
291 if (p > 0) {
292 sol[0] = 0.0 - b/3;
293 return 1; // only a single real root at t=0
294 }
295
296 // two additional real roots at t = sqrt(-p) and t = -sqrt(-p)
297 using std::sqrt;
298 sol[0] = -sqrt(-p) - b/3;
299 sol[1] = 0.0 - b/3;
300 sol[2] = sqrt(-p) - b/3;
301
302 return 3;
303 }
304
305 // At this point
306 //
307 // t^3 + p*t + q = 0
308 //
309 // with p != 0 and q != 0 holds. Introducing the variables u and v
310 // with the properties
311 //
312 // u + v = t and 3*u*v + p = 0
313 //
314 // leads to
315 //
316 // u^3 + v^3 + q = 0 .
317 //
318 // multiplying both sides with u^3 and taking advantage of the
319 // fact that u*v = -p/3 leads to
320 //
321 // u^6 + q*u^3 - p^3/27 = 0
322 //
323 // Now, substituting u^3 = w yields
324 //
325 // w^2 + q*w - p^3/27 = 0
326 //
327 // This is a quadratic equation with the solutions
328 //
329 // w = -q/2 + sqrt(q^2/4 + p^3/27) and
330 // w = -q/2 - sqrt(q^2/4 + p^3/27)
331 //
332 // Since w is equivalent to u^3 it is sufficient to only look at
333 // one of the two cases. Then, there are still 2 cases: positive
334 // and negative discriminant.
335 Scalar wDisc = q*q/4 + p*p*p/27;
336 if (wDisc >= 0) { // the positive discriminant case:
337 // calculate the cube root of - q/2 + sqrt(q^2/4 + p^3/27)
338 using std::cbrt;
339 using std::sqrt;
340
341 // Choose the root that is safe against the loss of precision
342 // that can cause wDisc to be equal to q*q/4 despite p != 0.
343 // We do not want u to be zero in that case. Mathematically,
344 // we are happy with either root.
345 const Scalar u = [&]{
346 return q > 0 ? cbrt(-0.5*q - sqrt(wDisc)) : cbrt(-0.5*q + sqrt(wDisc));
347 }();
348 // at this point, u != 0 since p^3 = 0 is necessary in order
349 // for u = 0 to hold, so
350 sol[0] = u - p/(3*u) - b/3;
351 // does not produce a division by zero. the remaining two
352 // roots of u are rotated by +- 2/3*pi in the complex plane
353 // and thus not considered here
354 invertCubicPolynomialPostProcess_(sol, 1, a, b, c, d, numPostProcessIterations);
355 return 1;
356 }
357 else { // the negative discriminant case:
358 Scalar uCubedRe = - q/2;
359 using std::sqrt;
360 Scalar uCubedIm = sqrt(-wDisc);
361 // calculate the cube root of - q/2 + sqrt(q^2/4 + p^3/27)
362 using std::cbrt;
363 Scalar uAbs = cbrt(sqrt(uCubedRe*uCubedRe + uCubedIm*uCubedIm));
364 using std::atan2;
365 Scalar phi = atan2(uCubedIm, uCubedRe)/3;
366
367 // calculate the length and the angle of the primitive root
368
369 // with the definitions from above it follows that
370 //
371 // x = u - p/(3*u) - b/3
372 //
373 // where x and u are complex numbers. Rewritten in polar form
374 // this is equivalent to
375 //
376 // x = |u|*e^(i*phi) - p*e^(-i*phi)/(3*|u|) - b/3 .
377 //
378 // Factoring out the e^ terms and subtracting the additional
379 // terms, yields
380 //
381 // x = (e^(i*phi) + e^(-i*phi))*(|u| - p/(3*|u|)) - y - b/3
382 //
383 // with
384 //
385 // y = - |u|*e^(-i*phi) + p*e^(i*phi)/(3*|u|) .
386 //
387 // The crucial observation is the fact that y is the conjugate
388 // of - x + b/3. This means that after taking advantage of the
389 // relation
390 //
391 // e^(i*phi) + e^(-i*phi) = 2*cos(phi)
392 //
393 // the equation
394 //
395 // x = 2*cos(phi)*(|u| - p / (3*|u|)) - conj(x) - 2*b/3
396 //
397 // holds. Since |u|, p, b and cos(phi) are real numbers, it
398 // follows that Im(x) = - Im(x) and thus Im(x) = 0. This
399 // implies
400 //
401 // Re(x) = x = cos(phi)*(|u| - p / (3*|u|)) - b/3 .
402 //
403 // Considering the fact that u is a cubic root, we have three
404 // values for phi which differ by 2/3*pi. This allows to
405 // calculate the three real roots of the polynomial:
406 for (int i = 0; i < 3; ++i) {
407 using std::cos;
408 sol[i] = cos(phi)*(uAbs - p/(3*uAbs)) - b/3;
409 phi += 2*M_PI/3;
410 }
411
412 // post process the obtained solution to increase numerical
413 // precision
414 invertCubicPolynomialPostProcess_(sol, 3, a, b, c, d, numPostProcessIterations);
415
416 // sort the result
417 using std::sort;
418 sort(sol, sol + 3);
419
420 return 3;
421 }
422
423 // NOT REACHABLE!
424 return 0;
425}
426
439template <class Scalar, int dim>
440bool isLarger(const Dune::FieldVector<Scalar, dim> &pos,
441 const Dune::FieldVector<Scalar, dim> &smallerVec)
442{
443 for (int i=0; i < dim; i++)
444 {
445 if (pos[i]<= smallerVec[i])
446 {
447 return false;
448 }
449 }
450 return true;
451}
452
465template <class Scalar, int dim>
466bool isSmaller(const Dune::FieldVector<Scalar, dim> &pos,
467 const Dune::FieldVector<Scalar, dim> &largerVec)
468{
469 for (int i=0; i < dim; i++)
470 {
471 if (pos[i]>= largerVec[i])
472 {
473 return false;
474 }
475 }
476 return true;
477}
478
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)
499{
500 if (isLarger(pos, smallerVec) && isSmaller(pos, largerVec))
501 {
502 return true;
503 }
504 else
505 return false;
506}
507
509namespace InterpolationPolicy { struct Linear; }
510
517template <class Policy = InterpolationPolicy::Linear, class Scalar, class ... Parameter>
518Scalar interpolate(Scalar ip, Parameter&& ... params)
519{ return Policy::interpolate(ip, std::forward<Parameter>(params) ...); }
520
525namespace InterpolationPolicy {
526
531struct Linear
532{
538 template<class Scalar>
539 static constexpr Scalar interpolate(Scalar ip, const std::array<Scalar, 2>& params)
540 {
541 return params[0]*(1.0 - ip) + params[1]*ip;
542 }
543};
544
550{
558 template<class Scalar, class RandomAccessContainer0, class RandomAccessContainer1>
559 static constexpr Scalar interpolate(Scalar ip, const RandomAccessContainer0& range, const RandomAccessContainer1& values)
560 {
561 // check bounds
562 if (ip > range.back()) return values.back();
563 if (ip < range[0]) return values[0];
564
565 // if we are within bounds find the index of the lower bound
566 const auto lookUpIndex = std::distance(range.begin(), std::lower_bound(range.begin(), range.end(), ip));
567 if (lookUpIndex == 0)
568 return values[0];
569
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]}});
572 }
573
574 template<class Scalar, class RandomAccessContainer>
575 static constexpr Scalar interpolate(Scalar ip, const std::pair<RandomAccessContainer, RandomAccessContainer>& table)
576 {
577 const auto& [range, values] = table;
578 return interpolate(ip, range, values);
579 }
580};
581
582} // end namespace InterpolationPolicy
583
593template <class Scalar>
594std::vector<Scalar> linspace(const Scalar begin, const Scalar end,
595 std::size_t samples,
596 bool endPoint = true)
597{
598 using std::max;
599 samples = max(std::size_t{2}, samples); // only makes sense for 2 or more 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;
605 return vec;
606}
607
608
621template <class Scalar>
622Scalar antoine(Scalar temperature,
623 Scalar A,
624 Scalar B,
625 Scalar C)
626{
627 const Scalar ln10 = 2.3025850929940459;
628 using std::exp;
629 return exp(ln10*(A - B/(C + temperature)));
630}
631
640template<class ValueType>
641constexpr int sign(const ValueType& value) noexcept
642{
643 return (ValueType(0) < value) - (value < ValueType(0));
644}
645
653template <class Scalar>
654Dune::FieldVector<Scalar, 3> crossProduct(const Dune::FieldVector<Scalar, 3> &vec1,
655 const Dune::FieldVector<Scalar, 3> &vec2)
656{
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]};
660}
661
669template <class Scalar>
670Scalar crossProduct(const Dune::FieldVector<Scalar, 2> &vec1,
671 const Dune::FieldVector<Scalar, 2> &vec2)
672{ return vec1[0]*vec2[1]-vec1[1]*vec2[0]; }
673
682template <class Scalar>
683Scalar tripleProduct(const Dune::FieldVector<Scalar, 3> &vec1,
684 const Dune::FieldVector<Scalar, 3> &vec2,
685 const Dune::FieldVector<Scalar, 3> &vec3)
686{ return crossProduct<Scalar>(vec1, vec2)*vec3; }
687
694template <class Scalar, int m, int n>
695Dune::FieldMatrix<Scalar, n, m> getTransposed(const Dune::FieldMatrix<Scalar, m, n>& M)
696{
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)
700 T[j][i] = M[i][j];
701
702 return T;
703}
704
711template <class Scalar>
712Dune::DynamicMatrix<Scalar> getTransposed(const Dune::DynamicMatrix<Scalar>& M)
713{
714 std::size_t rows_T = M.M();
715 std::size_t cols_T = M.N();
716
717 Dune::DynamicMatrix<Scalar> M_T(rows_T, cols_T, 0.0);
718
719 for (std::size_t i = 0; i < rows_T; ++i)
720 for (std::size_t j = 0; j < cols_T; ++j)
721 M_T[i][j] = M[j][i];
722
723 return M_T;
724}
725
733template <class Scalar>
734Dune::DynamicMatrix<Scalar> multiplyMatrices(const Dune::DynamicMatrix<Scalar> &M1,
735 const Dune::DynamicMatrix<Scalar> &M2)
736{
737 using size_type = typename Dune::DynamicMatrix<Scalar>::size_type;
738 const size_type rows = M1.N();
739 const size_type cols = M2.M();
740
741 DUNE_ASSERT_BOUNDS(M1.M() == M2.N());
742
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];
748
749 return result;
750}
751
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)
762{
763 using size_type = typename Dune::FieldMatrix<Scalar, rows1, cols2>::size_type;
764
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];
770
771 return result;
772}
773
774
781template <class MatrixType>
782typename Dune::DenseMatrix<MatrixType>::field_type
783trace(const Dune::DenseMatrix<MatrixType>& M)
784{
785 const auto rows = M.N();
786 DUNE_ASSERT_BOUNDS(rows == M.M()); // rows == cols
787
788 using MatType = Dune::DenseMatrix<MatrixType>;
789 typename MatType::field_type trace = 0.0;
790
791 for (typename MatType::size_type i = 0; i < rows; ++i)
792 trace += M[i][i];
793
794 return trace;
795}
796
810template <class MAT, class V>
811typename Dune::DenseVector<V>::derived_type
812mv(const Dune::DenseMatrix<MAT>& M,
813 const Dune::DenseVector<V>& v)
814{
815 typename Dune::DenseVector<V>::derived_type res(v);
816 M.mv(v, res);
817 return res;
818}
819
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)
841{
842 typename Dune::DenseVector<V>::derived_type res(v);
843 res *= m;
844 return res;
845}
846
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)
866{
867 return v1*mv(M, v2);
868}
869
889template <class V1, class FieldScalar, class V2>
890typename std::enable_if_t<Dune::IsNumber<FieldScalar>::value, FieldScalar>
891vtmv(const Dune::DenseVector<V1>& v1,
892 const FieldScalar m,
893 const Dune::DenseVector<V2>& v2)
894{
895 return m*(v1*v2);
896}
897
909template <class Scalar>
910std::array<Scalar, 2> linearRegression(const std::vector<Scalar>& x,
911 const std::vector<Scalar>& y)
912{
913 if (x.size() != y.size())
914 DUNE_THROW(Dune::InvalidStateException, "x and y array must have the same length.");
915
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();
918
919 // calculate temporary variables necessary for slope computation
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); }
923 );
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); }
927 );
928
929 // compute slope and intercept of the regression line
930 const Scalar slope = numerator / denominator;
931 const Scalar intercept = averageY - slope * averageX;
932
933 return {intercept, slope};
934}
935
936} // end namespace Dumux
937
938#endif
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:294
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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 > &params)
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