3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
lowrekepsilonncmodel.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 *****************************************************************************/
28#ifndef DUMUX_LOWREKEPSILON_NC_MODEL_HH
29#define DUMUX_LOWREKEPSILON_NC_MODEL_HH
30
35
36#include "iofields.hh"
37
38namespace Dumux {
39
41// properties for the single-phase, multi-component low-Re k-epsilon model
43namespace Properties {
44
46// Type tags
48
49// Create new type tags
50namespace TTag {
52struct LowReKEpsilonNC { using InheritsFrom = std::tuple<NavierStokesNC>; };
53} // end namespace TTag
54
56// default property values
58
68template<int dimension, int nComp, bool useMoles, int replaceCompEqIdx>
69struct LowReKEpsilonNCModelTraits : NavierStokesNCModelTraits<dimension, nComp, useMoles, replaceCompEqIdx>
70{
73 static constexpr int numEq() { return dimension+nComp+2; }
74
76 static constexpr bool usesTurbulenceModel() { return true; }
77
79 static constexpr auto turbulenceModel()
81
84};
85
87template<class TypeTag>
88struct ModelTraits<TypeTag, TTag::LowReKEpsilonNC>
89{
90private:
92 static constexpr int dimension = GridView::dimension;
94 static constexpr int numComponents = FluidSystem::numComponents;
95 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
96 static constexpr int replaceCompEqIdx = getPropValue<TypeTag, Properties::ReplaceCompEqIdx>();
97public:
99};
100
102template<class TypeTag>
103struct VolumeVariables<TypeTag, TTag::LowReKEpsilonNC>
104{
105private:
110 static_assert(FSY::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid system");
111 static_assert(FST::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid state");
112 static_assert(FSY::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid system");
113 static_assert(FST::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid state");
115
117 template<class BaseTraits, class DT>
118 struct NCTraits : public BaseTraits { using DiffusionType = DT; };
119
121public:
123};
124
126template<class TypeTag>
127struct LocalResidual<TypeTag, TTag::LowReKEpsilonNC>
128{
129private:
131public:
133};
134
136template<class TypeTag>
137struct FluxVariables<TypeTag, TTag::LowReKEpsilonNC>
138{
139private:
141public:
143};
144
146template<class TypeTag>
147struct IOFields<TypeTag, TTag::LowReKEpsilonNC> { using type = FreeflowNCIOFields<LowReKEpsilonIOFields, true/*turbulenceModel*/>; };
148
150// Property values for non-isothermal multi-component low-Re k-epsilon model
152
153// Create new type tags
154namespace TTag {
156struct LowReKEpsilonNCNI { using InheritsFrom = std::tuple<LowReKEpsilonNC, NavierStokesNCNI>; };
157} // end namespace TTag
158
160template<class TypeTag>
161struct ModelTraits<TypeTag, TTag::LowReKEpsilonNCNI>
162{
163private:
165 static constexpr int dim = GridView::dimension;
167 static constexpr int numComponents = FluidSystem::numComponents;
168 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
169 static constexpr int replaceCompEqIdx = getPropValue<TypeTag, Properties::ReplaceCompEqIdx>();
171public:
173};
174
176template<class TypeTag>
177struct VolumeVariables<TypeTag, TTag::LowReKEpsilonNCNI>
178{
179private:
184 static_assert(FSY::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid system");
185 static_assert(FST::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid state");
186 static_assert(FSY::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid system");
187 static_assert(FST::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid state");
189
191 template<class BaseTraits, class DT>
192 struct NCNITraits : public BaseTraits { using DiffusionType = DT; };
194public:
196};
197
199template<class TypeTag>
200struct LocalResidual<TypeTag, TTag::LowReKEpsilonNCNI>
201{
202private:
204public:
206};
207
209template<class TypeTag>
210struct FluxVariables<TypeTag, TTag::LowReKEpsilonNCNI>
211{
212private:
214public:
216};
217
219template<class TypeTag>
220struct IOFields<TypeTag, TTag::LowReKEpsilonNCNI>
221{
222private:
223 using IsothermalIOFields = FreeflowNCIOFields<LowReKEpsilonIOFields, true/*turbulenceModel*/>;
224public:
225 using type = FreeflowNonIsothermalIOFields<IsothermalIOFields, true/*turbulenceModel*/>;
226};
227
228} // end namespace Properties
229} // end namespace Dumux
230
231#endif
A single-phase, multi-component free-flow model.
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Traits class encapsulating model specifications.
Definition: common/properties.hh:53
A class helping models to define input and output fields.
Definition: common/properties.hh:63
Definition: common/properties.hh:74
The secondary variables within a sub-control volume.
Definition: common/properties.hh:107
Container storing the different types of flux variables.
Definition: common/properties.hh:113
The flux variables class for the multi-component free-flow model using the staggered grid discretizat...
Definition: freeflow/compositional/fluxvariables.hh:34
Adds I/O fields specific to the FreeflowNC model.
Definition: freeflow/compositional/iofields.hh:38
Element-wise calculation of the multi-component free-flow residual for models using the staggered dis...
Definition: freeflow/compositional/localresidual.hh:36
The type tags for the single-phase, multi-component isothermal low-Re k-epsilon model.
Definition: lowrekepsilonncmodel.hh:52
std::tuple< NavierStokesNC > InheritsFrom
Definition: lowrekepsilonncmodel.hh:52
Traits for the low-Reynolds k-epsilon multi-component modelstates some specifics of the isothermal mu...
Definition: lowrekepsilonncmodel.hh:70
static constexpr int numEq()
Definition: lowrekepsilonncmodel.hh:73
static constexpr bool usesTurbulenceModel()
The model does include a turbulence model.
Definition: lowrekepsilonncmodel.hh:76
static constexpr auto turbulenceModel()
return the type of turbulence model used
Definition: lowrekepsilonncmodel.hh:79
The type tags for the single-phase, multi-component non-isothermal low-Re k-epsilon models.
Definition: lowrekepsilonncmodel.hh:156
std::tuple< LowReKEpsilonNC, NavierStokesNCNI > InheritsFrom
Definition: lowrekepsilonncmodel.hh:156
Traits for the multi-component free-flow model.
Definition: navierstokesncmodel.hh:88
Volume variables for the single-phase, multi-component free-flow model.
Definition: freeflow/compositional/volumevariables.hh:40
Traits class for the volume variables of the Navier-Stokes model.
Definition: freeflow/navierstokes/model.hh:123
Adds I/O fields specific to non-isothermal free-flow models.
Definition: freeflow/nonisothermal/iofields.hh:38
Specifies a number properties of non-isothermal free-flow flow models based on the specifics of a giv...
Definition: freeflow/nonisothermal/model.hh:59
The common indices for isothermal two-equation RANS models.
Definition: freeflow/rans/twoeq/indices.hh:41
The flux variables class for the low-Reynolds k-epsilon model using the staggered grid discretization...
Definition: freeflow/rans/twoeq/lowrekepsilon/fluxvariables.hh:34
Adds I/O fields for the low-Re k-epsilon turbulence model.
Definition: freeflow/rans/twoeq/lowrekepsilon/iofields.hh:36
Element-wise calculation of the residual for low-Reynolds k-epsilon models using the staggered discre...
Definition: freeflow/rans/twoeq/lowrekepsilon/localresidual.hh:36
Volume variables for the isothermal single-phase low-Re k-epsilons model.
Definition: freeflow/rans/twoeq/lowrekepsilon/volumevariables.hh:40
A single-phase, isothermal low-Reynolds k-epsilon model.
Declares all properties used in Dumux.
Adds I/O fields specific to the tracer model.