3.1-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 <cmath>
29#include <utility>
30
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>
36
37namespace Dumux {
38
48template <class Scalar>
49constexpr Scalar arithmeticMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
50{
51 static_assert(Dune::IsNumber<Scalar>::value, "The arguments x, y, wx, and wy have to be numbers!");
52
53 if (x*y <= 0)
54 return 0;
55 return (x * wx + y * wy)/(wx + wy);
56}
57
67template <class Scalar>
68constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
69{
70 static_assert(Dune::IsNumber<Scalar>::value, "The arguments x, y, wx, and wy have to be numbers!");
71
72 if (x*y <= 0)
73 return 0;
74 return (wx + wy) * x * y / (wy * x + wx * y);
75}
76
85template <class Scalar>
86Scalar geometricMean(Scalar x, Scalar y) noexcept
87{
88 static_assert(Dune::IsNumber<Scalar>::value, "The arguments x and y have to be numbers!");
89
90 if (x*y <= 0)
91 return 0;
92 using std::sqrt;
93 return sqrt(x*y)*sign(x);
94}
95
107template <class Scalar, int m, int n>
108void harmonicMeanMatrix(Dune::FieldMatrix<Scalar, m, n> &K,
109 const Dune::FieldMatrix<Scalar, m, n> &Ki,
110 const Dune::FieldMatrix<Scalar, m, n> &Kj)
111{
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])) {
115 K[rowIdx][colIdx] =
116 harmonicMean(Ki[rowIdx][colIdx],
117 Kj[rowIdx][colIdx]);
118 }
119 else
120 K[rowIdx][colIdx] = Ki[rowIdx][colIdx];
121 }
122 }
123}
124
139template <class Scalar, class SolContainer>
140int invertLinearPolynomial(SolContainer &sol,
141 Scalar a,
142 Scalar b)
143{
144 if (a == 0.0)
145 return 0;
146
147 sol[0] = -b/a;
148 return 1;
149}
150
167template <class Scalar, class SolContainer>
168int invertQuadraticPolynomial(SolContainer &sol,
169 Scalar a,
170 Scalar b,
171 Scalar c)
172{
173 // check for a line
174 if (a == 0.0)
175 return invertLinearPolynomial(sol, b, c);
176
177 // discriminant
178 Scalar Delta = b*b - 4*a*c;
179 if (Delta < 0)
180 return 0; // no real roots
181
182 using std::sqrt;
183 Delta = sqrt(Delta);
184 sol[0] = (- b + Delta)/(2*a);
185 sol[1] = (- b - Delta)/(2*a);
186
187 // sort the result
188 if (sol[0] > sol[1])
189 {
190 using std::swap;
191 swap(sol[0], sol[1]);
192 }
193 return 2; // two real roots
194}
195
197template <class Scalar, class SolContainer>
198void invertCubicPolynomialPostProcess_(SolContainer &sol,
199 int numSol,
200 Scalar a,
201 Scalar b,
202 Scalar c,
203 Scalar d)
204{
205 // do one Newton iteration on the analytic solution if the
206 // precision is increased
207 for (int i = 0; i < numSol; ++i) {
208 Scalar x = sol[i];
209 Scalar fOld = d + x*(c + x*(b + x*a));
210
211 Scalar fPrime = c + x*(2*b + x*3*a);
212 if (fPrime == 0.0)
213 continue;
214 x -= fOld/fPrime;
215
216 Scalar fNew = d + x*(c + x*(b + x*a));
217 using std::abs;
218 if (abs(fNew) < abs(fOld))
219 sol[i] = x;
220 }
221}
223
241template <class Scalar, class SolContainer>
242int invertCubicPolynomial(SolContainer *sol,
243 Scalar a,
244 Scalar b,
245 Scalar c,
246 Scalar d)
247{
248 // reduces to a quadratic polynomial
249 if (a == 0)
250 return invertQuadraticPolynomial(sol, b, c, d);
251
252 // normalize the polynomial
253 b /= a;
254 c /= a;
255 d /= a;
256 a = 1;
257
258 // get rid of the quadratic term by subsituting x = t - b/3
259 Scalar p = c - b*b/3;
260 Scalar q = d + (2*b*b*b - 9*b*c)/27;
261
262 // now we are at the form t^3 + p*t + q = 0. First we handle some
263 // special cases to avoid divisions by zero later...
264 if (p == 0.0 && q == 0.0) {
265 // t^3 = 0, i.e. triple root at t = 0
266 sol[0] = sol[1] = sol[2] = 0.0 - b/3;
267 return 3;
268 }
269 else if (p == 0.0 && q != 0.0) {
270 // t^3 + q = 0,
271 //
272 // i. e. single real root at t=curt(q)
273 using std::cbrt;
274 Scalar t = cbrt(q);
275 sol[0] = t - b/3;
276
277 return 1;
278 }
279 else if (p != 0.0 && q == 0.0) {
280 // t^3 + p*t = 0 = t*(t^2 + p),
281 //
282 // i. e. roots at t = 0, t^2 + p = 0
283 if (p > 0) {
284 sol[0] = 0.0 - b/3;
285 return 1; // only a single real root at t=0
286 }
287
288 // two additional real roots at t = sqrt(-p) and t = -sqrt(-p)
289 using std::sqrt;
290 sol[0] = -sqrt(-p) - b/3;
291 sol[1] = 0.0 - b/3;
292 sol[2] = sqrt(-p) - b/3;
293
294 return 3;
295 }
296
297 // At this point
298 //
299 // t^3 + p*t + q = 0
300 //
301 // with p != 0 and q != 0 holds. Introducing the variables u and v
302 // with the properties
303 //
304 // u + v = t and 3*u*v + p = 0
305 //
306 // leads to
307 //
308 // u^3 + v^3 + q = 0 .
309 //
310 // multiplying both sides with u^3 and taking advantage of the
311 // fact that u*v = -p/3 leads to
312 //
313 // u^6 + q*u^3 - p^3/27 = 0
314 //
315 // Now, substituting u^3 = w yields
316 //
317 // w^2 + q*w - p^3/27 = 0
318 //
319 // This is a quadratic equation with the solutions
320 //
321 // w = -q/2 + sqrt(q^2/4 + p^3/27) and
322 // w = -q/2 - sqrt(q^2/4 + p^3/27)
323 //
324 // Since w is equivalent to u^3 it is sufficient to only look at
325 // one of the two cases. Then, there are still 2 cases: positive
326 // and negative discriminant.
327 Scalar wDisc = q*q/4 + p*p*p/27;
328 if (wDisc >= 0) { // the positive discriminant case:
329 // calculate the cube root of - q/2 + sqrt(q^2/4 + p^3/27)
330 using std::cbrt;
331 using std::sqrt;
332 Scalar u = cbrt(-q/2 + sqrt(wDisc));
333
334 // at this point, u != 0 since p^3 = 0 is necessary in order
335 // for u = 0 to hold, so
336 sol[0] = u - p/(3*u) - b/3;
337 // does not produce a division by zero. the remaining two
338 // roots of u are rotated by +- 2/3*pi in the complex plane
339 // and thus not considered here
340 invertCubicPolynomialPostProcess_(sol, 1, a, b, c, d);
341 return 1;
342 }
343 else { // the negative discriminant case:
344 Scalar uCubedRe = - q/2;
345 using std::sqrt;
346 Scalar uCubedIm = sqrt(-wDisc);
347 // calculate the cube root of - q/2 + sqrt(q^2/4 + p^3/27)
348 using std::cbrt;
349 Scalar uAbs = cbrt(sqrt(uCubedRe*uCubedRe + uCubedIm*uCubedIm));
350 using std::atan2;
351 Scalar phi = atan2(uCubedIm, uCubedRe)/3;
352
353 // calculate the length and the angle of the primitive root
354
355 // with the definitions from above it follows that
356 //
357 // x = u - p/(3*u) - b/3
358 //
359 // where x and u are complex numbers. Rewritten in polar form
360 // this is equivalent to
361 //
362 // x = |u|*e^(i*phi) - p*e^(-i*phi)/(3*|u|) - b/3 .
363 //
364 // Factoring out the e^ terms and subtracting the additional
365 // terms, yields
366 //
367 // x = (e^(i*phi) + e^(-i*phi))*(|u| - p/(3*|u|)) - y - b/3
368 //
369 // with
370 //
371 // y = - |u|*e^(-i*phi) + p*e^(i*phi)/(3*|u|) .
372 //
373 // The crucial observation is the fact that y is the conjugate
374 // of - x + b/3. This means that after taking advantage of the
375 // relation
376 //
377 // e^(i*phi) + e^(-i*phi) = 2*cos(phi)
378 //
379 // the equation
380 //
381 // x = 2*cos(phi)*(|u| - p / (3*|u|)) - conj(x) - 2*b/3
382 //
383 // holds. Since |u|, p, b and cos(phi) are real numbers, it
384 // follows that Im(x) = - Im(x) and thus Im(x) = 0. This
385 // implies
386 //
387 // Re(x) = x = cos(phi)*(|u| - p / (3*|u|)) - b/3 .
388 //
389 // Considering the fact that u is a cubic root, we have three
390 // values for phi which differ by 2/3*pi. This allows to
391 // calculate the three real roots of the polynomial:
392 for (int i = 0; i < 3; ++i) {
393 using std::cos;
394 sol[i] = cos(phi)*(uAbs - p/(3*uAbs)) - b/3;
395 phi += 2*M_PI/3;
396 }
397
398 // post process the obtained solution to increase numerical
399 // precision
400 invertCubicPolynomialPostProcess_(sol, 3, a, b, c, d);
401
402 // sort the result
403 using std::sort;
404 sort(sol, sol + 3);
405
406 return 3;
407 }
408
409 // NOT REACHABLE!
410 return 0;
411}
412
425template <class Scalar, int dim>
426bool isLarger(const Dune::FieldVector<Scalar, dim> &pos,
427 const Dune::FieldVector<Scalar, dim> &smallerVec)
428{
429 for (int i=0; i < dim; i++)
430 {
431 if (pos[i]<= smallerVec[i])
432 {
433 return false;
434 }
435 }
436 return true;
437}
438
451template <class Scalar, int dim>
452bool isSmaller(const Dune::FieldVector<Scalar, dim> &pos,
453 const Dune::FieldVector<Scalar, dim> &largerVec)
454{
455 for (int i=0; i < dim; i++)
456 {
457 if (pos[i]>= largerVec[i])
458 {
459 return false;
460 }
461 }
462 return true;
463}
464
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)
485{
486 if (isLarger(pos, smallerVec) && isSmaller(pos, largerVec))
487 {
488 return true;
489 }
490 else
491 return false;
492}
493
495namespace InterpolationPolicy { struct Linear; }
496
503template <class Policy = InterpolationPolicy::Linear, class Scalar, class Parameter>
504Scalar interpolate(Scalar ip, Parameter&& params)
505{ return Policy::interpolate(ip, std::forward<Parameter>(params)); }
506
511namespace InterpolationPolicy {
512
517struct Linear
518{
524 template<class Scalar>
525 static constexpr Scalar interpolate(Scalar ip, const std::array<Scalar, 2>& params)
526 {
527 return params[0]*(1.0 - ip) + params[1]*ip;
528 }
529};
530
536{
543 template<class Scalar, class RandomAccessContainer>
544 static constexpr Scalar interpolate(Scalar ip, const std::pair<RandomAccessContainer, RandomAccessContainer>& table)
545 {
546 const auto& range = table.first;
547 const auto& values = table.second;
548
549 // check bounds
550 if (ip > range.back()) return values.back();
551 if (ip < range[0]) return values[0];
552
553 // if we are within bounds find the index of the lower bound
554 const auto lookUpIndex = std::distance(range.begin(), std::lower_bound(range.begin(), range.end(), ip));
555 if (lookUpIndex == 0)
556 return values[0];
557
558 const auto ipLinear = (ip - range[lookUpIndex-1])/(range[lookUpIndex] - range[lookUpIndex-1]);
559 return Dumux::interpolate<Linear>(ipLinear, std::array<Scalar, 2>{{values[lookUpIndex-1], values[lookUpIndex]}});
560 }
561};
562
563} // end namespace InterpolationPolicy
564
573template <class Scalar>
574std::vector<Scalar> linspace(const Scalar begin, const Scalar end, std::size_t samples)
575{
576 using std::max;
577 samples = max(std::size_t{2}, samples); // only makes sense for 2 or more samples
578 const Scalar delta = (end-begin)/static_cast<Scalar>(samples-1);
579 std::vector<Scalar> vec(samples);
580 for (std::size_t i = 0; i < samples; ++i)
581 vec[i] = begin + i*delta;
582 return vec;
583}
584
585
598template <class Scalar>
599Scalar antoine(Scalar temperature,
600 Scalar A,
601 Scalar B,
602 Scalar C)
603{
604 const Scalar ln10 = 2.3025850929940459;
605 using std::exp;
606 return exp(ln10*(A - B/(C + temperature)));
607}
608
617template<class ValueType>
618constexpr int sign(const ValueType& value) noexcept
619{
620 return (ValueType(0) < value) - (value < ValueType(0));
621}
622
630template <class Scalar>
631Dune::FieldVector<Scalar, 3> crossProduct(const Dune::FieldVector<Scalar, 3> &vec1,
632 const Dune::FieldVector<Scalar, 3> &vec2)
633{
634 return {vec1[1]*vec2[2]-vec1[2]*vec2[1],
635 vec1[2]*vec2[0]-vec1[0]*vec2[2],
636 vec1[0]*vec2[1]-vec1[1]*vec2[0]};
637}
638
646template <class Scalar>
647Scalar crossProduct(const Dune::FieldVector<Scalar, 2> &vec1,
648 const Dune::FieldVector<Scalar, 2> &vec2)
649{ return vec1[0]*vec2[1]-vec1[1]*vec2[0]; }
650
659template <class Scalar>
660Scalar tripleProduct(const Dune::FieldVector<Scalar, 3> &vec1,
661 const Dune::FieldVector<Scalar, 3> &vec2,
662 const Dune::FieldVector<Scalar, 3> &vec3)
663{ return crossProduct<Scalar>(vec1, vec2)*vec3; }
664
671template <class Scalar, int m, int n>
672Dune::FieldMatrix<Scalar, n, m> getTransposed(const Dune::FieldMatrix<Scalar, m, n>& M)
673{
674 Dune::FieldMatrix<Scalar, n, m> T;
675 for (std::size_t i = 0; i < m; ++i)
676 for (std::size_t j = 0; j < n; ++j)
677 T[j][i] = M[i][j];
678
679 return T;
680}
681
688template <class Scalar>
689Dune::DynamicMatrix<Scalar> getTransposed(const Dune::DynamicMatrix<Scalar>& M)
690{
691 std::size_t rows_T = M.M();
692 std::size_t cols_T = M.N();
693
694 Dune::DynamicMatrix<Scalar> M_T(rows_T, cols_T, 0.0);
695
696 for (std::size_t i = 0; i < rows_T; ++i)
697 for (std::size_t j = 0; j < cols_T; ++j)
698 M_T[i][j] = M[j][i];
699
700 return M_T;
701}
702
710template <class Scalar>
711Dune::DynamicMatrix<Scalar> multiplyMatrices(const Dune::DynamicMatrix<Scalar> &M1,
712 const Dune::DynamicMatrix<Scalar> &M2)
713{
714 using size_type = typename Dune::DynamicMatrix<Scalar>::size_type;
715 const size_type rows = M1.N();
716 const size_type cols = M2.M();
717
718 DUNE_ASSERT_BOUNDS(M1.M() == M2.N());
719
720 Dune::DynamicMatrix<Scalar> result(rows, cols, 0.0);
721 for (size_type i = 0; i < rows; i++)
722 for (size_type j = 0; j < cols; j++)
723 for (size_type k = 0; k < M1.M(); k++)
724 result[i][j] += M1[i][k]*M2[k][j];
725
726 return result;
727}
728
736template <class Scalar, int rows1, int cols1, int cols2>
737Dune::FieldMatrix<Scalar, rows1, cols2> multiplyMatrices(const Dune::FieldMatrix<Scalar, rows1, cols1> &M1,
738 const Dune::FieldMatrix<Scalar, cols1, cols2> &M2)
739{
740 using size_type = typename Dune::FieldMatrix<Scalar, rows1, cols2>::size_type;
741
742 Dune::FieldMatrix<Scalar, rows1, cols2> result(0.0);
743 for (size_type i = 0; i < rows1; i++)
744 for (size_type j = 0; j < cols2; j++)
745 for (size_type k = 0; k < cols1; k++)
746 result[i][j] += M1[i][k]*M2[k][j];
747
748 return result;
749}
750
751
758template <class MatrixType>
759typename Dune::DenseMatrix<MatrixType>::field_type
760trace(const Dune::DenseMatrix<MatrixType>& M)
761{
762 const auto rows = M.N();
763 DUNE_ASSERT_BOUNDS(rows == M.M()); // rows == cols
764
765 using MatType = Dune::DenseMatrix<MatrixType>;
766 typename MatType::field_type trace = 0.0;
767
768 for (typename MatType::size_type i = 0; i < rows; ++i)
769 trace += M[i][i];
770
771 return trace;
772}
773
787template <class MAT, class V>
788typename Dune::DenseVector<V>::derived_type
789mv(const Dune::DenseMatrix<MAT>& M,
790 const Dune::DenseVector<V>& v)
791{
792 typename Dune::DenseVector<V>::derived_type res(v);
793 M.mv(v, res);
794 return res;
795}
796
814template <class FieldScalar, class V>
815typename std::enable_if_t<Dune::IsNumber<FieldScalar>::value,
816 typename Dune::DenseVector<V>::derived_type>
817mv(const FieldScalar m, const Dune::DenseVector<V>& v)
818{
819 typename Dune::DenseVector<V>::derived_type res(v);
820 res *= m;
821 return res;
822}
823
838template <class V1, class MAT, class V2>
839typename Dune::DenseMatrix<MAT>::value_type
840vtmv(const Dune::DenseVector<V1>& v1,
841 const Dune::DenseMatrix<MAT>& M,
842 const Dune::DenseVector<V2>& v2)
843{
844 return v1*mv(M, v2);
845}
846
866template <class V1, class FieldScalar, class V2>
867typename std::enable_if_t<Dune::IsNumber<FieldScalar>::value, FieldScalar>
868vtmv(const Dune::DenseVector<V1>& v1,
869 const FieldScalar m,
870 const Dune::DenseVector<V2>& v2)
871{
872 return m*(v1*v2);
873}
874
875} // end namespace Dumux
876
877#endif
std::vector< Scalar > linspace(const Scalar begin, const Scalar end, std::size_t samples)
Generates linearly spaced vectors.
Definition: math.hh:574
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:840
Dune::DynamicMatrix< Scalar > multiplyMatrices(const Dune::DynamicMatrix< Scalar > &M1, const Dune::DynamicMatrix< Scalar > &M2)
Multiply two dynamic matrices.
Definition: math.hh:711
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:631
int invertCubicPolynomial(SolContainer *sol, Scalar a, Scalar b, Scalar c, Scalar d)
Invert a cubic polynomial analytically.
Definition: math.hh:242
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:789
Dune::DenseMatrix< MatrixType >::field_type trace(const Dune::DenseMatrix< MatrixType > &M)
Trace of a dense matrix.
Definition: math.hh:760
Dune::FieldMatrix< Scalar, n, m > getTransposed(const Dune::FieldMatrix< Scalar, m, n > &M)
Transpose a FieldMatrix.
Definition: math.hh:672
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
Scalar interpolate(Scalar ip, Parameter &&params)
a generic function to interpolate given a set of parameters and an interpolation point
Definition: math.hh:504
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:660
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:618
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:599
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
make the local view function available whenever we use the grid geometry
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:518
static constexpr Scalar interpolate(Scalar ip, const std::array< Scalar, 2 > &params)
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)
interpolate linearly in a piecewise linear function (tabularized function)
Definition: math.hh:544