24#ifndef DUMUX_FREEFLOW_NC_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_FREEFLOW_NC_STAGGERED_FLUXVARIABLES_HH
42template<
class TypeTag, DiscretizationMethod discMethod>
45template<
class TypeTag>
52 using FVElementGeometry =
typename GridGeometry::LocalView;
53 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
54 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
60 static constexpr auto numComponents = ModelTraits::numFluidComponents();
61 static constexpr bool useMoles = ModelTraits::useMoles();
67 template<
class ElementVolumeVariables,
class ElementFaceVariables,
class FluxVariablesCache>
69 const Element& element,
70 const FVElementGeometry& fvGeometry,
71 const ElementVolumeVariables& elemVolVars,
72 const ElementFaceVariables& elemFaceVars,
73 const SubControlVolumeFace& scvf,
74 const FluxVariablesCache& fluxVarsCache)
76 CellCenterPrimaryVariables flux(0.0);
78 const auto diffusiveFluxes = MolecularDiffusionType::flux(problem, element, fvGeometry, elemVolVars, scvf);
80 static constexpr auto referenceSystemFormulation = MolecularDiffusionType::referenceSystemFormulation();
84 auto upwindTerm = [compIdx](
const auto& volVars)
86 const auto density =
useMoles ? volVars.molarDensity() : volVars.density();
87 const auto fraction =
useMoles ? volVars.moleFraction(compIdx) : volVars.massFraction(compIdx);
88 return density * fraction;
91 flux[compIdx] = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTerm);
96 flux[compIdx] +=
useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx) : diffusiveFluxes[compIdx];
99 flux[compIdx] +=
useMoles ? diffusiveFluxes[compIdx] : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
101 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
107 flux[ModelTraits::replaceCompEqIdx()] = std::accumulate(flux.begin(), flux.end(), 0.0);
The available discretization methods in Dumux.
Base class for the flux variables living on a sub control volume face.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
NavierStokesFluxVariablesImpl< TypeTag, GetPropType< TypeTag, Properties::GridGeometry >::discMethod > NavierStokesFluxVariables
The flux variables class for the Navier-Stokes model. This is a convenience alias for that actual,...
Definition freeflow/navierstokes/fluxvariables.hh:45
DiscretizationMethod
The available discretization methods in Dumux.
Definition method.hh:37
@ staggered
Definition method.hh:38
@ massAveraged
Definition referencesystemformulation.hh:46
@ molarAveraged
Definition referencesystemformulation.hh:46
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition propertysystem.hh:150
The flux variables class for the multi-component free-flow model using the staggered grid discretizat...
Definition freeflow/compositional/fluxvariables.hh:34
GetPropType< TypeTag, Properties::MolecularDiffusionType > MolecularDiffusionType
Definition freeflow/compositional/staggered/fluxvariables.hh:62
static constexpr auto numComponents
Definition freeflow/compositional/staggered/fluxvariables.hh:60
CellCenterPrimaryVariables computeMassFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Computes the flux for the cell center residual.
Definition freeflow/compositional/staggered/fluxvariables.hh:68
static constexpr bool useMoles
Definition freeflow/compositional/staggered/fluxvariables.hh:61
Declares all properties used in Dumux.
The flux variables class for the Navier-Stokes model. This is a convenience alias for that actual,...