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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 "laplaceFilter.H"
31 #include "fvm.H"
32 #include "fvc.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(laplaceFilter, 0);
39  addToRunTimeSelectionTable(LESfilter, laplaceFilter, dictionary);
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::laplaceFilter::laplaceFilter(const fvMesh& mesh, scalar widthCoeff)
46 :
47  LESfilter(mesh),
48  widthCoeff_(widthCoeff),
49  coeff_
50  (
51  IOobject
52  (
53  "laplaceFilterCoeff",
54  mesh.time().timeName(),
55  mesh
56  ),
57  mesh,
59  calculatedFvPatchScalarField::typeName
60  )
61 {
62  coeff_.ref() = pow(mesh.V(), 2.0/3.0)/widthCoeff_;
63 }
64 
65 
66 Foam::laplaceFilter::laplaceFilter(const fvMesh& mesh, const dictionary& bd)
67 :
68  LESfilter(mesh),
69  widthCoeff_
70  (
71  bd.optionalSubDict(type() + "Coeffs").get<scalar>("widthCoeff")
72  ),
73  coeff_
74  (
75  IOobject
76  (
77  "laplaceFilterCoeff",
78  mesh.time().timeName(),
79  mesh
80  ),
81  mesh,
83  calculatedFvPatchScalarField::typeName
84  )
85 {
86  coeff_.ref() = pow(mesh.V(), 2.0/3.0)/widthCoeff_;
87 }
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  bd.optionalSubDict(type() + "Coeffs").readEntry("widthCoeff", widthCoeff_);
95 }
96 
97 
98 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
99 
100 Foam::tmp<Foam::volScalarField> Foam::laplaceFilter::operator()
101 (
102  const tmp<volScalarField>& unFilteredField
103 ) const
104 {
105  correctBoundaryConditions(unFilteredField);
106 
107  tmp<volScalarField> filteredField =
108  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
109 
110  unFilteredField.clear();
111 
112  return filteredField;
113 }
114 
115 
116 Foam::tmp<Foam::volVectorField> Foam::laplaceFilter::operator()
117 (
118  const tmp<volVectorField>& unFilteredField
119 ) const
120 {
121  correctBoundaryConditions(unFilteredField);
122 
123  tmp<volVectorField> filteredField =
124  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
125 
126  unFilteredField.clear();
127 
128  return filteredField;
129 }
130 
131 
132 Foam::tmp<Foam::volSymmTensorField> Foam::laplaceFilter::operator()
133 (
134  const tmp<volSymmTensorField>& unFilteredField
135 ) const
136 {
137  correctBoundaryConditions(unFilteredField);
138 
139  tmp<volSymmTensorField> filteredField =
140  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
141 
142  unFilteredField.clear();
143 
144  return filteredField;
145 }
146 
147 
148 Foam::tmp<Foam::volTensorField> Foam::laplaceFilter::operator()
149 (
150  const tmp<volTensorField>& unFilteredField
151 ) const
152 {
153  correctBoundaryConditions(unFilteredField);
154 
155  tmp<volTensorField> filteredField =
156  unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
157 
158  unFilteredField.clear();
159 
160  return filteredField;
161 }
162 
163 
164 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::laplaceFilter::read
virtual void read(const dictionary &)
Read the LESfilter dictionary.
Definition: laplaceFilter.C:92
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
calculatedFvPatchFields.H
Foam::LESfilter::mesh
const fvMesh & mesh() const
Return mesh reference.
Definition: LESfilter.H:138
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:47
correctBoundaryConditions
cellMask correctBoundaryConditions()
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
timeName
word timeName
Definition: getTimeIndex.H:3
fvm.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::LESfilter
Abstract class for LES filters.
Definition: LESfilter.H:57
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
laplaceFilter.H
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:749
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::dictionary::optionalSubDict
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:577
fvc.H
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:179