version 3.8
freeflow/compositional/staggered/fluxvariables.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_FREEFLOW_NC_STAGGERED_FLUXVARIABLES_HH
13#define DUMUX_FREEFLOW_NC_STAGGERED_FLUXVARIABLES_HH
14
15#include <numeric>
21
22namespace Dumux {
23
29// forward declaration
30template<class TypeTag, class DiscretizationMethod>
31class FreeflowNCFluxVariablesImpl;
32
33template<class TypeTag>
34class FreeflowNCFluxVariablesImpl<TypeTag, DiscretizationMethods::Staggered>
35: public NavierStokesFluxVariables<TypeTag>
36{
40 using FVElementGeometry = typename GridGeometry::LocalView;
41 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
42 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
43 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
47
48public:
49 static constexpr auto numComponents = ModelTraits::numFluidComponents();
50 static constexpr bool useMoles = ModelTraits::useMoles();
52
56 template<class ElementVolumeVariables, class ElementFaceVariables, class FluxVariablesCache>
57 CellCenterPrimaryVariables computeMassFlux(const Problem& problem,
58 const Element& element,
59 const FVElementGeometry& fvGeometry,
60 const ElementVolumeVariables& elemVolVars,
61 const ElementFaceVariables& elemFaceVars,
62 const SubControlVolumeFace& scvf,
63 const FluxVariablesCache& fluxVarsCache)
64 {
65 CellCenterPrimaryVariables flux(0.0);
66
67 const auto diffusiveFluxes = MolecularDiffusionType::flux(problem, element, fvGeometry, elemVolVars, scvf);
68
69 static constexpr auto referenceSystemFormulation = MolecularDiffusionType::referenceSystemFormulation();
70
71 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
72 {
73 auto upwindTerm = [compIdx](const auto& volVars)
74 {
75 const auto density = useMoles ? volVars.molarDensity() : volVars.density();
76 const auto fraction = useMoles ? volVars.moleFraction(compIdx) : volVars.massFraction(compIdx);
77 return density * fraction;
78 };
79
80 flux[compIdx] = ParentType::advectiveFluxForCellCenter(problem, fvGeometry, elemVolVars, elemFaceVars, scvf, upwindTerm);
81
82 // check for the reference system and adapt units of the diffusive flux accordingly.
83 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
84 {
85 flux[compIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx) : diffusiveFluxes[compIdx];
86 }
87 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
88 flux[compIdx] += useMoles ? diffusiveFluxes[compIdx] : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
89 else
90 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
91 }
92
93
94 // in case one balance is substituted by the total mass balance
95 if (ModelTraits::replaceCompEqIdx() < numComponents)
96 {
97 // accumulate fluxes to a total mass based flux
98 Scalar totalMassFlux = 0.0;
99 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
100 totalMassFlux += useMoles ? flux[compIdx]*FluidSystem::molarMass(compIdx) : flux[compIdx];
101
102 flux[ModelTraits::replaceCompEqIdx()] = totalMassFlux;
103 }
104
105 return flux;
106 }
107};
108
109} // end namespace
110
111#endif
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:57
GetPropType< TypeTag, Properties::MolecularDiffusionType > MolecularDiffusionType
Definition: freeflow/compositional/staggered/fluxvariables.hh:51
The flux variables class for the multi-component free-flow model using the staggered grid discretizat...
Definition: freeflow/compositional/fluxvariables.hh:22
The flux variables class for the Navier-Stokes model using the staggered grid discretization.
Definition: freeflow/navierstokes/fluxvariables.hh:23
Defines all properties used in Dumux.
Base class for the flux variables living on a sub control volume face.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
Definition: adapt.hh:17
The reference frameworks and formulations available for splitting total fluxes into a advective and d...