version 3.10-dev
freeflow/rans/twoeq/lowrekepsilon/model.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//
66#ifndef DUMUX_LOWREKEPSILON_MODEL_HH
67#define DUMUX_LOWREKEPSILON_MODEL_HH
68
74
75#include "problem.hh"
76#include "fluxvariables.hh"
77#include "localresidual.hh"
78#include "volumevariables.hh"
79#include "iofields.hh"
80
81namespace Dumux {
82namespace Properties {
83
90template<int dimension>
92{
94 static constexpr int dim() { return dimension; }
95
98 static constexpr int numEq() { return dim()+1+2; }
99
101 static constexpr int numFluidComponents() { return 1; }
102
105
107 static constexpr auto turbulenceModel()
109};
110
112// default property values for the isothermal low-Reynolds k-epsilon model
114
115// Create new type tags
116namespace TTag {
118struct LowReKEpsilon { using InheritsFrom = std::tuple<RANS>; };
119} // end namespace TTag
120
122template<class TypeTag>
123struct ModelTraits<TypeTag, TTag::LowReKEpsilon>
124{
125private:
127 static constexpr int dim = GridView::dimension;
128public:
130};
131
133template<class TypeTag>
134struct FluxVariables<TypeTag, TTag::LowReKEpsilon>
135{
136private:
138public:
140};
141
143template<class TypeTag>
144struct LocalResidual<TypeTag, TTag::LowReKEpsilon>
145{
146private:
148public:
150};
151
153template<class TypeTag>
154struct VolumeVariables<TypeTag, TTag::LowReKEpsilon>
155{
156private:
161
162 static_assert(FSY::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid system");
163 static_assert(FST::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid state");
164 static_assert(!FSY::isMiscible(), "The Navier-Stokes model only works with immiscible fluid systems.");
165
168public:
170};
171
173template<class TypeTag>
174struct IOFields<TypeTag, TTag::LowReKEpsilon> { using type = LowReKEpsilonIOFields; };
175
177// default property values for the non-isothermal low-Reynolds k-epsilon model
179
180// Create new type tags
181namespace TTag {
183struct LowReKEpsilonNI { using InheritsFrom = std::tuple<LowReKEpsilon, RANSNI>; };
184} // end namespace TTag
185
187template<class TypeTag>
188struct ModelTraits<TypeTag, TTag::LowReKEpsilonNI>
189{
190private:
192 static constexpr int dim = GridView::dimension;
194public:
196};
197
199template<class TypeTag>
200struct VolumeVariables<TypeTag, TTag::LowReKEpsilonNI>
201{
202private:
207
208 static_assert(FSY::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid system");
209 static_assert(FST::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid state");
210 static_assert(!FSY::isMiscible(), "The Navier-Stokes model only works with immiscible fluid systems.");
211
214public:
216};
217
219template<class TypeTag>
220struct IOFields<TypeTag, TTag::LowReKEpsilonNI> { using type = FreeflowNonIsothermalIOFields<LowReKEpsilonIOFields, true/*turbulenceModel*/>; };
221
222} // end properties
223} // end namespace
224
225#endif
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
The flux variables class for the Navier-Stokes model using the staggered grid discretization.
Definition: freeflow/navierstokes/fluxvariables.hh:23
Element-wise calculation of the Navier-Stokes residual for models using the staggered discretization.
Definition: freeflow/navierstokes/localresidual.hh:23
Volume variables for the single-phase Navier-Stokes model.
Definition: freeflow/navierstokes/volumevariables.hh:26
Defines all properties used in Dumux.
Defines a type tag and some properties for free flow models.
A single-phase, isothermal Reynolds-Averaged Navier-Stokes model.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: adapt.hh:17
Base class for the flux variables in porous medium models.
Adds I/O fields specific to the tracer model.
Element-wise calculation of the local residual for problems using fully implicit tracer model.
Python wrapper for volume variables (finite volume schemes)
TODO: docme!
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
The common indices for the isothermal Navier-Stokes model.
Definition: freeflow/navierstokes/indices.hh:25
Traits class for the volume variables of the Navier-Stokes model.
Definition: freeflow/navierstokes/model.hh:111
Traits for the low-Reynolds k-epsilon model.
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:92
static constexpr int dim()
The dimension of the model.
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:94
static constexpr int numEq()
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:98
static constexpr int numFluidComponents()
The number of components.
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:101
static constexpr auto turbulenceModel()
return the type of turbulence model used
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:107
Traits for the Reynolds-averaged Navier-Stokes model.
Definition: freeflow/rans/model.hh:61
The type tag for the single-phase, isothermal low-Reynolds k-epsilon model.
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:118
std::tuple< RANS > InheritsFrom
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:118
The type tag for the single-phase, non-isothermal low-Reynolds k-epsilon model.
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:183
std::tuple< LowReKEpsilon, RANSNI > InheritsFrom
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:183
The common indices for isothermal two-equation RANS models.
Definition: freeflow/rans/twoeq/indices.hh:29
The available free flow turbulence models in Dumux.