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-2016 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 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 // turbulence constants - file-scope
46 static const scalar Cmu(0.09);
47 static const scalar kappa(0.41);
48 
49 
50 int main(int argc, char *argv[])
51 {
52  argList::addNote
53  (
54  "apply a simplified boundary-layer model to the velocity and\n"
55  "turbulence fields based on the 1/7th power-law."
56  );
57 
58  argList::addOption
59  (
60  "ybl",
61  "scalar",
62  "specify the boundary-layer thickness"
63  );
64  argList::addOption
65  (
66  "Cbl",
67  "scalar",
68  "boundary-layer thickness as Cbl * mean distance to wall"
69  );
70  argList::addBoolOption
71  (
72  "writenut",
73  "write nut field"
74  );
75 
76  #include "setRootCase.H"
77 
78  if (!args.optionFound("ybl") && !args.optionFound("Cbl"))
79  {
81  << "Neither option 'ybl' or 'Cbl' have been provided to calculate "
82  << "the boundary-layer thickness.\n"
83  << "Please choose either 'ybl' OR 'Cbl'."
84  << exit(FatalError);
85  }
86  else if (args.optionFound("ybl") && args.optionFound("Cbl"))
87  {
89  << "Both 'ybl' and 'Cbl' have been provided to calculate "
90  << "the boundary-layer thickness.\n"
91  << "Please choose either 'ybl' OR 'Cbl'."
92  << exit(FatalError);
93  }
94 
95  #include "createTime.H"
96  #include "createMesh.H"
97  #include "createFields.H"
98 
99  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
100 
101  // Modify velocity by applying a 1/7th power law boundary-layer
102  // u/U0 = (y/ybl)^(1/7)
103  // assumes U0 is the same as the current cell velocity
104 
105  Info<< "Setting boundary layer velocity" << nl << endl;
106  scalar yblv = ybl.value();
107  forAll(U, celli)
108  {
109  if (y[celli] <= yblv)
110  {
111  mask[celli] = 1;
112  U[celli] *= ::pow(y[celli]/yblv, (1.0/7.0));
113  }
114  }
115  mask.correctBoundaryConditions();
116 
117  Info<< "Writing U\n" << endl;
118  U.write();
119 
120  // Update/re-write phi
121  #include "createPhi.H"
122  phi.write();
123 
124  singlePhaseTransportModel laminarTransport(U, phi);
125 
126  autoPtr<incompressible::turbulenceModel> turbulence
127  (
129  );
130 
131  if (isA<incompressible::RASModel>(turbulence()))
132  {
133  // Calculate nut - reference nut is calculated by the turbulence model
134  // on its construction
135  tmp<volScalarField> tnut = turbulence->nut();
136  volScalarField& nut = tnut.ref();
138  nut = (1 - mask)*nut + mask*sqr(kappa*min(y, ybl))*::sqrt(2)*S;
139 
140  // do not correct BC - wall functions will 'undo' manipulation above
141  // by using nut from turbulence model
142 
143  if (args.optionFound("writenut"))
144  {
145  Info<< "Writing nut" << endl;
146  nut.write();
147  }
148 
149 
150  //--- Read and modify turbulence fields
151 
152  // Turbulence k
153  tmp<volScalarField> tk = turbulence->k();
154  volScalarField& k = tk.ref();
155  scalar ck0 = pow025(Cmu)*kappa;
156  k = (1 - mask)*k + mask*sqr(nut/(ck0*min(y, ybl)));
157 
158  // do not correct BC - operation may use inconsistent fields wrt these
159  // local manipulations
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 = tepsilon.ref();
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 - wall functions will use non-updated k from
173  // turbulence model
174  // epsilon.correctBoundaryConditions();
175 
176  Info<< "Writing epsilon\n" << endl;
177  epsilon.write();
178 
179  // Turbulence omega
180  IOobject omegaHeader
181  (
182  "omega",
183  runTime.timeName(),
184  mesh,
185  IOobject::MUST_READ,
186  IOobject::NO_WRITE,
187  false
188  );
189 
190  if (omegaHeader.headerOk())
191  {
192  volScalarField omega(omegaHeader, mesh);
193  dimensionedScalar k0("VSMALL", k.dimensions(), VSMALL);
194  omega = (1 - mask)*omega + mask*epsilon/(Cmu*k + k0);
195 
196  // do not correct BC - wall functions will use non-updated k from
197  // turbulence model
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.headerOk())
216  {
217  volScalarField nuTilda(nuTildaHeader, mesh);
218  nuTilda = nut;
219 
220  // do not correct BC
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 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
surfaceScalarField & phi
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:428
U
Definition: pEqn.H:83
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:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
dimensionedScalar pow025(const dimensionedScalar &ds)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
scalar y
singlePhaseTransportModel laminarTransport(U, phi)
dynamicFvMesh & mesh
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
static const char nl
Definition: Ostream.H:262
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::argList args(argc, argv)
scalar nut