openfoam/applications/utilities/postProcessing/turbulence/R/R.C

170 lines
4.3 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
R
Description
Calculates and writes the Reynolds stress R for the current time step.
Compressible modes is automatically selected based on the existence of the
"thermophysicalProperties" dictionary required to construct the
thermodynamics package.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "turbulentTransportModel.H"
#include "turbulentFluidThermoModel.H"
#include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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
)
{
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"
#include "setRootCase.H"
#include "createTime.H"
instantList timeDirs = timeSelector::select0(runTime, args);
#include "createNamedMesh.H"
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
(
IOobject
(
basicThermo::dictName,
runTime.constant(),
mesh
).headerOk()
)
{
calcCompressibleR(mesh, runTime, U);
}
else
{
calcIncompressibleR(mesh, runTime, U);
}
}
else
{
Info<< " no " << UHeader.name() << " field" << endl;
}
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //