Loading...
Searching...
No Matches
sigma.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2022 Upstream CFD GmbH
9 Copyright (C) 2022-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "sigma.H"
30#include "DESModel.H"
31#include "fvOptions.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
37namespace LESModels
38{
39
40// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41
42template<class BasicTurbulenceModel>
44{
45 this->nut_ =
46 sqr(this->delta()*Csigma_)
48 this->nut_.correctBoundaryConditions();
49 fv::options::New(this->mesh_).correct(this->nut_);
51 BasicTurbulenceModel::correctNut();
52}
53
54
55// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56
57template<class BasicTurbulenceModel> sigma<BasicTurbulenceModel>::sigma
58(
59 const alphaField& alpha,
60 const rhoField& rho,
61 const volVectorField& U,
62 const surfaceScalarField& alphaRhoPhi,
64 const transportModel& transport,
65 const word& propertiesName,
66 const word& type
67)
68:
69 LESeddyViscosity<BasicTurbulenceModel>
70 (
71 type,
72 alpha,
73 rho,
74 U,
75 alphaRhoPhi,
76 phi,
77 transport,
78 propertiesName
79 ),
80
81 Ck_
82 (
83 dimensioned<scalar>::getOrAddToDict
84 (
85 "Ck",
86 this->coeffDict_,
87 0.094
88 )
89 ),
90
91 Cw_
92 (
93 dimensioned<scalar>::getOrAddToDict
94 (
95 "Cw",
96 this->coeffDict_,
97 0.325
98 )
99 ),
100
101 Csigma_
102 (
103 dimensioned<scalar>::getOrAddToDict
104 (
105 "Csigma",
106 this->coeffDict_,
107 1.68
108 )
109 )
110{
111 if (type == typeName)
112 {
113 this->printCoeffs(type);
115}
116
117
118// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119
120template<class BasicTurbulenceModel>
122{
124 {
125 Ck_.readIfPresent(this->coeffDict());
126 Cw_.readIfPresent(this->coeffDict());
127 Csigma_.readIfPresent(this->coeffDict());
128
129 return true;
131
132 return false;
133}
134
135
136template<class BasicTurbulenceModel>
138{
140 (
141 IOobject::groupName("k", this->U_.group()),
142 (2.0*Ck_/this->Ce_)
143 *sqr(this->delta())
144 *magSqr(devSymm(fvc::grad(this->U_)))
145 );
146}
147
148
149template<class BasicTurbulenceModel>
151{
153 (
154 IOobject::groupName("epsilon", this->U_.group()),
155 this->Ce_*k()*sqrt(k())/this->delta()
156 );
157}
158
159
160template<class BasicTurbulenceModel>
162{
164 correctNut();
165}
166
167
168// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169
170} // End namespace LESModels
171} // End namespace Foam
172
173// ************************************************************************* //
label k
scalar delta
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
static tmp< volScalarField > Ssigma(const volTensorField &gradU, const dimensionedScalar &coeff)
Return modified strain rate.
Definition DESModel.C:82
Eddy viscosity LES SGS model base class.
virtual bool read()
Read model coefficients if they have changed.
BasicTurbulenceModel::alphaField alphaField
Definition sigma.H:113
BasicTurbulenceModel::rhoField rhoField
Definition sigma.H:114
dimensionedScalar Ck_
Definition sigma.H:98
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition sigma.C:130
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition sigma.C:154
virtual tmp< volScalarField > epsilon() const
Return SGS disipation rate.
Definition sigma.C:143
dimensionedScalar Csigma_
Definition sigma.H:100
virtual void correctNut()
Update the SGS eddy-viscosity.
Definition sigma.C:36
dimensionedScalar Cw_
Definition sigma.H:99
BasicTurbulenceModel::transportModel transportModel
Definition sigma.H:115
virtual bool read()
Read model coefficients if they have changed.
Definition sigma.C:114
Generic dimensioned Type class.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
Namespace for LES SGS models.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
volScalarField & alpha