version 3.10-dev
komegancmodel.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_KOMEGA_NC_MODEL_HH
17#define DUMUX_KOMEGA_NC_MODEL_HH
18
23
24#include "iofields.hh"
25
26namespace Dumux {
27
29// properties for the single-phase, multi-component k-omega model
31namespace Properties {
32
34// Type tags
36
37// Create new type tags
38namespace TTag {
40struct KOmegaNC { 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 KOmegaNCModelTraits : 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()
68 { return TurbulenceModel::komega; }
69
72};
73
75template<class TypeTag>
76struct ModelTraits<TypeTag, TTag::KOmegaNC>
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::KOmegaNC>
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; };
108public:
110};
111
113template<class TypeTag>
114struct LocalResidual<TypeTag, TTag::KOmegaNC>
115{
116private:
118public:
120};
121
123template<class TypeTag>
124struct FluxVariables<TypeTag, TTag::KOmegaNC>
125{
126private:
128public:
130};
131
133template<class TypeTag>
134struct IOFields<TypeTag, TTag::KOmegaNC> { using type = FreeflowNCIOFields<KOmegaIOFields, true/*turbulenceModel*/>; };
135
137// Property values for non-isothermal multi-component k-omega model
139
140// Create new type tags
141namespace TTag {
143struct KOmegaNCNI { using InheritsFrom = std::tuple<KOmegaNC, NavierStokesNCNI>; };
144} // end namespace TTag
145
147template<class TypeTag>
148struct ModelTraits<TypeTag, TTag::KOmegaNCNI>
149{
150private:
152 static constexpr int dim = GridView::dimension;
154 static constexpr int numComponents = FluidSystem::numComponents;
155 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
156 static constexpr int replaceCompEqIdx = getPropValue<TypeTag, Properties::ReplaceCompEqIdx>();
158
159public:
161};
162
164template<class TypeTag>
165struct VolumeVariables<TypeTag, TTag::KOmegaNCNI>
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::KOmegaNCNI>
189{
190private:
192public:
194};
195
197template<class TypeTag>
198struct FluxVariables<TypeTag, TTag::KOmegaNCNI>
199{
200private:
202public:
204};
205
207template<class TypeTag>
208struct IOFields<TypeTag, TTag::KOmegaNCNI>
209{
210private:
211 using IsothermalIOFields = FreeflowNCIOFields<KOmegaIOFields, 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 k-omega model using the staggered grid discretization.
Definition: freeflow/rans/twoeq/komega/fluxvariables.hh:22
Element-wise calculation of the residual for k-omega models using the staggered discretization.
Definition: freeflow/rans/twoeq/komega/localresidual.hh:24
Volume variables for the isothermal single-phase k-omega 2-Eq model.
Definition: freeflow/rans/twoeq/komega/volumevariables.hh:29
Defines all properties used in Dumux.
A single-phase, isothermal k-omega 2-Eq. 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 Reynolds-Averaged Navier-Stokes model.
Definition: freeflow/rans/twoeq/komega/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 k-omega multi-component modelstates some specifics of the isothermal multi-component k...
Definition: komegancmodel.hh:58
static constexpr int numEq()
Definition: komegancmodel.hh:61
static constexpr auto turbulenceModel()
return the type of turbulence model used
Definition: komegancmodel.hh:67
static constexpr bool usesTurbulenceModel()
The model does include a turbulence model.
Definition: komegancmodel.hh:64
The type tags for the single-phase, multi-component isothermal k-omega model.
Definition: komegancmodel.hh:40
std::tuple< NavierStokesNC > InheritsFrom
Definition: komegancmodel.hh:40
The type tags for the single-phase, multi-component non-isothermal k-omega models.
Definition: komegancmodel.hh:143
std::tuple< KOmegaNC, NavierStokesNCNI > InheritsFrom
Definition: komegancmodel.hh:143
The common indices for isothermal two-equation RANS models.
Definition: freeflow/rans/twoeq/indices.hh:29