24#ifndef DUMUX_FACE_TENSOR_AVERAGE_HH
25#define DUMUX_FACE_TENSOR_AVERAGE_HH
27#include <dune/common/fmatrix.hh>
40template <
class Scalar,
int dim>
43 const Dune::FieldVector<Scalar, dim>&
normal)
55template <
class Scalar,
int dim>
57 const Dune::FieldMatrix<Scalar, dim>& T2,
58 const Dune::FieldVector<Scalar, dim>&
normal)
61 Dune::FieldVector<Scalar, dim> tmp;
62 Dune::FieldVector<Scalar, dim> tmp2;
65 const Scalar alpha1 = tmp*
normal;
66 const Scalar alpha2 = tmp2*
normal;
69 const Scalar alphaAverage = 0.5*(alpha1 + alpha2);
71 Dune::FieldMatrix<Scalar, dim> T(0.0);
72 for (
int i = 0; i < dim; ++i)
74 for (
int j = 0; j < dim; ++j)
76 T[i][j] += 0.5*(T1[i][j] + T2[i][j]);
78 T[i][j] += alphaHarmonic - alphaAverage;
Define some often used mathematical functions.
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:36
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
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:41