applyBoundaryLayer.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
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 \*---------------------------------------------------------------------------*/
37 
38 #include "argList.H"
39 #include "viscosityModel.H"
41 #include "wallDist.H"
42 #include "bound.H"
43 #include "fvcFlux.H"
44 
45 using namespace Foam;
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 // Turbulence constants - file-scope
50 static const scalar Cmu(0.09);
51 static const scalar kappa(0.41);
52 
53 
54 int main(int argc, char *argv[])
55 {
57  (
58  "apply a simplified boundary-layer model to the velocity and\n"
59  "turbulence fields based on the 1/7th power-law."
60  );
61 
63  (
64  "ybl",
65  "scalar",
66  "specify the boundary-layer thickness"
67  );
69  (
70  "Cbl",
71  "scalar",
72  "boundary-layer thickness as Cbl * mean distance to wall"
73  );
75  (
76  "writenut",
77  "write nut field"
78  );
79 
80  #include "setRootCase.H"
81 
82  if (!args.optionFound("ybl") && !args.optionFound("Cbl"))
83  {
85  << "Neither option 'ybl' or 'Cbl' have been provided to calculate "
86  << "the boundary-layer thickness.\n"
87  << "Please choose either 'ybl' OR 'Cbl'."
88  << exit(FatalError);
89  }
90  else if (args.optionFound("ybl") && args.optionFound("Cbl"))
91  {
93  << "Both 'ybl' and 'Cbl' have been provided to calculate "
94  << "the boundary-layer thickness.\n"
95  << "Please choose either 'ybl' OR 'Cbl'."
96  << exit(FatalError);
97  }
98 
99  #include "createTime.H"
100  #include "createMeshNoChangers.H"
101  #include "createFields.H"
102 
103  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104 
105  // Modify velocity by applying a 1/7th power law boundary-layer
106  // u/U0 = (y/ybl)^(1/7)
107  // assumes U0 is the same as the current cell velocity
108 
109  Info<< "Setting boundary layer velocity" << nl << endl;
110  scalar yblv = ybl.value();
111  forAll(U, celli)
112  {
113  if (y[celli] <= yblv)
114  {
115  mask[celli] = 1;
116  U[celli] *= ::pow(y[celli]/yblv, (1.0/7.0));
117  }
118  }
119  mask.correctBoundaryConditions();
120 
121  Info<< "Writing U\n" << endl;
122  U.write();
123 
124  // Update/re-write phi
125  #include "createPhi.H"
126  phi.write();
127 
129 
131  (
133  );
134 
135  if (isA<incompressible::RASModel>(turbulence()))
136  {
137  // Calculate nut
138  turbulence->validate();
139  tmp<volScalarField> tnut = turbulence->nut();
140  volScalarField& nut = const_cast<volScalarField&>(tnut());
142  nut = (1 - mask)*nut + mask*sqr(kappa*min(y, ybl))*::sqrt(2)*S;
143 
144  // Do not correct BC - wall functions will 'undo' manipulation above
145  // by using nut from turbulence model
146 
147  if (args.optionFound("writenut"))
148  {
149  Info<< "Writing nut" << endl;
150  nut.write();
151  }
152 
153 
154  //--- Read and modify turbulence fields
155 
156  // Turbulence k
158  volScalarField& k = const_cast<volScalarField&>(tk());
159  scalar ck0 = pow025(Cmu)*kappa;
160  k = (1 - mask)*k + mask*sqr(nut/(ck0*min(y, ybl)));
161  k.correctBoundaryConditions();
162 
163  Info<< "Writing k\n" << endl;
164  k.write();
165 
166 
167  // Turbulence epsilon
168  tmp<volScalarField> tepsilon = turbulence->epsilon();
169  volScalarField& epsilon = const_cast<volScalarField&>(tepsilon());
170  scalar ce0 = ::pow(Cmu, 0.75)/kappa;
171  epsilon = (1 - mask)*epsilon + mask*ce0*k*sqrt(k)/min(y, ybl);
172 
173  // Do not correct BC - G set by the wall-functions is not available
174  // epsilon.correctBoundaryConditions();
175 
176  Info<< "Writing epsilon\n" << endl;
177  epsilon.write();
178 
179  // Turbulence omega
180  typeIOobject<volScalarField> omegaHeader
181  (
182  "omega",
183  runTime.name(),
184  mesh,
187  false
188  );
189 
190  if (omegaHeader.headerOk())
191  {
192  volScalarField omega(omegaHeader, mesh);
193 
194  const incompressible::RASModel& rasModel =
195  refCast<const incompressible::RASModel>(turbulence());
196 
197  omega = (1 - mask)*omega + mask*ce0*sqrt(k)/(Cmu*min(y, ybl));
198  bound(omega, rasModel.omegaMin());
199 
200  // Do not correct BC - G set by the wall-functions is not available
201  // omega.correctBoundaryConditions();
202 
203  Info<< "Writing omega\n" << endl;
204  omega.write();
205  }
206 
207  // Turbulence nuTilda
208  typeIOobject<volScalarField> nuTildaHeader
209  (
210  "nuTilda",
211  runTime.name(),
212  mesh,
215  false
216  );
217 
218  if (nuTildaHeader.headerOk())
219  {
220  volScalarField nuTilda(nuTildaHeader, mesh);
221  nuTilda = nut;
222 
223  // Do not correct BC - G set by the wall-functions is not available
224  // nuTilda.correctBoundaryConditions();
225 
226  Info<< "Writing nuTilda\n" << endl;
227  nuTilda.write();
228  }
229  }
230 
231  Info<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
232  << " ClockTime = " << runTime.elapsedClockTime() << " s"
233  << nl << endl;
234 
235  Info<< "End\n" << endl;
236 
237  return 0;
238 }
239 
240 
241 // ************************************************************************* //
scalar y
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Bound the given scalar field where it is below the specified minimum.
Generic GeometricField class.
virtual Ostream & write(const char)=0
Write character.
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
static autoPtr< incompressibleMomentumTransportModel > New(const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Return a reference to the selected turbulence model.
A class for managing temporary objects.
Definition: tmp.H:55
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
static autoPtr< viscosityModel > New(const fvMesh &mesh, const word &group=word::null)
Return a reference to the selected viscosity model.
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const scalar omega
const scalar nut
const scalar epsilon
int main(int argc, char *argv[])
Definition: financialFoam.C:44
Calculate the face-flux of the given field.
U
Definition: pEqn.H:72
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
RASModel< momentumTransportModel > RASModel
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
error FatalError
static const char nl
Definition: Ostream.H:260
dimensionedScalar pow025(const dimensionedScalar &ds)
Foam::argList args(argc, argv)
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.name(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating face flux\n"<< endl;surfaceScalarField phi(IOobject("phi", runTime.name(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar(mesh.Sf().dimensions() *U.dimensions(), 0));autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))