noiseFFT.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-2015 OpenFOAM Foundation
9 Copyright (C) 2016-2017 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "noiseFFT.H"
30#include "IFstream.H"
31#include "DynamicList.H"
33#include "HashSet.H"
34#include "fft.H"
35
36using namespace Foam::constant;
37
38// * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * //
39
40Foam::scalar Foam::noiseFFT::p0 = 2e-5;
41
42
44(
45 const label N,
46 const scalar deltaT
47)
48{
49 auto tf = tmp<scalarField>::New(N/2, Zero);
50 auto& f = tf.ref();
51
52 const scalar deltaf = 1.0/(N*deltaT);
53 forAll(f, i)
54 {
55 f[i] = i*deltaf;
56 }
57
58 return tf;
59}
60
61
63{
64 return 10*log10(PSDf/sqr(p0));
65}
66
67
69{
70 return 10*log10(Prms2/sqr(p0));
71}
72
73
75(
76 const scalarField& f,
77 const scalar fLower,
78 const scalar fUpper,
79 const scalar octave,
80 labelList& fBandIDs,
81 scalarField& fCentre
82)
83{
84 // Set the indices of to the lower frequency bands for the input frequency
85 // range. Ensure that the centre frequency passes though 1000 Hz
86
87 // Low frequency bound given by:
88 // fLow = f0*(2^(0.5*bandI/octave))
89
90 // Initial (lowest centre frequency)
91 scalar fTest = 15.625;
92
93 const scalar fRatio = pow(2, 1.0/octave);
94 const scalar fRatioL2C = pow(2, 0.5/octave);
95
96 // IDs of band IDs
97 labelHashSet bandIDs(f.size());
98
99 // Centre frequencies
101
102 // Convert to lower band limit
103 fTest /= fRatioL2C;
104
105 forAll(f, i)
106 {
107 if (f[i] >= fTest)
108 {
109 // Advance band if appropriate
110 while (f[i] > fTest)
111 {
112 fTest *= fRatio;
113 }
114 fTest /= fRatio;
115
116 if (bandIDs.insert(i))
117 {
118 // Also store (next) centre frequency
119 fc.append(fTest*fRatioL2C);
120 }
121 fTest *= fRatio;
122
123 if (fTest > fUpper)
124 {
125 break;
126 }
127 }
128 }
129
130 fBandIDs = bandIDs.sortedToc();
131
132 if (fc.size())
133 {
134 // Remove the last centre frequency (beyond upper frequency limit)
135 fc.remove();
136
137 fCentre.transfer(fc);
138 }
139}
140
141
142// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
143
144Foam::noiseFFT::noiseFFT(const scalar deltaT, const label windowSize)
145:
146 scalarField(),
147 deltaT_(deltaT)
148{
149 if (windowSize > 1)
150 {
151 planInfo_.active = true;
152 planInfo_.windowSize = windowSize;
153 planInfo_.in.setSize(windowSize);
154 planInfo_.out.setSize(windowSize);
155
156 // Using real to half-complex fftw 'kind'
157 planInfo_.plan =
158 fftw_plan_r2r_1d
159 (
160 windowSize,
161 planInfo_.in.data(),
162 planInfo_.out.data(),
163 FFTW_R2HC,
164 windowSize <= 8192 ? FFTW_MEASURE : FFTW_ESTIMATE
165 );
166 }
167 else
168 {
169 planInfo_.active = false;
170 }
171}
172
173
174// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
175
177{
178 if (planInfo_.active)
179 {
180 planInfo_.active = false;
181 fftw_destroy_plan(planInfo_.plan);
182 fftw_cleanup();
183 }
184}
185
186
187// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
188
190{
191 this->transfer(data);
192
193 scalarField& p = *this;
194 p -= average(p);
195}
196
197
198void Foam::noiseFFT::setData(const fileName& pFileName, const label skip)
199{
200 // Construct pressure data file
201 IFstream pFile(pFileName);
202
203 // Check pFile stream is OK
204 if (!pFile.good())
205 {
207 << "Cannot read file " << pFileName
208 << exit(FatalError);
209 }
210
211 if (skip)
212 {
213 scalar dummyt, dummyp;
214
215 for (label i = 0; i < skip; ++i)
216 {
217 pFile >> dummyt;
218
219 if (!pFile.good() || pFile.eof())
220 {
222 << "Number of points in file " << pFileName
223 << " is less than the number to be skipped = " << skip
224 << exit(FatalError);
225 }
226
227 pFile >> dummyp;
228 }
229 }
230
231 scalar t = 0, T0 = 0, T1 = 0;
232 DynamicList<scalar> pData(1024);
233 label i = 0;
234
235 while (!(pFile >> t).eof())
236 {
237 if (i == 0)
238 {
239 T0 = t;
240 }
241
242 T1 = t;
243 pFile >> pData(i++);
244 }
245
246 // Note: assumes fixed time step
247 deltaT_ = (T1 - T0)/pData.size();
248
249 this->transfer(pData);
250
251 scalarField& p = *this;
252 p -= average(p);
253}
254
255
257{
258 scalarField t(size());
259 forAll(t, i)
260 {
261 t[i] = i*deltaT_;
262 }
263
264 return graph
265 (
266 "p(t)",
267 "t [s]",
268 "p(t) [Pa]",
269 t,
270 *this
271 );
272}
273
274
276(
277 const tmp<scalarField>& tpn
278) const
279{
280 if (planInfo_.active)
281 {
282 const scalarField& pn = tpn();
283 if (pn.size() != planInfo_.windowSize)
284 {
286 << "Expected pressure data to have " << planInfo_.windowSize
287 << " values, but received " << pn.size() << " values"
288 << abort(FatalError);
289 }
290
291 List<double>& in = planInfo_.in;
292 const List<double>& out = planInfo_.out;
293 forAll(in, i)
294 {
295 in[i] = pn[i];
296 }
297 tpn.clear();
298
299 ::fftw_execute(planInfo_.plan);
300
301 const label n = planInfo_.windowSize;
302 const label nBy2 = n/2;
303 auto tresult = tmp<scalarField>::New(nBy2 + 1);
304 auto& result = tresult.ref();
305
306 // 0 th value = DC component
307 // nBy2 th value = real only if n is even
308 result[0] = out[0];
309 for (label i = 1; i <= nBy2; ++i)
310 {
311 const auto re = out[i];
312 const auto im = out[n - i];
313 result[i] = sqrt(re*re + im*im);
314 }
315
316 return tresult;
317 }
318
319 return mag(fft::realTransform1D(tpn));
320}
321
322
324{
325 const label N = window.nSamples();
326 const label nWindow = window.nWindow();
327
328 scalarField meanPf(N/2 + 1, Zero);
329
330 for (label windowI = 0; windowI < nWindow; ++windowI)
331 {
332 meanPf += Pf(window.apply<scalar>(*this, windowI));
333 }
334
335 meanPf /= scalar(nWindow);
336
337 scalar deltaf = 1.0/(N*deltaT_);
338 scalarField f(meanPf.size());
339 forAll(f, i)
340 {
341 f[i] = i*deltaf;
342 }
343
344 return graph
345 (
346 "P(f)",
347 "f [Hz]",
348 "P(f) [Pa]",
349 f,
350 meanPf
351 );
352}
353
354
356{
357 const label N = window.nSamples();
358 const label nWindow = window.nWindow();
359
360 scalarField RMSMeanPf(N/2 + 1, Zero);
361 for (label windowI = 0; windowI < nWindow; ++windowI)
362 {
363 RMSMeanPf += sqr(Pf(window.apply<scalar>(*this, windowI)));
364 }
365
366 RMSMeanPf = sqrt(RMSMeanPf/scalar(nWindow))/scalar(N);
367
368 scalar deltaf = 1.0/(N*deltaT_);
369 scalarField f(RMSMeanPf.size());
370 forAll(f, i)
371 {
372 f[i] = i*deltaf;
373 }
374
375 return graph
376 (
377 "Prms(f)",
378 "f [Hz]",
379 "Prms(f) [Pa]",
380 f,
381 RMSMeanPf
382 );
383}
384
385
387{
388 const label N = window.nSamples();
389 const label nWindow = window.nWindow();
390
391 scalarField psd(N/2 + 1, Zero);
392
393 for (label windowI = 0; windowI < nWindow; ++windowI)
394 {
395 psd += sqr(Pf(window.apply<scalar>(*this, windowI)));
396 }
397
398 scalar deltaf = 1.0/(N*deltaT_);
399 scalar fs = 1.0/deltaT_;
400 psd /= scalar(nWindow)*fs*N;
401
402 // Scaling due to use of 1-sided FFT
403 // Note: do not scale DC component
404 psd *= 2;
405 psd.first() /= 2;
406 psd.last() /= 2;
407
408 scalarField f(psd.size());
409
410 if (0) // if (debug)
411 {
412 Pout<< "Average PSD: " << average(psd) << endl;
413 }
414
415 forAll(f, i)
416 {
417 f[i] = i*deltaf;
418 }
419
420 return graph
421 (
422 "PSD(f)",
423 "f [Hz]",
424 "PSD(f) [PaPa_Hz]",
425 f,
426 psd
427 );
428}
429
430
432{
433 return graph
434 (
435 "PSD(f)",
436 "f [Hz]",
437 "PSD_dB(f) [dB_Hz]",
438 gPSDf.x(),
439 10*log10(gPSDf.y()/sqr(p0))
440 );
441}
442
443
445(
446 const graph& g,
447 const labelUList& freqBandIDs
448) const
449{
450 if (freqBandIDs.size() < 2)
451 {
453 << "Octave frequency bands are not defined "
454 << "- skipping octaves calculation"
455 << endl;
456
457 return graph
458 (
459 "octave",
460 "x",
461 "y",
462 scalarField(),
464 );
465 }
466
467 const scalarField& f = g.x();
468 const scalarField& data = g.y();
469
470 scalarField octData(freqBandIDs.size() - 1, Zero);
471 scalarField fm(freqBandIDs.size() -1, Zero);
472
473 for (label bandI = 0; bandI < freqBandIDs.size() - 1; ++bandI)
474 {
475 label fb0 = freqBandIDs[bandI];
476 label fb1 = freqBandIDs[bandI+1];
477 fm[bandI] = f[fb0];
478
479 if (fb0 == fb1) continue;
480
481 for (label freqI = fb0; freqI < fb1; ++freqI)
482 {
483 label f0 = f[freqI];
484 label f1 = f[freqI + 1];
485 scalar dataAve = 0.5*(data[freqI] + data[freqI + 1]);
486 octData[bandI] += dataAve*(f1 - f0);
487 }
488 }
489
490 return graph
491 (
492 "octaves(f)",
493 "fm [Hz]",
494 "octave data",
495 fm,
496 octData
497 );
498}
499
500
501Foam::scalar Foam::noiseFFT::dbToPa(const scalar db) const
502{
503 return p0*pow(10.0, db/20.0);
504}
505
506
508(
509 const tmp<scalarField>& db
510) const
511{
512 return p0*pow(10.0, db/20.0);
513}
514
515
516// ************************************************************************* //
label n
const uniformDimensionedVectorField & g
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
T remove()
Remove and return the last element. Fatal on an empty list.
Definition: DynamicListI.H:655
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
Input from file stream, using an ISstream.
Definition: IFstream.H:57
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:233
bool eof() const noexcept
True if end of input seen.
Definition: IOstream.H:239
void transfer(List< T > &list)
Definition: List.C:447
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
T & last()
Return the last element of the list.
Definition: UListI.H:216
Database for solution data, solver performance and other reduced data.
Definition: data.H:58
static tmp< complexField > realTransform1D(const scalarField &field)
Transform real-value data.
Definition: fft.C:120
A class for handling file names.
Definition: fileName.H:76
Class to create, store and output qgraph files.
Definition: graph.H:62
const scalarField & y() const
Definition: graph.C:144
const scalarField & x() const
Definition: graph.H:163
Performs FFT of pressure field to generate noise data.
Definition: noiseFFT.H:72
graph meanPf(const windowModel &window) const
Return the multi-window mean fft of the complete pressure data [Pa].
Definition: noiseFFT.C:323
static tmp< scalarField > frequencies(const label N, const scalar deltaT)
Return the FFT frequencies.
Definition: noiseFFT.C:44
static tmp< scalarField > SPL(const scalarField &Prms2)
Return the SPL [dB].
Definition: noiseFFT.C:68
static scalar p0
Reference pressure.
Definition: noiseFFT.H:109
void setData(scalarList &data)
Definition: noiseFFT.C:189
static void octaveBandInfo(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.
Definition: noiseFFT.C:75
graph pt() const
Return the graph of pressure as a function of time.
Definition: noiseFFT.C:256
graph octaves(const graph &g, const labelUList &freqBandIDs) const
Generate octave data.
Definition: noiseFFT.C:445
static tmp< scalarField > PSD(const scalarField &PSDf)
Return the PSD [dB/Hz].
Definition: noiseFFT.C:62
~noiseFFT()
Destructor. Cleanup/destroy plan.
Definition: noiseFFT.C:176
graph PSDf(const windowModel &window) const
Return the multi-window PSD (power spectral density) of the complete.
Definition: noiseFFT.C:386
scalar dbToPa(const scalar db) const
Convert the db into Pa.
Definition: noiseFFT.C:501
graph RMSmeanPf(const windowModel &window) const
Return the multi-window RMS mean fft of the complete pressure.
Definition: noiseFFT.C:355
tmp< scalarField > Pf(const tmp< scalarField > &pn) const
Return the fft of the given pressure data.
Definition: noiseFFT.C:276
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
Base class for windowing models.
Definition: windowModel.H:55
label nWindow() const
Return the number of windows.
Definition: windowModel.C:62
tmp< Field< Type > > apply(const Field< Type > &fld, const label windowI) const
Return the windowed data.
label nSamples() const
Return the number of samples in the window.
Definition: windowModel.C:56
volScalarField & p
const volScalarField & p0
Definition: EEqn.H:36
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define WarningInFunction
Report a warning using Foam::Warning.
Different types of constants.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar log10(const dimensionedScalar &ds)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
labelList f(nPoints)
scalar T0
Definition: createFields.H:22
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
const Vector< label > N(dict.get< Vector< label > >("N"))