math.hh 29.5 KB
Newer Older
1
2
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// 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
12
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
13
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
14
15
16
17
 *   GNU General Public License for more details.                            *
 *                                                                           *
 *   You should have received a copy of the GNU General Public License       *
 *   along with this program.  If not, see <http://www.gnu.org/licenses/>.   *
18
19
20
 *****************************************************************************/
/*!
 * \file
21
 * \ingroup Common
22
23
24
25
26
 * \brief Define some often used mathematical functions
 */
#ifndef DUMUX_MATH_HH
#define DUMUX_MATH_HH

27
#include <algorithm>
28
#include <numeric>
29
#include <cmath>
30
#include <utility>
31

Dennis Gläser's avatar
Dennis Gläser committed
32
#include <dune/common/typetraits.hh>
Andreas Lauser's avatar
Andreas Lauser committed
33
34
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
35
#include <dune/common/dynmatrix.hh>
36
#include <dune/common/float_cmp.hh>
37

38
39
namespace Dumux {

40
/*!
41
 * \ingroup Common
42
 * \brief Calculate the (weighted) arithmetic mean of two scalar values.
43
44
45
 *
 * \param x The first input value
 * \param y The second input value
46
47
 * \param wx The first weight
 * \param wy The second weight
48
49
 */
template <class Scalar>
50
constexpr Scalar arithmeticMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
51
{
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
    static_assert(Dune::IsNumber<Scalar>::value, "The arguments x, y, wx, and wy have to be numbers!");

    if (x*y <= 0)
        return 0;
    return (x * wx + y * wy)/(wx + wy);
}

/*!
 * \ingroup Common
 * \brief Calculate the (weighted) harmonic mean of two scalar values.
 *
 * \param x The first input value
 * \param y The second input value
 * \param wx The first weight
 * \param wy The second weight
 */
template <class Scalar>
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx = 1.0, Scalar wy = 1.0) noexcept
{
    static_assert(Dune::IsNumber<Scalar>::value, "The arguments x, y, wx, and wy have to be numbers!");
72

73
74
    if (x*y <= 0)
        return 0;
75
    return (wx + wy) * x * y / (wy * x + wx * y);
76
77
78
}

/*!
79
 * \ingroup Common
80
 * \brief Calculate the geometric mean of two scalar values.
81
82
83
 *
 * \param x The first input value
 * \param y The second input value
84
 * \note as std::sqrt is not constexpr this function is not constexpr
85
86
 */
template <class Scalar>
87
Scalar geometricMean(Scalar x, Scalar y) noexcept
88
{
89
90
    static_assert(Dune::IsNumber<Scalar>::value, "The arguments x and y have to be numbers!");

91
92
    if (x*y <= 0)
        return 0;
93
94
    using std::sqrt;
    return sqrt(x*y)*sign(x);
95
96
97
}

/*!
98
 * \ingroup Common
99
100
101
102
 * \brief Calculate the harmonic mean of a fixed-size matrix.
 *
 * This is done by calculating the harmonic mean for each entry
 * individually. The result is stored the first argument.
103
104
105
106
 *
 * \param K The matrix for storing the result
 * \param Ki The first input matrix
 * \param Kj The second input matrix
107
108
109
110
111
112
113
114
 */
template <class Scalar, int m, int n>
void harmonicMeanMatrix(Dune::FieldMatrix<Scalar, m, n> &K,
                        const Dune::FieldMatrix<Scalar, m, n> &Ki,
                        const Dune::FieldMatrix<Scalar, m, n> &Kj)
{
    for (int rowIdx=0; rowIdx < m; rowIdx++){
        for (int colIdx=0; colIdx< n; colIdx++){
115
            if (Dune::FloatCmp::ne<Scalar>(Ki[rowIdx][colIdx], Kj[rowIdx][colIdx])) {
116
117
118
119
120
121
122
123
124
125
                K[rowIdx][colIdx] =
                    harmonicMean(Ki[rowIdx][colIdx],
                                 Kj[rowIdx][colIdx]);
            }
            else
                K[rowIdx][colIdx] = Ki[rowIdx][colIdx];
        }
    }
}

