Skip to content
Snippets Groups Projects
Commit 7c3d2ec1 authored by andy's avatar andy
Browse files

ENH: Updated R utility to include compressible cases

parent 6adba987
Branches
Tags
No related merge requests found
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 "createPhi.H"
singlePhaseTransportModel laminarTransport(U, phi);
autoPtr<incompressible::turbulenceModel> model
(
incompressible::turbulenceModel::New(U, phi, laminarTransport)
);
Info<< "Writing R field" << nl << endl;
model->R()().write();
}
void calcCompressibleR
(
const fvMesh& mesh,
const Time& runTime,
const volVectorField& U
)
{
#include "createFields.H"
IOobject rhoHeader
(
"rho",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
);
if (!rhoHeader.headerOk())
{
Info<< " no " << rhoHeader.name() <<" field" << endl;
return;
}
Info<< "\nCalculating the Reynolds Stress R\n" << endl;
Info<< "Reading field rho\n" << endl;
volScalarField rho(rhoHeader, mesh);
volSymmTensorField R
#include "compressibleCreatePhi.H"
autoPtr<fluidThermo> pThermo(fluidThermo::New(mesh));
fluidThermo& thermo = pThermo();
autoPtr<compressible::turbulenceModel> model
(
IOobject
compressible::turbulenceModel::New
(
"R",
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::NO_READ,
IOobject::AUTO_WRITE
),
RASModel->R()
);
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;
}
}
R.write();
Info<< "End\n" << endl;
Info<< "End" << endl;
return 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)
);
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment