model.hh 8.86 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
 *   See the file COPYING for full copying permissions.                      *
 *                                                                           *
 *   This program is free software: you can redistribute it and/or modify    *
 *   it under the terms of the GNU General Public License as published by    *
 *   the Free Software Foundation, either version 2 of the License, or       *
 *   (at your option) any later version.                                     *
 *                                                                           *
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
 *   GNU General Public License for more details.                            *
 *                                                                           *
 *   You should have received a copy of the GNU General Public License       *
 *   along with this program.  If not, see <http://www.gnu.org/licenses/>.   *
 *****************************************************************************/
/*!
 * \file
 * \ingroup OneEqModel
 *
 * \brief A single-phase, isothermal one-equation turbulence model by Spalart-Allmaras
 *
 * \copydoc RANSModel
 *
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
 * This model, published by Spalart and Allmaras 1992 \cite Spalart1992a,
 * uses one additional PDE for a working variable \f$ \tilde{\nu} \f$.
 * This variable has the units of a viscosity and can be converted to the eddy
 * viscosity via a model function~(\f$ f_\text{v1} \f$):
 * \f[
 *  \nu_\text{t} = \tilde{\nu} f_\text{v1}
 * \f]
 *
 * Here, as proposed by Wilcox \cite Wilcox2008a and Versteeg \cite Versteeg2009a, the correction
 * term which account for the transition or trip, is dropped from the original equations,
 * such that the balance equation simplifies to:
 * \f[
 *   \frac{\partial \tilde{\nu}}{\partial t}
 *   + \nabla \cdot \left( \tilde{\nu} \textbf{v} \right)
 *   - c_\text{b1} \tilde{S} \tilde{\nu}
 *   - \frac{1}{\sigma_{\tilde{\nu}}} \nabla \cdot \left( \left[ \nu + \tilde{\nu} \right] \nabla \tilde{\nu} \right)
 *   - \frac{c_\text{b2}}{\sigma_{\tilde{\nu}}} \left| \nabla \tilde{\nu} \right|^2
 *   + c_\text{w1} f_\text{w} \frac{\tilde{\nu}^2}{y^2}
 *   = 0
 * \f]
 *
 * Here, a modified mean effective strain rate (\f$ \tilde{S} \f$) based on
 * the mean rotation rate tensor (\f$ \mathbf{\Omega} \f$) is used:
 * \f[
 *  \tilde{S} = \sqrt{2 \mathbf{\Omega} \cdot \mathbf{\Omega}} + \frac{\tilde{\nu}}{\kappa^2 y^2} f_\text{v2}
 * \f]
 * \f[
 *  \mathbf{\Omega} = \frac{1}{2} \left( \nabla \textbf{v}_\text{g}
 *                                  - \nabla \textbf{v}_\text{g}^\intercal \right)
 * \f]
 *
 * This balance equation is linked to the flow geometry by the distance to the closest wall ($y$).
 * Further, the model uses the following functions and expressions:
 * \f[ \chi = \frac{\tilde{\nu}}{\nu} \f]
 * \f[ f_\text{v1} = \frac{\chi^3}{\chi^3+c_\text{v1}^3} \f]
 * \f[ f_\text{v2} = 1 - \frac{\chi}{1+f_\text{v1}\chi} \f]
 * \f[ f_\text{w} = g_\text{w} \left( \frac{1+c_\text{w3}^6}{g^6_\text{w}+c_\text{w3}^6}
 *                             \right)^\frac{1}{6} \f]
 * \f[ g_\text{w} = r_\text{w} + c_\text{w2} (r_\text{w}^6 - r_\text{w}) \f]
 * \f[ r_\text{w} = \min \left[ \frac{\tilde{\nu}}{\tilde{S}\kappa^2 y^2},10\right] \f]
 * \f[ \sigma_{\tilde{\nu}} = \nicefrac{2}{3} \f]
 * \f[ c_\text{b1} = 0.1355 \f]
 * \f[ c_\text{b2} = 0.622 \f]
 * \f[ c_\text{v1} = 7.1 \f]
 * \f[ c_\text{w1} = \frac{c_\text{b1}}{\kappa^2}
 *                   + \frac{1+c_\text{b2}}{\sigma_{\tilde{\nu}}} \f]
 * \f[ c_\text{w2} = 0.3 \f]
 * \f[ c_\text{w3} = 2 \f]
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
 */

#ifndef DUMUX_ONEEQ_MODEL_HH
#define DUMUX_ONEEQ_MODEL_HH

#include <dumux/common/properties.hh>
#include <dumux/freeflow/properties.hh>
#include <dumux/freeflow/rans/model.hh>

#include "fluxvariables.hh"
#include "indices.hh"
#include "localresidual.hh"
#include "volumevariables.hh"
#include "vtkoutputfields.hh"

