version 3.10-dev
facetensoraverage.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
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
18namespace Dumux {
19
28template <class Scalar, int dim>
29Scalar faceTensorAverage(const Scalar T1,
30 const Scalar T2,
31 const Dune::FieldVector<Scalar, dim>& normal)
32{ return Dumux::harmonicMean(T1, T2); }
33
43template <class Scalar, int dim>
44Dune::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 Dune::FieldVector<Scalar, dim> tmp;
50 Dune::FieldVector<Scalar, dim> tmp2;
51 T1.mv(normal, tmp);
52 T2.mv(normal, tmp2);
53 const Scalar alpha1 = tmp*normal;
54 const Scalar alpha2 = tmp2*normal;
55
56 const Scalar alphaHarmonic = Dumux::harmonicMean(alpha1, alpha2);
57 const Scalar alphaAverage = 0.5*(alpha1 + alpha2);
58
59 Dune::FieldMatrix<Scalar, dim> T(0.0);
60 for (int i = 0; i < dim; ++i)
61 {
62 for (int j = 0; j < dim; ++j)
63 {
64 T[i][j] += 0.5*(T1[i][j] + T2[i][j]);
65 if (i == j)
66 T[i][j] += alphaHarmonic - alphaAverage;
67 }
68 }
69
70 return T;
71}
72
73} // end namespace Dumux
74
75#endif
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:57
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:26
Define some often used mathematical functions.
Definition: adapt.hh:17
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:29