126
/*!
127
 * \ingroup Common
128
129
 * \brief Invert a linear polynomial analytically
 *
130
131
132
133
134
135
136
137
138
 * The polynomial is defined as
 * \f[ p(x) = a\; x + b \f]
 *
 * This method Returns the number of solutions which are in the real
 * numbers, i.e. 1 except if the slope of the line is 0.
 *
 * \param sol Container into which the solutions are written
 * \param a The coefficiont for the linear term
 * \param b The coefficiont for the constant term
139
140
 */
template <class Scalar, class SolContainer>
141
int invertLinearPolynomial(SolContainer &sol,
142
143
144
145
146
147
148
149
150
151
152
                           Scalar a,
                           Scalar b)
{
    if (a == 0.0)
        return 0;

    sol[0] = -b/a;
    return 1;
}

/*!
153
 * \ingroup Common
154
155
 * \brief Invert a quadratic polynomial analytically
 *
156
157
158
 * The polynomial is defined as
 * \f[ p(x) = a\; x^2 + + b\;x + c \f]
 *
159
 * This method returns the number of solutions which are in the real
160
161
162
163
164
165
166
 * numbers. The "sol" argument contains the real roots of the parabola
 * in order with the smallest root first.
 *
 * \param sol Container into which the solutions are written
 * \param a The coefficiont for the quadratic term
 * \param b The coefficiont for the linear term
 * \param c The coefficiont for the constant term
167
168
 */
template <class Scalar, class SolContainer>
169
int invertQuadraticPolynomial(SolContainer &sol,
170
171
172
173
174
175
176
                              Scalar a,
                              Scalar b,
                              Scalar c)
{
    // check for a line
    if (a == 0.0)
        return invertLinearPolynomial(sol, b, c);
177

178
179
180
181
182
    // discriminant
    Scalar Delta = b*b - 4*a*c;
    if (Delta < 0)
        return 0; // no real roots

183
184
    using std::sqrt;
    Delta = sqrt(Delta);
185
186
    sol[0] = (- b + Delta)/(2*a);
    sol[1] = (- b - Delta)/(2*a);
187
188
189

    // sort the result
    if (sol[0] > sol[1])
190
191
192
193
    {
        using std::swap;
        swap(sol[0], sol[1]);
    }
194
195
196
    return 2; // two real roots
}

Thomas Fetzer's avatar
Thomas Fetzer committed
197
//! \cond false
198
199
200
template <class Scalar, class SolContainer>
void invertCubicPolynomialPostProcess_(SolContainer &sol,
                                       int numSol,
201
202
                                       Scalar a, Scalar b, Scalar c, Scalar d,
                                       std::size_t numIterations = 1)
203
{
204
205
    const auto eval = [&](auto x){ return d + x*(c + x*(b + x*a)); };
    const auto evalDeriv = [&](auto x){ return c + x*(2*b + x*3*a); };
206

207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
    // do numIterations Newton iterations on the analytic solution
    // and update result if the precision is increased
    for (int i = 0; i < numSol; ++i)
    {
        Scalar x = sol[i];
        Scalar fCurrent = eval(x);
        const Scalar fOld = fCurrent;
        for (int j = 0; j < numIterations; ++j)
        {
            const Scalar fPrime = evalDeriv(x);
            if (fPrime == 0.0)
                break;
            x -= fCurrent/fPrime;
            fCurrent = eval(x);
        }
222

223
        using std::abs;
224
        if (abs(fCurrent) < abs(fOld))
225
226
227
228
229
            sol[i] = x;
    }
}
//! \endcond

230
/*!
231
 * \ingroup Common
232
233
 * \brief Invert a cubic polynomial analytically
 *
234
235
236
 * The polynomial is defined as
 * \f[ p(x) = a\; x^3 + + b\;x^3 + c\;x + d \f]
 *
237
 * This method returns the number of solutions which are in the real
238
239
240
 * numbers. The "sol" argument contains the real roots of the cubic
 * polynomial in order with the smallest root first.
 *
241
242
243
 * \note The closer the roots are to each other the less
 *       precise the inversion becomes. Increase number of post-processing iterations for improved results.
 *
244
 * \param sol Container into which the solutions are written
245
246
247
248
249
 * \param a The coefficient for the cubic term
 * \param b The coefficient for the quadratic term
 * \param c The coefficient for the linear term
 * \param d The coefficient for the constant term
 * \param numPostProcessIterations The number of iterations to increase precision of the analytical result
250
251
 */
