All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2020 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"
28 #include "fvc.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  defineTypeNameAndDebug(simpleFilter, 0);
35  addToRunTimeSelectionTable(LESfilter, simpleFilter, dictionary);
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const fvMesh& mesh
44 )
45 :
46  LESfilter(mesh)
47 {}
48 
49 
51 :
52  LESfilter(mesh)
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
57 
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
63 
64 Foam::tmp<Foam::volScalarField> Foam::simpleFilter::operator()
65 (
66  const tmp<volScalarField>& unFilteredField
67 ) const
68 {
69  correctBoundaryConditions(unFilteredField);
70 
71  tmp<volScalarField> filteredField = fvc::surfaceSum
72  (
73  mesh().magSf()*fvc::interpolate(unFilteredField)
74  )/fvc::surfaceSum(mesh().magSf());
75 
76  unFilteredField.clear();
77 
78  return filteredField;
79 }
80 
81 
82 Foam::tmp<Foam::volVectorField> Foam::simpleFilter::operator()
83 (
84  const tmp<volVectorField>& unFilteredField
85 ) const
86 {
87  correctBoundaryConditions(unFilteredField);
88 
89  tmp<volVectorField> filteredField = fvc::surfaceSum
90  (
91  mesh().magSf()*fvc::interpolate(unFilteredField)
92  )/fvc::surfaceSum(mesh().magSf());
93 
94  unFilteredField.clear();
95 
96  return filteredField;
97 }
98 
99 
100 Foam::tmp<Foam::volSymmTensorField> Foam::simpleFilter::operator()
101 (
102  const tmp<volSymmTensorField>& unFilteredField
103 ) const
104 {
105  correctBoundaryConditions(unFilteredField);
106 
108  (
109  mesh().magSf()*fvc::interpolate(unFilteredField)
110  )/fvc::surfaceSum(mesh().magSf());
111 
112  unFilteredField.clear();
113 
114  return filteredField;
115 }
116 
117 
118 Foam::tmp<Foam::volTensorField> Foam::simpleFilter::operator()
119 (
120  const tmp<volTensorField>& unFilteredField
121 ) const
122 {
123  correctBoundaryConditions(unFilteredField);
124 
125  tmp<volTensorField> filteredField = fvc::surfaceSum
126  (
127  mesh().magSf()*fvc::interpolate(unFilteredField)
128  )/fvc::surfaceSum(mesh().magSf());
129 
130  unFilteredField.clear();
131 
132  return filteredField;
133 }
134 
135 
136 // ************************************************************************* //
Abstract class for LES filters.
Definition: LESfilter.H:54
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const fvMesh & mesh() const
Return mesh reference.
Definition: LESfilter.H:130
Macros for easy insertion into run-time selection tables.
dynamicFvMesh & mesh
simpleFilter(const fvMesh &mesh)
Construct from components.
Definition: simpleFilter.C:42
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
virtual void read(const dictionary &)
Read the LESfilter dictionary.
Definition: simpleFilter.C:58
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.
void correctBoundaryConditions(const tmp< GeoFieldType > &tgf) const
Temporary function to ensure the coupled boundary conditions of the.
Definition: LESfilter.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.