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-------------------------------------------------------------------------------
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 "laplaceFilter.H"
31#include "fvm.H"
32#include "fvc.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
46:
48 widthCoeff_(widthCoeff),
49 coeff_
50 (
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
67:
69 widthCoeff_
70 (
71 bd.optionalSubDict(type() + "Coeffs").get<scalar>("widthCoeff")
72 ),
73 coeff_
74 (
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
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
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
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
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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Abstract class for LES filters.
Definition: LESfilter.H:58
const fvMesh & mesh() const
Return mesh reference.
Definition: LESfilter.H:138
virtual bool read()
Re-read model coefficients if they have changed.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
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
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Laplace filter for LES.
Definition: laplaceFilter.H:61
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:47
Namespace for OpenFOAM.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
cellMask correctBoundaryConditions()