template <class Scalar, class SolContainer>
252
int invertCubicPolynomial(SolContainer *sol,
253
254
                          Scalar a, Scalar b, Scalar c, Scalar d,
                          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
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
    // normalize the polynomial
    b /= a;
    c /= a;
    d /= a;
    a = 1;

    // get rid of the quadratic term by subsituting x = t - b/3
    Scalar p = c - b*b/3;
    Scalar q = d + (2*b*b*b - 9*b*c)/27;

    // now we are at the form t^3 + p*t + q = 0. First we handle some
    // special cases to avoid divisions by zero later...
    if (p == 0.0 && q == 0.0) {
        // t^3 = 0, i.e. triple root at t = 0
        sol[0] = sol[1] = sol[2] = 0.0 - b/3;
        return 3;
    }
    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
282
        using std::cbrt;
        Scalar t = cbrt(q);
283
        sol[0] = t - b/3;
284

285
286
287
        return 1;
    }
    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
292
        if (p > 0) {
            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
298
        using std::sqrt;
        sol[0] = -sqrt(-p) - b/3;
299
        sol[1] = 0.0 - b/3;
300
        sol[2] = sqrt(-p) - b/3;
301

302
303
304
305
306
307
308
309
        return 3;
    }

    // At this point
    //
    // t^3 + p*t + q = 0
    //
    // 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
317
318
319
320
321
322
323
324
    //
    // u^3 + v^3 + q = 0 .
    //
    // multiplying both sides with u^3 and taking advantage of the
    // fact that u*v = -p/3 leads to
    //
    // u^6 + q*u^3 - p^3/27 = 0
    //
    // Now, substituting u^3 = w yields
    //
325
    // w^2 + q*w - p^3/27 = 0
326
    //
327
    // This is a quadratic equation with the solutions
328
    //
329
330
    // w = -q/2 + sqrt(q^2/4 + p^3/27) and
    // w = -q/2 - sqrt(q^2/4 + p^3/27)
331
    //
332
333
334
    // Since w is equivalent to u^3 it is sufficient to only look at
    // one of the two cases. Then, there are still 2 cases: positive
    // and negative discriminant.
335
    Scalar wDisc = q*q/4 + p*p*p/27;
336
337
    if (wDisc >= 0) { // the positive discriminant case:
        // calculate the cube root of - q/2 + sqrt(q^2/4 + p^3/27)
338
339
        using std::cbrt;
        using std::sqrt;
340

341
342
343
        // Choose the root that is safe against the loss of precision
        // that can cause wDisc to be equal to q*q/4 despite p != 0.
        // We do not want u to be zero in that case. Mathematically,
344
        // we are happy with either root.
345
346
347
        const Scalar u = [&]{
            return q > 0 ? cbrt(-0.5*q - sqrt(wDisc)) : cbrt(-0.5*q + sqrt(wDisc));
        }();
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
352
353
        // does not produce a division by zero. the remaining two
        // roots of u are rotated by +- 2/3*pi in the complex plane
        // and thus not considered here
354
        invertCubicPolynomialPostProcess_(sol, 1, a, b, c, d, numPostProcessIterations);
355
356
        return 1;
    }
