version 3.10-dev
kepsilonncmodel.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//
16#ifndef DUMUX_KEPSILON_NC_MODEL_HH
17#define DUMUX_KEPSILON_NC_MODEL_HH
18
23
24#include "iofields.hh"
25
26namespace Dumux {
27
29// properties for the single-phase, multi-component k-epsilon model
31namespace Properties {
32
34// Type tags
36
37// Create new type tags
38namespace TTag {
40struct KEpsilonNC { using InheritsFrom = std::tuple<NavierStokesNC>; };
41} // end namespace TTag
42
44// default property values
46
51template<int dimension, int nComp, bool useMoles, int replaceCompEqIdx>
52struct KEpsilonNCModelTraits : NavierStokesNCModelTraits<dimension, nComp, useMoles, replaceCompEqIdx>
53{
56 static constexpr int numEq() { return dimension+nComp+2; }
57
59 static constexpr bool usesTurbulenceModel() { return true; }
60
63
65 static constexpr auto turbulenceModel()
67};
68
70template<class TypeTag>
71struct ModelTraits<TypeTag, TTag::KEpsilonNC>
72{
73private:
75 static constexpr int dimension = GridView::dimension;
77 static constexpr int numComponents = FluidSystem::numComponents;
78 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
79 static constexpr int replaceCompEqIdx = getPropValue<TypeTag, Properties::ReplaceCompEqIdx>();
80public:
82};
83
85template<class TypeTag>
86struct VolumeVariables<TypeTag, TTag::KEpsilonNC>
87{
88private:
93 static_assert(FSY::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid system");
94 static_assert(FST::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid state");
95 static_assert(FSY::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid system");
96 static_assert(FST::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid state");
98
100 template<class BaseTraits, class DT>
101 struct NCTraits : public BaseTraits { using DiffusionType = DT; };
102
104public:
106};
107
109template<class TypeTag>
110struct LocalResidual<TypeTag, TTag::KEpsilonNC>
111{
112private:
114public:
116};
117
119template<class TypeTag>
120struct FluxVariables<TypeTag, TTag::KEpsilonNC>
121{
122private:
124public:
126};
127
129template<class TypeTag>
130struct IOFields<TypeTag, TTag::KEpsilonNC> { using type = FreeflowNCIOFields<KEpsilonIOFields, true/*turbulenceModel*/>; };
131
133// Property values for non-isothermal multi-component k-epsilon model
135
136// Create new type tags
137namespace TTag {
139struct KEpsilonNCNI { using InheritsFrom = std::tuple<KEpsilonNC, NavierStokesNCNI>; };
140} // end namespace TTag
141
143template<class TypeTag>
144struct ModelTraits<TypeTag, TTag::KEpsilonNCNI>
145{
146private:
148 static constexpr int dim = GridView::dimension;
150 static constexpr int numComponents = FluidSystem::numComponents;
151 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
152 static constexpr int replaceCompEqIdx = getPropValue<TypeTag, Properties::ReplaceCompEqIdx>();
154public:
156};
157
159template<class TypeTag>
160struct VolumeVariables<TypeTag, TTag::KEpsilonNCNI>
161{
162private:
167 static_assert(FSY::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid system");
168 static_assert(FST::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid state");
169 static_assert(FSY::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid system");
170 static_assert(FST::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid state");
172
174 template<class BaseTraits, class DT>
175 struct NCNITraits : public BaseTraits { using DiffusionType = DT; };
177public:
179};
180
182template<class TypeTag>
183struct LocalResidual<TypeTag, TTag::KEpsilonNCNI>
184{
185private:
187public:
189};
190
192template<class TypeTag>
193struct FluxVariables<TypeTag, TTag::KEpsilonNCNI>
194{
195private:
197public:
199};
200
202template<class TypeTag>
203struct IOFields<TypeTag, TTag::KEpsilonNCNI>
204{
205private:
206 using IsothermalIOFields = FreeflowNCIOFields<KEpsilonIOFields, true/*turbulenceModel*/>;
207public:
208 using type = FreeflowNonIsothermalIOFields<IsothermalIOFields, true/*turbulenceModel*/>;
209};
210
211} // end namespace Properties
212} // end namespace Dumux
213
214#endif
The flux variables class for the multi-component free-flow model using the staggered grid discretizat...
Definition: freeflow/compositional/fluxvariables.hh:22
Element-wise calculation of the multi-component free-flow residual for models using the staggered dis...
Definition: freeflow/compositional/localresidual.hh:24
Volume variables for the single-phase, multi-component free-flow model.
Definition: freeflow/compositional/volumevariables.hh:28
The flux variables class for the k-epsilon model using the staggered grid discretization.
Definition: freeflow/rans/twoeq/kepsilon/fluxvariables.hh:22
Element-wise calculation of the residual for k-epsilon models using the staggered discretization.
Definition: freeflow/rans/twoeq/kepsilon/localresidual.hh:24
Volume variables for the isothermal single-phase k-epsilon model.
Definition: freeflow/rans/twoeq/kepsilon/volumevariables.hh:27
Defines all properties used in Dumux.
A single-phase, isothermal k-epsilon model.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: adapt.hh:17
A single-phase, multi-component free-flow model.
Adds I/O fields specific to the tracer model.
Adds I/O fields specific to the FreeflowNC model.
Definition: freeflow/compositional/iofields.hh:26
Specifies a number properties of non-isothermal free-flow flow models based on the specifics of a giv...
Definition: freeflow/nonisothermal/model.hh:47
Adds I/O fields specific to non-isothermal free-flow models.
Definition: freeflow/nonisothermal/iofields.hh:26
Adds I/O fields for the k-epsilon turbulence model.
Definition: freeflow/rans/twoeq/kepsilon/iofields.hh:24
Traits for the multi-component free-flow model.
Definition: navierstokesncmodel.hh:76
Traits class for the volume variables of the Navier-Stokes model.
Definition: freeflow/navierstokes/model.hh:111
Traits for the low-Reynolds k-epsilon multi-component model states some specifics of the isothermal m...
Definition: kepsilonncmodel.hh:53
static constexpr int numEq()
Definition: kepsilonncmodel.hh:56
static constexpr bool usesTurbulenceModel()
The model does include a turbulence model.
Definition: kepsilonncmodel.hh:59
static constexpr auto turbulenceModel()
return the type of turbulence model used
Definition: kepsilonncmodel.hh:65
The type tags for the single-phase, multi-component isothermal k-epsilon model.
Definition: kepsilonncmodel.hh:40
std::tuple< NavierStokesNC > InheritsFrom
Definition: kepsilonncmodel.hh:40
The type tags for the single-phase, multi-component non-isothermal k-epsilon models.
Definition: kepsilonncmodel.hh:139
std::tuple< KEpsilonNC, NavierStokesNCNI > InheritsFrom
Definition: kepsilonncmodel.hh:139
The common indices for isothermal two-equation RANS models.
Definition: freeflow/rans/twoeq/indices.hh:29