3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
mpfa/computetransmissibility.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 *****************************************************************************/
26#ifndef DUMUX_DISCRETIZATION_CC_MPFA_COMPUTE_TRANSMISSIBILITY_HH
27#define DUMUX_DISCRETIZATION_CC_MPFA_COMPUTE_TRANSMISSIBILITY_HH
28
29#include <dune/common/typetraits.hh>
30#include <dune/common/fvector.hh>
31
32#include <dumux/common/math.hh>
34
35namespace Dumux {
36
49template<class EG, class IVSubControlVolume, class Tensor>
50Dune::FieldVector<typename Tensor::field_type, IVSubControlVolume::myDimension>
51computeMpfaTransmissibility(const IVSubControlVolume& scv,
52 const typename EG::SubControlVolumeFace& scvf,
53 const Tensor& t,
54 typename IVSubControlVolume::ctype extrusionFactor)
55{
56 Dune::FieldVector<typename Tensor::field_type, IVSubControlVolume::myDimension> wijk;
57 for (unsigned int dir = 0; dir < IVSubControlVolume::myDimension; ++dir)
58 wijk[dir] = vtmv(scvf.unitOuterNormal(), t, scv.nu(dir));
59
61 wijk *= Extrusion::area(scvf)*extrusionFactor;
62 wijk /= scv.detX();
63
64 return wijk;
65}
66
79template< class EG, class IVSubControlVolume, class Tensor, std::enable_if_t< Dune::IsNumber<Tensor>::value, int > = 1 >
80Dune::FieldVector<Tensor, IVSubControlVolume::myDimension>
81computeMpfaTransmissibility(const IVSubControlVolume& scv,
82 const typename EG::SubControlVolumeFace& scvf,
83 const Tensor& t,
84 typename IVSubControlVolume::ctype extrusionFactor)
85{
86 Dune::FieldVector<Tensor, IVSubControlVolume::myDimension> wijk;
87 for (unsigned int dir = 0; dir < IVSubControlVolume::myDimension; ++dir)
88 wijk[dir] = vtmv(scvf.unitOuterNormal(), t, scv.nu(dir));
89
91 wijk *= Extrusion::area(scvf)*extrusionFactor;
92 wijk /= scv.detX();
93
94 return wijk;
95}
96
97} // end namespace Dumux
98
99#endif
Define some often used mathematical functions.
Helper classes to compute the integration elements.
Dune::FieldVector< typename Tensor::field_type, IVSubControlVolume::myDimension > computeMpfaTransmissibility(const IVSubControlVolume &scv, const typename EG::SubControlVolumeFace &scvf, const Tensor &t, typename IVSubControlVolume::ctype extrusionFactor)
Free function to evaluate the Mpfa transmissibility associated with the flux (in the form of flux = t...
Definition: mpfa/computetransmissibility.hh:51
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:863
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
Traits extracting the public Extrusion type from T Defaults to NoExtrusion if no such type is found.
Definition: extrusion.hh:166