anisotropicFilter.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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"
27 #include "fvcSurfaceIntegrate.H"
28 #include "fvcSnGrad.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const fvMesh& mesh,
45  scalar widthCoeff
46 )
47 :
48  LESfilter(mesh),
49  widthCoeff_(widthCoeff),
50  coeff_
51  (
52  IOobject
53  (
54  "anisotropicFilterCoeff",
55  mesh.time().name(),
56  mesh
57  ),
58  mesh,
60  calculatedFvPatchVectorField::typeName
61  )
62 {
63  for (direction d=0; d<vector::nComponents; d++)
64  {
65  coeff_.primitiveFieldRef().replace
66  (
67  d,
68  (1/widthCoeff_)*
69  sqr
70  (
71  2.0*mesh.V()
72  /fvc::surfaceSum(mag(mesh.Sf().component(d)))().primitiveField()
73  )
74  );
75  }
76 }
77 
78 
80 (
81  const fvMesh& mesh,
82  const dictionary& bd
83 )
84 :
85  LESfilter(mesh),
86  widthCoeff_
87  (
88  bd.optionalSubDict(type() + "Coeffs").lookup<scalar>("widthCoeff")
89  ),
90  coeff_
91  (
92  IOobject
93  (
94  "anisotropicFilterCoeff",
95  mesh.time().name(),
96  mesh
97  ),
98  mesh,
100  calculatedFvPatchScalarField::typeName
101  )
102 {
103  for (direction d=0; d<vector::nComponents; d++)
104  {
105  coeff_.primitiveFieldRef().replace
106  (
107  d,
108  (1/widthCoeff_)*
109  sqr
110  (
111  2.0*mesh.V()
112  /fvc::surfaceSum(mag(mesh.Sf().component(d)))().primitiveField()
113  )
114  );
115  }
116 }
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 {
123  bd.optionalSubDict(type() + "Coeffs").lookup("widthCoeff") >> widthCoeff_;
124 }
125 
126 
127 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
128 
129 Foam::tmp<Foam::volScalarField> Foam::anisotropicFilter::operator()
130 (
131  const tmp<volScalarField>& unFilteredField
132 ) const
133 {
134  correctBoundaryConditions(unFilteredField);
135 
136  tmp<volScalarField> tmpFilteredField =
137  unFilteredField
138  + (
139  coeff_
141  (
142  mesh().Sf()
143  *fvc::snGrad(unFilteredField())
144  )
145  );
146 
147  unFilteredField.clear();
148 
149  return tmpFilteredField;
150 }
151 
152 
153 Foam::tmp<Foam::volVectorField> Foam::anisotropicFilter::operator()
154 (
155  const tmp<volVectorField>& unFilteredField
156 ) const
157 {
158  correctBoundaryConditions(unFilteredField);
159 
160  tmp<volVectorField> tmpFilteredField =
161  unFilteredField
162  + (
163  coeff_
165  (
166  mesh().Sf()
167  *fvc::snGrad(unFilteredField())
168  )
169  );
170 
171  unFilteredField.clear();
172 
173  return tmpFilteredField;
174 }
175 
176 
177 Foam::tmp<Foam::volSymmTensorField> Foam::anisotropicFilter::operator()
178 (
179  const tmp<volSymmTensorField>& unFilteredField
180 ) const
181 {
182  correctBoundaryConditions(unFilteredField);
183 
184  tmp<volSymmTensorField> tmpFilteredField
185  (
187  (
188  "anisotropicFilteredSymmTensorField",
189  mesh(),
190  unFilteredField().dimensions()
191  )
192  );
193 
194  for (direction d=0; d<symmTensor::nComponents; d++)
195  {
196  tmpFilteredField.ref().replace
197  (
198  d, anisotropicFilter::operator()(unFilteredField().component(d))
199  );
200  }
201 
202  unFilteredField.clear();
203 
204  return tmpFilteredField;
205 }
206 
207 
208 Foam::tmp<Foam::volTensorField> Foam::anisotropicFilter::operator()
209 (
210  const tmp<volTensorField>& unFilteredField
211 ) const
212 {
213  correctBoundaryConditions(unFilteredField);
214 
215  tmp<volTensorField> tmpFilteredField
216  (
218  (
219  "anisotropicFilteredTensorField",
220  mesh(),
221  unFilteredField().dimensions()
222  )
223  );
224 
225  for (direction d=0; d<tensor::nComponents; d++)
226  {
227  tmpFilteredField.ref().replace
228  (
229  d, anisotropicFilter::operator()(unFilteredField().component(d))
230  );
231  }
232 
233  unFilteredField.clear();
234 
235  return tmpFilteredField;
236 }
237 
238 
239 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:467
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Abstract class for LES filters.
Definition: LESfilter.H:55
const fvMesh & mesh() const
Return mesh reference.
Definition: LESfilter.H:130
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
anisotropic filter
anisotropicFilter(const fvMesh &mesh, scalar widthCoeff)
Construct from components.
virtual void read(const dictionary &)
Read the LESfilter dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1076
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const surfaceVectorField & Sf() const
Return cell face area vectors.
A class for managing temporary objects.
Definition: tmp.H:55
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Calculate the snGrad of the given volField.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
U correctBoundaryConditions()
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
void surfaceIntegrate(Field< Type > &ivf, const SurfaceField< Type > &ssf)
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
const dimensionSet dimLength
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
uint8_t direction
Definition: direction.H:45
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)