Loading...
Searching...
No Matches
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-2023 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 minFreq 25;
38 maxFreq 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 minFreq | Lower frequency bounds (fl) | no | 25
62 maxFreq | Upper frequency bounds (fu) | no | 10000
63 sampleFreq | Sample frequency | no | 0
64 startTime | Start time | no | 0
65 outputPrefix | Prefix applied to output files| no | ''
66 SPLweighting | Weighting: dBA, dBB, dBC, DBD | no | none
67 dBRef | Reference for dB calculation | no | 2e-5
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
75 If provided, the sampleFreq is used to define the deltaT between samples.
76 Otherwise the deltaT is determined from the time samples themselves.
77
78SourceFiles
79 noiseModel.C
80
81\*---------------------------------------------------------------------------*/
82
83#ifndef Foam_noiseModel_H
84#define Foam_noiseModel_H
85
86#include "writeFile.H"
87#include "dictionary.H"
88#include "scalarList.H"
89#include "instantList.H"
90#include "windowModel.H"
91#include "Enum.H"
92#include "Tuple2.H"
94#include <fftw3.h>
95
96// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97
98namespace Foam
99{
100
101/*---------------------------------------------------------------------------*\
102 Class noiseModel Declaration
103\*---------------------------------------------------------------------------*/
104
105class noiseModel
106:
108{
109public:
110
111 enum class weightingType
112 {
113 none,
114 dBA,
115 dBB,
116 dBC,
117 dBD
118 };
119
120 static const Enum<weightingType> weightingTypeNames_;
121
122 //- FFTW planner information
123 // Note: storage uses double for use directly with FFTW
124 struct planInfo
125 {
126 bool active;
127 label windowSize;
128 List<double> in;
129 List<double> out;
130 fftw_plan plan;
131 };
132
133 //- Octave band information
134 struct octaveBandInfo
135 {
136 label octave;
137
138 // IDs of bin boundaries in pressure data
140
141 // Centre frequencies for each bin
143 };
144
145
146protected:
147
148 // Protected Data
149
150 //- Copy of dictionary used for construction
151 const dictionary dict_;
152
153 //- Reference density (to convert from kinematic to static pressure)
154 scalar rhoRef_;
155
156 //- Number of samples in sampling window, default = 2^16
157 label nSamples_;
158
159 //- Lower frequency limit, default = 25Hz
160 scalar fLower_;
161
162 //- Upper frequency limit, default = 10kHz
163 scalar fUpper_;
164
165 //- Prescribed sample frequency
166 scalar sampleFreq_;
167
168 //- Start time, default = 0s
169 scalar startTime_;
170
171 //- Window model
172 autoPtr<windowModel> windowModelPtr_;
173
174 //- Weighting
176
177 //- Reference for dB calculation, default = 2e-5
178 scalar dBRef_;
179
180
181 // Data validation
182
183 //- Min pressure value
186 //- Min pressure value
188
189
190 // Write options
191
192 //- Output file prefix, default = ''
194
195 //- Write Prmsf; default = yes
196 bool writePrmsf_;
198 //- Write SPL; default = yes
201 //- Write PSD; default = yes
204 //- Write PSDf; default = yes
205 bool writePSDf_;
206
207 //- Write writeOctaves; default = yes
208 bool writeOctaves_;
210
211 // FFT
212
213 //- Plan information for FFTW
215
216
217 // Protected Member Functions
218
219 //- Check and return uniform time step
221 (
222 const scalarList& times
223 ) const;
224
225 //- Return true if all pressure data is within min/max bounds
226 bool validateBounds(const scalarList& p) const;
227
228 //- Find and return start time index
229 FOAM_DEPRECATED_FOR(2022-07, "instant::findStart() static method")
230 static label findStartTimeIndex
231 (
232 const instantList& allTimes,
233 const scalar startTime
234 )
235 {
236 return instant::findStart(allTimes, startTime);
237 }
239 //- Return the base output directory
240 fileName baseFileDir(const label dataseti) const;
241
242 //- Write output file header
244 (
245 Ostream& os,
246 const string& x,
247 const string& y,
248 const UList<Tuple2<string, token>>& headerValues
250 ) const;
251
252 // Write frequency-based data to file
254 (
255 Ostream& os,
256 const scalarField& f,
257 const scalarField& fx
258 ) const;
259
260 //- Create a field of equally spaced frequencies for the current set of
261 //- data - assumes a constant time step
264 const scalar deltaT,
265 const bool check
266 ) const;
267
268 //- Return a list of the frequency indices wrt f field that correspond
269 //- to the bands limits for a given octave
270 static void setOctaveBands
271 (
272 const scalarField& f,
273 const scalar fLower,
274 const scalar fUpper,
275 const scalar octave,
276 labelList& fBandIDs,
277 scalarField& fCentre
278 );
279
280 //- Generate octave data
282 (
283 const scalarField& data,
284 const scalarField& f,
285 const labelUList& freqBandIDs
286 ) const;
287
288 //- Return the fft of the given pressure data
289 tmp<scalarField> Pf(const scalarField& p) const;
290
291 //- Return the multi-window mean fft of the complete pressure data [Pa]
292 tmp<scalarField> meanPf(const scalarField& p) const;
293
294 //- Return the multi-window RMS mean fft of the complete pressure
295 //- data [Pa]
297
298 //- Return the multi-window Power Spectral Density (PSD) of the complete
299 //- pressure data [Pa^2/Hz]
300 tmp<scalarField> PSDf(const scalarField& p, const scalar deltaT) const;
301
302
303 // Weightings
305 //- A weighting function
306 scalar RAf(const scalar f) const;
307
308 //- A weighting as gain in dB
309 scalar gainA(const scalar f) const;
310
311 //- B weighting function
312 scalar RBf(const scalar f) const;
313
314 //- B weighting as gain in dB
315 scalar gainB(const scalar f) const;
316
317 //- C weighting function
318 scalar RCf(const scalar f) const;
320 //- C weighting as gain in dB
321 scalar gainC(const scalar f) const;
322
323 //- D weighting function
324 scalar RDf(const scalar f) const;
325
326 //- D weighting as gain in dB
327 scalar gainD(const scalar f) const;
328
329
330public:
331
332 //- Runtime type information
333 TypeName("noiseModel");
334
335 //- Run time selection table
337 (
338 autoPtr,
341 (
342 const dictionary& dict,
343 const objectRegistry& obr
344 ),
345 (dict, obr)
346 );
347
348
349 // Generated Methods
350
351 //- No copy construct
352 noiseModel(const noiseModel&) = delete;
353
354 //- No copy assignment
355 void operator=(const noiseModel&) = delete;
356
357
358 // Constructors
359
360 //- Constructor
362 (
363 const dictionary& dict,
364 const objectRegistry& obr,
365 const word& name,
366 const bool readFields = true
367 );
368
369 //- Selector
371 (
372 const dictionary& dict,
373 const objectRegistry& obr
374 );
375
376
377 //- Destructor
378 virtual ~noiseModel() = default;
379
380
381 // Public Member Functions
382
383 //- Read from dictionary
384 virtual bool read(const dictionary& dict);
385
386 //- Abstract call to calculate
387 virtual void calculate() = 0;
388
389 //- PSD [dB/Hz]
391
392 //- SPL [dB]
394 (
395 const scalarField& Prms2,
396 const scalar f
397 ) const;
398
399 //- SPL [dB]
401 (
402 const scalarField& Prms2,
403 const scalarField& f
404 ) const;
405
406 //- Clean up the FFTW
407 void cleanFFTW();
408
409 //- Helper function to check weightings
410 void writeWeightings() const;
411};
412
413
414// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
415
416} // End namespace Foam
417
418// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
419
420#endif
421
422// ************************************************************************* //
scalar y
Foam::label startTime
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition Tuple2.H:51
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A class for handling file names.
Definition fileName.H:75
Base class for writing single files from the function objects.
Definition writeFile.H:113
fileName baseFileDir() const
Return the base directory for output.
Definition writeFile.C:43
static label findStart(const UList< instant > &times, const scalar timeVal)
Find and return index of given start time (linear search).
Definition instant.C:37
Base class for noise models.
Definition noiseModel.H:178
tmp< scalarField > uniformFrequencies(const scalar deltaT, const bool check) const
Create a field of equally spaced frequencies for the current set of data - assumes a constant time st...
Definition noiseModel.C:295
scalar RAf(const scalar f) const
A weighting function.
Definition noiseModel.C:501
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:238
declareRunTimeSelectionTable(autoPtr, noiseModel, dictionary,(const dictionary &dict, const objectRegistry &obr),(dict, obr))
Run time selection table.
tmp< scalarField > octaves(const scalarField &data, const scalarField &f, const labelUList &freqBandIDs) const
Generate octave data.
Definition noiseModel.C:332
tmp< scalarField > RMSmeanPf(const scalarField &p) const
Return the multi-window RMS mean fft of the complete pressure data [Pa].
Definition noiseModel.C:447
virtual void calculate()=0
Abstract call to calculate.
scalar fUpper_
Upper frequency limit, default = 10kHz.
Definition noiseModel.H:248
scalar gainC(const scalar f) const
C weighting as gain in dB.
Definition noiseModel.C:567
tmp< Foam::scalarField > PSD(const scalarField &PSDf) const
PSD [dB/Hz].
Definition noiseModel.C:758
static autoPtr< noiseModel > New(const dictionary &dict, const objectRegistry &obr)
Selector.
tmp< scalarField > meanPf(const scalarField &p) const
Return the multi-window mean fft of the complete pressure data [Pa].
Definition noiseModel.C:426
scalar maxPressure_
Min pressure value.
Definition noiseModel.H:286
bool writePSD_
Write PSD; default = yes.
Definition noiseModel.H:309
scalar gainA(const scalar f) const
A weighting as gain in dB.
Definition noiseModel.C:518
scalar minPressure_
Min pressure value.
Definition noiseModel.H:281
static void setOctaveBands(const scalarField &f, const scalar fLower, const scalar fUpper, const scalar octave, labelList &fBandIDs, scalarField &fCentre)
Return a list of the frequency indices wrt f field that correspond to the bands limits for a given oc...
Definition noiseModel.C:48
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition noiseModel.C:649
bool writeOctaves_
Write writeOctaves; default = yes.
Definition noiseModel.H:319
fileName outputPrefix_
Output file prefix, default = ''.
Definition noiseModel.H:294
scalar RDf(const scalar f) const
D weighting function.
Definition noiseModel.C:578
tmp< scalarField > Pf(const scalarField &p) const
Return the fft of the given pressure data.
Definition noiseModel.C:383
scalar RCf(const scalar f) const
C weighting function.
Definition noiseModel.C:556
noiseModel(const noiseModel &)=delete
No copy construct.
tmp< scalarField > SPL(const scalarField &Prms2, const scalar f) const
SPL [dB].
Definition noiseModel.C:767
planInfo planInfo_
Plan information for FFTW.
Definition noiseModel.H:327
static label findStartTimeIndex(const instantList &allTimes, const scalar startTime)
Find and return start time index.
Definition noiseModel.H:350
void operator=(const noiseModel &)=delete
No copy assignment.
autoPtr< windowModel > windowModelPtr_
Window model.
Definition noiseModel.H:263
static const Enum< weightingType > weightingTypeNames_
Definition noiseModel.H:190
bool writeSPL_
Write SPL; default = yes.
Definition noiseModel.H:304
TypeName("noiseModel")
Runtime type information.
const dictionary dict_
Copy of dictionary used for construction.
Definition noiseModel.H:228
void writeFileHeader(Ostream &os, const string &x, const string &y, const UList< Tuple2< string, token > > &headerValues=UList< Tuple2< string, token > >::null()) const
Write output file header.
Definition noiseModel.C:250
bool writePrmsf_
Write Prmsf; default = yes.
Definition noiseModel.H:299
virtual ~noiseModel()=default
Destructor.
weightingType SPLweighting_
Weighting.
Definition noiseModel.H:268
scalar checkUniformTimeStep(const scalarList &times) const
Check and return uniform time step.
Definition noiseModel.C:177
bool writePSDf_
Write PSDf; default = yes.
Definition noiseModel.H:314
scalar gainD(const scalar f) const
D weighting as gain in dB.
Definition noiseModel.C:591
tmp< scalarField > PSDf(const scalarField &p, const scalar deltaT) const
Return the multi-window Power Spectral Density (PSD) of the complete pressure data [Pa^2/Hz].
Definition noiseModel.C:466
void writeWeightings() const
Helper function to check weightings.
Definition noiseModel.C:883
scalar startTime_
Start time, default = 0s.
Definition noiseModel.H:258
scalar sampleFreq_
Prescribed sample frequency.
Definition noiseModel.H:253
void cleanFFTW()
Clean up the FFTW.
Definition noiseModel.C:872
scalar dBRef_
Reference for dB calculation, default = 2e-5.
Definition noiseModel.H:273
void writeFreqDataToFile(Ostream &os, const scalarField &f, const scalarField &fx) const
Definition noiseModel.C:278
bool validateBounds(const scalarList &p) const
Return true if all pressure data is within min/max bounds.
Definition noiseModel.C:216
scalar fLower_
Lower frequency limit, default = 25Hz.
Definition noiseModel.H:243
scalar gainB(const scalar f) const
B weighting as gain in dB.
Definition noiseModel.C:545
scalar RBf(const scalar f) const
B weighting function.
Definition noiseModel.C:529
Registry of regIOobjects.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
List< label > labelList
A List of labels.
Definition List.H:62
static void check(const int retVal, const char *what)
List< instant > instantList
List of instants.
Definition instantList.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
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
#define FOAM_DEPRECATED_FOR(since, replacement)
Definition stdFoam.H:43
Octave band information.
Definition noiseModel.H:210
FFTW planner information.
Definition noiseModel.H:198
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68