applyBoundaryLayer.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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  applyBoundaryLayer
26 
27 Description
28  Apply a simplified boundary-layer model to the velocity and
29  turbulence fields based on the 1/7th power-law.
30 
31  The uniform boundary-layer thickness is either provided via the -ybl option
32  or calculated as the average of the distance to the wall scaled with
33  the thickness coefficient supplied via the option -Cbl. If both options
34  are provided -ybl is used.
35 
36  Compressible modes is automatically selected based on the existence of the
37  "thermophysicalProperties" dictionary required to construct the
38  thermodynamics package.
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #include "fvCFD.H"
46 #include "wallDist.H"
47 #include "processorFvPatchField.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 // turbulence constants - file-scope
52 static const scalar Cmu(0.09);
53 static const scalar kappa(0.41);
54 
55 void correctProcessorPatches(volScalarField& vf)
56 {
57  if (!Pstream::parRun())
58  {
59  return;
60  }
61 
62  // Not possible to use correctBoundaryConditions on fields as they may
63  // use local info as opposed to the constraint values employed here,
64  // but still need to update processor patches
65 
66  volScalarField::GeometricBoundaryField& bf = vf.boundaryField();
67 
68  forAll(bf, patchI)
69  {
70  if (isA<processorFvPatchField<scalar> >(bf[patchI]))
71  {
72  bf[patchI].initEvaluate();
73  }
74  }
75 
76  forAll(bf, patchI)
77  {
78  if (isA<processorFvPatchField<scalar> >(bf[patchI]))
79  {
80  bf[patchI].evaluate();
81  }
82  }
83 }
84 
85 
86 template<class TurbulenceModel>
88 (
89  TurbulenceModel& turbulence,
90  const volScalarField& mask,
91  const volScalarField& nut,
92  const volScalarField& y,
93  const dimensionedScalar& ybl,
94  const scalar Cmu,
95  const scalar kappa
96 )
97 {
98  // Turbulence k
99  tmp<volScalarField> tk = turbulence->k();
100  volScalarField& k = tk();
101  scalar ck0 = pow025(Cmu)*kappa;
102  k = (1 - mask)*k + mask*sqr(nut/(ck0*min(y, ybl)));
103 
104  // Do not correct BC
105  // - operation may use inconsistent fields wrt these local manipulations
106  // k.correctBoundaryConditions();
107  correctProcessorPatches(k);
108 
109  Info<< "Writing k\n" << endl;
110  k.write();
111 
112  return tk;
113 }
114 
115 
116 template<class TurbulenceModel>
118 (
119  TurbulenceModel& turbulence,
120  const volScalarField& mask,
121  const volScalarField& k,
122  const volScalarField& y,
123  const dimensionedScalar& ybl,
124  const scalar Cmu,
125  const scalar kappa
126 )
127 {
128  // Turbulence epsilon
129  tmp<volScalarField> tepsilon = turbulence->epsilon();
130  volScalarField& epsilon = tepsilon();
131  scalar ce0 = ::pow(Cmu, 0.75)/kappa;
132  epsilon = (1 - mask)*epsilon + mask*ce0*k*sqrt(k)/min(y, ybl);
133  epsilon.max(SMALL);
134 
135  // Do not correct BC
136  // - operation may use inconsistent fields wrt these local manipulations
137  // epsilon.correctBoundaryConditions();
138  correctProcessorPatches(epsilon);
139 
140  Info<< "Writing epsilon\n" << endl;
141  epsilon.write();
142 
143  return tepsilon;
144 }
145 
146 
147 void calcOmega
148 (
149  const fvMesh& mesh,
150  const volScalarField& mask,
151  const volScalarField& k,
152  const volScalarField& epsilon
153 )
154 {
155  // Turbulence omega
156  IOobject omegaHeader
157  (
158  "omega",
159  mesh.time().timeName(),
160  mesh,
161  IOobject::MUST_READ,
162  IOobject::NO_WRITE,
163  false
164  );
165 
166  if (omegaHeader.headerOk())
167  {
168  volScalarField omega(omegaHeader, mesh);
169  dimensionedScalar k0("SMALL", k.dimensions(), SMALL);
170 
171  omega = (1 - mask)*omega + mask*epsilon/(Cmu*k + k0);
172  omega.max(SMALL);
173 
174  // Do not correct BC
175  // - operation may use inconsistent fields wrt these local
176  // manipulations
177  // omega.correctBoundaryConditions();
178  correctProcessorPatches(omega);
179 
180  Info<< "Writing omega\n" << endl;
181  omega.write();
182  }
183 }
184 
185 
186 void setField
187 (
188  const fvMesh& mesh,
189  const word& fieldName,
190  const volScalarField& value
191 )
192 {
193  IOobject fldHeader
194  (
195  fieldName,
196  mesh.time().timeName(),
197  mesh,
198  IOobject::MUST_READ,
199  IOobject::NO_WRITE,
200  false
201  );
202 
203  if (fldHeader.headerOk())
204  {
205  volScalarField fld(fldHeader, mesh);
206  fld = value;
207 
208  // Do not correct BC
209  // - operation may use inconsistent fields wrt these local
210  // manipulations
211  // fld.correctBoundaryConditions();
212  correctProcessorPatches(fld);
213 
214  Info<< "Writing " << fieldName << nl << endl;
215  fld.write();
216  }
217 }
218 
219 
220 void calcCompressible
221 (
222  const fvMesh& mesh,
223  const volScalarField& mask,
224  const volVectorField& U,
225  const volScalarField& y,
226  const dimensionedScalar& ybl
227 )
228 {
229  const Time& runTime = mesh.time();
230 
231  autoPtr<fluidThermo> pThermo(fluidThermo::New(mesh));
232  fluidThermo& thermo = pThermo();
233  volScalarField rho(thermo.rho());
234 
235  // Update/re-write phi
236  #include "compressibleCreatePhi.H"
237  phi.write();
238 
239  autoPtr<compressible::turbulenceModel> turbulence
240  (
242  (
243  rho,
244  U,
245  phi,
246  thermo
247  )
248  );
249 
250  // Calculate nut - reference nut is calculated by the turbulence model
251  // on its construction
252  tmp<volScalarField> tnut = turbulence->nut();
253 
254  volScalarField& nut = tnut();
256  nut = (1 - mask)*nut + mask*sqr(kappa*min(y, ybl))*::sqrt(2)*S;
257 
258  // Do not correct BC - wall functions will 'undo' manipulation above
259  // by using nut from turbulence model
260  correctProcessorPatches(nut);
261  nut.write();
262 
263  tmp<volScalarField> k =
264  calcK(turbulence, mask, nut, y, ybl, Cmu, kappa);
265  tmp<volScalarField> epsilon =
266  calcEpsilon(turbulence, mask, k, y, ybl, Cmu, kappa);
267  calcOmega(mesh, mask, k, epsilon);
268  setField(mesh, "nuTilda", nut);
269 }
270 
271 
272 void calcIncompressible
273 (
274  const fvMesh& mesh,
275  const volScalarField& mask,
276  const volVectorField& U,
277  const volScalarField& y,
278  const dimensionedScalar& ybl
279 )
280 {
281  const Time& runTime = mesh.time();
282 
283  // Update/re-write phi
284  #include "createPhi.H"
285  phi.write();
286 
287  singlePhaseTransportModel laminarTransport(U, phi);
288 
289  autoPtr<incompressible::turbulenceModel> turbulence
290  (
292  );
293 
294  tmp<volScalarField> tnut = turbulence->nut();
295 
296  // Calculate nut - reference nut is calculated by the turbulence model
297  // on its construction
298  volScalarField& nut = tnut();
300  nut = (1 - mask)*nut + mask*sqr(kappa*min(y, ybl))*::sqrt(2)*S;
301 
302  // Do not correct BC - wall functions will 'undo' manipulation above
303  // by using nut from turbulence model
304  correctProcessorPatches(nut);
305  nut.write();
306 
307  tmp<volScalarField> k =
308  calcK(turbulence, mask, nut, y, ybl, Cmu, kappa);
309  tmp<volScalarField> epsilon =
310  calcEpsilon(turbulence, mask, k, y, ybl, Cmu, kappa);
311  calcOmega(mesh, mask, k, epsilon);
312  setField(mesh, "nuTilda", nut);
313 }
314 
315 
316 int main(int argc, char *argv[])
317 {
318  argList::addNote
319  (
320  "apply a simplified boundary-layer model to the velocity and\n"
321  "turbulence fields based on the 1/7th power-law."
322  );
323 
324  #include "addRegionOption.H"
325 
326  argList::addOption
327  (
328  "ybl",
329  "scalar",
330  "specify the boundary-layer thickness"
331  );
332  argList::addOption
333  (
334  "Cbl",
335  "scalar",
336  "boundary-layer thickness as Cbl * mean distance to wall"
337  );
338 
339  #include "setRootCase.H"
340 
341  if (!args.optionFound("ybl") && !args.optionFound("Cbl"))
342  {
344  << "Neither option 'ybl' or 'Cbl' have been provided to calculate "
345  << "the boundary-layer thickness.\n"
346  << "Please choose either 'ybl' OR 'Cbl'."
347  << exit(FatalError);
348  }
349  else if (args.optionFound("ybl") && args.optionFound("Cbl"))
350  {
352  << "Both 'ybl' and 'Cbl' have been provided to calculate "
353  << "the boundary-layer thickness.\n"
354  << "Please choose either 'ybl' OR 'Cbl'."
355  << exit(FatalError);
356  }
357 
358  #include "createTime.H"
359  #include "createNamedMesh.H"
360  #include "createFields.H"
361 
362  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
363 
364  // Modify velocity by applying a 1/7th power law boundary-layer
365  // u/U0 = (y/ybl)^(1/7)
366  // assumes U0 is the same as the current cell velocity
367 
368  Info<< "Setting boundary layer velocity" << nl << endl;
369  scalar yblv = ybl.value();
370  forAll(U, cellI)
371  {
372  if (y[cellI] <= yblv)
373  {
374  mask[cellI] = 1;
375  U[cellI] *= ::pow(y[cellI]/yblv, (1.0/7.0));
376  }
377  }
378  mask.correctBoundaryConditions();
379 
380  Info<< "Writing U\n" << endl;
381  U.write();
382 
383 
384  if
385  (
386  IOobject
387  (
389  runTime.constant(),
390  mesh
391  ).headerOk()
392  )
393  {
394  calcCompressible(mesh, mask, U, y, ybl);
395  }
396  else
397  {
398  calcIncompressible(mesh, mask, U, y, ybl);
399  }
400 
401 
402  Info<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
403  << " ClockTime = " << runTime.elapsedClockTime() << " s"
404  << nl << endl;
405 
406  Info<< "End\n" << endl;
407 
408  return 0;
409 }
410 
411 
412 // ************************************************************************* //
wallDist.H
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:82
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
pThermo
Info<< "Reading thermophysical properties\n"<< endl;autoPtr< psiThermo > pThermo(psiThermo::New(mesh))
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
processorFvPatchField.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
singlePhaseTransportModel.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
turbulentTransportModel.H
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
U
U
Definition: pEqn.H:46
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:131
dictName
const word dictName("particleTrackDict")
nut
nut
Definition: createTDFields.H:71
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
addRegionOption.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
setField
surfacesMesh setField(triSurfaceToAgglom)
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::FatalError
error FatalError
createNamedMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
compressibleCreatePhi.H
Creates and initialises the face-flux field phi.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
rho
rho
Definition: pEqn.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
setRootCase.H
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
epsilon
epsilon
Definition: createTDFields.H:56
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::isA
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
fvCFD.H
args
Foam::argList args(argc, argv)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
turbulentFluidThermoModel.H
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104