version 3.8
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// 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_LOWREKEPSILON_NC_MODEL_HH
17#define DUMUX_LOWREKEPSILON_NC_MODEL_HH
18
23
24#include "iofields.hh"
25
26namespace Dumux {
27
29// properties for the single-phase, multi-component low-Re k-epsilon model
31namespace Properties {
32
34// Type tags
36
37// Create new type tags
38namespace TTag {
40struct LowReKEpsilonNC { using InheritsFrom = std::tuple<NavierStokesNC>; };
41} // end namespace TTag
42
44// default property values
46
56template<int dimension, int nComp, bool useMoles, int replaceCompEqIdx>
57struct LowReKEpsilonNCModelTraits : NavierStokesNCModelTraits<dimension, nComp, useMoles, replaceCompEqIdx>
58{
61 static constexpr int numEq() { return dimension+nComp+2; }
62
64 static constexpr bool usesTurbulenceModel() { return true; }
65
67 static constexpr auto turbulenceModel()
69
72};
73
75template<class TypeTag>
76struct ModelTraits<TypeTag, TTag::LowReKEpsilonNC>
77{
78private:
80 static constexpr int dimension = GridView::dimension;
82 static constexpr int numComponents = FluidSystem::numComponents;
83 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
84 static constexpr int replaceCompEqIdx = getPropValue<TypeTag, Properties::ReplaceCompEqIdx>();
85public:
87};
88
90template<class TypeTag>
91struct VolumeVariables<TypeTag, TTag::LowReKEpsilonNC>
92{
93private:
98 static_assert(FSY::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid system");
99 static_assert(FST::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid state");
100 static_assert(FSY::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid system");
101 static_assert(FST::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid state");
103
105 template<class BaseTraits, class DT>
106 struct NCTraits : public BaseTraits { using DiffusionType = DT; };
107
109public:
111};
112
114template<class TypeTag>
115struct LocalResidual<TypeTag, TTag::LowReKEpsilonNC>
116{
117private:
119public:
121};
122
124template<class TypeTag>
125struct FluxVariables<TypeTag, TTag::LowReKEpsilonNC>
126{
127private:
129public:
131};
132
134template<class TypeTag>
135struct IOFields<TypeTag, TTag::LowReKEpsilonNC> { using type = FreeflowNCIOFields<LowReKEpsilonIOFields, true/*turbulenceModel*/>; };
136
138// Property values for non-isothermal multi-component low-Re k-epsilon model
140
141// Create new type tags
142namespace TTag {
144struct LowReKEpsilonNCNI { using InheritsFrom = std::tuple<LowReKEpsilonNC, NavierStokesNCNI>; };
145} // end namespace TTag
146
148template<class TypeTag>
149struct ModelTraits<TypeTag, TTag::LowReKEpsilonNCNI>
150{
151private:
153 static constexpr int dim = GridView::dimension;
155 static constexpr int numComponents = FluidSystem::numComponents;
156 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
157 static constexpr int replaceCompEqIdx = getPropValue<TypeTag, Properties::ReplaceCompEqIdx>();
159public:
161};
162
164template<class TypeTag>
165struct VolumeVariables<TypeTag, TTag::LowReKEpsilonNCNI>
166{
167private:
172 static_assert(FSY::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid system");
173 static_assert(FST::numComponents == MT::numFluidComponents(), "Number of components mismatch between model and fluid state");
174 static_assert(FSY::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid system");
175 static_assert(FST::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid state");
177
179 template<class BaseTraits, class DT>
180 struct NCNITraits : public BaseTraits { using DiffusionType = DT; };
182public:
184};
185
187template<class TypeTag>
188struct LocalResidual<TypeTag, TTag::LowReKEpsilonNCNI>
189{
190private:
192public:
194};
195
197template<class TypeTag>
198struct FluxVariables<TypeTag, TTag::LowReKEpsilonNCNI>
199{
200private:
202public:
204};
205
207template<class TypeTag>
208struct IOFields<TypeTag, TTag::LowReKEpsilonNCNI>
209{
210private:
211 using IsothermalIOFields = FreeflowNCIOFields<LowReKEpsilonIOFields, true/*turbulenceModel*/>;
212public:
213 using type = FreeflowNonIsothermalIOFields<IsothermalIOFields, true/*turbulenceModel*/>;
214};
215
216} // end namespace Properties
217} // end namespace Dumux
218
219#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 low-Reynolds k-epsilon model using the staggered grid discretization...
Definition: freeflow/rans/twoeq/lowrekepsilon/fluxvariables.hh:22
Element-wise calculation of the residual for low-Reynolds k-epsilon models using the staggered discre...
Definition: freeflow/rans/twoeq/lowrekepsilon/localresidual.hh:24
Volume variables for the isothermal single-phase low-Re k-epsilons model.
Definition: freeflow/rans/twoeq/lowrekepsilon/volumevariables.hh:28
Defines all properties used in Dumux.
A single-phase, isothermal low-Reynolds 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 low-Re k-epsilon turbulence model.
Definition: freeflow/rans/twoeq/lowrekepsilon/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 modelstates some specifics of the isothermal mu...
Definition: lowrekepsilonncmodel.hh:58
static constexpr int numEq()
Definition: lowrekepsilonncmodel.hh:61
static constexpr bool usesTurbulenceModel()
The model does include a turbulence model.
Definition: lowrekepsilonncmodel.hh:64
static constexpr auto turbulenceModel()
return the type of turbulence model used
Definition: lowrekepsilonncmodel.hh:67
The type tags for the single-phase, multi-component isothermal low-Re k-epsilon model.
Definition: lowrekepsilonncmodel.hh:40
std::tuple< NavierStokesNC > InheritsFrom
Definition: lowrekepsilonncmodel.hh:40
The type tags for the single-phase, multi-component non-isothermal low-Re k-epsilon models.
Definition: lowrekepsilonncmodel.hh:144
std::tuple< LowReKEpsilonNC, NavierStokesNCNI > InheritsFrom
Definition: lowrekepsilonncmodel.hh:144
The common indices for isothermal two-equation RANS models.
Definition: freeflow/rans/twoeq/indices.hh:29