evaluateNearWall.H
Go to the documentation of this file.
1 {
2  // Evaluate near-wall behaviour
3 
4  scalar nu = turbulence->nu()().boundaryField()[patchId][faceId];
5  scalar nut = turbulence->nut()().boundaryField()[patchId][faceId];
6  symmTensor R = turbulence->devReff()().boundaryField()[patchId][faceId];
7  scalar epsilon = turbulence->epsilon()()[cellId];
8 // scalar omega = turbulence->omega()()[cellId];
9  scalar k = turbulence->k()()[cellId];
10  scalar magUp = mag(U[cellId] - U.boundaryField()[patchId][faceId]);
11 
13 
14  scalar uTau = ::sqrt(mag(tauw));
15 
16  scalar yPlus = uTau*y[cellId]/(nu + ROOTVSMALL);
17 
18  scalar uPlus = magUp/(uTau + ROOTVSMALL);
19 
20  scalar nutPlus = nut/nu;
21 
22  scalar kPlus = k/(sqr(uTau) + ROOTVSMALL);
23 
24  scalar epsilonPlus = epsilon*nu/(pow4(uTau) + ROOTVSMALL);
25 
26 // scalar omegaPlus = omega*nu/(sqr(uTau) + ROOTVSMALL);
27 
28  scalar Rey = magUp*y[cellId]/nu;
29 
30  Info<< "Rey = " << Rey << ", uTau = " << uTau << ", nut+ = " << nutPlus
31  << ", y+ = " << yPlus << ", u+ = " << uPlus
32  << ", k+ = " << kPlus << ", epsilon+ = " << epsilonPlus
33  << endl;
34 }
nutPlus
scalar nutPlus
Definition: evaluateNearWall.H:20
turbulence
Info<< "Reading field U\n"<< 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("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
epsilonPlus
scalar epsilonPlus
Definition: evaluateNearWall.H:24
uPlus
scalar uPlus
Definition: evaluateNearWall.H:18
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
magUp
scalar magUp
Definition: evaluateNearWall.H:10
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
R
symmTensor R
Universal gas constant: default SI units: [J/mol/K].
Definition: evaluateNearWall.H:6
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
faceId
label faceId(-1)
kPlus
scalar kPlus
Definition: evaluateNearWall.H:22
flowDirection
vector flowDirection
Definition: createFields.H:41
Foam::symmTensor
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:59
tauw
scalar tauw
Definition: evaluateNearWall.H:12
uTau
scalar uTau
Definition: evaluateNearWall.H:14
cellId
label cellId
Definition: interrogateWallPatches.H:67
U
U
Definition: pEqn.H:72
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Rey
scalar Rey
Definition: evaluateNearWall.H:28
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
patchId
label patchId(-1)
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7
nut
scalar nut
Definition: evaluateNearWall.H:5
k
scalar k
Boltzmann constant.
Definition: evaluateNearWall.H:9
wallNormal
vector wallNormal(Zero)
y
scalar y
Definition: LISASMDCalcMethod1.H:14