357
358
    else { // the negative discriminant case:
        Scalar uCubedRe = - q/2;
359
360
        using std::sqrt;
        Scalar uCubedIm = sqrt(-wDisc);
361
        // calculate the cube root of - q/2 + sqrt(q^2/4 + p^3/27)
362
363
364
365
        using std::cbrt;
        Scalar uAbs  = cbrt(sqrt(uCubedRe*uCubedRe + uCubedIm*uCubedIm));
        using std::atan2;
        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
372
        //
        // x = u - p/(3*u) - b/3
        //
373
        // where x and u are complex numbers. Rewritten in polar form
374
375
        // this is equivalent to
        //
376
        // x = |u|*e^(i*phi) - p*e^(-i*phi)/(3*|u|) - b/3 .
377
        //
378
379
        // Factoring out the e^ terms and subtracting the additional
        // 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
388
389
        // The crucial observation is the fact that y is the conjugate
        // of - x + b/3. This means that after taking advantage of the
        // 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
398
399
        // holds. Since |u|, p, b and cos(phi) are real numbers, it
        // follows that Im(x) = - Im(x) and thus Im(x) = 0. This
        // implies
400
        //
401
        // Re(x) = x = cos(phi)*(|u| - p / (3*|u|)) - b/3 .
402
403
404
405
406
        //
        // Considering the fact that u is a cubic root, we have three
        // values for phi which differ by 2/3*pi. This allows to
        // calculate the three real roots of the polynomial:
        for (int i = 0; i < 3; ++i) {
407
408
            using std::cos;
            sol[i] = cos(phi)*(uAbs - p/(3*uAbs)) - b/3;
409
410
            phi += 2*M_PI/3;
        }
411

412
413
        // post process the obtained solution to increase numerical
        // precision
414
        invertCubicPolynomialPostProcess_(sol, 3, a, b, c, d, numPostProcessIterations);
415

416
        // sort the result
417
418
        using std::sort;
        sort(sol, sol + 3);
419

420
421
422
423
424
425
        return 3;
    }

    // NOT REACHABLE!
    return 0;
}
426

427
/*!
428
 * \ingroup Common
429
430
431
432
433
434
435
 * \brief Comparison of two position vectors
 *
 * Compares an current position vector with a reference vector, and returns true
 * if the position vector is larger.
 * "Larger" in this case means that all the entries of each spacial dimension are
 * larger compared to the reference vector.
 *
436
437
 * \param pos Vektor holding the current Position that is to be checked
 * \param smallerVec Reference vector, holding the minimum values for comparison.
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
 */
template <class Scalar, int dim>
bool isLarger(const Dune::FieldVector<Scalar, dim> &pos,
              const Dune::FieldVector<Scalar, dim> &smallerVec)
{
    for (int i=0; i < dim; i++)
    {
        if (pos[i]<= smallerVec[i])
        {
            return false;
        }
    }
    return true;
}

/*!
454
 * \ingroup Common
455
456
457
458
459
460
461
 * \brief Comparison of two position vectors
 *
 * Compares an current position vector with a reference vector, and returns true
 * if the position vector is smaller.
 * "Smaller" in this case means that all the entries of each spacial dimension are
 * smaller in comparison with the reference vector.
 *
462
463
 * \param pos Vektor holding the current Position that is to be checked
 * \param largerVec Reference vector, holding the maximum values for comparison.
464
465
466
 */
template <class Scalar, int dim>
bool isSmaller(const Dune::FieldVector<Scalar, dim> &pos,
Andreas Lauser's avatar
Andreas Lauser committed
467
               const Dune::FieldVector<Scalar, dim> &largerVec)
468
469
470
471
472
473
474
475
476
477
478
479
{
    for (int i=0; i < dim; i++)
    {
        if (pos[i]>= largerVec[i])
        {
            return false;
        }
    }
    return true;
}

/*!
480
 * \ingroup Common
481
482
483
484
485
486
487
488
489
490
 * \brief Comparison of three position vectors
 *
 * Compares an current position vector with two reference vector, and returns true
 * if the position vector lies in between them.
 * "Between" in this case means that all the entries of each spacial dimension are
 * smaller in comparison with the larger reference vector as well as larger campared
 * to the smaller reference.
 * This is comfortable to cheack weather the current position is located inside or
 * outside of a lense with different properties.
 *
491
492
493
 * \param pos Vektor holding the current Position that is to be checked
 * \param smallerVec Reference vector, holding the minimum values for comparison.
 * \param largerVec Reference vector, holding the maximum values for comparison.
494
495
496
497
498
499
500
501
502
503
504
505
506
507
 */
