anisotropicFilter.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-2017 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 "anisotropicFilter.H"
31#include "wallFvPatch.H"
32#include "fvc.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
46(
47 const fvMesh& mesh,
48 scalar widthCoeff
49)
50:
52 widthCoeff_(widthCoeff),
53 coeff_
54 (
56 (
57 "anisotropicFilterCoeff",
58 mesh.time().timeName(),
59 mesh
60 ),
61 mesh,
63 calculatedFvPatchVectorField::typeName
64 )
65{
66 for (direction d=0; d<vector::nComponents; d++)
67 {
69 (
70 d,
71 (1/widthCoeff_)*
72 sqr
73 (
74 2.0*mesh.V()
75 /fvc::surfaceSum(mag(mesh.Sf().component(d)))().primitiveField()
76 )
77 );
78 }
79}
80
81
83(
84 const fvMesh& mesh,
85 const dictionary& bd
86)
87:
89 widthCoeff_
90 (
91 bd.optionalSubDict(type() + "Coeffs").get<scalar>("widthCoeff")
92 ),
93 coeff_
94 (
96 (
97 "anisotropicFilterCoeff",
98 mesh.time().timeName(),
99 mesh
100 ),
101 mesh,
103 calculatedFvPatchScalarField::typeName
104 )
105{
106 for (direction d=0; d<vector::nComponents; d++)
107 {
109 (
110 d,
111 (1/widthCoeff_)*
112 sqr
113 (
114 2.0*mesh.V()
115 /fvc::surfaceSum(mag(mesh.Sf().component(d)))().primitiveField()
116 )
117 );
118 }
119}
120
121
122// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123
125{
126 bd.optionalSubDict(type() + "Coeffs").readEntry("widthCoeff", widthCoeff_);
127}
128
129
130// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
131
133(
134 const tmp<volScalarField>& unFilteredField
135) const
136{
137 correctBoundaryConditions(unFilteredField);
138
139 tmp<volScalarField> tmpFilteredField =
140 unFilteredField
141 + (
142 coeff_
144 (
145 mesh().Sf()
146 *fvc::snGrad(unFilteredField())
147 )
148 );
149
150 unFilteredField.clear();
151
152 return tmpFilteredField;
153}
154
155
157(
158 const tmp<volVectorField>& unFilteredField
159) const
160{
161 correctBoundaryConditions(unFilteredField);
162
163 tmp<volVectorField> tmpFilteredField =
164 unFilteredField
165 + (
166 coeff_
168 (
169 mesh().Sf()
170 *fvc::snGrad(unFilteredField())
171 )
172 );
173
174 unFilteredField.clear();
175
176 return tmpFilteredField;
177}
178
179
181(
182 const tmp<volSymmTensorField>& unFilteredField
183) const
184{
185 correctBoundaryConditions(unFilteredField);
186
187 tmp<volSymmTensorField> tmpFilteredField
188 (
190 (
192 (
193 "anisotropicFilteredSymmTensorField",
194 mesh().time().timeName(),
195 mesh()
196 ),
197 mesh(),
198 unFilteredField().dimensions()
199 )
200 );
201
202 for (direction d=0; d<symmTensor::nComponents; d++)
203 {
204 tmpFilteredField.ref().replace
205 (
206 d, anisotropicFilter::operator()(unFilteredField().component(d))
207 );
208 }
209
210 unFilteredField.clear();
211
212 return tmpFilteredField;
213}
214
215
217(
218 const tmp<volTensorField>& unFilteredField
219) const
220{
221 correctBoundaryConditions(unFilteredField);
222
223 tmp<volTensorField> tmpFilteredField
224 (
226 (
228 (
229 "anisotropicFilteredTensorField",
230 mesh().time().timeName(),
231 mesh()
232 ),
233 mesh(),
234 unFilteredField().dimensions()
235 )
236 );
237
238 for (direction d=0; d<tensor::nComponents; d++)
239 {
240 tmpFilteredField.ref().replace
241 (
242 d, anisotropicFilter::operator()(unFilteredField().component(d))
243 );
244 }
245
246 unFilteredField.clear();
247
248 return tmpFilteredField;
249}
250
251
252// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:557
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Abstract class for LES filters.
Definition: LESfilter.H:58
const fvMesh & mesh() const
Return mesh reference.
Definition: LESfilter.H:138
virtual bool read()
Re-read model coefficients if they have changed.
anisotropic filter
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:577
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const surfaceVectorField & Sf() const
Return cell face area vectors.
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
T & ref() const
Definition: tmpI.H:227
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
uint8_t direction
Definition: direction.H:56
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
cellMask correctBoundaryConditions()