3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_FACE_TENSOR_AVERAGE_HH
25#define DUMUX_FACE_TENSOR_AVERAGE_HH
26
27#include <dune/common/fmatrix.hh>
28#include <dumux/common/math.hh>
29
30namespace Dumux {
31
40template <class Scalar, int dim>
41Scalar faceTensorAverage(const Scalar T1,
42 const Scalar T2,
43 const Dune::FieldVector<Scalar, dim>& normal)
44{ return Dumux::harmonicMean(T1, T2); }
45
55template <class Scalar, int dim>
56Dune::FieldMatrix<Scalar, dim> faceTensorAverage(const Dune::FieldMatrix<Scalar, dim>& T1,
57 const Dune::FieldMatrix<Scalar, dim>& T2,
58 const Dune::FieldVector<Scalar, dim>& normal)
59{
60 // determine nT*k*n
61 Dune::FieldVector<Scalar, dim> tmp;
62 Dune::FieldVector<Scalar, dim> tmp2;
63 T1.mv(normal, tmp);
64 T2.mv(normal, tmp2);
65 const Scalar alpha1 = tmp*normal;
66 const Scalar alpha2 = tmp2*normal;
67
68 const Scalar alphaHarmonic = Dumux::harmonicMean(alpha1, alpha2);
69 const Scalar alphaAverage = 0.5*(alpha1 + alpha2);
70
71 Dune::FieldMatrix<Scalar, dim> T(0.0);
72 for (int i = 0; i < dim; ++i)
73 {
74 for (int j = 0; j < dim; ++j)
75 {
76 T[i][j] += 0.5*(T1[i][j] + T2[i][j]);
77 if (i == j)
78 T[i][j] += alphaHarmonic - alphaAverage;
79 }
80 }
81
82 return T;
83}
84
85} // end namespace Dumux
86
87#endif
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:38
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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