template <class Scalar, int dim>
bool isBetween(const Dune::FieldVector<Scalar, dim> &pos,
              const Dune::FieldVector<Scalar, dim> &smallerVec,
              const Dune::FieldVector<Scalar, dim> &largerVec)
{
   if (isLarger(pos, smallerVec) && isSmaller(pos, largerVec))
       {
           return true;
       }
   else
       return false;
}

508
509
510
511
512
513
514
515
516
//! forward declaration of the linear interpolation policy (default)
namespace InterpolationPolicy { struct Linear; }

/*!
 * \ingroup Common
 * \brief a generic function to interpolate given a set of parameters and an interpolation point
 * \param params the parameters used for interpolation (depends on the policy used)
 * \param ip the interpolation point
 */
517
518
519
template <class Policy = InterpolationPolicy::Linear, class Scalar, class ... Parameter>
Scalar interpolate(Scalar ip, Parameter&& ... params)
{ return Policy::interpolate(ip, std::forward<Parameter>(params) ...); }
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535

/*!
 * \ingroup Common
 * \brief Interpolation policies
 */
namespace InterpolationPolicy  {

/*!
 * \ingroup Common
 * \brief interpolate linearly between two given values
 */
struct Linear
{
    /*!
     * \brief interpolate linearly between two given values
     * \param ip the interpolation point in [0,1]
536
     * \param params array with the lower and upper bound
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
     */
    template<class Scalar>
    static constexpr Scalar interpolate(Scalar ip, const std::array<Scalar, 2>& params)
    {
        return params[0]*(1.0 - ip) + params[1]*ip;
    }
};

/*!
 * \ingroup Common
 * \brief interpolate linearly in a piecewise linear function (tabularized function)
 */
struct LinearTable
{
    /*!
     * \brief interpolate linearly in a piecewise linear function (tabularized function)
     * \param ip the interpolation point
554
555
     * \param range positions of values
     * \param values values to interpolate from
556
557
     * \note if the interpolation point is out of bounds this will return the bounds
     */
558
559
    template<class Scalar, class RandomAccessContainer0, class RandomAccessContainer1>
    static constexpr Scalar interpolate(Scalar ip, const RandomAccessContainer0& range, const RandomAccessContainer1& values)
560
561
562
563
564
565
566
567
568
569
570
571
572
    {
        // check bounds
        if (ip > range.back()) return values.back();
        if (ip < range[0]) return values[0];

        // if we are within bounds find the index of the lower bound
        const auto lookUpIndex = std::distance(range.begin(), std::lower_bound(range.begin(), range.end(), ip));
        if (lookUpIndex == 0)
            return values[0];

        const auto ipLinear = (ip - range[lookUpIndex-1])/(range[lookUpIndex] - range[lookUpIndex-1]);
        return Dumux::interpolate<Linear>(ipLinear, std::array<Scalar, 2>{{values[lookUpIndex-1], values[lookUpIndex]}});
    }
573
574
575
576
577
578
579

    template<class Scalar, class RandomAccessContainer>
    static constexpr Scalar interpolate(Scalar ip, const std::pair<RandomAccessContainer, RandomAccessContainer>& table)
    {
        const auto& [range, values] = table;
        return interpolate(ip, range, values);
    }
580
581
582
583
};

} // end namespace InterpolationPolicy

584
585
586
587
588
589
590
/*!
 * \ingroup Common
 * \brief Generates linearly spaced vectors
 *
 * \param begin The first value in the vector
 * \param end The last value in the vector
 * \param samples The size of the vector
591
 * \param endPoint if the range is including the interval's end point or not
592
593
 */
template <class Scalar>
594
595
596
std::vector<Scalar> linspace(const Scalar begin, const Scalar end,
                             std::size_t samples,
                             bool endPoint = true)
597
598
599
{
    using std::max;
    samples = max(std::size_t{2}, samples); // only makes sense for 2 or more samples
600
601
    const Scalar divisor = endPoint ? samples-1 : samples;
    const Scalar delta = (end-begin)/divisor;
602
603
604
605
606
607
    std::vector<Scalar> vec(samples);
    for (std::size_t i = 0; i < samples; ++i)
        vec[i] = begin + i*delta;
    return vec;
}

608
609

