laplaceFilter.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 "laplaceFilter.H"
27 #include "fvcLaplacian.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 Foam::laplaceFilter::laplaceFilter(const fvMesh& mesh, scalar widthCoeff)
42 :
43  LESfilter(mesh),
44  widthCoeff_(widthCoeff),
45  coeff_
46  (
47  IOobject
48  (
49  "laplaceFilterCoeff",
50  mesh.time().name(),
51  mesh
52  ),
53  mesh,
55  calculatedFvPatchScalarField::typeName
56  )
57 {
58  coeff_.ref() = pow(mesh.V(), 2.0/3.0)/widthCoeff_;
59 }
60 
61 
63 :
64  LESfilter(mesh),
65  widthCoeff_
66  (
67  bd.optionalSubDict(type() + "Coeffs").lookup<scalar>("widthCoeff")
68  ),
69  coeff_
70  (
71  IOobject
72  (
73  "laplaceFilterCoeff",
74  mesh.time().name(),
75  mesh
76  ),
77  mesh,
79  calculatedFvPatchScalarField::typeName
80  )
81 {
82  coeff_.ref() = pow(mesh.V(), 2.0/3.0)/widthCoeff_;
83 }
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 {
90  bd.optionalSubDict(type() + "Coeffs").lookup("widthCoeff") >> widthCoeff_;
91 }
92 
93 
94 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
95 
96 Foam::tmp<Foam::volScalarField> Foam::laplaceFilter::operator()
97 (
98  const tmp<volScalarField>& unFilteredField
99 ) const
100 {
101  correctBoundaryConditions(unFilteredField);
102 
103  tmp<volScalarField> filteredField =
104  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
105 
106  unFilteredField.clear();
107 
108  return filteredField;
109 }
110 
111 
112 Foam::tmp<Foam::volVectorField> Foam::laplaceFilter::operator()
113 (
114  const tmp<volVectorField>& unFilteredField
115 ) const
116 {
117  correctBoundaryConditions(unFilteredField);
118 
119  tmp<volVectorField> filteredField =
120  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
121 
122  unFilteredField.clear();
123 
124  return filteredField;
125 }
126 
127 
128 Foam::tmp<Foam::volSymmTensorField> Foam::laplaceFilter::operator()
129 (
130  const tmp<volSymmTensorField>& unFilteredField
131 ) const
132 {
133  correctBoundaryConditions(unFilteredField);
134 
135  tmp<volSymmTensorField> filteredField =
136  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
137 
138  unFilteredField.clear();
139 
140  return filteredField;
141 }
142 
143 
144 Foam::tmp<Foam::volTensorField> Foam::laplaceFilter::operator()
145 (
146  const tmp<volTensorField>& unFilteredField
147 ) const
148 {
149  correctBoundaryConditions(unFilteredField);
150 
151  tmp<volTensorField> filteredField =
152  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
153 
154  unFilteredField.clear();
155 
156  return filteredField;
157 }
158 
159 
160 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Internal & ref()
Return a reference to the dimensioned internal 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
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.
Laplace filter for LES.
Definition: laplaceFilter.H:59
laplaceFilter(const fvMesh &mesh, scalar widthCoeff)
Construct from components.
Definition: laplaceFilter.C:41
virtual void read(const dictionary &)
Read the LESfilter dictionary.
Definition: laplaceFilter.C:88
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the laplacian of the given field.
U correctBoundaryConditions()
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
Namespace for OpenFOAM.
const dimensionSet dimLength
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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)