3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 *****************************************************************************/
78#ifndef DUMUX_LOWREKEPSILON_MODEL_HH
79#define DUMUX_LOWREKEPSILON_MODEL_HH
80
86
87#include "problem.hh"
88#include "fluxvariables.hh"
89#include "localresidual.hh"
90#include "volumevariables.hh"
91#include "iofields.hh"
92
93namespace Dumux {
94namespace Properties {
95
102template<int dimension>
104{
106 static constexpr int dim() { return dimension; }
107
110 static constexpr int numEq() { return dim()+1+2; }
111
113 static constexpr int numFluidComponents() { return 1; }
114
117
119 static constexpr auto turbulenceModel()
121};
122
124// default property values for the isothermal low-Reynolds k-epsilon model
126
127// Create new type tags
128namespace TTag {
130struct LowReKEpsilon { using InheritsFrom = std::tuple<RANS>; };
131} // end namespace TTag
132
134template<class TypeTag>
135struct ModelTraits<TypeTag, TTag::LowReKEpsilon>
136{
137private:
139 static constexpr int dim = GridView::dimension;
140public:
142};
143
145template<class TypeTag>
146struct FluxVariables<TypeTag, TTag::LowReKEpsilon>
147{
148private:
150public:
152};
153
155template<class TypeTag>
156struct LocalResidual<TypeTag, TTag::LowReKEpsilon>
157{
158private:
160public:
162};
163
165template<class TypeTag>
166struct VolumeVariables<TypeTag, TTag::LowReKEpsilon>
167{
168private:
173
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");
176 static_assert(!FSY::isMiscible(), "The Navier-Stokes model only works with immiscible fluid systems.");
177
180public:
182};
183
185template<class TypeTag>
186struct IOFields<TypeTag, TTag::LowReKEpsilon> { using type = LowReKEpsilonIOFields; };
187
189// default property values for the non-isothermal low-Reynolds k-epsilon model
191
192// Create new type tags
193namespace TTag {
195struct LowReKEpsilonNI { using InheritsFrom = std::tuple<LowReKEpsilon, RANSNI>; };
196} // end namespace TTag
197
199template<class TypeTag>
200struct ModelTraits<TypeTag, TTag::LowReKEpsilonNI>
201{
202private:
204 static constexpr int dim = GridView::dimension;
206public:
208};
209
211template<class TypeTag>
212struct VolumeVariables<TypeTag, TTag::LowReKEpsilonNI>
213{
214private:
219
220 static_assert(FSY::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid system");
221 static_assert(FST::numPhases == MT::numFluidPhases(), "Number of phases mismatch between model and fluid state");
222 static_assert(!FSY::isMiscible(), "The Navier-Stokes model only works with immiscible fluid systems.");
223
226public:
228};
229
231template<class TypeTag>
232struct IOFields<TypeTag, TTag::LowReKEpsilonNI> { using type = FreeflowNonIsothermalIOFields<LowReKEpsilonIOFields, true/*turbulenceModel*/>; };
233
234} // end properties
235} // end namespace
236
237#endif
The available free flow turbulence models in Dumux.
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:103
Container storing the different types of flux variables.
Definition: common/properties.hh:109
The flux variables class for the Navier-Stokes model using the staggered grid discretization.
Definition: freeflow/navierstokes/fluxvariables.hh:35
The common indices for the isothermal Navier-Stokes model.
Definition: freeflow/navierstokes/indices.hh:37
Element-wise calculation of the Navier-Stokes residual for models using the staggered discretization.
Definition: freeflow/navierstokes/localresidual.hh:35
Traits class for the volume variables of the Navier-Stokes model.
Definition: freeflow/navierstokes/model.hh:125
Volume variables for the single-phase Navier-Stokes model.
Definition: freeflow/navierstokes/volumevariables.hh:38
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
Traits for the Reynolds-averaged Navier-Stokes model.
Definition: freeflow/rans/model.hh:73
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
Traits for the low-Reynolds k-epsilon model.
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:104
static constexpr int dim()
The dimension of the model.
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:106
static constexpr int numEq()
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:110
static constexpr int numFluidComponents()
The number of components.
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:113
static constexpr auto turbulenceModel()
return the type of turbulence model used
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:119
The type tag for the single-phase, isothermal low-Reynolds k-epsilon model.
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:130
std::tuple< RANS > InheritsFrom
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:130
The type tag for the single-phase, non-isothermal low-Reynolds k-epsilon model.
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:195
std::tuple< LowReKEpsilon, RANSNI > InheritsFrom
Definition: freeflow/rans/twoeq/lowrekepsilon/model.hh:195
Volume variables for the isothermal single-phase low-Re k-epsilons model.
Definition: freeflow/rans/twoeq/lowrekepsilon/volumevariables.hh:40
Declares 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.
Base class for all porous media problems.
Adds I/O fields specific to the tracer model.
Base class for the flux variables in porous medium models.
Element-wise calculation of the local residual for problems using fully implicit tracer model.
Base class for the model specific class which provides access to all volume averaged quantities.