Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
Andrew Heather
committed
\\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM 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 3 of the License, or
(at your option) any later version.
OpenFOAM 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 OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
applyBoundaryLayer
Description
Apply a simplified boundary-layer model to the velocity and
turbulence fields based on the 1/7th power-law.
The uniform boundary-layer thickness is either provided via the -ybl option
or calculated as the average of the distance to the wall scaled with
the thickness coefficient supplied via the option -Cbl. If both options
are provided -ybl is used.
Andrew Heather
committed
Compressible modes is automatically selected based on the existence of the
"thermophysicalProperties" dictionary required to construct the
thermodynamics package.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "singlePhaseTransportModel.H"
#include "turbulentTransportModel.H"
mattijs
committed
#include "turbulentFluidThermoModel.H"
mattijs
committed
#include "processorFvPatchField.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// turbulence constants - file-scope
static const scalar Cmu(0.09);
static const scalar kappa(0.41);
mattijs
committed
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
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
161
162
163
164
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
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
void correctProcessorPatches(volScalarField& vf)
{
if (!Pstream::parRun())
{
return;
}
// Not possible to use correctBoundaryConditions on fields as they may
// use local info as opposed to the constraint values employed here,
// but still need to update processor patches
volScalarField::GeometricBoundaryField& bf = vf.boundaryField();
forAll(bf, patchI)
{
if (isA<processorFvPatchField<scalar> >(bf[patchI]))
{
bf[patchI].initEvaluate();
}
}
forAll(bf, patchI)
{
if (isA<processorFvPatchField<scalar> >(bf[patchI]))
{
bf[patchI].evaluate();
}
}
}
template<class TurbulenceModel>
Foam::tmp<Foam::volScalarField> calcK
(
TurbulenceModel& turbulence,
const volScalarField& mask,
const volScalarField& nut,
const volScalarField& y,
const dimensionedScalar& ybl,
const scalar Cmu,
const scalar kappa
)
{
// Turbulence k
tmp<volScalarField> tk = turbulence->k();
volScalarField& k = tk();
scalar ck0 = pow025(Cmu)*kappa;
k = (1 - mask)*k + mask*sqr(nut/(ck0*min(y, ybl)));
// Do not correct BC
// - operation may use inconsistent fields wrt these local manipulations
// k.correctBoundaryConditions();
correctProcessorPatches(k);
Info<< "Writing k\n" << endl;
k.write();
return tk;
}
template<class TurbulenceModel>
Foam::tmp<Foam::volScalarField> calcEpsilon
(
TurbulenceModel& turbulence,
const volScalarField& mask,
const volScalarField& k,
const volScalarField& y,
const dimensionedScalar& ybl,
const scalar Cmu,
const scalar kappa
)
{
// Turbulence epsilon
tmp<volScalarField> tepsilon = turbulence->epsilon();
volScalarField& epsilon = tepsilon();
scalar ce0 = ::pow(Cmu, 0.75)/kappa;
epsilon = (1 - mask)*epsilon + mask*ce0*k*sqrt(k)/min(y, ybl);
epsilon.max(SMALL);
// Do not correct BC
// - operation may use inconsistent fields wrt these local manipulations
// epsilon.correctBoundaryConditions();
correctProcessorPatches(epsilon);
Info<< "Writing epsilon\n" << endl;
epsilon.write();
return tepsilon;
}
void calcOmega
(
const fvMesh& mesh,
const volScalarField& mask,
const volScalarField& k,
const volScalarField& epsilon
)
{
// Turbulence omega
IOobject omegaHeader
(
"omega",
mesh.time().timeName(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
);
if (omegaHeader.headerOk())
{
volScalarField omega(omegaHeader, mesh);
dimensionedScalar k0("SMALL", k.dimensions(), SMALL);
omega = (1 - mask)*omega + mask*epsilon/(Cmu*k + k0);
omega.max(SMALL);
// Do not correct BC
// - operation may use inconsistent fields wrt these local
// manipulations
// omega.correctBoundaryConditions();
correctProcessorPatches(omega);
Info<< "Writing omega\n" << endl;
omega.write();
}
}
void setField
(
const fvMesh& mesh,
const word& fieldName,
const volScalarField& value
)
{
IOobject fldHeader
(
fieldName,
mesh.time().timeName(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
);
if (fldHeader.headerOk())
{
volScalarField fld(fldHeader, mesh);
fld = value;
// Do not correct BC
// - operation may use inconsistent fields wrt these local
// manipulations
// fld.correctBoundaryConditions();
correctProcessorPatches(fld);
Info<< "Writing " << fieldName << nl << endl;
fld.write();
}
}
void calcCompressible
(
const fvMesh& mesh,
const volScalarField& mask,
const volVectorField& U,
const volScalarField& y,
const dimensionedScalar& ybl
)
{
const Time& runTime = mesh.time();
autoPtr<fluidThermo> pThermo(fluidThermo::New(mesh));
fluidThermo& thermo = pThermo();
volScalarField rho(thermo.rho());
// Update/re-write phi
#include "compressibleCreatePhi.H"
phi.write();
autoPtr<compressible::turbulenceModel> turbulence
(
compressible::turbulenceModel::New
(
rho,
U,
phi,
thermo
)
);
// Calculate nut - reference nut is calculated by the turbulence model
// on its construction
tmp<volScalarField> tnut = turbulence->nut();
volScalarField& nut = tnut();
volScalarField S(mag(dev(symm(fvc::grad(U)))));
nut = (1 - mask)*nut + mask*sqr(kappa*min(y, ybl))*::sqrt(2)*S;
// Do not correct BC - wall functions will 'undo' manipulation above
// by using nut from turbulence model
correctProcessorPatches(nut);
nut.write();
tmp<volScalarField> k =
calcK(turbulence, mask, nut, y, ybl, Cmu, kappa);
tmp<volScalarField> epsilon =
calcEpsilon(turbulence, mask, k, y, ybl, Cmu, kappa);
calcOmega(mesh, mask, k, epsilon);
setField(mesh, "nuTilda", nut);
}
void calcIncompressible
(
const fvMesh& mesh,
const volScalarField& mask,
const volVectorField& U,
const volScalarField& y,
const dimensionedScalar& ybl
)
{
const Time& runTime = mesh.time();
// Update/re-write phi
#include "createPhi.H"
phi.write();
singlePhaseTransportModel laminarTransport(U, phi);
autoPtr<incompressible::turbulenceModel> turbulence
(
incompressible::turbulenceModel::New(U, phi, laminarTransport)
);
tmp<volScalarField> tnut = turbulence->nut();
// Calculate nut - reference nut is calculated by the turbulence model
// on its construction
volScalarField& nut = tnut();
volScalarField S(mag(dev(symm(fvc::grad(U)))));
nut = (1 - mask)*nut + mask*sqr(kappa*min(y, ybl))*::sqrt(2)*S;
// Do not correct BC - wall functions will 'undo' manipulation above
// by using nut from turbulence model
correctProcessorPatches(nut);
nut.write();
tmp<volScalarField> k =
calcK(turbulence, mask, nut, y, ybl, Cmu, kappa);
tmp<volScalarField> epsilon =
calcEpsilon(turbulence, mask, k, y, ybl, Cmu, kappa);
calcOmega(mesh, mask, k, epsilon);
setField(mesh, "nuTilda", nut);
}
int main(int argc, char *argv[])
{
argList::addNote
(
"apply a simplified boundary-layer model to the velocity and\n"
"turbulence fields based on the 1/7th power-law."
);
mattijs
committed
#include "addRegionOption.H"
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
argList::addOption
(
"ybl",
"scalar",
"specify the boundary-layer thickness"
);
argList::addOption
(
"Cbl",
"scalar",
"boundary-layer thickness as Cbl * mean distance to wall"
);
#include "setRootCase.H"
if (!args.optionFound("ybl") && !args.optionFound("Cbl"))
{
FatalErrorIn(args.executable())
<< "Neither option 'ybl' or 'Cbl' have been provided to calculate "
<< "the boundary-layer thickness.\n"
<< "Please choose either 'ybl' OR 'Cbl'."
<< exit(FatalError);
}
else if (args.optionFound("ybl") && args.optionFound("Cbl"))
{
FatalErrorIn(args.executable())
<< "Both 'ybl' and 'Cbl' have been provided to calculate "
<< "the boundary-layer thickness.\n"
<< "Please choose either 'ybl' OR 'Cbl'."
<< exit(FatalError);
}
#include "createTime.H"
mattijs
committed
#include "createNamedMesh.H"
#include "createFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// Modify velocity by applying a 1/7th power law boundary-layer
// u/U0 = (y/ybl)^(1/7)
// assumes U0 is the same as the current cell velocity
Info<< "Setting boundary layer velocity" << nl << endl;
scalar yblv = ybl.value();
forAll(U, cellI)
{
if (y[cellI] <= yblv)
{
mask[cellI] = 1;
U[cellI] *= ::pow(y[cellI]/yblv, (1.0/7.0));
}
}
mask.correctBoundaryConditions();
Info<< "Writing U\n" << endl;
U.write();
Andrew Heather
committed
if
(
IOobject
(
basicThermo::dictName,
runTime.constant(),
mesh
).headerOk()
)
mattijs
committed
calcCompressible(mesh, mask, U, y, ybl);
mattijs
committed
else
{
calcIncompressible(mesh, mask, U, y, ybl);
}
Info<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //