/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2017 OpenFOAM Foundation ------------------------------------------------------------------------------- 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 "laplaceFilter.H" #include "addToRunTimeSelectionTable.H" #include "calculatedFvPatchFields.H" #include "fvm.H" #include "fvc.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // namespace Foam { defineTypeNameAndDebug(laplaceFilter, 0); addToRunTimeSelectionTable(LESfilter, laplaceFilter, dictionary); } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::laplaceFilter::laplaceFilter(const fvMesh& mesh, scalar widthCoeff) : LESfilter(mesh), widthCoeff_(widthCoeff), coeff_ ( IOobject ( "laplaceFilterCoeff", mesh.time().timeName(), mesh ), mesh, dimensionedScalar(dimLength*dimLength, Zero), calculatedFvPatchScalarField::typeName ) { coeff_.ref() = pow(mesh.V(), 2.0/3.0)/widthCoeff_; } Foam::laplaceFilter::laplaceFilter(const fvMesh& mesh, const dictionary& bd) : LESfilter(mesh), widthCoeff_ ( bd.optionalSubDict(type() + "Coeffs").get("widthCoeff") ), coeff_ ( IOobject ( "laplaceFilterCoeff", mesh.time().timeName(), mesh ), mesh, dimensionedScalar(dimLength*dimLength, Zero), calculatedFvPatchScalarField::typeName ) { coeff_.ref() = pow(mesh.V(), 2.0/3.0)/widthCoeff_; } // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // void Foam::laplaceFilter::read(const dictionary& bd) { bd.optionalSubDict(type() + "Coeffs").readEntry("widthCoeff", widthCoeff_); } // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * // Foam::tmp Foam::laplaceFilter::operator() ( const tmp& unFilteredField ) const { correctBoundaryConditions(unFilteredField); tmp filteredField = unFilteredField() + fvc::laplacian(coeff_, unFilteredField()); unFilteredField.clear(); return filteredField; } Foam::tmp Foam::laplaceFilter::operator() ( const tmp& unFilteredField ) const { correctBoundaryConditions(unFilteredField); tmp filteredField = unFilteredField() + fvc::laplacian(coeff_, unFilteredField()); unFilteredField.clear(); return filteredField; } Foam::tmp Foam::laplaceFilter::operator() ( const tmp& unFilteredField ) const { correctBoundaryConditions(unFilteredField); tmp filteredField = unFilteredField() + fvc::laplacian(coeff_, unFilteredField()); unFilteredField.clear(); return filteredField; } Foam::tmp Foam::laplaceFilter::operator() ( const tmp& unFilteredField ) const { correctBoundaryConditions(unFilteredField); tmp filteredField = unFilteredField() + fvc::laplacian(coeff_, unFilteredField()); unFilteredField.clear(); return filteredField; } // ************************************************************************* //