/*!
610
 * \ingroup Common
611
612
613
614
 * \brief Evaluates the Antoine equation used to calculate the vapour
 *        pressure of various liquids.
 *
 * See http://en.wikipedia.org/wiki/Antoine_equation
615
616
617
618
619
 *
 * \param temperature The temperature [K] of the fluid
 * \param A The first coefficient for the Antoine equation
 * \param B The first coefficient for the Antoine equation
 * \param C The first coefficient for the Antoine equation
620
621
622
623
624
625
626
627
 */
template <class Scalar>
Scalar antoine(Scalar temperature,
               Scalar A,
               Scalar B,
               Scalar C)
{
    const Scalar ln10 = 2.3025850929940459;
628
629
630
631
632
    using std::exp;
    return exp(ln10*(A - B/(C + temperature)));
}

/*!
633
 * \ingroup Common
634
635
636
637
638
639
640
 * \brief Sign or signum function.
 *
 * Returns 1 for a positive argument.
 * Returns -1 for a negative argument.
 * Returns 0 if the argument is zero.
 */
template<class ValueType>
641
constexpr int sign(const ValueType& value) noexcept
642
643
{
    return (ValueType(0) < value) - (value < ValueType(0));
644
645
}

646
/*!
647
 * \ingroup Common
648
 * \brief Cross product of two vectors in three-dimensional Euclidean space
649
650
651
652
 *
 * \param vec1 The first vector
 * \param vec2 The second vector
 */
653
654
655
template <class Scalar>
Dune::FieldVector<Scalar, 3> crossProduct(const Dune::FieldVector<Scalar, 3> &vec1,
                                          const Dune::FieldVector<Scalar, 3> &vec2)
656
{
657
658
659
    return {vec1[1]*vec2[2]-vec1[2]*vec2[1],
            vec1[2]*vec2[0]-vec1[0]*vec2[2],
            vec1[0]*vec2[1]-vec1[1]*vec2[0]};
660
661
}

662
/*!
663
 * \ingroup Common
664
665
666
667
668
669
670
671
672
 * \brief Cross product of two vectors in two-dimensional Euclidean space retuning scalar
 *
 * \param vec1 The first vector
 * \param vec2 The second vector
 */
template <class Scalar>
Scalar crossProduct(const Dune::FieldVector<Scalar, 2> &vec1,
                    const Dune::FieldVector<Scalar, 2> &vec2)
{   return vec1[0]*vec2[1]-vec1[1]*vec2[0]; }
673

Dennis Gläser's avatar
Dennis Gläser committed
674
/*!
675
 * \ingroup Common
Dennis Gläser's avatar
Dennis Gläser committed
676
677
678
679
680
681
682
683
 * \brief Triple product of three vectors in three-dimensional Euclidean space retuning scalar
 *
 * \param vec1 The first vector
 * \param vec2 The second vector
 * \param vec3 The third vector
 */
template <class Scalar>
Scalar tripleProduct(const Dune::FieldVector<Scalar, 3> &vec1,
Dennis Gläser's avatar
Dennis Gläser committed
684
685
                     const Dune::FieldVector<Scalar, 3> &vec2,
                     const Dune::FieldVector<Scalar, 3> &vec3)
Dennis Gläser's avatar
Dennis Gläser committed
686
{   return crossProduct<Scalar>(vec1, vec2)*vec3; }
687

688
/*!
689
 * \ingroup Common
690
691
692
693
694
 * \brief Transpose a FieldMatrix
 *
 * \param M The matrix to be transposed
 */
template <class Scalar, int m, int n>
695
Dune::FieldMatrix<Scalar, n, m> getTransposed(const Dune::FieldMatrix<Scalar, m, n>& M)
696
{
697
698
699
    Dune::FieldMatrix<Scalar, n, m> T;
    for (std::size_t i = 0; i < m; ++i)
        for (std::size_t j = 0; j < n; ++j)
700
701
702
703
704
            T[j][i] = M[i][j];

    return T;
}

705
/*!
706
 * \ingroup Common
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
 * \brief Transpose a DynamicMatrix
 *
 * \param M The matrix to be transposed
 */
