Backport of wallShearStress to allow compresible flow

This commit is contained in:
Dominik Christ 2014-10-02 17:23:22 +01:00
parent 3a22faaa90
commit 494d070d56
2 changed files with 141 additions and 15 deletions

View file

@ -1,10 +1,15 @@
EXE_INC = \ EXE_INC = \
-I$(LIB_SRC)/transportModels \ -I$(LIB_SRC)/transportModels \
-I$(LIB_SRC)/turbulenceModels \ -I$(LIB_SRC)/turbulenceModels \
-I$(LIB_SRC)/turbulenceModels/incompressible/RAS/RASModel \ -I$(LIB_SRC)/turbulenceModels/incompressible/RAS/derivedFvPatchFields/wallFunctions/nutWallFunctions \
-I$(LIB_SRC)/turbulenceModels/compressible/RAS/derivedFvPatchFields/wallFunctions/mutWallFunctions \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
-I$(LIB_SRC)/finiteVolume/lnInclude -I$(LIB_SRC)/finiteVolume/lnInclude
EXE_LIBS = \ EXE_LIBS = \
-lincompressibleRASModels \ -lincompressibleRASModels \
-lincompressibleTransportModels \ -lincompressibleTransportModels \
-lbasicThermophysicalModels \
-lspecie \
-lcompressibleRASModels \
-lfiniteVolume -lfiniteVolume

View file

@ -25,23 +25,126 @@ Application
wallShearStress wallShearStress
Description Description
Calculates and writes the wall shear stress, for the specified times. Calculates and reports wall shear stress for all patches, for the
specified times when using RAS turbulence models.
Default behaviour assumes operating in incompressible mode.
Use the -compressible option for compressible RAS cases.
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "fvCFD.H" #include "fvCFD.H"
#include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H" #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
#include "RASModel.H" #include "incompressible/RAS/RASModel/RASModel.H"
#include "basicPsiThermo.H"
#include "compressible/RAS/RASModel/RASModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
void calcIncompressible
(
const fvMesh& mesh,
const Time& runTime,
const volVectorField& U,
volVectorField& wallShearStress
)
{
#include "createPhi.H"
singlePhaseTransportModel laminarTransport(U, phi);
autoPtr<incompressible::RASModel> model
(
incompressible::RASModel::New(U, phi, laminarTransport)
);
const volSymmTensorField Reff(model->devReff());
forAll(wallShearStress.boundaryField(), patchI)
{
wallShearStress.boundaryField()[patchI] =
(
-mesh.Sf().boundaryField()[patchI]
/mesh.magSf().boundaryField()[patchI]
) & Reff.boundaryField()[patchI];
}
}
void calcCompressible
(
const fvMesh& mesh,
const Time& runTime,
const volVectorField& U,
volVectorField& wallShearStress
)
{
IOobject rhoHeader
(
"rho",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
);
if (!rhoHeader.headerOk())
{
Info<< " no rho field" << endl;
return;
}
Info<< "Reading field rho\n" << endl;
volScalarField rho(rhoHeader, mesh);
#include "compressibleCreatePhi.H"
autoPtr<basicPsiThermo> pThermo
(
basicPsiThermo::New(mesh)
);
basicPsiThermo& thermo = pThermo();
autoPtr<compressible::RASModel> model
(
compressible::RASModel::New
(
rho,
U,
phi,
thermo
)
);
const volSymmTensorField Reff(model->devRhoReff());
forAll(wallShearStress.boundaryField(), patchI)
{
wallShearStress.boundaryField()[patchI] =
(
-mesh.Sf().boundaryField()[patchI]
/mesh.magSf().boundaryField()[patchI]
) & Reff.boundaryField()[patchI];
}
}
int main(int argc, char *argv[]) int main(int argc, char *argv[])
{ {
timeSelector::addOptions(); timeSelector::addOptions();
#include "addRegionOption.H"
argList::validOptions.insert("compressible","");
#include "setRootCase.H" #include "setRootCase.H"
#include "createTime.H" #include "createTime.H"
instantList timeDirs = timeSelector::select0(runTime, args); instantList timeDirs = timeSelector::select0(runTime, args);
#include "createMesh.H" #include "createNamedMesh.H"
bool compressible = args.optionFound("compressible");
forAll(timeDirs, timeI) forAll(timeDirs, timeI)
{ {
@ -49,10 +152,6 @@ int main(int argc, char *argv[])
Info<< "Time = " << runTime.timeName() << endl; Info<< "Time = " << runTime.timeName() << endl;
mesh.readUpdate(); mesh.readUpdate();
#include "createFields.H"
volSymmTensorField Reff(RASModel->devReff());
volVectorField wallShearStress volVectorField wallShearStress
( (
IOobject IOobject
@ -67,19 +166,41 @@ int main(int argc, char *argv[])
dimensionedVector dimensionedVector
( (
"wallShearStress", "wallShearStress",
Reff.dimensions(), sqr(dimLength)/sqr(dimTime),
vector::zero vector::zero
) )
); );
forAll(wallShearStress.boundaryField(), patchi) IOobject UHeader
(
"U",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
);
if (UHeader.headerOk())
{ {
wallShearStress.boundaryField()[patchi] = Info<< "Reading field U\n" << endl;
( volVectorField U(UHeader, mesh);
-mesh.Sf().boundaryField()[patchi]
/mesh.magSf().boundaryField()[patchi] if (compressible)
) & Reff.boundaryField()[patchi]; {
calcCompressible(mesh, runTime, U, wallShearStress);
}
else
{
calcIncompressible(mesh, runTime, U, wallShearStress);
}
} }
else
{
Info<< " no U field" << endl;
}
Info<< "Writing wall shear stress to field " << wallShearStress.name()
<< nl << endl;
wallShearStress.write(); wallShearStress.write();
} }