DeardorffDiffStress.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2021 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 "DeardorffDiffStress.H"
30#include "fvOptions.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace LESModels
38{
39
40// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41
42template<class BasicTurbulenceModel>
44{
45 this->nut_ = Ck_*sqrt(this->k())*this->delta();
46 this->nut_.correctBoundaryConditions();
47 fv::options::New(this->mesh_).correct(this->nut_);
48
49 BasicTurbulenceModel::correctNut();
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
55template<class BasicTurbulenceModel>
57(
58 const alphaField& alpha,
59 const rhoField& rho,
60 const volVectorField& U,
61 const surfaceScalarField& alphaRhoPhi,
63 const transportModel& transport,
64 const word& propertiesName,
65 const word& type
66)
67:
68 ReynoldsStress<LESModel<BasicTurbulenceModel>>
69 (
70 type,
71 alpha,
72 rho,
73 U,
74 alphaRhoPhi,
75 phi,
76 transport,
77 propertiesName
78 ),
79
80 Ck_
81 (
82 dimensioned<scalar>::getOrAddToDict
83 (
84 "Ck",
85 this->coeffDict_,
86 0.094
87 )
88 ),
89 Cm_
90 (
91 dimensioned<scalar>::getOrAddToDict
92 (
93 "Cm",
94 this->coeffDict_,
95 4.13
96 )
97 ),
98 Ce_
99 (
100 dimensioned<scalar>::getOrAddToDict
101 (
102 "Ce",
103 this->coeffDict_,
104 1.05
105 )
106 ),
107 Cs_
108 (
109 dimensioned<scalar>::getOrAddToDict
110 (
111 "Cs",
112 this->coeffDict_,
113 0.25
114 )
115 )
116{
117 if (type == typeName)
118 {
119 this->printCoeffs(type);
120 this->boundNormalStress(this->R_);
121 }
122}
123
124
125// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
126
127template<class BasicTurbulenceModel>
129{
131 {
132 Ck_.readIfPresent(this->coeffDict());
133 Cm_.readIfPresent(this->coeffDict());
134 Ce_.readIfPresent(this->coeffDict());
135 Cs_.readIfPresent(this->coeffDict());
136
137 return true;
138 }
139
140 return false;
141}
142
143
144template<class BasicTurbulenceModel>
146{
147 if (!this->turbulence_)
148 {
149 return;
150 }
151
152 // Local references
153 const alphaField& alpha = this->alpha_;
154 const rhoField& rho = this->rho_;
155 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
156 const volVectorField& U = this->U_;
157 volSymmTensorField& R = this->R_;
159
161
163 const volTensorField& gradU = tgradU();
164
165 volSymmTensorField D(symm(gradU));
166
167 volSymmTensorField P(-twoSymm(R & gradU));
168
169 volScalarField k(this->k());
170
172 (
174 + fvm::div(alphaRhoPhi, R)
175 - fvm::laplacian(I*this->nu() + Cs_*(k/this->epsilon())*R, R)
176 + fvm::Sp(Cm_*alpha*rho*sqrt(k)/this->delta(), R)
177 ==
178 alpha*rho*P
179 + (4.0/5.0)*alpha*rho*k*D
180 - ((2.0/3.0)*(1.0 - Cm_/this->Ce_)*I)*(alpha*rho*this->epsilon())
181 + fvOptions(alpha, rho, R)
182 );
183
184 REqn.ref().relax();
185 fvOptions.constrain(REqn.ref());
186 REqn.ref().solve();
187 fvOptions.correct(R);
188 this->boundNormalStress(R);
189
190 correctNut();
191}
192
193
194// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195
196} // End namespace LESModels
197} // End namespace Foam
198
199// ************************************************************************* //
label k
scalar delta
#define R(A, B, C, D, E, F, K, M)
fv::options & fvOptions
surfaceScalarField & phi
Templated abstract base class for LES SGS models.
Definition: LESModel.H:62
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:34
Differential SGS Stress Equation Model for incompressible and compressible flows.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual void correct()
Correct sub-grid stress, eddy-Viscosity and related properties.
virtual void correctNut()
Update the eddy-viscosity.
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Read model coefficients if they have changed.
Reynolds-stress turbulence model base class.
void boundNormalStress(volSymmTensorField &R) const
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Generic dimensioned Type class.
Finite-volume options.
Definition: fvOptions.H:59
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
thermo correct()
scalar epsilon
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:94
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar sqrt(const dimensionedScalar &ds)
volScalarField & nu
volScalarField & alpha
const dimensionedScalar & D