noiseModel.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) 2015-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::noiseModel
28
29Description
30 Base class for noise models.
31
32 Data is read from a dictionary, e.g.
33
34 \verbatim
35 rhoRef 1;
36 N 4096;
37 fl 25;
38 fu 10000;
39 startTime 0;
40
41 outputPrefix "test1";
42
43 SPLweighting dBA;
44
45 // Optional write options dictionary
46 writeOptions
47 {
48 writePrmsf no;
49 writeSPL yes;
50 writePSD yes;
51 writePSDf no;
52 writeOctaves yes;
53 }
54 \endverbatim
55
56 where
57 \table
58 Property | Description | Required | Default value
59 rhoRef | Reference density | no | 1
60 N | Number of samples in sampling window | no | 65536 (2^16)
61 fl | Lower frequency bounds | no | 25
62 fu | Upper frequency bounds | no | 10000
63 startTime | Start time | no | 0
64 outputPrefix | Prefix applied to output files| no | ''
65 SPLweighting | Weighting: dBA, dBB, dBC, DBD | no | none
66 dBRef | Reference for dB calculation | no | 2e-5
67 graphFormat | Graph format | no | raw
68 writePrmsf | Write Prmsf data | no | yes
69 writeSPL | Write SPL data | no | yes
70 writePSD | Write PSD data | no | yes
71 writePSDf | Write PSDf data | no | yes
72 writeOctaves | Write octaves data | no | yes
73 \endtable
74
75SourceFiles
76 noiseModel.C
77
78\*---------------------------------------------------------------------------*/
79
80#ifndef noiseModel_H
81#define noiseModel_H
82
83#include "dictionary.H"
84#include "scalarList.H"
85#include "instantList.H"
86#include "windowModel.H"
87#include "Enum.H"
89#include "graph.H"
90#include <fftw3.h>
91
92// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
93
94namespace Foam
95{
96
97/*---------------------------------------------------------------------------*\
98 Class noiseModel Declaration
99\*---------------------------------------------------------------------------*/
100
101class noiseModel
102{
103public:
104
105 enum class weightingType
106 {
107 none,
108 dBA,
109 dBB,
110 dBC,
111 dBD
112 };
113
114 static const Enum<weightingType> weightingTypeNames_;
115
116 //- FFTW planner information
117 // Note: storage uses double for use directly with FFTW
118 struct planInfo
119 {
120 bool active;
121 label windowSize;
122 List<double> in;
123 List<double> out;
124 fftw_plan plan;
125 };
126
127 //- Octave band information
128 struct octaveBandInfo
129 {
130 label octave;
131
132 // IDs of bin boundaries in pressure data
134
135 // Centre frequencies for each bin
137 };
138
139
140private:
141
142 // Private Member Functions
143
144 //- No copy construct
145 noiseModel(const noiseModel&) = delete;
146
147 //- No copy assignment
148 void operator=(const noiseModel&) = delete;
149
150
151protected:
152
153 // Protected Data
154
155 //- Copy of dictionary used for construction
156 const dictionary dict_;
157
158 //- Reference density (to convert from kinematic to static pressure)
159 scalar rhoRef_;
160
161 //- Number of samples in sampling window, default = 2^16
162 label nSamples_;
163
164 //- Lower frequency limit, default = 25Hz
165 scalar fLower_;
166
167 //- Upper frequency limit, default = 10kHz
168 scalar fUpper_;
169
170 //- Start time, default = 0s
171 scalar startTime_;
172
173 //- Window model
174 autoPtr<windowModel> windowModelPtr_;
176 //- Graph format
178
179 //- Weighting
181
182 //- Reference for dB calculation, default = 2e-5
183 scalar dBRef_;
184
185
186 // Data validation
187
188 //- Min pressure value
189 scalar minPressure_;
190
191 //- Min pressure value
193
195 // Write options
197 //- Output file prefix, default = ''
199
200 //- Write Prmsf; default = yes
201 bool writePrmsf_;
203 //- Write SPL; default = yes
205
206 //- Write PSD; default = yes
208
209 //- Write PSDf; default = yes
211
212 //- Write writeOctaves; default = yes
213 bool writeOctaves_;
214
215
216 // FFT
217
218 //- Plan information for FFTW
219 mutable planInfo planInfo_;
220
221
222 // Protected Member Functions
223
224 //- Helper function to read write options and provide info feedback
225 void readWriteOption
226 (
227 const dictionary& dict,
228 const word& lookup,
229 bool& option
230 ) const;
231
232 //- Check and return uniform time step
234 (
235 const scalarList& times
236 ) const;
237
238 //- Return true if all pressure data is within min/max bounds
239 bool validateBounds(const scalarList& p) const;
240
241 //- Find and return start time index
243 (
244 const instantList& allTimes,
245 const scalar startTime
246 ) const;
247
248 //- Return the base output directory
249 fileName baseFileDir(const label dataseti) const;
250
251 //- Create a field of equally spaced frequencies for the current set of
252 //- data - assumes a constant time step
255 const scalar deltaT,
256 const bool check
257 ) const;
258
259 //- Return a list of the frequency indices wrt f field that correspond
260 //- to the bands limits for a given octave
261 static void setOctaveBands
262 (
264 const scalar fLower,
265 const scalar fUpper,
266 const scalar octave,
267 labelList& fBandIDs,
268 scalarField& fCentre
269 );
270
271 //- Generate octave data
273 (
274 const scalarField& data,
276 const labelUList& freqBandIDs
277 ) const;
279 //- Return the fft of the given pressure data
280 tmp<scalarField> Pf(const scalarField& p) const;
282 //- Return the multi-window mean fft of the complete pressure data [Pa]
283 tmp<scalarField> meanPf(const scalarField& p) const;
285 //- Return the multi-window RMS mean fft of the complete pressure
286 //- data [Pa]
288
289 //- Return the multi-window Power Spectral Density (PSD) of the complete
290 //- pressure data [Pa^2/Hz]
291 tmp<scalarField> PSDf(const scalarField& p, const scalar deltaT) const;
292
294 // Weightings
295
296 //- A weighting function
297 scalar RAf(const scalar f) const;
298
299 //- A weighting as gain in dB
300 scalar gainA(const scalar f) const;
301
302 //- B weighting function
303 scalar RBf(const scalar f) const;
304
305 //- B weighting as gain in dB
306 scalar gainB(const scalar f) const;
307
308 //- C weighting function
309 scalar RCf(const scalar f) const;
310
311 //- C weighting as gain in dB
312 scalar gainC(const scalar f) const;
313
314 //- D weighting function
315 scalar RDf(const scalar f) const;
316
317 //- D weighting as gain in dB
318 scalar gainD(const scalar f) const;
319
320
321public:
322
323 //- Runtime type information
324 TypeName("noiseModel");
325
326 //- Run time selection table
328 (
329 autoPtr,
332 (
333 const dictionary& dict
334 ),
335 (dict)
336 );
337
338 //- Selector
339 static autoPtr<noiseModel> New(const dictionary& dict);
340
341 //- Constructor
342 noiseModel(const dictionary& dict, const bool readFields = true);
343
344 //- Destructor
345 virtual ~noiseModel() = default;
346
347
348 // Public Member Functions
349
350 //- Read from dictionary
351 virtual bool read(const dictionary& dict);
352
353 //- Abstract call to calculate
354 virtual void calculate() = 0;
355
356 //- PSD [dB/Hz]
358
359 //- SPL [dB]
361 (
362 const scalarField& Prms2,
363 const scalar f
364 ) const;
365
366 //- SPL [dB]
368 (
369 const scalarField& Prms2,
370 const scalarField& f
371 ) const;
372
373 //- Clean up the FFTW
374 void cleanFFTW();
375
376 //- Helper function to check weightings
377 void writeWeightings() const;
378};
379
380
381// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
382
383} // End namespace Foam
384
385// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
386
387#endif
388
389// ************************************************************************* //
Foam::label startTime
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Database for solution data, solver performance and other reduced data.
Definition: data.H:58
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A class for handling file names.
Definition: fileName.H:76
Base class for noise models.
Definition: noiseModel.H:176
tmp< scalarField > uniformFrequencies(const scalar deltaT, const bool check) const
Definition: noiseModel.C:277
label findStartTimeIndex(const instantList &allTimes, const scalar startTime) const
Find and return start time index.
Definition: noiseModel.C:243
scalar RAf(const scalar f) const
A weighting function.
Definition: noiseModel.C:483
scalar rhoRef_
Reference density (to convert from kinematic to static pressure)
Definition: noiseModel.H:233
label nSamples_
Number of samples in sampling window, default = 2^16.
Definition: noiseModel.H:236
tmp< scalarField > octaves(const scalarField &data, const scalarField &f, const labelUList &freqBandIDs) const
Generate octave data.
Definition: noiseModel.C:314
tmp< scalarField > RMSmeanPf(const scalarField &p) const
Definition: noiseModel.C:429
virtual void calculate()=0
Abstract call to calculate.
scalar fUpper_
Upper frequency limit, default = 10kHz.
Definition: noiseModel.H:242
scalar gainC(const scalar f) const
C weighting as gain in dB.
Definition: noiseModel.C:549
tmp< Foam::scalarField > PSD(const scalarField &PSDf) const
PSD [dB/Hz].
Definition: noiseModel.C:713
static autoPtr< noiseModel > New(const dictionary &dict)
Selector.
Definition: noiseModelNew.C:32
tmp< scalarField > meanPf(const scalarField &p) const
Return the multi-window mean fft of the complete pressure data [Pa].
Definition: noiseModel.C:408
scalar maxPressure_
Min pressure value.
Definition: noiseModel.H:266
bool writePSD_
Write PSD; default = yes.
Definition: noiseModel.H:281
scalar gainA(const scalar f) const
A weighting as gain in dB.
Definition: noiseModel.C:500
scalar minPressure_
Min pressure value.
Definition: noiseModel.H:263
word graphFormat_
Graph format.
Definition: noiseModel.H:251
static void setOctaveBands(const scalarField &f, const scalar fLower, const scalar fUpper, const scalar octave, labelList &fBandIDs, scalarField &fCentre)
Definition: noiseModel.C:55
fileName baseFileDir(const label dataseti) const
Return the base output directory.
Definition: noiseModel.C:262
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: noiseModel.C:624
bool writeOctaves_
Write writeOctaves; default = yes.
Definition: noiseModel.H:287
fileName outputPrefix_
Output file prefix, default = ''.
Definition: noiseModel.H:272
scalar RDf(const scalar f) const
D weighting function.
Definition: noiseModel.C:560
tmp< scalarField > Pf(const scalarField &p) const
Return the fft of the given pressure data.
Definition: noiseModel.C:365
scalar RCf(const scalar f) const
C weighting function.
Definition: noiseModel.C:538
tmp< scalarField > SPL(const scalarField &Prms2, const scalar f) const
SPL [dB].
Definition: noiseModel.C:722
planInfo planInfo_
Plan information for FFTW.
Definition: noiseModel.H:293
autoPtr< windowModel > windowModelPtr_
Window model.
Definition: noiseModel.H:248
static const Enum< weightingType > weightingTypeNames_
Definition: noiseModel.H:188
bool writeSPL_
Write SPL; default = yes.
Definition: noiseModel.H:278
TypeName("noiseModel")
Runtime type information.
const dictionary dict_
Copy of dictionary used for construction.
Definition: noiseModel.H:230
declareRunTimeSelectionTable(autoPtr, noiseModel, dictionary,(const dictionary &dict),(dict))
Run time selection table.
bool writePrmsf_
Write Prmsf; default = yes.
Definition: noiseModel.H:275
virtual ~noiseModel()=default
Destructor.
weightingType SPLweighting_
Weighting.
Definition: noiseModel.H:254
scalar checkUniformTimeStep(const scalarList &times) const
Check and return uniform time step.
Definition: noiseModel.C:186
void readWriteOption(const dictionary &dict, const word &lookup, bool &option) const
Helper function to read write options and provide info feedback.
Definition: noiseModel.C:163
bool writePSDf_
Write PSDf; default = yes.
Definition: noiseModel.H:284
scalar gainD(const scalar f) const
D weighting as gain in dB.
Definition: noiseModel.C:573
tmp< scalarField > PSDf(const scalarField &p, const scalar deltaT) const
Definition: noiseModel.C:448
void writeWeightings() const
Helper function to check weightings.
Definition: noiseModel.C:838
scalar startTime_
Start time, default = 0s.
Definition: noiseModel.H:245
void cleanFFTW()
Clean up the FFTW.
Definition: noiseModel.C:827
scalar dBRef_
Reference for dB calculation, default = 2e-5.
Definition: noiseModel.H:257
bool validateBounds(const scalarList &p) const
Return true if all pressure data is within min/max bounds.
Definition: noiseModel.C:220
scalar fLower_
Lower frequency limit, default = 25Hz.
Definition: noiseModel.H:239
scalar gainB(const scalar f) const
B weighting as gain in dB.
Definition: noiseModel.C:527
scalar RBf(const scalar f) const
B weighting function.
Definition: noiseModel.C:511
Lookup type of boundary radiation properties.
Definition: lookup.H:66
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
static void check(const int retVal, const char *what)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
labelList f(nPoints)
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
dictionary dict
Octave band information.
Definition: noiseModel.H:203
FFTW planner information.
Definition: noiseModel.H:193
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73