laplaceFilter.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 "laplaceFilter.H"
29 #include "fvm.H"
30 #include "fvc.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(laplaceFilter, 0);
37  addToRunTimeSelectionTable(LESfilter, laplaceFilter, dictionary);
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 Foam::laplaceFilter::laplaceFilter(const fvMesh& mesh, scalar widthCoeff)
44 :
45  LESfilter(mesh),
46  widthCoeff_(widthCoeff),
47  coeff_
48  (
49  IOobject
50  (
51  "laplaceFilterCoeff",
52  mesh.time().timeName(),
53  mesh
54  ),
55  mesh,
57  calculatedFvPatchScalarField::typeName
58  )
59 {
60  coeff_.ref() = pow(mesh.V(), 2.0/3.0)/widthCoeff_;
61 }
62 
63 
64 Foam::laplaceFilter::laplaceFilter(const fvMesh& mesh, const dictionary& bd)
65 :
66  LESfilter(mesh),
67  widthCoeff_(readScalar(bd.subDict(type() + "Coeffs").lookup("widthCoeff"))),
68  coeff_
69  (
70  IOobject
71  (
72  "laplaceFilterCoeff",
73  mesh.time().timeName(),
74  mesh
75  ),
76  mesh,
78  calculatedFvPatchScalarField::typeName
79  )
80 {
81  coeff_.ref() = pow(mesh.V(), 2.0/3.0)/widthCoeff_;
82 }
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
88 {
89  bd.subDict(type() + "Coeffs").lookup("widthCoeff") >> widthCoeff_;
90 }
91 
92 
93 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
94 
95 Foam::tmp<Foam::volScalarField> Foam::laplaceFilter::operator()
96 (
97  const tmp<volScalarField>& unFilteredField
98 ) const
99 {
100  correctBoundaryConditions(unFilteredField);
101 
102  tmp<volScalarField> filteredField =
103  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
104 
105  unFilteredField.clear();
106 
107  return filteredField;
108 }
109 
110 
111 Foam::tmp<Foam::volVectorField> Foam::laplaceFilter::operator()
112 (
113  const tmp<volVectorField>& unFilteredField
114 ) const
115 {
116  correctBoundaryConditions(unFilteredField);
117 
118  tmp<volVectorField> filteredField =
119  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
120 
121  unFilteredField.clear();
122 
123  return filteredField;
124 }
125 
126 
127 Foam::tmp<Foam::volSymmTensorField> Foam::laplaceFilter::operator()
128 (
129  const tmp<volSymmTensorField>& unFilteredField
130 ) const
131 {
132  correctBoundaryConditions(unFilteredField);
133 
134  tmp<volSymmTensorField> filteredField =
135  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
136 
137  unFilteredField.clear();
138 
139  return filteredField;
140 }
141 
142 
143 Foam::tmp<Foam::volTensorField> Foam::laplaceFilter::operator()
144 (
145  const tmp<volTensorField>& unFilteredField
146 ) const
147 {
148  correctBoundaryConditions(unFilteredField);
149 
150  tmp<volTensorField> filteredField =
151  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
152 
153  unFilteredField.clear();
154 
155  return filteredField;
156 }
157 
158 
159 // ************************************************************************* //
Abstract class for LES filters.
Definition: LESfilter.H:54
void correctBoundaryConditions(const tmp< GeoFieldType > &tgf) const
Temporary function to ensure the coupled boundary conditions of the.
Definition: LESfilter.H:79
virtual void read(const dictionary &)
Read the LESfilter dictionary.
Definition: laplaceFilter.C:87
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
const fvMesh & mesh() const
Return mesh reference.
Definition: LESfilter.H:134
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
Macros for easy insertion into run-time selection tables.
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
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)
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
A class for managing temporary objects.
Definition: PtrList.H:54
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451