All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2021 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 "fvCFD.H"
41 #include "wallDist.H"
42 #include "bound.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 // Turbulence constants - file-scope
47 static const scalar Cmu(0.09);
48 static const scalar kappa(0.41);
49 
50 
51 int main(int argc, char *argv[])
52 {
53  argList::addNote
54  (
55  "apply a simplified boundary-layer model to the velocity and\n"
56  "turbulence fields based on the 1/7th power-law."
57  );
58 
59  argList::addOption
60  (
61  "ybl",
62  "scalar",
63  "specify the boundary-layer thickness"
64  );
65  argList::addOption
66  (
67  "Cbl",
68  "scalar",
69  "boundary-layer thickness as Cbl * mean distance to wall"
70  );
71  argList::addBoolOption
72  (
73  "writenut",
74  "write nut field"
75  );
76 
77  #include "setRootCase.H"
78 
79  if (!args.optionFound("ybl") && !args.optionFound("Cbl"))
80  {
82  << "Neither option 'ybl' or 'Cbl' have been provided to calculate "
83  << "the boundary-layer thickness.\n"
84  << "Please choose either 'ybl' OR 'Cbl'."
85  << exit(FatalError);
86  }
87  else if (args.optionFound("ybl") && args.optionFound("Cbl"))
88  {
90  << "Both 'ybl' and 'Cbl' have been provided to calculate "
91  << "the boundary-layer thickness.\n"
92  << "Please choose either 'ybl' OR 'Cbl'."
93  << exit(FatalError);
94  }
95 
96  #include "createTime.H"
97  #include "createMesh.H"
98  #include "createFields.H"
99 
100  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
101 
102  // Modify velocity by applying a 1/7th power law boundary-layer
103  // u/U0 = (y/ybl)^(1/7)
104  // assumes U0 is the same as the current cell velocity
105 
106  Info<< "Setting boundary layer velocity" << nl << endl;
107  scalar yblv = ybl.value();
108  forAll(U, celli)
109  {
110  if (y[celli] <= yblv)
111  {
112  mask[celli] = 1;
113  U[celli] *= ::pow(y[celli]/yblv, (1.0/7.0));
114  }
115  }
116  mask.correctBoundaryConditions();
117 
118  Info<< "Writing U\n" << endl;
119  U.write();
120 
121  // Update/re-write phi
122  #include "createPhi.H"
123  phi.write();
124 
125  singlePhaseTransportModel laminarTransport(U, phi);
126 
127  autoPtr<incompressible::momentumTransportModel> turbulence
128  (
130  );
131 
132  if (isA<incompressible::RASModel>(turbulence()))
133  {
134  // Calculate nut
135  turbulence->validate();
136  tmp<volScalarField> tnut = turbulence->nut();
137  volScalarField& nut = const_cast<volScalarField&>(tnut());
139  nut = (1 - mask)*nut + mask*sqr(kappa*min(y, ybl))*::sqrt(2)*S;
140 
141  // Do not correct BC - wall functions will 'undo' manipulation above
142  // by using nut from turbulence model
143 
144  if (args.optionFound("writenut"))
145  {
146  Info<< "Writing nut" << endl;
147  nut.write();
148  }
149 
150 
151  //--- Read and modify turbulence fields
152 
153  // Turbulence k
154  tmp<volScalarField> tk = turbulence->k();
155  volScalarField& k = const_cast<volScalarField&>(tk());
156  scalar ck0 = pow025(Cmu)*kappa;
157  k = (1 - mask)*k + mask*sqr(nut/(ck0*min(y, ybl)));
158  k.correctBoundaryConditions();
159 
160  Info<< "Writing k\n" << endl;
161  k.write();
162 
163 
164  // Turbulence epsilon
165  tmp<volScalarField> tepsilon = turbulence->epsilon();
166  volScalarField& epsilon = const_cast<volScalarField&>(tepsilon());
167  scalar ce0 = ::pow(Cmu, 0.75)/kappa;
168  epsilon = (1 - mask)*epsilon + mask*ce0*k*sqrt(k)/min(y, ybl);
169 
170  // Do not correct BC - G set by the wall-functions is not available
171  // epsilon.correctBoundaryConditions();
172 
173  Info<< "Writing epsilon\n" << endl;
174  epsilon.write();
175 
176  // Turbulence omega
177  IOobject omegaHeader
178  (
179  "omega",
180  runTime.timeName(),
181  mesh,
182  IOobject::MUST_READ,
183  IOobject::NO_WRITE,
184  false
185  );
186 
187  if (omegaHeader.typeHeaderOk<volScalarField>(true))
188  {
189  volScalarField omega(omegaHeader, mesh);
190 
191  const incompressible::RASModel& rasModel =
192  refCast<const incompressible::RASModel>(turbulence());
193 
194  omega = (1 - mask)*omega + mask*ce0*sqrt(k)/(Cmu*min(y, ybl));
195  bound(omega, rasModel.omegaMin());
196 
197  // Do not correct BC - G set by the wall-functions is not available
198  // omega.correctBoundaryConditions();
199 
200  Info<< "Writing omega\n" << endl;
201  omega.write();
202  }
203 
204  // Turbulence nuTilda
205  IOobject nuTildaHeader
206  (
207  "nuTilda",
208  runTime.timeName(),
209  mesh,
210  IOobject::MUST_READ,
211  IOobject::NO_WRITE,
212  false
213  );
214 
215  if (nuTildaHeader.typeHeaderOk<volScalarField>(true))
216  {
217  volScalarField nuTilda(nuTildaHeader, mesh);
218  nuTilda = nut;
219 
220  // Do not correct BC - G set by the wall-functions is not available
221  // nuTilda.correctBoundaryConditions();
222 
223  Info<< "Writing nuTilda\n" << endl;
224  nuTilda.write();
225  }
226  }
227 
228  Info<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
229  << " ClockTime = " << runTime.elapsedClockTime() << " s"
230  << nl << endl;
231 
232  Info<< "End\n" << endl;
233 
234  return 0;
235 }
236 
237 
238 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
dimensionedSymmTensor sqr(const dimensionedVector &dv)
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensionedScalar pow025(const dimensionedScalar &ds)
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
scalar y
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
dynamicFvMesh & mesh
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const GeometricField< Type, fvPatchField, volMesh > &)
phi
Definition: correctPhi.H:3
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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
static const char nl
Definition: Ostream.H:260
Bound the given scalar field if it has gone unbounded.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::argList args(argc, argv)
singlePhaseTransportModel laminarTransport(U, phi)
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
scalar nut