proudmanAcousticPower.H
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) 2019-2021 OpenCFD Ltd.
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
26Class
27 Foam::functionObjects::proudmanAcousticPower
28
29Group
30 grpFieldFunctionObjects
31
32Description
33 Computes the acoustic power due to the volume of isotropic turbulence
34 using Proudman's formula.
35
36 The acoustic power, i.e. \f$ P_A \f$ [\f$W/m^3\f$], in terms of turbulent
37 kinetic energy, i.e. \f$ k \f$, and turbulent kinetic energy dissipation
38 rate, i.e. \f$ \epsilon \f$, is given as:
39
40 \f[
41 P_A = \alpha_\epsilon \rho \epsilon M_t^5
42 \f]
43
44 where \f$ \alpha_\epsilon = 0.1 \f$ is a constant and
45
46 \f[
47 M_t = \frac{\sqrt{2 k}}{a_0}
48 \f]
49
50 with \f$ a_0 \f$ the speed of sound. The acoustic power is also output in
51 dB using:
52
53 \f[
54 L_P = 10 \log \frac{P_A}{P_{ref}}
55 \f]
56
57 where \f$ P_{ref} = 1e^{-12} \f$ [\f$W/m^3\f$] is a constant.
58
59 Operands:
60 \table
61 Operand | Type | Location
62 input | volScalarField | $FOAM_CASE/<time>/<inpField>
63 output file | - | -
64 output field | volScalarField | $FOAM_CASE/<time>/<outField>
65 \endtable
66
67Usage
68 Minimal example by using \c system/controlDict.functions:
69 \verbatim
70 proudmanAcousticPower1
71 {
72 // Mandatory entries (unmodifiable)
73 type proudmanAcousticPower;
74 libs (fieldFunctionObjects);
75
76 // Optional entries (runtime modifiable)
77 alphaEps 0.1;
78 // For incompressible flow simulations
79 rhoInf 1.225;
80 aRef 340;
81 // Turbulence field names (if not retrieved from the turb model)
82 k kMean;
83 epsilon epsilonMean;
84 omega none; // omegaMean
85
86 // Optional (inherited) entries
87 ...
88 }
89 \endverbatim
90
91 where the entries mean:
92 \table
93 Property | Description | Type | Req'd | Dflt
94 type | Type name: proudmanAcousticPower | word | yes | -
95 libs | Library name: fieldFunctionObjects | word | yes | -
96 rhoInf | Freestream density (for incompressible) | scalar <!--
97 --> | conditional | -
98 aRef | Speed of sound (incompressible) | scalar <!--
99 --> | conditional | -
100 alphaEps | Empirical model coefficient | scalar | no | 0.1
101 k | Turbulence k field name | word | no | none
102 epsilon | Turbulence epsilon field name | word | no | none
103 omega | Turbulence omega field name | word | no | none
104 \endtable
105
106 The inherited entries are elaborated in:
107 - \link functionObject.H \endlink
108
109 Usage by the \c postProcess utility is not available.
110
111Note
112 The freestream density and reference speed of sound are only necessary
113 when a thermodynamics package is unavailable, typically for incompressible
114 cases.
115
116See also
117 - Foam::functionObject
118 - Foam::functionObjects::fvMeshFunctionObject
119 - ExtendedCodeGuide::functionObjects::field::proudmanAcousticPower
120
121SourceFiles
122 proudmanAcousticPower.C
123
124\*---------------------------------------------------------------------------*/
125
126#ifndef functionObjects_proudmanAcousticPower_H
127#define functionObjects_proudmanAcousticPower_H
128
129#include "fvMeshFunctionObject.H"
130#include "volFieldsFwd.H"
131
132// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133
134namespace Foam
135{
136namespace functionObjects
137{
138
139/*---------------------------------------------------------------------------*\
140 Class proudmanAcousticPower Declaration
141\*---------------------------------------------------------------------------*/
142
143class proudmanAcousticPower
144:
145 public fvMeshFunctionObject
146{
147 // Private Data
148
149 //- Empirical model coefficient
150 scalar alphaEps_;
151
152 //- Freestream density (incompressible calcs only)
153 dimensionedScalar rhoInf_;
154
155 //- Reference speed of sound (incompressible calcs only)
156 dimensionedScalar aRef_;
157
158 //- Name of turbulence k field; default = none
159 word kName_;
160
161 //- Name of turbulence epsilon field; default = none
162 word epsilonName_;
163
164 //- Name of turbulence omega field; default = none
165 word omegaName_;
166
167
168 // Private Member Functions
169
170 //- Multiply the field by density and return
171 tmp<volScalarField> rhoScale(const tmp<volScalarField>& fld) const;
172
173 //- Speed of sound
174 tmp<volScalarField> a() const;
175
176 //- Turbulence kinetic energy dissipation rate
177 tmp<volScalarField> k() const;
178
179 //- Turbulence dissipation
180 tmp<volScalarField> epsilon() const;
181
182
183public:
184
185 //- Runtime type information
186 TypeName("proudmanAcousticPower");
187
188
189 // Constructors
190
191 //- Construct from Time and dictionary
193 (
194 const word& name,
195 const Time& runTime,
196 const dictionary&
197 );
198
199 //- No copy construct
200 proudmanAcousticPower(const proudmanAcousticPower&) = delete;
201
202 //- No copy assignment
203 void operator=(const proudmanAcousticPower&) = delete;
204
205
206 //- Destructor
207 virtual ~proudmanAcousticPower() = default;
208
209
210 // Member Functions
211
212 //- Read the Proudman acoustic power data
213 virtual bool read(const dictionary&);
214
215 //- Calculate the Proudman acoustic power
216 virtual bool execute();
217
218 //- Write the Proudman acoustic power
219 virtual bool write();
220};
221
222
223// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224
225} // End namespace functionObjects
226} // End namespace Foam
227
228// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229
230#endif
231
232// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const word & name() const noexcept
Return the name of this functionObject.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Computes the acoustic power due to the volume of isotropic turbulence using Proudman's formula.
proudmanAcousticPower(const proudmanAcousticPower &)=delete
No copy construct.
virtual ~proudmanAcousticPower()=default
Destructor.
void operator=(const proudmanAcousticPower &)=delete
No copy assignment.
proudmanAcousticPower(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
virtual bool execute()
Calculate the Proudman acoustic power.
TypeName("proudmanAcousticPower")
Runtime type information.
virtual bool write()
Write the Proudman acoustic power.
virtual bool read(const dictionary &)
Read the Proudman acoustic power data.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
engineTime & runTime
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73