mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: Updated R utility to include compressible cases
This commit is contained in:
@ -1,13 +1,19 @@
|
||||
EXE_INC = \
|
||||
-I$(LIB_SRC)/postProcessing/postCalc \
|
||||
-I$(LIB_SRC)/turbulenceModels \
|
||||
-I$(LIB_SRC)/turbulenceModels/incompressible/RAS/RASModel \
|
||||
-I$(LIB_SRC)/meshTools/lnInclude \
|
||||
-I$(LIB_SRC)/transportModels \
|
||||
-I$(LIB_SRC)/turbulenceModels \
|
||||
-I$(LIB_SRC)/turbulenceModels/compressible/RAS/derivedFvPatchFields/wallFunctions/mutWallFunctions \
|
||||
-I$(LIB_SRC)/turbulenceModels/incompressible/RAS/derivedFvPatchFields/wallFunctions/nutWallFunctions \
|
||||
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
|
||||
-I$(LIB_SRC)/finiteVolume/lnInclude
|
||||
|
||||
EXE_LIBS = \
|
||||
$(FOAM_LIBBIN)/postCalc.o \
|
||||
-lincompressibleRASModels \
|
||||
-lincompressibleTransportModels \
|
||||
-lincompressibleRASModels \
|
||||
-lfluidThermophysicalModels \
|
||||
-lspecie \
|
||||
-lcompressibleRASModels \
|
||||
-lfiniteVolume \
|
||||
-lgenericPatchFields
|
||||
-lgenericPatchFields \
|
||||
-lmeshTools \
|
||||
-lsampling
|
||||
|
||||
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -29,35 +29,140 @@ Description
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "calc.H"
|
||||
#include "fvCFD.H"
|
||||
#include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
|
||||
#include "RASModel.H"
|
||||
#include "incompressible/turbulenceModel/turbulenceModel.H"
|
||||
|
||||
#include "fluidThermo.H"
|
||||
#include "compressible/turbulenceModel/turbulenceModel.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
|
||||
void calcIncompressibleR
|
||||
(
|
||||
const fvMesh& mesh,
|
||||
const Time& runTime,
|
||||
const volVectorField& U
|
||||
)
|
||||
{
|
||||
#include "createFields.H"
|
||||
#include "createPhi.H"
|
||||
|
||||
Info<< "\nCalculating the Reynolds Stress R\n" << endl;
|
||||
singlePhaseTransportModel laminarTransport(U, phi);
|
||||
|
||||
volSymmTensorField R
|
||||
autoPtr<incompressible::turbulenceModel> model
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"R",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::NO_READ,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
RASModel->R()
|
||||
incompressible::turbulenceModel::New(U, phi, laminarTransport)
|
||||
);
|
||||
|
||||
R.write();
|
||||
Info<< "Writing R field" << nl << endl;
|
||||
|
||||
Info<< "End" << endl;
|
||||
model->R()().write();
|
||||
}
|
||||
|
||||
|
||||
void calcCompressibleR
|
||||
(
|
||||
const fvMesh& mesh,
|
||||
const Time& runTime,
|
||||
const volVectorField& U
|
||||
)
|
||||
{
|
||||
IOobject rhoHeader
|
||||
(
|
||||
"rho",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::MUST_READ,
|
||||
IOobject::NO_WRITE
|
||||
);
|
||||
|
||||
if (!rhoHeader.headerOk())
|
||||
{
|
||||
Info<< " no " << rhoHeader.name() <<" field" << endl;
|
||||
return;
|
||||
}
|
||||
|
||||
Info<< "Reading field rho\n" << endl;
|
||||
volScalarField rho(rhoHeader, mesh);
|
||||
|
||||
#include "compressibleCreatePhi.H"
|
||||
|
||||
autoPtr<fluidThermo> pThermo(fluidThermo::New(mesh));
|
||||
fluidThermo& thermo = pThermo();
|
||||
|
||||
autoPtr<compressible::turbulenceModel> model
|
||||
(
|
||||
compressible::turbulenceModel::New
|
||||
(
|
||||
rho,
|
||||
U,
|
||||
phi,
|
||||
thermo
|
||||
)
|
||||
);
|
||||
|
||||
Info<< "Writing R field" << nl << endl;
|
||||
|
||||
model->R()().write();
|
||||
}
|
||||
|
||||
|
||||
int main(int argc, char *argv[])
|
||||
{
|
||||
timeSelector::addOptions();
|
||||
|
||||
#include "addRegionOption.H"
|
||||
|
||||
argList::addBoolOption
|
||||
(
|
||||
"compressible",
|
||||
"calculate compressible R"
|
||||
);
|
||||
|
||||
#include "setRootCase.H"
|
||||
#include "createTime.H"
|
||||
instantList timeDirs = timeSelector::select0(runTime, args);
|
||||
#include "createNamedMesh.H"
|
||||
|
||||
const bool compressible = args.optionFound("compressible");
|
||||
|
||||
forAll(timeDirs, timeI)
|
||||
{
|
||||
runTime.setTime(timeDirs[timeI], timeI);
|
||||
Info<< "Time = " << runTime.timeName() << endl;
|
||||
|
||||
IOobject UHeader
|
||||
(
|
||||
"U",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::MUST_READ,
|
||||
IOobject::NO_WRITE
|
||||
);
|
||||
|
||||
if (UHeader.headerOk())
|
||||
{
|
||||
Info<< "Reading field " << UHeader.name() << nl << endl;
|
||||
volVectorField U(UHeader, mesh);
|
||||
|
||||
if (compressible)
|
||||
{
|
||||
calcCompressibleR(mesh, runTime, U);
|
||||
}
|
||||
else
|
||||
{
|
||||
calcIncompressibleR(mesh, runTime, U);
|
||||
}
|
||||
}
|
||||
else
|
||||
{
|
||||
Info<< " no " << UHeader.name() << " field" << endl;
|
||||
}
|
||||
}
|
||||
|
||||
Info<< "End\n" << endl;
|
||||
|
||||
return 0;
|
||||
}
|
||||
|
||||
|
||||
|
||||
@ -1,22 +0,0 @@
|
||||
Info<< "Reading field U\n" << endl;
|
||||
volVectorField U
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"U",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::MUST_READ,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
mesh
|
||||
);
|
||||
|
||||
#include "createPhi.H"
|
||||
|
||||
singlePhaseTransportModel laminarTransport(U, phi);
|
||||
|
||||
autoPtr<incompressible::RASModel> RASModel
|
||||
(
|
||||
incompressible::RASModel::New(U, phi, laminarTransport)
|
||||
);
|
||||
Reference in New Issue
Block a user