template <class Scalar>
Dune::DynamicMatrix<Scalar> getTransposed(const Dune::DynamicMatrix<Scalar>& M)
{
    std::size_t rows_T = M.M();
    std::size_t cols_T = M.N();

    Dune::DynamicMatrix<Scalar> M_T(rows_T, cols_T, 0.0);

    for (std::size_t i = 0; i < rows_T; ++i)
        for (std::size_t j = 0; j < cols_T; ++j)
            M_T[i][j] = M[j][i];

    return M_T;
}

726
/*!
727
 * \ingroup Common
728
729
730
731
732
733
734
735
736
 * \brief Multiply two dynamic matrices
 *
 * \param M1 The first dynamic matrix
 * \param M2 The second dynamic matrix (to be multiplied to M1 from the right side)
 */
template <class Scalar>
Dune::DynamicMatrix<Scalar> multiplyMatrices(const Dune::DynamicMatrix<Scalar> &M1,
                                             const Dune::DynamicMatrix<Scalar> &M2)
{
737
738
739
    using size_type = typename Dune::DynamicMatrix<Scalar>::size_type;
    const size_type rows = M1.N();
    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
745
746
    for (size_type i = 0; i < rows; i++)
        for (size_type j = 0; j < cols; j++)
            for (size_type k = 0; k < M1.M(); k++)
747
748
749
750
751
                result[i][j] += M1[i][k]*M2[k][j];

    return result;
}

752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
/*!
 * \ingroup Common
 * \brief Multiply two field matrices
 *
 * \param M1 The first field matrix
 * \param M2 The second field matrix (to be multiplied to M1 from the right side)
 */
template <class Scalar, int rows1, int cols1, int cols2>
Dune::FieldMatrix<Scalar, rows1, cols2> multiplyMatrices(const Dune::FieldMatrix<Scalar, rows1, cols1> &M1,
                                                         const Dune::FieldMatrix<Scalar, cols1, cols2> &M2)
{
    using size_type = typename Dune::FieldMatrix<Scalar, rows1, cols2>::size_type;

    Dune::FieldMatrix<Scalar, rows1, cols2> result(0.0);
    for (size_type i = 0; i < rows1; i++)
        for (size_type j = 0; j < cols2; j++)
            for (size_type k = 0; k < cols1; k++)
                result[i][j] += M1[i][k]*M2[k][j];

    return result;
}


775
/*!
776
 * \ingroup Common
777
 * \brief Trace of a dense matrix
778
 *
779
 * \param M The dense matrix
780
 */
781
782
783
template <class MatrixType>
typename Dune::DenseMatrix<MatrixType>::field_type
trace(const Dune::DenseMatrix<MatrixType>& M)
784
{
785
786
    const auto rows = M.N();
    DUNE_ASSERT_BOUNDS(rows == M.M()); // rows == cols
787

788
789
    using MatType = Dune::DenseMatrix<MatrixType>;
    typename MatType::field_type trace = 0.0;
790

791
    for (typename MatType::size_type i = 0; i < rows; ++i)
792
793
794
795
796
        trace += M[i][i];

    return trace;
}

Dennis Gläser's avatar
Dennis Gläser committed
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
/*!
 * \ingroup Common
 * \brief Returns the result of the projection of
 *        a vector v with a Matrix M.
 *
 *        Note: We use DenseVector and DenseMatrix here so that
 *        it can be used with the statically and dynamically
 *        allocated Dune Vectors/Matrices. Size mismatch
 *        assertions are done in the respective Dune classes.
 *
 * \param M The matrix
 * \param v The vector
 */
template <class MAT, class V>
typename Dune::DenseVector<V>::derived_type
mv(const Dune::DenseMatrix<MAT>& M,
   const Dune::DenseVector<V>& v)
{
    typename Dune::DenseVector<V>::derived_type res(v);
    M.mv(v, res);
    return res;
}

/*!
 * \ingroup Common
 * \brief Returns the result of a vector v multiplied by a scalar m.
 *
 *        Note: We use DenseVector and DenseMatrix here so that
 *        it can be used with the statically and dynamically
 *        allocated Dune Vectors/Matrices. Size mismatch
 *        assertions are done in the respective Dune classes.
 *
 * \note We need the enable_if to make sure that only Scalars
 *       fit here. Matrix types are forced to use the above
 *       mv(DenseMatrix, DenseVector) instead.
 *
 * \param m The scale factor
 * \param v The vector
 */

