26#ifndef DUMUX_DISCRETIZATION_CC_MPFA_COMPUTE_TRANSMISSIBILITY_HH
27#define DUMUX_DISCRETIZATION_CC_MPFA_COMPUTE_TRANSMISSIBILITY_HH
29#include <dune/common/typetraits.hh>
30#include <dune/common/fvector.hh>
49template<
class EG,
class IVSubControlVolume,
class Tensor>
50[[deprecated(
"Will be removed after release 3.6. Use interface with additional fvGeometry parameter instead.")]]
51Dune::FieldVector<typename Tensor::field_type, IVSubControlVolume::myDimension>
53 const typename EG::SubControlVolumeFace& scvf,
55 typename IVSubControlVolume::ctype extrusionFactor)
57 Dune::FieldVector<typename Tensor::field_type, IVSubControlVolume::myDimension> wijk;
58 for (
unsigned int dir = 0; dir < IVSubControlVolume::myDimension; ++dir)
59 wijk[dir] =
vtmv(scvf.unitOuterNormal(), t, scv.nu(dir));
62 wijk *= Extrusion::area(scvf)*extrusionFactor;
81template<
class EG,
class IVSubControlVolume,
class Tensor>
82Dune::FieldVector<typename Tensor::field_type, IVSubControlVolume::myDimension>
84 const IVSubControlVolume& scv,
85 const typename EG::SubControlVolumeFace& scvf,
87 typename IVSubControlVolume::ctype extrusionFactor)
89 Dune::FieldVector<typename Tensor::field_type, IVSubControlVolume::myDimension> wijk;
90 for (
unsigned int dir = 0; dir < IVSubControlVolume::myDimension; ++dir)
91 wijk[dir] =
vtmv(scvf.unitOuterNormal(), t, scv.nu(dir));
94 wijk *= Extrusion::area(fvGeometry, scvf)*extrusionFactor;
112template< class EG, class IVSubControlVolume, class Tensor, std::enable_if_t< Dune::IsNumber<Tensor>::value,
int > = 1 >
113[[deprecated(
"Will be removed after release 3.6. Use interface with additional fvGeometry parameter instead.")]]
114Dune::FieldVector<Tensor, IVSubControlVolume::myDimension>
116 const typename EG::SubControlVolumeFace& scvf,
118 typename IVSubControlVolume::ctype extrusionFactor)
120 Dune::FieldVector<Tensor, IVSubControlVolume::myDimension> wijk;
121 for (
unsigned int dir = 0; dir < IVSubControlVolume::myDimension; ++dir)
122 wijk[dir] =
vtmv(scvf.unitOuterNormal(), t, scv.nu(dir));
125 wijk *= Extrusion::area(scvf)*extrusionFactor;
144template< class EG, class IVSubControlVolume, class Tensor, std::enable_if_t< Dune::IsNumber<Tensor>::value,
int > = 1 >
145Dune::FieldVector<Tensor, IVSubControlVolume::myDimension>
147 const IVSubControlVolume& scv,
148 const typename EG::SubControlVolumeFace& scvf,
150 typename IVSubControlVolume::ctype extrusionFactor)
152 Dune::FieldVector<Tensor, IVSubControlVolume::myDimension> wijk;
153 for (
unsigned int dir = 0; dir < IVSubControlVolume::myDimension; ++dir)
154 wijk[dir] =
vtmv(scvf.unitOuterNormal(), t, scv.nu(dir));
157 wijk *= Extrusion::area(fvGeometry, scvf)*extrusionFactor;
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:52
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
Traits extracting the public Extrusion type from T Defaults to NoExtrusion if no such type is found.
Definition: extrusion.hh:240