Loading...
Searching...
No Matches
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{
38 defineTypeNameAndDebug(laplaceFilter, 0);
39 addToRunTimeSelectionTable(LESfilter, laplaceFilter, dictionary);
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45Foam::laplaceFilter::laplaceFilter(const fvMesh& mesh, scalar widthCoeff)
46:
48 widthCoeff_(widthCoeff),
49 coeff_
50 (
52 (
53 "laplaceFilterCoeff",
54 mesh.time().timeName(),
55 mesh
56 ),
57 mesh,
59 )
60{
61 coeff_.ref() = pow(mesh.V(), 2.0/3.0)/widthCoeff_;
62}
63
64
65Foam::laplaceFilter::laplaceFilter(const fvMesh& mesh, const dictionary& bd)
66:
68 widthCoeff_
69 (
70 bd.optionalSubDict(type() + "Coeffs").get<scalar>("widthCoeff")
71 ),
72 coeff_
73 (
75 (
76 "laplaceFilterCoeff",
77 mesh.time().timeName(),
78 mesh
79 ),
80 mesh,
82 )
84 coeff_.ref() = pow(mesh.V(), 2.0/3.0)/widthCoeff_;
85}
86
87
88// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89
92 bd.optionalSubDict(type() + "Coeffs").readEntry("widthCoeff", widthCoeff_);
93}
94
95
96// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
97
98Foam::tmp<Foam::volScalarField> Foam::laplaceFilter::operator()
99(
100 const tmp<volScalarField>& unFilteredField
101) const
102{
103 correctBoundaryConditions(unFilteredField);
104
105 tmp<volScalarField> filteredField =
106 unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
108 unFilteredField.clear();
109
110 return filteredField;
111}
112
113
114Foam::tmp<Foam::volVectorField> Foam::laplaceFilter::operator()
115(
116 const tmp<volVectorField>& unFilteredField
117) const
118{
119 correctBoundaryConditions(unFilteredField);
120
121 tmp<volVectorField> filteredField =
122 unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
123
124 unFilteredField.clear();
125
126 return filteredField;
127}
128
129
130Foam::tmp<Foam::volSymmTensorField> Foam::laplaceFilter::operator()
131(
132 const tmp<volSymmTensorField>& unFilteredField
133) const
134{
135 correctBoundaryConditions(unFilteredField);
136
137 tmp<volSymmTensorField> filteredField =
138 unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
140 unFilteredField.clear();
141
142 return filteredField;
143}
144
145
146Foam::tmp<Foam::volTensorField> Foam::laplaceFilter::operator()
147(
148 const tmp<volTensorField>& unFilteredField
149) const
150{
151 correctBoundaryConditions(unFilteredField);
152
153 tmp<volTensorField> filteredField =
154 unFilteredField() + fvc::laplacian(coeff_, unFilteredField());
155
156 unFilteredField.clear();
157
158 return filteredField;
159}
160
161
162// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Abstract class for LES filters.
Definition LESfilter.H:54
const fvMesh & mesh() const
Return mesh reference.
Definition LESfilter.H:151
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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:560
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Laplace filter for LES.
virtual void read(const dictionary &)
Read the LESfilter dictionary.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Namespace for OpenFOAM.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
cellMask correctBoundaryConditions()