/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2013-2014 OpenFOAM Foundation \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd. ------------------------------------------------------------------------------- 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 "Lambda2.H" #include "volFields.H" #include "dictionary.H" #include "zeroGradientFvPatchFields.H" #include "fvcGrad.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // namespace Foam { defineTypeNameAndDebug(Lambda2, 0); } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::Lambda2::Lambda2 ( const word& name, const objectRegistry& obr, const dictionary& dict, const bool loadFromFiles ) : name_(name), obr_(obr), active_(true), UName_("U"), resultName_(name), log_(true) { // Check if the available mesh is an fvMesh, otherwise deactivate if (!isA(obr_)) { active_ = false; WarningIn ( "Lambda2::Lambda2" "(" "const word&, " "const objectRegistry&, " "const dictionary&, " "const bool" ")" ) << "No fvMesh available, deactivating " << name_ << nl << endl; } read(dict); if (active_) { const fvMesh& mesh = refCast(obr_); volScalarField* Lambda2Ptr ( new volScalarField ( IOobject ( resultName_, mesh.time().timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE ), mesh, dimensionedScalar("0", dimless/sqr(dimTime), 0.0) ) ); mesh.objectRegistry::store(Lambda2Ptr); } } // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // Foam::Lambda2::~Lambda2() {} // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // void Foam::Lambda2::read(const dictionary& dict) { if (active_) { log_.readIfPresent("log", dict); dict.readIfPresent("UName", UName_); if (!dict.readIfPresent("resultName", resultName_)) { resultName_ = name_; if (UName_ != "U") { resultName_ = resultName_ + "(" + UName_ + ")"; } } } } void Foam::Lambda2::execute() { if (active_) { const fvMesh& mesh = refCast(obr_); const volVectorField& U = mesh.lookupObject(UName_); const volTensorField gradU(fvc::grad(U)); const volTensorField SSplusWW ( (symm(gradU) & symm(gradU)) + (skew(gradU) & skew(gradU)) ); volScalarField& Lambda2 = const_cast ( mesh.lookupObject(resultName_) ); Lambda2 = -eigenValues(SSplusWW)().component(vector::Y); } } void Foam::Lambda2::end() { if (active_) { execute(); } } void Foam::Lambda2::timeSet() { // Do nothing } void Foam::Lambda2::write() { if (active_) { const volScalarField& Lambda2 = obr_.lookupObject(resultName_); if (log_) Info << type() << " " << name_ << " output:" << nl << " writing field " << Lambda2.name() << nl << endl; Lambda2.write(); } } // ************************************************************************* //