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-2025 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 template<class Type>
43 inline Foam::tmp<Foam::VolField<Type>> Foam::simpleFilter::filter
44 (
45  const tmp<VolField<Type>>& unFilteredField
46 ) const
47 {
48  correctBoundaryConditions(unFilteredField);
49 
50  tmp<VolField<Type>> filteredField
51  (
53  (
54  "simpleFilteredField",
56  (
57  mesh().magSf()*fvc::interpolate(unFilteredField)
58  )/fvc::surfaceSum(mesh().magSf()),
60  )
61  );
62 
63  unFilteredField.clear();
64 
65  return filteredField;
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
72 (
73  const fvMesh& mesh
74 )
75 :
77 {}
78 
79 
81 :
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
93 
94 Foam::tmp<Foam::volScalarField> Foam::simpleFilter::operator()
95 (
96  const tmp<volScalarField>& unFilteredField
97 ) const
98 {
99  return filter(unFilteredField);
100 }
101 
102 
103 Foam::tmp<Foam::volVectorField> Foam::simpleFilter::operator()
104 (
105  const tmp<volVectorField>& unFilteredField
106 ) const
107 {
108  return filter(unFilteredField);
109 }
110 
111 
112 Foam::tmp<Foam::volSymmTensorField> Foam::simpleFilter::operator()
113 (
114  const tmp<volSymmTensorField>& unFilteredField
115 ) const
116 {
117  return filter(unFilteredField);
118 }
119 
120 
121 Foam::tmp<Foam::volTensorField> Foam::simpleFilter::operator()
122 (
123  const tmp<volTensorField>& unFilteredField
124 ) const
125 {
126  return filter(unFilteredField);
127 }
128 
129 
130 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Abstract class for LES filters.
Definition: LESfilter.H:55
void correctBoundaryConditions(const tmp< GeoFieldType > &tgf) const
Temporary function to ensure the coupled boundary conditions of the.
Definition: LESfilter.H:72
const fvMesh & mesh() const
Return mesh reference.
Definition: LESfilter.H:130
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Simple top-hat filter used in dynamic LES models.
Definition: simpleFilter.H:53
simpleFilter(const fvMesh &mesh)
Construct from components.
Definition: simpleFilter.C:72
virtual void read(const dictionary &)
Read the LESfilter dictionary.
Definition: simpleFilter.C:88
A class for managing temporary objects.
Definition: tmp.H:55
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolInternalField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
Namespace for OpenFOAM.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
defineTypeNameAndDebug(combustionModel, 0)