ENH: updated buoyantBoussinesqPimpleFoam solver to include thermal wall functions

on (kinematic) turbulent thermal conductivity, kappat and updated tutorial
This commit is contained in:
andy
2010-10-01 17:22:59 +01:00
parent 6ab81bf32a
commit 8f03a3258b
3 changed files with 64 additions and 5 deletions

View File

@ -1,9 +1,8 @@
{
volScalarField kappaEff
(
"kappaEff",
turbulence->nu()/Pr + turbulence->nut()/Prt
);
kappat = turbulence->nut()/Prt;
kappat.correctBoundaryConditions();
volScalarField kappaEff("kappaEff", turbulence->nu()/Pr + kappat);
fvScalarMatrix TEqn
(

View File

@ -52,6 +52,21 @@
incompressible::RASModel::New(U, phi, laminarTransport)
);
// kinematic turbulent thermal thermal conductivity m2/s
Info<< "Reading field kappat\n" << endl;
volScalarField kappat
(
IOobject
(
"kappat",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
label pRefCell = 0;
scalar pRefValue = 0.0;
setRefCell

View File

@ -0,0 +1,45 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class volScalarField;
location "0";
object kappat;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
dimensions [0 2 -1 0 0 0 0];
internalField uniform 0;
boundaryField
{
floor
{
type kappatJayatillekeWallFunction;
Prt 0.85;
value uniform 0;
}
ceiling
{
type kappatJayatillekeWallFunction;
Prt 0.85;
value uniform 0;
}
fixedWalls
{
type kappatJayatillekeWallFunction;
Prt 0.85;
value uniform 0;
}
}
// ************************************************************************* //