25#ifndef DUMUX_DISCRETIZATION_BOX_DISPERSION_FLUX_HH
26#define DUMUX_DISCRETIZATION_BOX_DISPERSION_FLUX_HH
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
41template<
class TypeTag,
class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
48template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
56 using FVElementGeometry =
typename GridGeometry::LocalView;
57 using SubControlVolume =
typename GridGeometry::SubControlVolume;
58 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
62 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
63 using FluxVarCache =
typename GridFluxVariablesCache::FluxVariablesCache;
68 using Element =
typename GridView::template Codim<0>::Entity;
70 using Indices =
typename ModelTraits::Indices;
72 enum { dim = GridView::dimension} ;
73 enum { dimWorld = GridView::dimensionworld} ;
76 numPhases = ModelTraits::numFluidPhases(),
77 numComponents = ModelTraits::numFluidComponents()
80 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
81 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
82 using HeatFluxScalar = Scalar;
90 {
return referenceSystem; }
99 const Element& element,
100 const FVElementGeometry& fvGeometry,
101 const ElementVolumeVariables& elemVolVars,
102 const SubControlVolumeFace& scvf,
104 const ElementFluxVariablesCache& elemFluxVarsCache)
106 ComponentFluxVector componentFlux(0.0);
108 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
109 const auto& shapeValues = fluxVarsCache.shapeValues();
112 Scalar rhoMassOrMole(0.0);
113 for (
auto&& scv : scvs(fvGeometry))
116 rhoMassOrMole += rho * shapeValues[scv.indexInElement()][0];
119 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
122 const auto& dispersionTensor =
123 ModelTraits::CompositionalDispersionModel::compositionalDispersionTensor(problem, scvf, fvGeometry,
124 elemVolVars, elemFluxVarsCache,
128 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
129 for (
auto&& scv : scvs(fvGeometry))
131 const auto x =
massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx);
132 gradX.axpy(x, fluxVarsCache.gradN(scv.indexInElement()));
136 componentFlux[compIdx] = -1.0 * rhoMassOrMole *
vtmv(scvf.unitOuterNormal(), dispersionTensor, gradX) * scvf.area();
137 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
138 componentFlux[phaseIdx] -= componentFlux[compIdx];
140 return componentFlux;
148 const Element& element,
149 const FVElementGeometry& fvGeometry,
150 const ElementVolumeVariables& elemVolVars,
151 const SubControlVolumeFace& scvf,
153 const ElementFluxVariablesCache& elemFluxVarsCache)
156 const auto& dispersionTensor =
157 ModelTraits::ThermalDispersionModel::thermalDispersionTensor(problem, scvf, fvGeometry,
158 elemVolVars, elemFluxVarsCache,
161 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
162 Dune::FieldVector<Scalar, GridView::dimensionworld> gradTemp(0.0);
163 for (
auto&& scv : scvs(fvGeometry))
164 gradTemp.axpy(elemVolVars[scv].
temperature(), fluxVarsCache.gradN(scv.indexInElement()));
167 return -1.0*
vtmv(scvf.unitOuterNormal(), dispersionTensor, gradTemp)*Extrusion::area(fvGeometry, scvf);
Define some often used mathematical functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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:66
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:55
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:45
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
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
Definition: box/dispersionflux.hh:42
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: box/dispersionflux.hh:89
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:98
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:147
Traits of a flux variables type.
Definition: flux/traits.hh:44
static constexpr bool hasStationaryVelocityField()
Definition: flux/traits.hh:45
Declares all properties used in Dumux.