Go to the documentation of this file.
52 static const scalar Cmu(0.09);
53 static const scalar
kappa(0.41);
57 if (!Pstream::parRun())
66 volScalarField::GeometricBoundaryField& bf = vf.boundaryField();
70 if (
isA<processorFvPatchField<scalar> >(bf[patchI]))
72 bf[patchI].initEvaluate();
78 if (
isA<processorFvPatchField<scalar> >(bf[patchI]))
80 bf[patchI].evaluate();
86 template<
class TurbulenceModel>
107 correctProcessorPatches(
k);
116 template<
class TurbulenceModel>
129 tmp<volScalarField> tepsilon =
turbulence->epsilon();
138 correctProcessorPatches(
epsilon);
159 mesh.time().timeName(),
166 if (omegaHeader.headerOk())
171 omega = (1 - mask)*omega + mask*
epsilon/(Cmu*
k + k0);
178 correctProcessorPatches(omega);
189 const word& fieldName,
196 mesh.time().timeName(),
203 if (fldHeader.headerOk())
212 correctProcessorPatches(
fld);
220 void calcCompressible
229 const Time& runTime =
mesh.time();
239 autoPtr<compressible::turbulenceModel>
turbulence
260 correctProcessorPatches(
nut);
263 tmp<volScalarField>
k =
272 void calcIncompressible
281 const Time& runTime =
mesh.time();
284 #include "createPhi.H"
289 autoPtr<incompressible::turbulenceModel>
turbulence
304 correctProcessorPatches(
nut);
307 tmp<volScalarField>
k =
316 int main(
int argc,
char *argv[])
320 "apply a simplified boundary-layer model to the velocity and\n"
321 "turbulence fields based on the 1/7th power-law."
330 "specify the boundary-layer thickness"
336 "boundary-layer thickness as Cbl * mean distance to wall"
344 <<
"Neither option 'ybl' or 'Cbl' have been provided to calculate "
345 <<
"the boundary-layer thickness.\n"
346 <<
"Please choose either 'ybl' OR 'Cbl'."
352 <<
"Both 'ybl' and 'Cbl' have been provided to calculate "
353 <<
"the boundary-layer thickness.\n"
354 <<
"Please choose either 'ybl' OR 'Cbl'."
360 #include "createFields.H"
368 Info<<
"Setting boundary layer velocity" <<
nl <<
endl;
369 scalar yblv = ybl.value();
372 if (
y[cellI] <= yblv)
375 U[cellI] *=
::pow(
y[cellI]/yblv, (1.0/7.0));
378 mask.correctBoundaryConditions();
394 calcCompressible(
mesh, mask,
U,
y, ybl);
398 calcIncompressible(
mesh, mask,
U,
y, ybl);
402 Info<<
nl <<
"ExecutionTime = " << runTime.elapsedCpuTime() <<
" s"
403 <<
" ClockTime = " << runTime.elapsedClockTime() <<
" s"
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Info<< "Reading thermophysical properties\n"<< endl;autoPtr< psiThermo > pThermo(psiThermo::New(mesh))
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
dimensionedScalar pow025(const dimensionedScalar &ds)
const word dictName("particleTrackDict")
autoPtr< compressible::turbulenceModel > turbulence
virtual Ostream & write(const token &)=0
Write next token to stream.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
surfacesMesh setField(triSurfaceToAgglom)
int main(int argc, char *argv[])
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Creates and initialises the face-flux field phi.
errorManipArg< error, int > exit(error &err, const int errNo=1)
GeometricField< vector, fvPatchField, volMesh > volVectorField
singlePhaseTransportModel laminarTransport(U, phi)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
bool optionFound(const word &opt) const
Return true if the named option is found.
Foam::argList args(argc, argv)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)