namespace Dumux
{
namespace Properties {

/*!
 * \ingroup OneEqModel
 * \brief Traits for the Spalart-Allmaras model
97
98
 *
 * \tparam dimension The dimension of the problem
99
 */
100
101
template<int dimension>
struct OneEqModelTraits : RANSModelTraits<dimension>
102
103
104
105
106
107
108
109
110
111
112
113
{
    //! The dimension of the model
    static constexpr int dim() { return dimension; }

    //! There are as many momentum balance equations as dimensions,
    //! one mass balance equation and two turbulent transport equations
    static constexpr int numEq() { return dim()+1+1; }

    //! The number of components
    static constexpr int numComponents() { return 1; }

    //! the indices
114
    using Indices = OneEqIndices<dim(), numComponents()>;
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
};

///////////////////////////////////////////////////////////////////////////
// default property values for the isothermal Spalart-Allmaras model
///////////////////////////////////////////////////////////////////////////

//! The type tag for the single-phase, isothermal Spalart-Allmaras model
NEW_TYPE_TAG(OneEq, INHERITS_FROM(RANS));

//!< states some specifics of the isothermal Spalart-Allmaras model
SET_PROP(OneEq, ModelTraits)
{
private:
    using GridView = typename GET_PROP_TYPE(TypeTag, FVGridGeometry)::GridView;
    static constexpr int dim = GridView::dimension;
public:
131
    using type = OneEqModelTraits<dim>;
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
};

//! The flux variables
SET_PROP(OneEq, FluxVariables)
{
private:
    using BaseFluxVariables = NavierStokesFluxVariables<TypeTag>;
public:
    using type = OneEqFluxVariables<TypeTag, BaseFluxVariables>;
};

//! The local residual
SET_PROP(OneEq, LocalResidual)
{
private:
    using BaseLocalResidual = NavierStokesResidual<TypeTag>;
public:
    using type = OneEqResidual<TypeTag, BaseLocalResidual>;
};

//! Set the volume variables property
SET_PROP(OneEq, VolumeVariables)
{
private:
    using PV = typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
    using FSY = typename GET_PROP_TYPE(TypeTag, FluidSystem);
    using FST = typename GET_PROP_TYPE(TypeTag, FluidState);
    using MT = typename GET_PROP_TYPE(TypeTag, ModelTraits);

161
162
163
164
    static_assert(FSY::numPhases == MT::numPhases(), "Number of phases mismatch between model and fluid system");
    static_assert(FST::numPhases == MT::numPhases(), "Number of phases mismatch between model and fluid state");
    static_assert(!FSY::isMiscible(), "The Navier-Stokes model only works with immiscible fluid systems.");

165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
    using Traits = NavierStokesVolumeVariablesTraits<PV, FSY, FST, MT>;
    using NSVolVars = NavierStokesVolumeVariables<Traits>;
public:
    using type = OneEqVolumeVariables<Traits, NSVolVars>;
};

//! The specific vtk output fields
SET_PROP(OneEq, VtkOutputFields)
{
private:
    using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
public:
    using type = OneEqVtkOutputFields<FVGridGeometry>;
};

//////////////////////////////////////////////////////////////////
// default property values for the non-isothermal Spalart-Allmaras model
//////////////////////////////////////////////////////////////////

//! The type tag for the single-phase, isothermal Spalart-Allmaras model
NEW_TYPE_TAG(OneEqNI, INHERITS_FROM(RANSNI));

//! The model traits of the non-isothermal model
SET_PROP(OneEqNI, ModelTraits)
{
private:
    using GridView = typename GET_PROP_TYPE(TypeTag, FVGridGeometry)::GridView;
    static constexpr int dim = GridView::dimension;
193
    using IsothermalTraits = OneEqModelTraits<dim>;
194
195
196
197
198
199
200
201
202
203
204
205
206
public:
    using type = FreeflowNIModelTraits<IsothermalTraits>;
};

//! Set the volume variables property
SET_PROP(OneEqNI, VolumeVariables)
{
private:
    using PV = typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
    using FSY = typename GET_PROP_TYPE(TypeTag, FluidSystem);
    using FST = typename GET_PROP_TYPE(TypeTag, FluidState);
    using MT = typename GET_PROP_TYPE(TypeTag, ModelTraits);

207
208
209
210
    static_assert(FSY::numPhases == MT::numPhases(), "Number of phases mismatch between model and fluid system");
    static_assert(FST::numPhases == MT::numPhases(), "Number of phases mismatch between model and fluid state");
    static_assert(!FSY::isMiscible(), "The Navier-Stokes model only works with immiscible fluid systems.");

211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
    using Traits = NavierStokesVolumeVariablesTraits<PV, FSY, FST, MT>;
    using NSVolVars = NavierStokesVolumeVariables<Traits>;
public:
    using type = OneEqVolumeVariables<Traits, NSVolVars>;
};

//! The specific non-isothermal vtk output fields
SET_PROP(OneEqNI, VtkOutputFields)
{
private:
    using ModelTraits = typename GET_PROP_TYPE(TypeTag, ModelTraits);
    using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
    using IsothermalFields = OneEqVtkOutputFields<FVGridGeometry>;
public:
    using type = FreeflowNonIsothermalVtkOutputFields<IsothermalFields, ModelTraits>;
};

// \}
}

} // end namespace

#endif