initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / turbulenceModels / LES / LESfilters / simpleFilter / simpleFilter.C
blob0b8463efeb3cd4a2252402f9f9c967ff0decc8ef
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 "simpleFilter.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvc.H"
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 namespace Foam
35     defineTypeNameAndDebug(simpleFilter, 0);
36     addToRunTimeSelectionTable(LESfilter, simpleFilter, dictionary);
40 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
42 Foam::simpleFilter::simpleFilter
44     const fvMesh& mesh
47     LESfilter(mesh)
51 Foam::simpleFilter::simpleFilter(const fvMesh& mesh, const dictionary&)
53     LESfilter(mesh)
57 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
59 void Foam::simpleFilter::read(const dictionary&)
63 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
65 Foam::tmp<Foam::volScalarField> Foam::simpleFilter::operator()
67     const tmp<volScalarField>& unFilteredField
68 ) const
70     tmp<volScalarField> filteredField = fvc::surfaceSum
71     (
72         mesh().magSf()*fvc::interpolate(unFilteredField)
73     )/fvc::surfaceSum(mesh().magSf());
75     unFilteredField.clear();
77     return filteredField;
81 Foam::tmp<Foam::volVectorField> Foam::simpleFilter::operator()
83     const tmp<volVectorField>& unFilteredField
84 ) const
86     tmp<volVectorField> filteredField = fvc::surfaceSum
87     (
88         mesh().magSf()*fvc::interpolate(unFilteredField)
89     )/fvc::surfaceSum(mesh().magSf());
91     unFilteredField.clear();
93     return filteredField;
97 Foam::tmp<Foam::volSymmTensorField> Foam::simpleFilter::operator()
99     const tmp<volSymmTensorField>& unFilteredField
100 ) const
102     tmp<volSymmTensorField> filteredField = fvc::surfaceSum
103     (
104         mesh().magSf()*fvc::interpolate(unFilteredField)
105     )/fvc::surfaceSum(mesh().magSf());
107     unFilteredField.clear();
109     return filteredField;
113 Foam::tmp<Foam::volTensorField> Foam::simpleFilter::operator()
115     const tmp<volTensorField>& unFilteredField
116 ) const
118     tmp<volTensorField> filteredField = fvc::surfaceSum
119     (
120         mesh().magSf()*fvc::interpolate(unFilteredField)
121     )/fvc::surfaceSum(mesh().magSf());
123     unFilteredField.clear();
125     return filteredField;
129 // ************************************************************************* //