Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
34namespace Foam
35{
36 defineTypeNameAndDebug(simpleFilter, 0);
37 addToRunTimeSelectionTable(LESfilter, simpleFilter, dictionary);
38}
39
40
41// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42
43Foam::simpleFilter::simpleFilter
44(
45 const fvMesh& mesh
46)
47:
48 LESfilter(mesh)
49{}
50
51
52Foam::simpleFilter::simpleFilter(const fvMesh& mesh, const dictionary&)
54 LESfilter(mesh)
55{}
56
57
58// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61{}
62
63
64// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
65
66Foam::tmp<Foam::volScalarField> Foam::simpleFilter::operator()
67(
68 const tmp<volScalarField>& unFilteredField
69) const
70{
71 correctBoundaryConditions(unFilteredField);
72
74 (
75 mesh().magSf()*fvc::interpolate(unFilteredField)
76 )/fvc::surfaceSum(mesh().magSf());
78 unFilteredField.clear();
79
80 return filteredField;
81}
82
83
84Foam::tmp<Foam::volVectorField> Foam::simpleFilter::operator()
85(
86 const tmp<volVectorField>& unFilteredField
87) const
88{
89 correctBoundaryConditions(unFilteredField);
90
92 (
93 mesh().magSf()*fvc::interpolate(unFilteredField)
94 )/fvc::surfaceSum(mesh().magSf());
95
96 unFilteredField.clear();
97
98 return filteredField;
99}
100
101
102Foam::tmp<Foam::volSymmTensorField> Foam::simpleFilter::operator()
103(
104 const tmp<volSymmTensorField>& unFilteredField
105) const
106{
107 correctBoundaryConditions(unFilteredField);
108
109 tmp<volSymmTensorField> filteredField = fvc::surfaceSum
110 (
111 mesh().magSf()*fvc::interpolate(unFilteredField)
112 )/fvc::surfaceSum(mesh().magSf());
114 unFilteredField.clear();
115
116 return filteredField;
117}
118
119
120Foam::tmp<Foam::volTensorField> Foam::simpleFilter::operator()
121(
122 const tmp<volTensorField>& unFilteredField
123) const
124{
125 correctBoundaryConditions(unFilteredField);
126
128 (
129 mesh().magSf()*fvc::interpolate(unFilteredField)
130 )/fvc::surfaceSum(mesh().magSf());
131
132 unFilteredField.clear();
133
134 return filteredField;
135}
136
137
138// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Abstract class for LES filters.
Definition LESfilter.H:54
const fvMesh & mesh() const
Return mesh reference.
Definition LESfilter.H:151
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Simple top-hat filter used in dynamic LES models.
virtual void read(const dictionary &)
Read the LESfilter dictionary.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
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.
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Namespace for OpenFOAM.
cellMask correctBoundaryConditions()