3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_FREEFLOW_NC_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_FREEFLOW_NC_STAGGERED_FLUXVARIABLES_HH
26
27#include <numeric>
33
34namespace Dumux {
35
36// forward declaration
37template<class TypeTag, DiscretizationMethod discMethod>
38class FreeflowNCFluxVariablesImpl;
39
44template<class TypeTag>
46: public NavierStokesFluxVariables<TypeTag>
47{
51 using FVElementGeometry = typename GridGeometry::LocalView;
52 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
53 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
54 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
57
58public:
59 static constexpr auto numComponents = ModelTraits::numFluidComponents();
60 static constexpr bool useMoles = ModelTraits::useMoles();
62
66 template<class ElementVolumeVariables, class ElementFaceVariables, class FluxVariablesCache>
67 CellCenterPrimaryVariables computeMassFlux(const Problem& problem,
68 const Element& element,
69 const FVElementGeometry& fvGeometry,
70 const ElementVolumeVariables& elemVolVars,
71 const ElementFaceVariables& elemFaceVars,
72 const SubControlVolumeFace& scvf,
73 const FluxVariablesCache& fluxVarsCache)
74 {
75 CellCenterPrimaryVariables flux(0.0);
76
77 const auto diffusiveFluxes = MolecularDiffusionType::flux(problem, element, fvGeometry, elemVolVars, scvf);
78
79 static constexpr auto referenceSystemFormulation = MolecularDiffusionType::referenceSystemFormulation();
80
81 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
82 {
83 auto upwindTerm = [compIdx](const auto& volVars)
84 {
85 const auto density = useMoles ? volVars.molarDensity() : volVars.density();
86 const auto fraction = useMoles ? volVars.moleFraction(compIdx) : volVars.massFraction(compIdx);
87 return density * fraction;
88 };
89
90 flux[compIdx] = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTerm);
91
92 //check for the reference system and adapt units of the diffusive flux accordingly.
93 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
94 {
95 flux[compIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx) : diffusiveFluxes[compIdx];
96 }
97 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
98 flux[compIdx] += useMoles ? diffusiveFluxes[compIdx] : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
99 else
100 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
101 }
102
103 // in case one balance is substituted by the total mass balance
104 if (ModelTraits::replaceCompEqIdx() < numComponents)
105 {
106 flux[ModelTraits::replaceCompEqIdx()] = std::accumulate(flux.begin(), flux.end(), 0.0);
107 }
108
109 return flux;
110 }
111};
112
113} // end namespace
114
115#endif
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...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Definition: freeflow/compositional/fluxvariables.hh:34
GetPropType< TypeTag, Properties::MolecularDiffusionType > MolecularDiffusionType
Definition: freeflow/compositional/staggered/fluxvariables.hh:61
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:67
Definition: freeflow/navierstokes/fluxvariables.hh:35
Declares all properties used in Dumux.