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 -------------------------------------------------------------------------------
10 License
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 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(anisotropicFilter, 0);
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::anisotropicFilter::anisotropicFilter
46 (
47  const fvMesh& mesh,
48  scalar widthCoeff
49 )
50 :
51  LESfilter(mesh),
52  widthCoeff_(widthCoeff),
53  coeff_
54  (
55  IOobject
56  (
57  "anisotropicFilterCoeff",
58  mesh.time().timeName(),
59  mesh
60  ),
61  mesh,
63  )
64 {
65  for (direction d=0; d<vector::nComponents; d++)
66  {
67  coeff_.primitiveFieldRef().replace
68  (
69  d,
70  (1/widthCoeff_)*
71  sqr
72  (
73  2.0*mesh.V()
74  /fvc::surfaceSum(mag(mesh.Sf().component(d)))().primitiveField()
75  )
76  );
77  }
78 }
79 
80 
81 Foam::anisotropicFilter::anisotropicFilter
82 (
83  const fvMesh& mesh,
84  const dictionary& bd
85 )
86 :
87  LESfilter(mesh),
88  widthCoeff_
89  (
90  bd.optionalSubDict(type() + "Coeffs").get<scalar>("widthCoeff")
91  ),
92  coeff_
93  (
94  IOobject
95  (
96  "anisotropicFilterCoeff",
97  mesh.time().timeName(),
98  mesh
99  ),
100  mesh,
102  )
103 {
104  for (direction d=0; d<vector::nComponents; d++)
105  {
106  coeff_.primitiveFieldRef().replace
107  (
108  d,
109  (1/widthCoeff_)*
110  sqr
111  (
112  2.0*mesh.V()
113  /fvc::surfaceSum(mag(mesh.Sf().component(d)))().primitiveField()
114  )
115  );
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 
122 void Foam::anisotropicFilter::read(const dictionary& bd)
123 {
124  bd.optionalSubDict(type() + "Coeffs").readEntry("widthCoeff", widthCoeff_);
125 }
126 
127 
128 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
129 
130 Foam::tmp<Foam::volScalarField> Foam::anisotropicFilter::operator()
131 (
132  const tmp<volScalarField>& unFilteredField
133 ) const
134 {
135  correctBoundaryConditions(unFilteredField);
136 
137  tmp<volScalarField> tmpFilteredField =
138  unFilteredField
139  + (
140  coeff_
142  (
143  mesh().Sf()
144  *fvc::snGrad(unFilteredField())
145  )
146  );
147 
148  unFilteredField.clear();
149 
150  return tmpFilteredField;
151 }
152 
153 
154 Foam::tmp<Foam::volVectorField> Foam::anisotropicFilter::operator()
155 (
156  const tmp<volVectorField>& unFilteredField
157 ) const
158 {
159  correctBoundaryConditions(unFilteredField);
160 
161  tmp<volVectorField> tmpFilteredField =
162  unFilteredField
163  + (
164  coeff_
166  (
167  mesh().Sf()
168  *fvc::snGrad(unFilteredField())
169  )
170  );
171 
172  unFilteredField.clear();
173 
174  return tmpFilteredField;
175 }
176 
177 
178 Foam::tmp<Foam::volSymmTensorField> Foam::anisotropicFilter::operator()
179 (
180  const tmp<volSymmTensorField>& unFilteredField
181 ) const
182 {
183  correctBoundaryConditions(unFilteredField);
184 
185  tmp<volSymmTensorField> tmpFilteredField
186  (
188  (
189  IOobject
190  (
191  "anisotropicFilteredSymmTensorField",
192  mesh().time().timeName(),
193  mesh()
194  ),
195  mesh(),
196  unFilteredField().dimensions()
197  )
198  );
199 
200  for (direction d=0; d<symmTensor::nComponents; d++)
201  {
202  tmpFilteredField.ref().replace
203  (
204  d, anisotropicFilter::operator()(unFilteredField().component(d))
205  );
206  }
207 
208  unFilteredField.clear();
209 
210  return tmpFilteredField;
211 }
212 
213 
214 Foam::tmp<Foam::volTensorField> Foam::anisotropicFilter::operator()
215 (
216  const tmp<volTensorField>& unFilteredField
217 ) const
218 {
219  correctBoundaryConditions(unFilteredField);
220 
221  tmp<volTensorField> tmpFilteredField
222  (
223  new volTensorField
224  (
225  IOobject
226  (
227  "anisotropicFilteredTensorField",
228  mesh().time().timeName(),
229  mesh()
230  ),
231  mesh(),
232  unFilteredField().dimensions()
233  )
234  );
235 
236  for (direction d=0; d<tensor::nComponents; d++)
237  {
238  tmpFilteredField.ref().replace
239  (
240  d, anisotropicFilter::operator()(unFilteredField().component(d))
241  );
242  }
243 
244  unFilteredField.clear();
245 
246  return tmpFilteredField;
247 }
248 
249 
250 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Abstract class for LES filters.
Definition: LESfilter.H:53
const surfaceVectorField & Sf() const
Return cell face area vectors.
uint8_t direction
Definition: direction.H:46
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const fvMesh & mesh() const
Return mesh reference.
Definition: LESfilter.H:151
cellMask correctBoundaryConditions()
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:620
Macros for easy insertion into run-time selection tables.
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
anisotropic filter
dynamicFvMesh & mesh
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:109
defineTypeNameAndDebug(combustionModel, 0)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void read(const dictionary &)
Read the LESfilter dictionary.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition: tmpI.H:289
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
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:40
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Replace specified field component with content from another field.