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-2013 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 
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 {
61 }
62 
63 
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 {
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  tmp<volScalarField> filteredField =
101  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
102 
103  unFilteredField.clear();
104 
105  return filteredField;
106 }
107 
108 
109 Foam::tmp<Foam::volVectorField> Foam::laplaceFilter::operator()
110 (
111  const tmp<volVectorField>& unFilteredField
112 ) const
113 {
114  tmp<volVectorField> filteredField =
115  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
116 
117  unFilteredField.clear();
118 
119  return filteredField;
120 }
121 
122 
123 Foam::tmp<Foam::volSymmTensorField> Foam::laplaceFilter::operator()
124 (
125  const tmp<volSymmTensorField>& unFilteredField
126 ) const
127 {
128  tmp<volSymmTensorField> filteredField =
129  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
130 
131  unFilteredField.clear();
132 
133  return filteredField;
134 }
135 
136 
137 Foam::tmp<Foam::volTensorField> Foam::laplaceFilter::operator()
138 (
139  const tmp<volTensorField>& unFilteredField
140 ) const
141 {
142  tmp<volTensorField> filteredField =
143  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
144 
145  unFilteredField.clear();
146 
147  return filteredField;
148 }
149 
150 
151 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::laplaceFilter::read
virtual void read(const dictionary &)
Read the LESfilter dictionary.
Definition: laplaceFilter.C:87
Foam::GeometricField::dimensionedInternalField
DimensionedInternalField & dimensionedInternalField()
Return dimensioned internal field.
Definition: GeometricField.C:713
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
calculatedFvPatchFields.H
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::laplaceFilter::widthCoeff_
scalar widthCoeff_
Definition: laplaceFilter.H:62
Foam::LESfilter::mesh
const fvMesh & mesh() const
Return mesh reference.
Definition: LESfilter.H:117
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:199
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
fvm.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::LESfilter
Abstract class for LES filters.
Definition: LESfilter.H:54
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
laplaceFilter.H
Foam::laplaceFilter::laplaceFilter
laplaceFilter(const laplaceFilter &)
Disallow default bitwise copy construct and assignment.
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::laplaceFilter::coeff_
volScalarField coeff_
Definition: laplaceFilter.H:63
fvc.H
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress