13#ifndef DUMUX_DISCRETIZATION_BOX_DISPERSION_FLUX_HH
14#define DUMUX_DISCRETIZATION_BOX_DISPERSION_FLUX_HH
16#include <dune/common/fvector.hh>
17#include <dune/common/fmatrix.hh>
29template<
class TypeTag,
class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
36template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
44 using FVElementGeometry =
typename GridGeometry::LocalView;
45 using SubControlVolume =
typename GridGeometry::SubControlVolume;
46 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
50 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
51 using FluxVarCache =
typename GridFluxVariablesCache::FluxVariablesCache;
56 using Element =
typename GridView::template Codim<0>::Entity;
58 using Indices =
typename ModelTraits::Indices;
60 enum { dim = GridView::dimension} ;
61 enum { dimWorld = GridView::dimensionworld} ;
64 numPhases = ModelTraits::numFluidPhases(),
65 numComponents = ModelTraits::numFluidComponents()
68 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
69 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
70 using HeatFluxScalar = Scalar;
78 {
return referenceSystem; }
87 const Element& element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const SubControlVolumeFace& scvf,
92 const ElementFluxVariablesCache& elemFluxVarsCache)
94 ComponentFluxVector componentFlux(0.0);
96 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
97 const auto& shapeValues = fluxVarsCache.shapeValues();
100 Scalar rhoMassOrMole(0.0);
101 for (
auto&& scv : scvs(fvGeometry))
104 rhoMassOrMole += rho * shapeValues[scv.indexInElement()][0];
107 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
110 const auto& dispersionTensor =
111 ModelTraits::CompositionalDispersionModel::compositionalDispersionTensor(problem, scvf, fvGeometry,
112 elemVolVars, elemFluxVarsCache,
116 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
117 for (
auto&& scv : scvs(fvGeometry))
119 const auto x =
massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx);
120 gradX.axpy(x, fluxVarsCache.gradN(scv.indexInElement()));
124 componentFlux[compIdx] = -1.0 * rhoMassOrMole *
vtmv(scvf.unitOuterNormal(), dispersionTensor, gradX) * scvf.area();
125 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
126 componentFlux[phaseIdx] -= componentFlux[compIdx];
128 return componentFlux;
136 const Element& element,
137 const FVElementGeometry& fvGeometry,
138 const ElementVolumeVariables& elemVolVars,
139 const SubControlVolumeFace& scvf,
141 const ElementFluxVariablesCache& elemFluxVarsCache)
144 const auto& dispersionTensor =
145 ModelTraits::ThermalDispersionModel::thermalDispersionTensor(problem, scvf, fvGeometry,
146 elemVolVars, elemFluxVarsCache,
149 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
150 Dune::FieldVector<Scalar, GridView::dimensionworld> gradTemp(0.0);
151 for (
auto&& scv : scvs(fvGeometry))
152 gradTemp.axpy(elemVolVars[scv].
temperature(), fluxVarsCache.gradN(scv.indexInElement()));
155 return -1.0*
vtmv(scvf.unitOuterNormal(), dispersionTensor, gradTemp)*Extrusion::area(fvGeometry, scvf);
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: box/dispersionflux.hh:77
static ComponentFluxVector compositionalDispersionFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the dispersive fluxes of all components within a fluid phase across the given sub-control vol...
Definition: box/dispersionflux.hh:86
static HeatFluxScalar thermalDispersionFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the thermal dispersive flux across the given sub-control volume face.
Definition: box/dispersionflux.hh:135
Definition: box/dispersionflux.hh:30
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
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:851
VolumeVariables::PrimaryVariables::value_type massOrMoleFraction(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx, const int compIdx)
returns the mass or mole fraction to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:54
VolumeVariables::PrimaryVariables::value_type massOrMolarDensity(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx)
evaluates the density to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:43
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:33
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Define some often used mathematical functions.
The available discretization methods in Dumux.
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:39
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Traits of a flux variables type.
Definition: flux/traits.hh:32
static constexpr bool hasStationaryVelocityField()
Definition: flux/traits.hh:33