initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / LES / LESfilters / laplaceFilter / laplaceFilter.C
blob9b21fe9cee5187f0e3138376d18b6e897bae3e0a
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "laplaceFilter.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "calculatedFvPatchFields.H"
30 #include "fvm.H"
31 #include "fvc.H"
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 namespace Foam
37     defineTypeNameAndDebug(laplaceFilter, 0);
38     addToRunTimeSelectionTable(LESfilter, laplaceFilter, dictionary);
42 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
44 Foam::laplaceFilter::laplaceFilter(const fvMesh& mesh, scalar widthCoeff)
46     LESfilter(mesh),
47     widthCoeff_(widthCoeff),
48     coeff_
49     (
50         IOobject
51         (
52             "laplaceFilterCoeff",
53             mesh.time().timeName(),
54             mesh
55         ),
56         mesh,
57         dimensionedScalar("zero", dimLength*dimLength, 0),
58         calculatedFvPatchScalarField::typeName
59     )
61     coeff_.internalField() = pow(mesh.V(), 2.0/3.0)/widthCoeff_;
65 Foam::laplaceFilter::laplaceFilter(const fvMesh& mesh, const dictionary& bd)
67     LESfilter(mesh),
68     widthCoeff_(readScalar(bd.subDict(type() + "Coeffs").lookup("widthCoeff"))),
69     coeff_
70     (
71         IOobject
72         (
73             "laplaceFilterCoeff",
74             mesh.time().timeName(),
75             mesh
76         ),
77         mesh,
78         dimensionedScalar("zero", dimLength*dimLength, 0),
79         calculatedFvPatchScalarField::typeName
80     )
82     coeff_.internalField() = pow(mesh.V(), 2.0/3.0)/widthCoeff_;
86 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
88 void Foam::laplaceFilter::read(const dictionary& bd)
90     bd.subDict(type() + "Coeffs").lookup("widthCoeff") >> widthCoeff_;
94 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
96 Foam::tmp<Foam::volScalarField> Foam::laplaceFilter::operator()
98     const tmp<volScalarField>& unFilteredField
99 ) const
101     tmp<volScalarField> filteredField =
102         unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
104     unFilteredField.clear();
106     return filteredField;
110 Foam::tmp<Foam::volVectorField> Foam::laplaceFilter::operator()
112     const tmp<volVectorField>& unFilteredField
113 ) const
115     tmp<volVectorField> filteredField =
116         unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
118     unFilteredField.clear();
120     return filteredField;
124 Foam::tmp<Foam::volSymmTensorField> Foam::laplaceFilter::operator()
126     const tmp<volSymmTensorField>& unFilteredField
127 ) const
129     tmp<volSymmTensorField> filteredField =
130         unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
132     unFilteredField.clear();
134     return filteredField;
138 Foam::tmp<Foam::volTensorField> Foam::laplaceFilter::operator()
140     const tmp<volTensorField>& unFilteredField
141 ) const
143     tmp<volTensorField> filteredField =
144         unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
146     unFilteredField.clear();
148     return filteredField;
152 // ************************************************************************* //