/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2013-2016 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 . \*---------------------------------------------------------------------------*/ #include "Peclet.H" #include "volFields.H" #include "dictionary.H" #include "surfaceFields.H" #include "turbulentTransportModel.H" #include "turbulentFluidThermoModel.H" #include "surfaceInterpolate.H" #include "addToRunTimeSelectionTable.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // namespace Foam { namespace functionObjects { defineTypeNameAndDebug(Peclet, 0); addToRunTimeSelectionTable ( functionObject, Peclet, dictionary ); } } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::functionObjects::Peclet::Peclet ( const word& name, const Time& runTime, const dictionary& dict ) : fvMeshFunctionObject(name, runTime, dict) { read(dict); surfaceScalarField* PecletPtr ( new surfaceScalarField ( IOobject ( type(), mesh_.time().timeName(), mesh_, IOobject::NO_READ, IOobject::NO_WRITE ), mesh_, dimensionedScalar("0", dimless, 0.0) ) ); mesh_.objectRegistry::store(PecletPtr); } // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // Foam::functionObjects::Peclet::~Peclet() {} // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // bool Foam::functionObjects::Peclet::read(const dictionary& dict) { phiName_ = dict.lookupOrDefault("phiName", "phi"); rhoName_ = dict.lookupOrDefault("rhoName", "rho"); return true; } bool Foam::functionObjects::Peclet::execute(const bool postProcess) { typedef compressible::turbulenceModel cmpTurbModel; typedef incompressible::turbulenceModel icoTurbModel; tmp nuEff; if (mesh_.foundObject(turbulenceModel::propertiesName)) { const cmpTurbModel& model = mesh_.lookupObject ( turbulenceModel::propertiesName ); const volScalarField& rho = mesh_.lookupObject(rhoName_); nuEff = model.muEff()/rho; } else if ( mesh_.foundObject(turbulenceModel::propertiesName) ) { const icoTurbModel& model = mesh_.lookupObject ( turbulenceModel::propertiesName ); nuEff = model.nuEff(); } else if (mesh_.foundObject("transportProperties")) { const dictionary& model = mesh_.lookupObject("transportProperties"); nuEff = tmp ( new volScalarField ( IOobject ( "nuEff", mesh_.time().timeName(), mesh_, IOobject::NO_READ, IOobject::NO_WRITE ), mesh_, dimensionedScalar(model.lookup("nu")) ) ); } else { FatalErrorInFunction << "Unable to determine the viscosity" << exit(FatalError); } const surfaceScalarField& phi = mesh_.lookupObject(phiName_); surfaceScalarField& Peclet = const_cast ( mesh_.lookupObject(type()) ); Peclet = mag(phi) /( mesh_.magSf() *mesh_.surfaceInterpolation::deltaCoeffs() *fvc::interpolate(nuEff) ); return true; } bool Foam::functionObjects::Peclet::write(const bool postProcess) { const surfaceScalarField& Peclet = mesh_.lookupObject(type()); Info<< type() << " " << name() << " output:" << nl << " writing field " << Peclet.name() << nl << endl; Peclet.write(); return true; } // ************************************************************************* //