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 | Copyright (C) 2011-2016 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 "anisotropicFilter.H"
29 #include "wallFvPatch.H"
30 #include "fvc.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(anisotropicFilter, 0);
37  addToRunTimeSelectionTable(LESfilter, anisotropicFilter, dictionary);
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 Foam::anisotropicFilter::anisotropicFilter
44 (
45  const fvMesh& mesh,
46  scalar widthCoeff
47 )
48 :
49  LESfilter(mesh),
50  widthCoeff_(widthCoeff),
51  coeff_
52  (
53  IOobject
54  (
55  "anisotropicFilterCoeff",
56  mesh.time().timeName(),
57  mesh
58  ),
59  mesh,
61  calculatedFvPatchVectorField::typeName
62  )
63 {
64  for (direction d=0; d<vector::nComponents; d++)
65  {
66  coeff_.primitiveFieldRef().replace
67  (
68  d,
69  (1/widthCoeff_)*
70  sqr
71  (
72  2.0*mesh.V()
73  /fvc::surfaceSum(mag(mesh.Sf().component(d)))().primitiveField()
74  )
75  );
76  }
77 }
78 
79 
80 Foam::anisotropicFilter::anisotropicFilter
81 (
82  const fvMesh& mesh,
83  const dictionary& bd
84 )
85 :
86  LESfilter(mesh),
87  widthCoeff_(readScalar(bd.subDict(type() + "Coeffs").lookup("widthCoeff"))),
88  coeff_
89  (
90  IOobject
91  (
92  "anisotropicFilterCoeff",
93  mesh.time().timeName(),
94  mesh
95  ),
96  mesh,
98  calculatedFvPatchScalarField::typeName
99  )
100 {
101  for (direction d=0; d<vector::nComponents; d++)
102  {
103  coeff_.primitiveFieldRef().replace
104  (
105  d,
106  (1/widthCoeff_)*
107  sqr
108  (
109  2.0*mesh.V()
110  /fvc::surfaceSum(mag(mesh.Sf().component(d)))().primitiveField()
111  )
112  );
113  }
114 }
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
120 {
121  bd.subDict(type() + "Coeffs").lookup("widthCoeff") >> widthCoeff_;
122 }
123 
124 
125 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
126 
127 Foam::tmp<Foam::volScalarField> Foam::anisotropicFilter::operator()
128 (
129  const tmp<volScalarField>& unFilteredField
130 ) const
131 {
132  correctBoundaryConditions(unFilteredField);
133 
134  tmp<volScalarField> tmpFilteredField =
135  unFilteredField
136  + (
137  coeff_
139  (
140  mesh().Sf()
141  *fvc::snGrad(unFilteredField())
142  )
143  );
144 
145  unFilteredField.clear();
146 
147  return tmpFilteredField;
148 }
149 
150 
151 Foam::tmp<Foam::volVectorField> Foam::anisotropicFilter::operator()
152 (
153  const tmp<volVectorField>& unFilteredField
154 ) const
155 {
156  correctBoundaryConditions(unFilteredField);
157 
158  tmp<volVectorField> tmpFilteredField =
159  unFilteredField
160  + (
161  coeff_
163  (
164  mesh().Sf()
165  *fvc::snGrad(unFilteredField())
166  )
167  );
168 
169  unFilteredField.clear();
170 
171  return tmpFilteredField;
172 }
173 
174 
175 Foam::tmp<Foam::volSymmTensorField> Foam::anisotropicFilter::operator()
176 (
177  const tmp<volSymmTensorField>& unFilteredField
178 ) const
179 {
180  correctBoundaryConditions(unFilteredField);
181 
182  tmp<volSymmTensorField> tmpFilteredField
183  (
185  (
186  IOobject
187  (
188  "anisotropicFilteredSymmTensorField",
189  mesh().time().timeName(),
190  mesh()
191  ),
192  mesh(),
193  unFilteredField().dimensions()
194  )
195  );
196 
197  for (direction d=0; d<symmTensor::nComponents; d++)
198  {
199  tmpFilteredField.ref().replace
200  (
201  d, anisotropicFilter::operator()(unFilteredField().component(d))
202  );
203  }
204 
205  unFilteredField.clear();
206 
207  return tmpFilteredField;
208 }
209 
210 
211 Foam::tmp<Foam::volTensorField> Foam::anisotropicFilter::operator()
212 (
213  const tmp<volTensorField>& unFilteredField
214 ) const
215 {
216  correctBoundaryConditions(unFilteredField);
217 
218  tmp<volTensorField> tmpFilteredField
219  (
220  new volTensorField
221  (
222  IOobject
223  (
224  "anisotropicFilteredTensorField",
225  mesh().time().timeName(),
226  mesh()
227  ),
228  mesh(),
229  unFilteredField().dimensions()
230  )
231  );
232 
233  for (direction d=0; d<tensor::nComponents; d++)
234  {
235  tmpFilteredField.ref().replace
236  (
237  d, anisotropicFilter::operator()(unFilteredField().component(d))
238  );
239  }
240 
241  unFilteredField.clear();
242 
243  return tmpFilteredField;
244 }
245 
246 
247 // ************************************************************************* //
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
Abstract class for LES filters.
Definition: LESfilter.H:54
uint8_t direction
Definition: direction.H:46
const surfaceVectorField & Sf() const
Return cell face area vectors.
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
Macros for easy insertion into run-time selection tables.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:96
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:91
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
U correctBoundaryConditions()
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void read(const dictionary &)
Read the LESfilter dictionary.
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for managing temporary objects.
Definition: PtrList.H:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451