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