/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2012 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 "DESModelRegions.H" #include "volFields.H" #include "compressible/turbulenceModel/turbulenceModel.H" #include "compressible/LES/DESModel/DESModel.H" #include "incompressible/turbulenceModel/turbulenceModel.H" #include "incompressible/LES/DESModel/DESModel.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // defineTypeNameAndDebug(Foam::DESModelRegions, 0); // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // void Foam::DESModelRegions::makeFile() { // Create the output file if not already created if (outputFilePtr_.empty()) { if (debug) { Info<< "Creating output file." << endl; } // File update if (Pstream::master()) { fileName outputDir; word startTimeName = obr_.time().timeName(obr_.time().startTime().value()); if (Pstream::parRun()) { // Put in undecomposed case (Note: gives problems for // distributed data running) outputDir = obr_.time().path()/".."/name_/startTimeName; } else { outputDir = obr_.time().path()/name_/startTimeName; } // Create directory if does not exist mkDir(outputDir); // Open new file at start up outputFilePtr_.reset(new OFstream(outputDir/(type() + ".dat"))); // Add headers to output data outputFilePtr_() << "# DES model region coverage (% volume)" << nl << "# time " << token::TAB << "LES" << token::TAB << "RAS" << endl; } } } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::DESModelRegions::DESModelRegions ( const word& name, const objectRegistry& obr, const dictionary& dict, const bool loadFromFiles ) : name_(name), obr_(obr), active_(true), log_(false), outputFilePtr_(NULL) { // Check if the available mesh is an fvMesh, otherwise deactivate if (!isA(obr_)) { active_ = false; WarningIn ( "DESModelRegions::DESModelRegions" "(" "const word&, " "const objectRegistry&, " "const dictionary&, " "const bool" ")" ) << "No fvMesh available, deactivating." << nl << endl; } makeFile(); read(dict); } // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // Foam::DESModelRegions::~DESModelRegions() {} // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // void Foam::DESModelRegions::read(const dictionary& dict) { if (active_) { log_ = dict.lookupOrDefault("log", false); } } void Foam::DESModelRegions::execute() { // Do nothing - only valid on write } void Foam::DESModelRegions::end() { // Do nothing - only valid on write } void Foam::DESModelRegions::write() { typedef incompressible::turbulenceModel icoModel; typedef incompressible::DESModel icoDESModel; typedef compressible::turbulenceModel cmpModel; typedef compressible::DESModel cmpDESModel; if (active_) { const fvMesh& mesh = refCast(obr_); if (log_) { Info<< type() << " output:" << nl; } tmp result; label DESpresent = false; if (mesh.foundObject("turbulenceModel")) { const icoModel& model = mesh.lookupObject("turbulenceModel"); if (isA(model)) { const icoDESModel& des = dynamic_cast(model); result = des.LESRegion(); DESpresent = true; } } else if (mesh.foundObject("turbulenceModel")) { const cmpModel& model = mesh.lookupObject("turbulenceModel"); if (isA(model)) { const cmpDESModel& des = dynamic_cast(model); result = des.LESRegion(); DESpresent = true; } } if (DESpresent) { scalar prc = gSum(result().internalField()*mesh.V())/gSum(mesh.V())*100.0; if (Pstream::master()) { outputFilePtr_() << obr_.time().timeName() << token::TAB << prc << token::TAB << 100.0 - prc << endl; } if (log_) { Info<< " LES = " << prc << " % (volume)" << nl << " RES = " << 100.0 - prc << " % (volume)" << nl << endl; } } else { if (log_) { Info<< " No DES turbulence model found in database" << nl << endl; } } } } // ************************************************************************* //