simpleFilter.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "simpleFilter.H"
30 #include "fvc.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(simpleFilter, 0);
37  addToRunTimeSelectionTable(LESfilter, simpleFilter, dictionary);
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 Foam::simpleFilter::simpleFilter
44 (
45  const fvMesh& mesh
46 )
47 :
49 {}
50 
51 
52 Foam::simpleFilter::simpleFilter(const fvMesh& mesh, const dictionary&)
53 :
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
59 
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
65 
66 Foam::tmp<Foam::volScalarField> Foam::simpleFilter::operator()
67 (
68  const tmp<volScalarField>& unFilteredField
69 ) const
70 {
71  correctBoundaryConditions(unFilteredField);
72 
73  tmp<volScalarField> filteredField = fvc::surfaceSum
74  (
75  mesh().magSf()*fvc::interpolate(unFilteredField)
76  )/fvc::surfaceSum(mesh().magSf());
77 
78  unFilteredField.clear();
79 
80  return filteredField;
81 }
82 
83 
84 Foam::tmp<Foam::volVectorField> Foam::simpleFilter::operator()
85 (
86  const tmp<volVectorField>& unFilteredField
87 ) const
88 {
89  correctBoundaryConditions(unFilteredField);
90 
91  tmp<volVectorField> filteredField = fvc::surfaceSum
92  (
93  mesh().magSf()*fvc::interpolate(unFilteredField)
94  )/fvc::surfaceSum(mesh().magSf());
95 
96  unFilteredField.clear();
97 
98  return filteredField;
99 }
100 
101 
102 Foam::tmp<Foam::volSymmTensorField> Foam::simpleFilter::operator()
103 (
104  const tmp<volSymmTensorField>& unFilteredField
105 ) const
106 {
107  correctBoundaryConditions(unFilteredField);
108 
110  (
111  mesh().magSf()*fvc::interpolate(unFilteredField)
112  )/fvc::surfaceSum(mesh().magSf());
113 
114  unFilteredField.clear();
115 
116  return filteredField;
117 }
118 
119 
120 Foam::tmp<Foam::volTensorField> Foam::simpleFilter::operator()
121 (
122  const tmp<volTensorField>& unFilteredField
123 ) const
124 {
125  correctBoundaryConditions(unFilteredField);
126 
127  tmp<volTensorField> filteredField = fvc::surfaceSum
128  (
129  mesh().magSf()*fvc::interpolate(unFilteredField)
130  )/fvc::surfaceSum(mesh().magSf());
131 
132  unFilteredField.clear();
133 
134  return filteredField;
135 }
136 
137 
138 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::fvc::surfaceSum
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcSurfaceIntegrate.C:135
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
correctBoundaryConditions
cellMask correctBoundaryConditions()
Foam::simpleFilter::read
virtual void read(const dictionary &)
Read the LESfilter dictionary.
Definition: simpleFilter.C:60
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::LESfilter
Abstract class for LES filters.
Definition: LESfilter.H:57
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvc::interpolate
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
simpleFilter.H
fvc.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)