Pe: Create and write volScalarField in addition to the surfaceScalarField for ease of post-processing

Update handling of turbulence models
This commit is contained in:
Henry
2015-04-25 16:29:54 +01:00
parent 209b3740db
commit 992195c6bc
3 changed files with 49 additions and 121 deletions

View File

@ -1,9 +1,10 @@
EXE_INC = \
-I$(LIB_SRC)/postProcessing/postCalc \
-I$(LIB_SRC)/transportModels \
-I$(LIB_SRC)/transportModels/incompressible/singlePhaseTransportModel \
-I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/incompressible/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/compressible/lnInclude \
-I$(LIB_SRC)/transportModels \
-I$(LIB_SRC)/transportModels/compressible/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
-I$(LIB_SRC)/finiteVolume/lnInclude \
@ -20,4 +21,5 @@ EXE_LIBS = \
-lspecie \
-lfiniteVolume \
-lgenericPatchFields \
-lmeshTools
-lmeshTools \
-lsampling

View File

@ -25,20 +25,21 @@ Application
Pe
Description
Calculates and writes the Pe number as a surfaceScalarField obtained from
field phi.
Calculates the Peclet number Pe from the flux phi and writes the maximum
value, the surfaceScalarField Pef and volScalarField Pe.
The -noWrite option just outputs the max/min values without writing
the field.
With the -noWrite option just outputs the max value without writing
the fields.
\*---------------------------------------------------------------------------*/
#include "calc.H"
#include "fvc.H"
#include "surfaceInterpolate.H"
#include "fvcAverage.H"
#include "singlePhaseTransportModel.H"
#include "turbulentTransportModel.H"
#include "turbulentFluidThermoModel.H"
#include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -73,18 +74,9 @@ void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
mesh
);
IOobject RASPropertiesHeader
IOobject turbulencePropertiesHeader
(
"RASProperties",
runTime.constant(),
mesh,
IOobject::MUST_READ_IF_MODIFIED,
IOobject::NO_WRITE
);
IOobject LESPropertiesHeader
(
"LESProperties",
"turbulenceProperties",
runTime.constant(),
mesh,
IOobject::MUST_READ_IF_MODIFIED,
@ -95,15 +87,13 @@ void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
{
if (RASPropertiesHeader.headerOk())
if (turbulencePropertiesHeader.headerOk())
{
IOdictionary RASProperties(RASPropertiesHeader);
singlePhaseTransportModel laminarTransport(U, phi);
autoPtr<incompressible::RASModel> RASModel
autoPtr<incompressible::turbulenceModel> turbulenceModel
(
incompressible::RASModel::New
incompressible::turbulenceModel::New
(
U,
phi,
@ -117,47 +107,15 @@ void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
(
IOobject
(
"Pe",
"Pef",
runTime.timeName(),
mesh,
IOobject::NO_READ
mesh
),
mag(phi)
/(
mesh.magSf()
* mesh.surfaceInterpolation::deltaCoeffs()
* fvc::interpolate(RASModel->nuEff())
)
)
);
}
else if (LESPropertiesHeader.headerOk())
{
IOdictionary LESProperties(LESPropertiesHeader);
singlePhaseTransportModel laminarTransport(U, phi);
autoPtr<incompressible::LESModel> sgsModel
(
incompressible::LESModel::New(U, phi, laminarTransport)
);
PePtr.set
(
new surfaceScalarField
(
IOobject
(
"Pe",
runTime.timeName(),
mesh,
IOobject::NO_READ
),
mag(phi)
/(
mesh.magSf()
* mesh.surfaceInterpolation::deltaCoeffs()
* fvc::interpolate(sgsModel->nuEff())
* fvc::interpolate(turbulenceModel->nuEff())
)
)
);
@ -184,10 +142,9 @@ void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
(
IOobject
(
"Pe",
"Pef",
runTime.timeName(),
mesh,
IOobject::NO_READ
mesh
),
mag(phi)
/(
@ -201,10 +158,8 @@ void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
}
else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
{
if (RASPropertiesHeader.headerOk())
if (turbulencePropertiesHeader.headerOk())
{
IOdictionary RASProperties(RASPropertiesHeader);
autoPtr<fluidThermo> thermo(fluidThermo::New(mesh));
volScalarField rho
@ -218,9 +173,9 @@ void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
thermo->rho()
);
autoPtr<compressible::RASModel> RASModel
autoPtr<compressible::turbulenceModel> turbulenceModel
(
compressible::RASModel::New
compressible::turbulenceModel::New
(
rho,
U,
@ -235,58 +190,15 @@ void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
(
IOobject
(
"Pe",
runTime.timeName(),
mesh,
IOobject::NO_READ
),
mag(phi)
/(
mesh.magSf()
* mesh.surfaceInterpolation::deltaCoeffs()
* fvc::interpolate(RASModel->muEff())
)
)
);
}
else if (LESPropertiesHeader.headerOk())
{
IOdictionary LESProperties(LESPropertiesHeader);
autoPtr<fluidThermo> thermo(fluidThermo::New(mesh));
volScalarField rho
(
IOobject
(
"rho",
"Pef",
runTime.timeName(),
mesh
),
thermo->rho()
);
autoPtr<compressible::LESModel> sgsModel
(
compressible::LESModel::New(rho, U, phi, thermo())
);
PePtr.set
(
new surfaceScalarField
(
IOobject
(
"Pe",
runTime.timeName(),
mesh,
IOobject::NO_READ
),
mag(phi)
/(
mesh.magSf()
* mesh.surfaceInterpolation::deltaCoeffs()
* fvc::interpolate(sgsModel->muEff())
* fvc::interpolate(turbulenceModel->muEff())
)
)
);
@ -313,10 +225,9 @@ void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
(
IOobject
(
"Pe",
"Pef",
runTime.timeName(),
mesh,
IOobject::NO_READ
mesh
),
mag(phi)
/(
@ -335,19 +246,34 @@ void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
<< abort(FatalError);
}
Info<< "Pe max : " << max(PePtr()).value() << endl;
Info<< " Pe max : " << max(PePtr()).value() << endl;
if (writeResults)
{
Info<< " Writing surfaceScalarField : "
<< PePtr().name() << endl;
PePtr().write();
volScalarField Pe
(
IOobject
(
"Pe",
runTime.timeName(),
mesh
),
fvc::average(PePtr())
);
Info<< " Writing volScalarField : "
<< Pe.name() << endl;
Pe.write();
}
}
else
{
Info<< " No phi" << endl;
}
Info<< "\nEnd\n" << endl;
}

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -57,7 +57,7 @@ namespace Foam
int main(int argc, char *argv[])
{
Foam::timeSelector::addOptions();
# include "addRegionOption.H"
#include "addRegionOption.H"
Foam::argList::addBoolOption
(
"noWrite",