Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
src
turbulenceModels
LES
LESfilters
simpleFilter
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 | Copyright (C) 1991-2010 OpenCFD Ltd.
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 <
OpenFOAM/addToRunTimeSelectionTable.H
>
28
#include <
finiteVolume/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
41
Foam::simpleFilter::simpleFilter
42
(
43
const
fvMesh
&
mesh
44
)
45
:
46
LESfilter
(mesh)
47
{}
48
49
50
Foam::simpleFilter::simpleFilter(
const
fvMesh
& mesh,
const
dictionary
&)
51
:
52
LESfilter
(mesh)
53
{}
54
55
56
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
57
58
void
Foam::simpleFilter::read
(
const
dictionary
&)
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
tmp<volScalarField>
filteredField =
fvc::surfaceSum
70
(
71
mesh
().magSf()*
fvc::interpolate
(unFilteredField)
72
)/
fvc::surfaceSum
(
mesh
().magSf());
73
74
unFilteredField.clear();
75
76
return
filteredField;
77
}
78
79
80
Foam::tmp<Foam::volVectorField>
Foam::simpleFilter::operator()
81
(
82
const
tmp<volVectorField>
& unFilteredField
83
)
const
84
{
85
tmp<volVectorField>
filteredField =
fvc::surfaceSum
86
(
87
mesh
().magSf()*
fvc::interpolate
(unFilteredField)
88
)/
fvc::surfaceSum
(
mesh
().magSf());
89
90
unFilteredField.clear();
91
92
return
filteredField;
93
}
94
95
96
Foam::tmp<Foam::volSymmTensorField>
Foam::simpleFilter::operator()
97
(
98
const
tmp<volSymmTensorField>
& unFilteredField
99
)
const
100
{
101
tmp<volSymmTensorField>
filteredField =
fvc::surfaceSum
102
(
103
mesh
().magSf()*
fvc::interpolate
(unFilteredField)
104
)/
fvc::surfaceSum
(
mesh
().magSf());
105
106
unFilteredField.clear();
107
108
return
filteredField;
109
}
110
111
112
Foam::tmp<Foam::volTensorField>
Foam::simpleFilter::operator()
113
(
114
const
tmp<volTensorField>
& unFilteredField
115
)
const
116
{
117
tmp<volTensorField>
filteredField =
fvc::surfaceSum
118
(
119
mesh
().magSf()*
fvc::interpolate
(unFilteredField)
120
)/
fvc::surfaceSum
(
mesh
().magSf());
121
122
unFilteredField.clear();
123
124
return
filteredField;
125
}
126
127
128
// ************************ vim: set sw=4 sts=4 et: ************************ //