Loading...
Searching...
No Matches
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{
38 defineTypeNameAndDebug(anisotropicFilter, 0);
39 addToRunTimeSelectionTable(LESfilter, anisotropicFilter, dictionary);
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45Foam::anisotropicFilter::anisotropicFilter
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 )
64{
65 for (direction d=0; d<vector::nComponents; d++)
66 {
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
81Foam::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 {
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
124 bd.optionalSubDict(type() + "Coeffs").readEntry("widthCoeff", widthCoeff_);
125}
126
127
128// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
129
130Foam::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 );
148 unFilteredField.clear();
149
150 return tmpFilteredField;
151}
152
153
154Foam::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
178Foam::tmp<Foam::volSymmTensorField> Foam::anisotropicFilter::operator()
179(
180 const tmp<volSymmTensorField>& unFilteredField
181) const
182{
183 correctBoundaryConditions(unFilteredField);
184
185 auto tmpFilteredField = volSymmTensorField::New
186 (
187 "anisotropicFilteredSymmTensorField",
189 mesh(),
190 unFilteredField().dimensions()
191 );
192
193 for (direction d=0; d<symmTensor::nComponents; d++)
194 {
195 tmpFilteredField.ref().replace
196 (
197 d, anisotropicFilter::operator()(unFilteredField().component(d))
198 );
199 }
201 unFilteredField.clear();
202
203 return tmpFilteredField;
204}
205
206
207Foam::tmp<Foam::volTensorField> Foam::anisotropicFilter::operator()
208(
209 const tmp<volTensorField>& unFilteredField
210) const
211{
212 correctBoundaryConditions(unFilteredField);
213
214 auto tmpFilteredField = volTensorField::New
215 (
216 "anisotropicFilteredTensorField",
218 mesh(),
219 unFilteredField().dimensions()
220 );
221
222 for (direction d=0; d<tensor::nComponents; d++)
223 {
224 tmpFilteredField.ref().replace
225 (
226 d, anisotropicFilter::operator()(unFilteredField().component(d))
227 );
228 }
229
230 unFilteredField.clear();
231
232 return tmpFilteredField;
233}
234
235
236// ************************************************************************* //
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:619
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< symmTensor >::calculatedType())
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Abstract class for LES filters.
Definition LESfilter.H:54
const fvMesh & mesh() const
Return mesh reference.
Definition LESfilter.H:151
static constexpr direction nComponents
Number of components in this vector space.
virtual void read(const dictionary &)
Read the LESfilter dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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:560
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
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:40
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)
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 POSIX.C:801
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
uint8_t direction
Definition direction.H:49
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
cellMask correctBoundaryConditions()