Line | Branch | Exec | Source |
---|---|---|---|
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 | // SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup Flux | ||
10 | * \brief A free function to average a Tensor at an interface. | ||
11 | */ | ||
12 | #ifndef DUMUX_FACE_TENSOR_AVERAGE_HH | ||
13 | #define DUMUX_FACE_TENSOR_AVERAGE_HH | ||
14 | |||
15 | #include <dune/common/fmatrix.hh> | ||
16 | #include <dumux/common/math.hh> | ||
17 | |||
18 | namespace Dumux { | ||
19 | |||
20 | /*! | ||
21 | * \brief Average of a discontinuous scalar field at discontinuity interface | ||
22 | * (for compatibility reasons with the function below) | ||
23 | * \return the harmonic average of the scalars | ||
24 | * \param T1 first scalar parameter | ||
25 | * \param T2 second scalar parameter | ||
26 | * \param normal The unit normal vector of the interface | ||
27 | */ | ||
28 | template <class Scalar, int dim> | ||
29 |
8/12✓ Branch 0 taken 507995251 times.
✓ Branch 1 taken 88996997 times.
✓ Branch 2 taken 280074103 times.
✓ Branch 3 taken 92608639 times.
✓ Branch 4 taken 271690584 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2146886 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 4239872 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 10976 times.
✗ Branch 11 not taken.
|
1247763308 | Scalar faceTensorAverage(const Scalar T1, |
30 | const Scalar T2, | ||
31 | const Dune::FieldVector<Scalar, dim>& normal) | ||
32 |
6/6✓ Branch 0 taken 1528946 times.
✓ Branch 1 taken 619999756 times.
✓ Branch 2 taken 1513314 times.
✓ Branch 3 taken 4890838 times.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 47 times.
|
801908862 | { return Dumux::harmonicMean(T1, T2); } |
33 | |||
34 | /*! | ||
35 | * \brief Average of a discontinuous tensorial field at discontinuity interface | ||
36 | * \note We do a harmonic average of the part normal to the interface (alpha*I) and | ||
37 | * an arithmetic average of the tangential part (T - alpha*I). | ||
38 | * \return the averaged tensor | ||
39 | * \param T1 first tensor | ||
40 | * \param T2 second tensor | ||
41 | * \param normal The unit normal vector of the interface | ||
42 | */ | ||
43 | template <class Scalar, int dim> | ||
44 | 11490912 | Dune::FieldMatrix<Scalar, dim> faceTensorAverage(const Dune::FieldMatrix<Scalar, dim>& T1, | |
45 | const Dune::FieldMatrix<Scalar, dim>& T2, | ||
46 | const Dune::FieldVector<Scalar, dim>& normal) | ||
47 | { | ||
48 | // determine nT*k*n | ||
49 | 11490912 | Dune::FieldVector<Scalar, dim> tmp; | |
50 | 11490912 | Dune::FieldVector<Scalar, dim> tmp2; | |
51 |
2/2✓ Branch 0 taken 14523024 times.
✓ Branch 1 taken 11490912 times.
|
26013936 | T1.mv(normal, tmp); |
52 |
2/2✓ Branch 0 taken 14523024 times.
✓ Branch 1 taken 11490912 times.
|
26013936 | T2.mv(normal, tmp2); |
53 |
2/2✓ Branch 0 taken 14523024 times.
✓ Branch 1 taken 3032112 times.
|
17555136 | const Scalar alpha1 = tmp*normal; |
54 |
1/2✓ Branch 0 taken 11490912 times.
✗ Branch 1 not taken.
|
11490912 | const Scalar alpha2 = tmp2*normal; |
55 | |||
56 | 11490912 | const Scalar alphaHarmonic = Dumux::harmonicMean(alpha1, alpha2); | |
57 | 11490912 | const Scalar alphaAverage = 0.5*(alpha1 + alpha2); | |
58 | |||
59 | 11490912 | Dune::FieldMatrix<Scalar, dim> T(0.0); | |
60 |
2/2✓ Branch 0 taken 14523024 times.
✓ Branch 1 taken 11490912 times.
|
26013936 | for (int i = 0; i < dim; ++i) |
61 | { | ||
62 |
2/2✓ Branch 0 taken 20587248 times.
✓ Branch 1 taken 14523024 times.
|
35110272 | for (int j = 0; j < dim; ++j) |
63 | { | ||
64 |
2/2✓ Branch 0 taken 6064224 times.
✓ Branch 1 taken 6064224 times.
|
20587248 | T[i][j] += 0.5*(T1[i][j] + T2[i][j]); |
65 |
2/2✓ Branch 0 taken 6064224 times.
✓ Branch 1 taken 6064224 times.
|
12128448 | if (i == j) |
66 | 14523024 | T[i][j] += alphaHarmonic - alphaAverage; | |
67 | } | ||
68 | } | ||
69 | |||
70 | 11490912 | return T; | |
71 | } | ||
72 | |||
73 | } // end namespace Dumux | ||
74 | |||
75 | #endif | ||
76 |