simpleFilter.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "simpleFilter.H"
27 #include "fvcSurfaceIntegrate.H"
28 #include "surfaceInterpolate.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const fvMesh& mesh
45 )
46 :
47  LESfilter(mesh)
48 {}
49 
50 
52 :
53  LESfilter(mesh)
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
58 
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
64 
65 Foam::tmp<Foam::volScalarField> Foam::simpleFilter::operator()
66 (
67  const tmp<volScalarField>& unFilteredField
68 ) const
69 {
70  correctBoundaryConditions(unFilteredField);
71 
72  tmp<volScalarField> filteredField = fvc::surfaceSum
73  (
74  mesh().magSf()*fvc::interpolate(unFilteredField)
75  )/fvc::surfaceSum(mesh().magSf());
76 
77  unFilteredField.clear();
78 
79  return filteredField;
80 }
81 
82 
83 Foam::tmp<Foam::volVectorField> Foam::simpleFilter::operator()
84 (
85  const tmp<volVectorField>& unFilteredField
86 ) const
87 {
88  correctBoundaryConditions(unFilteredField);
89 
90  tmp<volVectorField> filteredField = fvc::surfaceSum
91  (
92  mesh().magSf()*fvc::interpolate(unFilteredField)
93  )/fvc::surfaceSum(mesh().magSf());
94 
95  unFilteredField.clear();
96 
97  return filteredField;
98 }
99 
100 
101 Foam::tmp<Foam::volSymmTensorField> Foam::simpleFilter::operator()
102 (
103  const tmp<volSymmTensorField>& unFilteredField
104 ) const
105 {
106  correctBoundaryConditions(unFilteredField);
107 
109  (
110  mesh().magSf()*fvc::interpolate(unFilteredField)
111  )/fvc::surfaceSum(mesh().magSf());
112 
113  unFilteredField.clear();
114 
115  return filteredField;
116 }
117 
118 
119 Foam::tmp<Foam::volTensorField> Foam::simpleFilter::operator()
120 (
121  const tmp<volTensorField>& unFilteredField
122 ) const
123 {
124  correctBoundaryConditions(unFilteredField);
125 
126  tmp<volTensorField> filteredField = fvc::surfaceSum
127  (
128  mesh().magSf()*fvc::interpolate(unFilteredField)
129  )/fvc::surfaceSum(mesh().magSf());
130 
131  unFilteredField.clear();
132 
133  return filteredField;
134 }
135 
136 
137 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Abstract class for LES filters.
Definition: LESfilter.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Simple top-hat filter used in dynamic LES models.
Definition: simpleFilter.H:54
simpleFilter(const fvMesh &mesh)
Construct from components.
Definition: simpleFilter.C:43
virtual void read(const dictionary &)
Read the LESfilter dictionary.
Definition: simpleFilter.C:59
A class for managing temporary objects.
Definition: tmp.H:55
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
U correctBoundaryConditions()
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
defineTypeNameAndDebug(combustionModel, 0)