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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 Application
28  applyBoundaryLayer
29 
30 Group
31  grpPreProcessingUtilities
32 
33 Description
34  Apply a simplified boundary-layer model to the velocity and
35  turbulence fields based on the 1/7th power-law.
36 
37  The uniform boundary-layer thickness is either provided via the -ybl option
38  or calculated as the average of the distance to the wall scaled with
39  the thickness coefficient supplied via the option -Cbl. If both options
40  are provided -ybl is used.
41 
42  Compressible modes is automatically selected based on the existence of the
43  "thermophysicalProperties" dictionary required to construct the
44  thermodynamics package.
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #include "fvCFD.H"
52 #include "wallDist.H"
53 #include "processorFvPatch.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 // Turbulence constants - file-scope
58 static const scalar Cmu(0.09);
59 static const scalar kappa(0.41);
60 
61 
62 template<class Type>
63 void correctProcessorPatches(GeometricField<Type, fvPatchField, volMesh>& fld)
64 {
65  if (UPstream::parRun())
66  {
67  fld.boundaryFieldRef().template evaluateCoupled<processorFvPatch>();
68  }
69 }
70 
71 
72 void blendField
73 (
74  const word& fieldName,
75  const fvMesh& mesh,
76  const scalarField& mask,
77  const scalarField& boundaryLayerField
78 )
79 {
80  IOobject fieldHeader
81  (
82  fieldName,
83  mesh.time().timeName(),
84  mesh,
85  IOobject::MUST_READ,
86  IOobject::NO_WRITE,
87  IOobject::NO_REGISTER
88  );
89 
90  if (fieldHeader.typeHeaderOk<volScalarField>(true))
91  {
92  volScalarField fld(fieldHeader, mesh);
93  scalarField& pf = fld.primitiveFieldRef();
94  pf = (1 - mask)*pf + mask*boundaryLayerField;
95  fld.clamp_min(SMALL);
96 
97  // Correct the processor patches only.
98  // Do not correct BC
99  // - operation may use inconsistent fields wrt these local
100  // manipulations
101  correctProcessorPatches(fld);
102 
103  Info<< "Writing " << fieldName << nl << endl;
104  fld.write();
105  }
106 }
107 
108 
109 void calcOmegaField
110 (
111  const fvMesh& mesh,
112  const scalarField& mask,
113  const scalarField& kBL,
114  const scalarField& epsilonBL
115 )
116 {
117  // Turbulence omega
118  IOobject omegaHeader
119  (
120  "omega",
121  mesh.time().timeName(),
122  mesh,
123  IOobject::MUST_READ,
124  IOobject::NO_WRITE,
125  IOobject::NO_REGISTER
126  );
127 
128  if (omegaHeader.typeHeaderOk<volScalarField>(true))
129  {
130  volScalarField omega(omegaHeader, mesh);
131  scalarField& pf = omega.primitiveFieldRef();
132 
133  pf = (1 - mask)*pf + mask*epsilonBL/(Cmu*kBL + SMALL);
134  omega.clamp_min(SMALL);
135 
136  // Correct the processor patches only.
137  // Do not correct BC
138  // - operation may use inconsistent fields wrt these local
139  // manipulations
140  correctProcessorPatches(omega);
141 
142  Info<< "Writing omega\n" << endl;
143  omega.write();
144  }
145 }
146 
147 
148 void setField
149 (
150  const fvMesh& mesh,
151  const word& fieldName,
152  const volScalarField& value
153 )
154 {
155  IOobject fldHeader
156  (
157  fieldName,
158  mesh.time().timeName(),
159  mesh,
160  IOobject::MUST_READ,
161  IOobject::NO_WRITE,
162  IOobject::NO_REGISTER
163  );
164 
165  if (fldHeader.typeHeaderOk<volScalarField>(true))
166  {
167  volScalarField fld(fldHeader, mesh);
168  fld = value;
169 
170  // Correct the processor patches only.
171  // Do not correct BC
172  // - operation may use inconsistent fields wrt these local
173  // manipulations
174  correctProcessorPatches(fld);
175 
176  Info<< "Writing " << fieldName << nl << endl;
177  fld.write();
178  }
179 }
180 
181 
182 tmp<volScalarField> calcNut
183 (
184  const fvMesh& mesh,
185  const volVectorField& U
186 )
187 {
188  const Time& runTime = mesh.time();
189 
190  if
191  (
192  IOobject
193  (
195  runTime.constant(),
196  mesh
197  ).typeHeaderOk<IOdictionary>(true)
198  )
199  {
200  // Compressible
201  autoPtr<fluidThermo> pThermo(fluidThermo::New(mesh));
202  fluidThermo& thermo = pThermo();
203  volScalarField rho(thermo.rho());
204 
205  // Update/re-write phi
206  #include "compressibleCreatePhi.H"
207  phi.write();
208 
209  autoPtr<compressible::turbulenceModel> turbulence
210  (
212  (
213  rho,
214  U,
215  phi,
216  thermo
217  )
218  );
219 
220  // Correct nut
221  turbulence->validate();
222 
223  return tmp<volScalarField>::New(turbulence->nut());
224  }
225  else
226  {
227  // Incompressible
228 
229  // Update/re-write phi
230  #include "createPhi.H"
231  phi.write();
232 
233  singlePhaseTransportModel laminarTransport(U, phi);
234 
235  autoPtr<incompressible::turbulenceModel> turbulence
236  (
238  );
239 
240  // Correct nut
241  turbulence->validate();
242 
243  return tmp<volScalarField>::New(turbulence->nut());
244  }
245 }
246 
247 
248 int main(int argc, char *argv[])
249 {
250  argList::addNote
251  (
252  "Apply a simplified boundary-layer model to the velocity and"
253  " turbulence fields based on the 1/7th power-law."
254  );
255 
256  #include "addRegionOption.H"
257 
258  argList::addOption
259  (
260  "ybl",
261  "scalar",
262  "Specify the boundary-layer thickness"
263  );
264  argList::addOption
265  (
266  "Cbl",
267  "scalar",
268  "Boundary-layer thickness as Cbl * mean distance to wall"
269  );
270  argList::addBoolOption
271  (
272  "writeTurbulenceFields", // (until 1906 was write-nut)
273  "Write the turbulence fields"
274  );
275  argList::addOptionCompat
276  (
277  "writeTurbulenceFields", {"write-nut", 1906}
278  );
279 
280  #include "setRootCase.H"
281 
282  if (!args.found("ybl") && !args.found("Cbl"))
283  {
285  << "Neither option 'ybl' or 'Cbl' have been provided to calculate "
286  << "the boundary-layer thickness.\n"
287  << "Please choose either 'ybl' OR 'Cbl'."
288  << exit(FatalError);
289  }
290  else if (args.found("ybl") && args.found("Cbl"))
291  {
293  << "Both 'ybl' and 'Cbl' have been provided to calculate "
294  << "the boundary-layer thickness.\n"
295  << "Please choose either 'ybl' OR 'Cbl'."
296  << exit(FatalError);
297  }
298 
299  const bool writeTurbulenceFields = args.found("writeTurbulenceFields");
300 
301  #include "createTime.H"
302  #include "createNamedMesh.H"
303  #include "createFields.H"
304 
305  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306 
307  // Modify velocity by applying a 1/7th power law boundary-layer
308  // u/U0 = (y/ybl)^(1/7)
309  // assumes U0 is the same as the current cell velocity
310  Info<< "Setting boundary layer velocity" << nl << endl;
311  const scalar yblv = ybl.value();
312  forAll(U, celli)
313  {
314  if ((y[celli] > 0) && (y[celli] <= yblv))
315  {
316  mask[celli] = 1;
317  U[celli] *= ::pow(y[celli]/yblv, (1.0/7.0));
318  }
319  }
320  mask.correctBoundaryConditions();
321  correctProcessorPatches(U);
322 
323  if (writeTurbulenceFields)
324  {
325  // Retrieve nut from turbulence model
326  volScalarField nut(calcNut(mesh, U));
327 
328  // Blend nut using boundary layer profile
329  volScalarField S("S", mag(devSymm(fvc::grad(U))));
330  nut = (1 - mask)*nut + mask*sqr(kappa*min(y, ybl))*::sqrt(2)*S;
331 
332  // Do not correct BC - wall functions will 'undo' manipulation above
333  // by using nut from turbulence model
334  correctProcessorPatches(nut);
335 
336  Info<< "Writing nut\n" << endl;
337  nut.write();
338 
339  // Boundary layer turbulence kinetic energy
340  scalar ck0 = pow025(Cmu)*kappa;
341  scalarField kBL(sqr(nut/(ck0*min(y, ybl))));
342 
343  // Boundary layer turbulence dissipation
344  scalar ce0 = ::pow(Cmu, 0.75)/kappa;
345  scalarField epsilonBL(ce0*kBL*sqrt(kBL)/min(y, ybl));
346 
347  // Process fields if they are present
348  blendField("k", mesh, mask, kBL);
349  blendField("epsilon", mesh, mask, epsilonBL);
350  calcOmegaField(mesh, mask, kBL, epsilonBL);
351  setField(mesh, "nuTilda", nut);
352  }
353 
354  // Write the updated U field
355  Info<< "Writing U\n" << endl;
356  U.write();
357 
358  Info<< nl;
359  runTime.printExecutionTime(Info);
360 
361  Info<< "End\n" << endl;
362 
363  return 0;
364 }
365 
366 
367 // ************************************************************************* //
Info<< "Reading thermophysical properties\"<< endl;autoPtr< psiReactionThermo > pThermo(psiReactionThermo::New(mesh))
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void clamp_min(const Type &lower)
Impose lower (floor) clamp on the field values (in-place)
Definition: Field.C:655
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
const word dictName("faMeshDefinition")
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Required Variables.
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:481
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
scalar y
singlePhaseTransportModel laminarTransport(U, phi)
dynamicFvMesh & mesh
surfacesMesh setField(triSurfaceToAgglom)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::argList args(argc, argv)
Creates and initialises the face-flux field phi.
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
scalar nut
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:51