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 "fvcFlux.H"
43 
44 using namespace Foam;
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 // Turbulence constants - file-scope
49 static const scalar Cmu(0.09);
50 static const scalar kappa(0.41);
51 
52 
53 int main(int argc, char *argv[])
54 {
56  (
57  "apply a simplified boundary-layer model to the velocity and\n"
58  "turbulence fields based on the 1/7th power-law."
59  );
60 
62  (
63  "ybl",
64  "scalar",
65  "specify the boundary-layer thickness"
66  );
68  (
69  "Cbl",
70  "scalar",
71  "boundary-layer thickness as Cbl * mean distance to wall"
72  );
74  (
75  "writenut",
76  "write nut field"
77  );
78 
79  #include "setRootCase.H"
80 
81  if (!args.optionFound("ybl") && !args.optionFound("Cbl"))
82  {
84  << "Neither option 'ybl' or 'Cbl' have been provided to calculate "
85  << "the boundary-layer thickness.\n"
86  << "Please choose either 'ybl' OR 'Cbl'."
87  << exit(FatalError);
88  }
89  else if (args.optionFound("ybl") && args.optionFound("Cbl"))
90  {
92  << "Both 'ybl' and 'Cbl' have been provided to calculate "
93  << "the boundary-layer thickness.\n"
94  << "Please choose either 'ybl' OR 'Cbl'."
95  << exit(FatalError);
96  }
97 
98  #include "createTime.H"
99  #include "createMeshNoChangers.H"
100  #include "createFields.H"
101 
102  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
103 
104  // Modify velocity by applying a 1/7th power law boundary-layer
105  // u/U0 = (y/ybl)^(1/7)
106  // assumes U0 is the same as the current cell velocity
107 
108  Info<< "Setting boundary layer velocity" << nl << endl;
109  scalar yblv = ybl.value();
110  forAll(U, celli)
111  {
112  if (y[celli] <= yblv)
113  {
114  mask[celli] = 1;
115  U[celli] *= ::pow(y[celli]/yblv, (1.0/7.0));
116  }
117  }
118  mask.correctBoundaryConditions();
119 
120  Info<< "Writing U\n" << endl;
121  U.write();
122 
123  // Update/re-write phi
124  #include "createPhi.H"
125  phi.write();
126 
128 
130  (
132  );
133 
134  if (isA<incompressible::RASModel>(turbulence()))
135  {
136  // Calculate nut
137  turbulence->validate();
138  tmp<volScalarField> tnut = turbulence->nut();
139  volScalarField& nut = const_cast<volScalarField&>(tnut());
141  nut = (1 - mask)*nut + mask*sqr(kappa*min(y, ybl))*::sqrt(2)*S;
142 
143  // Do not correct BC - wall functions will 'undo' manipulation above
144  // by using nut from turbulence model
145 
146  if (args.optionFound("writenut"))
147  {
148  Info<< "Writing nut" << endl;
149  nut.write();
150  }
151 
152 
153  //--- Read and modify turbulence fields
154 
155  // Turbulence k
157  volScalarField& k = const_cast<volScalarField&>(tk());
158  scalar ck0 = pow025(Cmu)*kappa;
159  k = (1 - mask)*k + mask*sqr(nut/(ck0*min(y, ybl)));
160  k.correctBoundaryConditions();
161 
162  Info<< "Writing k\n" << endl;
163  k.write();
164 
165 
166  // Turbulence epsilon
167  tmp<volScalarField> tepsilon = turbulence->epsilon();
168  volScalarField& epsilon = const_cast<volScalarField&>(tepsilon());
169  scalar ce0 = ::pow(Cmu, 0.75)/kappa;
170  epsilon = (1 - mask)*epsilon + mask*ce0*k*sqrt(k)/min(y, ybl);
171 
172  // Do not correct BC - G set by the wall-functions is not available
173  // epsilon.correctBoundaryConditions();
174 
175  Info<< "Writing epsilon\n" << endl;
176  epsilon.write();
177 
178  // Turbulence omega
179  typeIOobject<volScalarField> omegaHeader
180  (
181  "omega",
182  runTime.name(),
183  mesh,
186  false
187  );
188 
189  if (omegaHeader.headerOk())
190  {
191  volScalarField omega(omegaHeader, mesh);
192 
193  omega = (1 - mask)*omega + mask*ce0*sqrt(k)/(Cmu*min(y, ybl));
194 
195  // Do not correct BC - G set by the wall-functions is not available
196  // omega.correctBoundaryConditions();
197 
198  Info<< "Writing omega\n" << endl;
199  omega.write();
200  }
201 
202  // Turbulence nuTilda
203  typeIOobject<volScalarField> nuTildaHeader
204  (
205  "nuTilda",
206  runTime.name(),
207  mesh,
210  false
211  );
212 
213  if (nuTildaHeader.headerOk())
214  {
215  volScalarField nuTilda(nuTildaHeader, mesh);
216  nuTilda = nut;
217 
218  // Do not correct BC - G set by the wall-functions is not available
219  // nuTilda.correctBoundaryConditions();
220 
221  Info<< "Writing nuTilda\n" << endl;
222  nuTilda.write();
223  }
224  }
225 
226  Info<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
227  << " ClockTime = " << runTime.elapsedClockTime() << " s"
228  << nl << endl;
229 
230  Info<< "End\n" << endl;
231 
232  return 0;
233 }
234 
235 
236 // ************************************************************************* //
scalar y
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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.
Convenience class to handle the input of constant rotational speed. Reads an omega entry with default...
Definition: omega.H:54
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:334
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 > &)
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:257
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 > &)
error FatalError
static const char nl
Definition: Ostream.H:266
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))