template <class FieldScalar, class V>
typename std::enable_if_t<Dune::IsNumber<FieldScalar>::value,
                          typename Dune::DenseVector<V>::derived_type>
mv(const FieldScalar m, const Dune::DenseVector<V>& v)
{
    typename Dune::DenseVector<V>::derived_type res(v);
    res *= m;
    return res;
}

847
/*!
848
 * \ingroup Common
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
 * \brief Evaluates the scalar product of a vector v2, projected by
 *        a matrix M, with a vector v1.
 *
 *        Note: We use DenseVector and DenseMatrix here so that
 *        it can be used with the statically and dynamically
 *        allocated Dune Vectors/Matrices. Size mismatch
 *        assertions are done in the respective Dune classes.
 *
 * \param v1 The first vector
 * \param M The matrix
 * \param v2 The second vector
 */
template <class V1, class MAT, class V2>
typename Dune::DenseMatrix<MAT>::value_type
vtmv(const Dune::DenseVector<V1>& v1,
     const Dune::DenseMatrix<MAT>& M,
     const Dune::DenseVector<V2>& v2)
{
Dennis Gläser's avatar
Dennis Gläser committed
867
    return v1*mv(M, v2);
868
869
870
}

/*!
871
 * \ingroup Common
872
873
874
875
876
877
878
879
 * \brief Evaluates the scalar product of a vector v2, scaled by
 *        a scalar m, with a vector v1.
 *
 *        Note: We use DenseVector and DenseMatrix here so that
 *        it can be used with the statically and dynamically
 *        allocated Dune Vectors/Matrices. Size mismatch
 *        assertions are done in the respective Dune classes.
 *
880
881
882
883
 * \note We need the enable_if to make sure that only Scalars
 *       fit here. Matrix types are forced to use the above
 *       vtmv(DenseVector, DenseMatrix, DenseVector) instead.
 *
884
885
886
887
 * \param v1 The first vector
 * \param m The scale factor
 * \param v2 The second vector
 */
888

889
template <class V1, class FieldScalar, class V2>
890
891
typename std::enable_if_t<Dune::IsNumber<FieldScalar>::value, FieldScalar>
vtmv(const Dune::DenseVector<V1>& v1,
Dennis Gläser's avatar
Dennis Gläser committed
892
893
     const FieldScalar m,
     const Dune::DenseVector<V2>& v2)
894
895
896
897
{
    return m*(v1*v2);
}

898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
/*!
 * \ingroup Common
 * \brief Returns the [intercept, slope] of the regression line
 *        fitted to a set of (x, y) data points.
 *
 *        Note: We use least-square regression method to find
 *        the regression line.
 *
 * \param x x-values of the data set
 * \param y y-values of the data set
 */
template <class Scalar>
std::array<Scalar, 2> linearRegression(const std::vector<Scalar>& x,
                                       const std::vector<Scalar>& y)
{
    if (x.size() != y.size())
        DUNE_THROW(Dune::InvalidStateException, "x and y array must have the same length.");

    const Scalar averageX = std::accumulate(x.begin(), x.end(), 0.0)/x.size();
    const Scalar averageY = std::accumulate(y.begin(), y.end(), 0.0)/y.size();

    // calculate temporary variables necessary for slope computation
    const Scalar numerator = std::inner_product(
        x.begin(), x.end(), y.begin(), 0.0, std::plus<Scalar>(),
        [&](auto xx, auto yy) { return (xx - averageX) * (yy - averageY); }
    );
    const Scalar denominator = std::inner_product(
        x.begin(), x.end(), x.begin(), 0.0, std::plus<Scalar>(),
        [&](auto xx, auto yy) { return (xx - averageX) * (yy - averageX); }
    );

    // compute slope and intercept of the regression line
    const Scalar slope = numerator / denominator;
    const Scalar intercept = averageY - slope * averageX;

    return {intercept, slope};
}

936
} // end namespace Dumux
937
938

#endif