/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2013 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 . Application R Description Calculates and writes the Reynolds stress R for the current time step. \*---------------------------------------------------------------------------*/ #include "fvCFD.H" #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H" #include "incompressible/turbulenceModel/turbulenceModel.H" #include "fluidThermo.H" #include "compressible/turbulenceModel/turbulenceModel.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // void calcIncompressibleR ( const fvMesh& mesh, const Time& runTime, const volVectorField& U ) { #include "createPhi.H" singlePhaseTransportModel laminarTransport(U, phi); autoPtr 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 pThermo(fluidThermo::New(mesh)); fluidThermo& thermo = pThermo(); autoPtr 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; } // ************************************************************************* //