/
applyBoundaryLayer.C
236 lines (187 loc) · 7 KB
/
applyBoundaryLayer.C
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
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
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
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
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.
\*---------------------------------------------------------------------------*/
#include "argList.H"
#include "viscosityModel.H"
#include "incompressibleMomentumTransportModels.H"
#include "wallDist.H"
#include "fvcFlux.H"
using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// Turbulence constants - file-scope
static const scalar Cmu(0.09);
static const scalar kappa(0.41);
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."
);
argList::addOption
(
"ybl",
"scalar",
"specify the boundary-layer thickness"
);
argList::addOption
(
"Cbl",
"scalar",
"boundary-layer thickness as Cbl * mean distance to wall"
);
argList::addBoolOption
(
"writenut",
"write nut field"
);
#include "setRootCase.H"
if (!args.optionFound("ybl") && !args.optionFound("Cbl"))
{
FatalErrorInFunction
<< "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"))
{
FatalErrorInFunction
<< "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"
#include "createMeshNoChangers.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();
// Update/re-write phi
#include "createPhi.H"
phi.write();
autoPtr<viscosityModel> viscosity(viscosityModel::New(mesh));
autoPtr<incompressible::momentumTransportModel> turbulence
(
incompressible::momentumTransportModel::New(U, phi, viscosity)
);
if (isA<incompressible::RASModel>(turbulence()))
{
// Calculate nut
turbulence->validate();
tmp<volScalarField> tnut = turbulence->nut();
volScalarField& nut = const_cast<volScalarField&>(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
if (args.optionFound("writenut"))
{
Info<< "Writing nut" << endl;
nut.write();
}
//--- Read and modify turbulence fields
// Turbulence k
tmp<volScalarField> tk = turbulence->k();
volScalarField& k = const_cast<volScalarField&>(tk());
scalar ck0 = pow025(Cmu)*kappa;
k = (1 - mask)*k + mask*sqr(nut/(ck0*min(y, ybl)));
k.correctBoundaryConditions();
Info<< "Writing k\n" << endl;
k.write();
// Turbulence epsilon
tmp<volScalarField> tepsilon = turbulence->epsilon();
volScalarField& epsilon = const_cast<volScalarField&>(tepsilon());
scalar ce0 = ::pow(Cmu, 0.75)/kappa;
epsilon = (1 - mask)*epsilon + mask*ce0*k*sqrt(k)/min(y, ybl);
// Do not correct BC - G set by the wall-functions is not available
// epsilon.correctBoundaryConditions();
Info<< "Writing epsilon\n" << endl;
epsilon.write();
// Turbulence omega
typeIOobject<volScalarField> omegaHeader
(
"omega",
runTime.name(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
);
if (omegaHeader.headerOk())
{
volScalarField omega(omegaHeader, mesh);
omega = (1 - mask)*omega + mask*ce0*sqrt(k)/(Cmu*min(y, ybl));
// Do not correct BC - G set by the wall-functions is not available
// omega.correctBoundaryConditions();
Info<< "Writing omega\n" << endl;
omega.write();
}
// Turbulence nuTilda
typeIOobject<volScalarField> nuTildaHeader
(
"nuTilda",
runTime.name(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
);
if (nuTildaHeader.headerOk())
{
volScalarField nuTilda(nuTildaHeader, mesh);
nuTilda = nut;
// Do not correct BC - G set by the wall-functions is not available
// nuTilda.correctBoundaryConditions();
Info<< "Writing nuTilda\n" << endl;
nuTilda.write();
}
}
Info<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //