13#ifndef DUMUX_DISCRETIZATION_CC_TPFA_DISPERSION_FLUX_HH
14#define DUMUX_DISCRETIZATION_CC_TPFA_DISPERSION_FLUX_HH
16#include <dune/common/fvector.hh>
17#include <dune/common/fmatrix.hh>
30template<
class TypeTag,
class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
31class DispersionFluxImplementation;
37template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
45 using FVElementGeometry =
typename GridGeometry::LocalView;
46 using SubControlVolume =
typename GridGeometry::SubControlVolume;
47 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
51 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
52 using FluxVarCache =
typename GridFluxVariablesCache::FluxVariablesCache;
57 using Element =
typename GridView::template Codim<0>::Entity;
59 using Indices =
typename ModelTraits::Indices;
61 enum { dim = GridView::dimension} ;
62 enum { dimWorld = GridView::dimensionworld} ;
65 numPhases = ModelTraits::numFluidPhases(),
66 numComponents = ModelTraits::numFluidComponents()
69 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
70 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
71 using HeatFluxScalar = Scalar;
79 {
return referenceSystem; }
88 const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const SubControlVolumeFace& scvf,
93 const ElementFluxVariablesCache& elemFluxVarsCache)
95 if (scvf.numOutsideScvs() > 1 )
96 DUNE_THROW(Dune::NotImplemented,
"\n Dispersion using ccTPFA is only implemented for conforming grids.");
97 if (!stationaryVelocityField)
98 DUNE_THROW(Dune::NotImplemented,
"\n Dispersion using ccTPFA is only implemented for problems with stationary velocity fields");
100 ComponentFluxVector componentFlux(0.0);
101 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
102 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
105 const auto rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
106 const Scalar rho = 0.5*(rhoInside + rhoOutside);
108 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
110 const auto& dispersionTensor =
111 ModelTraits::CompositionalDispersionModel::compositionalDispersionTensor(problem, scvf, fvGeometry,
112 elemVolVars, elemFluxVarsCache,
114 const auto dij =
computeTpfaTransmissibility(fvGeometry, scvf, fvGeometry.scv(scvf.insideScvIdx()), dispersionTensor, insideVolVars.extrusionFactor());
116 const auto xInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
117 const auto xOutide =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
119 componentFlux[compIdx] = (rho * (xInside-xOutide) * dij) * scvf.area();
121 return componentFlux;
129 const Element& element,
130 const FVElementGeometry& fvGeometry,
131 const ElementVolumeVariables& elemVolVars,
132 const SubControlVolumeFace& scvf,
134 const ElementFluxVariablesCache& elemFluxVarsCache)
136 if (scvf.numOutsideScvs() > 1 )
137 DUNE_THROW(Dune::NotImplemented,
"\n Dispersion using ccTPFA is only implemented for conforming grids.");
138 if (!stationaryVelocityField)
139 DUNE_THROW(Dune::NotImplemented,
"\n Dispersion using ccTPFA is only implemented for problems with stationary velocity fields");
141 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
142 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
144 const auto& dispersionTensor =
145 ModelTraits::ThermalDispersionModel::thermalDispersionTensor(problem, scvf, fvGeometry,
146 elemVolVars, elemFluxVarsCache,
148 const auto dij =
computeTpfaTransmissibility(scvf, fvGeometry.scv(scvf.insideScvIdx()), dispersionTensor, insideVolVars.extrusionFactor());
151 const auto tInside = insideVolVars.temperature();
152 const auto tOutside = outsideVolVars.temperature();
155 return (tInside-tOutside) * dij * scvf.area();
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: cctpfa/dispersionflux.hh:87
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: cctpfa/dispersionflux.hh:78
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: cctpfa/dispersionflux.hh:128
Definition: box/dispersionflux.hh:30
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Tensor::field_type computeTpfaTransmissibility(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const typename FVElementGeometry::SubControlVolume &scv, const Tensor &T, typename FVElementGeometry::SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition: tpfa/computetransmissibility.hh:36
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.
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
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...