noiseModel.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) 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 \*---------------------------------------------------------------------------*/
27 
28 #include "noiseModel.H"
29 #include "functionObject.H"
30 #include "argList.H"
31 #include "fft.H"
32 #include "OFstream.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(noiseModel, 0);
39  defineRunTimeSelectionTable(noiseModel, dictionary);
40 }
41 
42 
45 ({
46  {weightingType::none, "dB"},
47  {weightingType::dBA, "dBA"},
48  {weightingType::dBB, "dBB"},
49  {weightingType::dBC, "dBC"},
50  {weightingType::dBD, "dBD"},
51 });
52 
53 
55 (
56  const scalarField& f,
57  const scalar fLower,
58  const scalar fUpper,
59  const scalar octave,
60  labelList& fBandIDs,
61  scalarField& fCentre
62 )
63 {
64  // Set the indices of to the lower frequency bands for the input frequency
65  // range. Ensure that the centre frequency passes though 1000 Hz
66 
67  // Low frequency bound given by:
68  // fLow = f0*(2^(0.5*bandI/octave))
69 
70  // Initial (lowest centre frequency)
71  scalar fTest = 15.625;
72 
73  const scalar fRatio = pow(2, 1.0/octave);
74  const scalar fRatioL2C = pow(2, 0.5/octave);
75 
76  // IDs of band IDs
77  labelHashSet bandIDs(f.size());
78 
79  // Centre frequencies
81 
82  // Convert to lower band limit
83  fTest /= fRatioL2C;
84 
85  forAll(f, i)
86  {
87  if (f[i] >= fTest)
88  {
89  // Advance band if appropriate
90  while (f[i] > fTest)
91  {
92  fTest *= fRatio;
93  }
94  fTest /= fRatio;
95 
96  if (bandIDs.insert(i))
97  {
98  // Also store (next) centre frequency
99  fc.append(fTest*fRatioL2C);
100  }
101  fTest *= fRatio;
102 
103  if (fTest > fUpper)
104  {
105  break;
106  }
107  }
108  }
109 
110  fBandIDs = bandIDs.sortedToc();
111 
112  if (fc.size())
113  {
114  // Remove the last centre frequency (beyond upper frequency limit)
115  fc.remove();
116 
117  fCentre.transfer(fc);
118  }
119 }
120 
121 
122 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
123 
125 (
126  const dictionary& dict,
127  const word& lookup,
128  bool& option
129 ) const
130 {
131  dict.readIfPresent(lookup, option);
132 
133  // Only writing on the master process
134  option = option && Pstream::master();
135 
136  if (option)
137  {
138  Info<< " " << lookup << ": " << "yes" << endl;
139  }
140  else
141  {
142  Info<< " " << lookup << ": " << "no" << endl;
143  }
144 }
145 
146 
148 (
149  const scalarList& times
150 ) const
151 {
152  scalar deltaT = -1.0;
153 
154  if (times.size() > 1)
155  {
156  // Uniform time step
157  deltaT = (times.last() - times.first())/scalar(times.size() - 1);
158 
159  for (label i = 1; i < times.size(); ++i)
160  {
161  scalar dT = times[i] - times[i-1];
162 
163  if (mag(dT/deltaT - 1) > 1e-8)
164  {
166  << "Unable to process data with a variable time step"
167  << exit(FatalError);
168  }
169  }
170  }
171  else
172  {
174  << "Unable to create FFT with a single value"
175  << exit(FatalError);
176  }
177 
178  return deltaT;
179 }
180 
181 
183 {
184  forAll(p, i)
185  {
186  if ((p[i] < minPressure_) || (p[i] > maxPressure_))
187  {
189  << "Pressure data at position " << i
190  << " is outside of permitted bounds:" << nl
191  << " pressure: " << p[i] << nl
192  << " minimum pressure: " << minPressure_ << nl
193  << " maximum pressure: " << maxPressure_ << nl
194  << endl;
195 
196  return false;
197  }
198  }
199 
200  return true;
201 }
202 
203 
205 (
206  const instantList& allTimes,
207  const scalar startTime
208 ) const
209 {
210  forAll(allTimes, timeI)
211  {
212  const instant& i = allTimes[timeI];
213 
214  if (i.value() >= startTime)
215  {
216  return timeI;
217  }
218  }
219 
220  return 0;
221 }
222 
223 
225 {
226  return
227  (
230  / "noise"
231  / outputPrefix_
232  / type()
233  / ("input" + Foam::name(dataseti))
234  );
235 }
236 
237 
239 (
240  const scalar deltaT
241 ) const
242 {
243  const auto& window = windowModelPtr_();
244  const label N = window.nSamples();
245 
246  auto tf = tmp<scalarField>::New(N/2 + 1, Zero);
247  auto& f = tf.ref();
248 
249  const scalar deltaf = 1.0/(N*deltaT);
250  forAll(f, i)
251  {
252  f[i] = i*deltaf;
253  }
254 
255  return tf;
256 }
257 
258 
260 (
261  const scalarField& data,
262  const scalarField& f,
263  const labelUList& freqBandIDs
264 ) const
265 {
266  auto toctData = tmp<scalarField>::New(freqBandIDs.size() - 1, Zero);
267 
268  if (freqBandIDs.size() < 2)
269  {
271  << "Octave frequency bands are not defined "
272  << "- skipping octaves calculation"
273  << endl;
274 
275  return toctData;
276  }
277 
278  auto& octData = toctData.ref();
279 
280  for (label bandI = 0; bandI < freqBandIDs.size() - 1; ++bandI)
281  {
282  label fb0 = freqBandIDs[bandI];
283  label fb1 = freqBandIDs[bandI+1];
284 
285  if (fb0 == fb1) continue;
286 
287  for (label freqI = fb0; freqI < fb1; ++freqI)
288  {
289  label f0 = f[freqI];
290  label f1 = f[freqI + 1];
291  scalar dataAve = 0.5*(data[freqI] + data[freqI + 1]);
292  octData[bandI] += dataAve*(f1 - f0);
293  }
294  }
295 
296  return toctData;
297 }
298 
299 
301 {
302  if (planInfo_.active)
303  {
304  if (p.size() != planInfo_.windowSize)
305  {
307  << "Expected pressure data to have " << planInfo_.windowSize
308  << " values, but received " << p.size() << " values"
309  << abort(FatalError);
310  }
311 
312  List<double>& in = planInfo_.in;
313  const List<double>& out = planInfo_.out;
314  forAll(in, i)
315  {
316  in[i] = p[i];
317  }
318 
319  ::fftw_execute(planInfo_.plan);
320 
321  const label n = planInfo_.windowSize;
322  const label nBy2 = n/2;
323  auto tresult = tmp<scalarField>::New(nBy2 + 1);
324  auto& result = tresult.ref();
325 
326  // 0 th value = DC component
327  // nBy2 th value = real only if n is even
328  result[0] = out[0];
329  for (label i = 1; i <= nBy2; ++i)
330  {
331  const auto re = out[i];
332  const auto im = out[n - i];
333  result[i] = sqrt(re*re + im*im);
334  }
335 
336  return tresult;
337  }
338 
339  return mag(fft::realTransform1D(p));
340 }
341 
342 
344 {
345  const auto& window = windowModelPtr_();
346  const label N = window.nSamples();
347  const label nWindow = window.nWindow();
348 
349  auto tmeanPf = tmp<scalarField>::New(N/2 + 1, Zero);
350  auto& meanPf = tmeanPf.ref();
351 
352  for (label windowI = 0; windowI < nWindow; ++windowI)
353  {
354  meanPf += Pf(window.apply<scalar>(p, windowI));
355  }
356 
357  meanPf /= scalar(nWindow);
358 
359  return tmeanPf;
360 }
361 
362 
364 (
365  const scalarField& p
366 ) const
367 {
368  const auto& window = windowModelPtr_();
369  const label N = window.nSamples();
370  const label nWindow = window.nWindow();
371 
372  scalarField RMSMeanPf(N/2 + 1, Zero);
373  for (label windowI = 0; windowI < nWindow; ++windowI)
374  {
375  RMSMeanPf += sqr(Pf(window.apply<scalar>(p, windowI)));
376  }
377 
378  return sqrt(RMSMeanPf/scalar(nWindow))/scalar(N);
379 }
380 
381 
383 (
384  const scalarField& p,
385  const scalar deltaT
386 ) const
387 {
388  const auto& window = windowModelPtr_();
389  const label N = window.nSamples();
390  const label nWindow = window.nWindow();
391 
392  auto tpsd = tmp<scalarField>::New(N/2 + 1, Zero);
393  auto& psd = tpsd.ref();
394 
395  for (label windowI = 0; windowI < nWindow; ++windowI)
396  {
397  psd += sqr(Pf(window.apply<scalar>(p, windowI)));
398  }
399 
400  scalar fs = 1.0/deltaT;
401  psd /= scalar(nWindow)*fs*N;
402 
403  // Scaling due to use of 1-sided FFT
404  // Note: do not scale DC component
405  psd *= 2;
406  psd.first() /= 2;
407  psd.last() /= 2;
408 
409  if (debug)
410  {
411  Pout<< "Average PSD: " << average(psd) << endl;
412  }
413 
414  return tpsd;
415 }
416 
417 
418 Foam::scalar Foam::noiseModel::RAf(const scalar f) const
419 {
420  const scalar c1 = sqr(12194.0);
421  const scalar c2 = sqr(20.6);
422  const scalar c3 = sqr(107.7);
423  const scalar c4 = sqr(739.9);
424 
425  const scalar f2 = f*f;
426 
427  return
428  c1*sqr(f2)
429  /(
430  (f2 + c2)*sqrt((f2 + c3)*(f2 + c4))*(f2 + c1)
431  );
432 }
433 
434 
435 Foam::scalar Foam::noiseModel::gainA(const scalar f) const
436 {
437  return 20*log10(RAf(f)) - 20*log10(RAf(1000));
438 }
439 
440 
441 Foam::scalar Foam::noiseModel::RBf(const scalar f) const
442 {
443  const scalar c1 = sqr(12194.0);
444  const scalar c2 = sqr(20.6);
445  const scalar c3 = sqr(158.5);
446 
447  const scalar f2 = f*f;
448 
449  return
450  c1*f2*f
451  /(
452  (f2 + c2)*sqrt(f2 + c3)*(f2 + c1)
453  );
454 }
455 
456 
457 Foam::scalar Foam::noiseModel::gainB(const scalar f) const
458 {
459  return 20*log10(RBf(f)) - 20*log10(RBf(1000));
460 }
461 
462 
463 Foam::scalar Foam::noiseModel::RCf(const scalar f) const
464 {
465  const scalar c1 = sqr(12194.0);
466  const scalar c2 = sqr(20.6);
467 
468  const scalar f2 = f*f;
469 
470  return c1*f2/((f2 + c2)*(f2 + c1));
471 }
472 
473 
474 Foam::scalar Foam::noiseModel::gainC(const scalar f) const
475 {
476  return 20*log10(RCf(f)) - 20*log10(RCf(1000));
477 }
478 
479 
480 Foam::scalar Foam::noiseModel::RDf(const scalar f) const
481 {
482  const scalar f2 = f*f;
483 
484  const scalar hf =
485  (sqr(1037918.48 - f2) + 1080768.16*f2)
486  /(sqr(9837328 - f2) + 11723776*f2);
487 
488  return
489  f/6.8966888496476e-5*Foam::sqrt(hf/((f2 + 79919.29)*(f2 + 1345600)));
490 }
491 
492 
493 Foam::scalar Foam::noiseModel::gainD(const scalar f) const
494 {
495  return 20*log10(RDf(f));
496 }
497 
498 
499 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
500 
501 Foam::noiseModel::noiseModel(const dictionary& dict, const bool readFields)
502 :
503  dict_(dict),
504  rhoRef_(1),
505  nSamples_(65536),
506  fLower_(25),
507  fUpper_(10000),
508  startTime_(0),
509  windowModelPtr_(),
510  graphFormat_("raw"),
511  SPLweighting_(weightingType::none),
512  minPressure_(-0.5*VGREAT),
513  maxPressure_(0.5*VGREAT),
514  outputPrefix_(),
515  writePrmsf_(true),
516  writeSPL_(true),
517  writePSD_(true),
518  writePSDf_(true),
519  writeOctaves_(true),
520  planInfo_()
521 {
522  planInfo_.active = false;
523 
524  if (readFields)
525  {
526  read(dict);
527  }
528 
529  if (debug)
530  {
531  writeWeightings();
532  }
533 }
534 
535 
536 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
537 
539 {
540  dict.readIfPresent("rhoRef", rhoRef_);
541  dict.readIfPresent("N", nSamples_);
542  dict.readIfPresent("fl", fLower_);
543  dict.readIfPresent("fu", fUpper_);
544  dict.readIfPresent("startTime", startTime_);
545  dict.readIfPresent("graphFormat", graphFormat_);
546  dict.readIfPresent("minPressure", minPressure_);
547  dict.readIfPresent("maxPressure", maxPressure_);
548  dict.readIfPresent("outputPrefix", outputPrefix_);
549 
550  if (fLower_ < 0)
551  {
553  << "fl: lower frequency bound must be greater than zero"
554  << exit(FatalIOError);
555  }
556 
557  if (fUpper_ < 0)
558  {
560  << "fu: upper frequency bound must be greater than zero"
561  << exit(FatalIOError);
562  }
563 
564  if (fUpper_ < fLower_)
565  {
567  << "fu: upper frequency bound must be greater than lower "
568  << "frequency bound (fl)"
569  << exit(FatalIOError);
570 
571  }
572 
573  Info<< " Frequency bounds:" << nl
574  << " lower: " << fLower_ << nl
575  << " upper: " << fUpper_ << endl;
576 
577  weightingTypeNames_.readIfPresent("SPLweighting", dict, SPLweighting_);
578 
579  Info<< " Weighting: " << weightingTypeNames_[SPLweighting_] << endl;
580 
581  Info<< " Write options:" << endl;
582  dictionary optDict(dict.subOrEmptyDict("writeOptions"));
583  readWriteOption(optDict, "writePrmsf", writePrmsf_);
584  readWriteOption(optDict, "writeSPL", writeSPL_);
585  readWriteOption(optDict, "writePSD", writePSD_);
586  readWriteOption(optDict, "writePSDf", writePSDf_);
587  readWriteOption(optDict, "writeOctaves", writeOctaves_);
588 
589 
590  windowModelPtr_ = windowModel::New(dict, nSamples_);
591 
592  cleanFFTW();
593 
594  const label windowSize = windowModelPtr_->nSamples();
595 
596  if (windowSize > 1)
597  {
598  planInfo_.active = true;
599  planInfo_.windowSize = windowSize;
600  planInfo_.in.setSize(windowSize);
601  planInfo_.out.setSize(windowSize);
602 
603  // Using real to half-complex fftw 'kind'
604  planInfo_.plan =
605  fftw_plan_r2r_1d
606  (
607  windowSize,
608  planInfo_.in.data(),
609  planInfo_.out.data(),
610  FFTW_R2HC,
611  windowSize <= 8192 ? FFTW_MEASURE : FFTW_ESTIMATE
612  );
613  }
614 
615  Info<< nl << endl;
616 
617  return true;
618 }
619 
620 
622 (
623  const scalarField& PSDf
624 ) const
625 {
626  return 10*log10(PSDf/sqr(2e-5));
627 }
628 
629 
631 (
632  const scalarField& Prms2,
633  const scalar f
634 ) const
635 {
636  tmp<scalarField> tspl(10*log10(Prms2/sqr(2e-5)));
637  scalarField& spl = tspl.ref();
638 
639  switch (SPLweighting_)
640  {
641  case weightingType::none:
642  {
643  break;
644  }
645  case weightingType::dBA:
646  {
647  spl += gainA(f);
648  break;
649  }
650  case weightingType::dBB:
651  {
652  spl += gainB(f);
653  break;
654  }
655  case weightingType::dBC:
656  {
657  spl += gainC(f);
658  break;
659  }
660  case weightingType::dBD:
661  {
662  spl += gainD(f);
663  break;
664  }
665  default:
666  {
668  << "Unknown weighting " << weightingTypeNames_[SPLweighting_]
669  << abort(FatalError);
670  }
671  }
672 
673  return tspl;
674 }
675 
676 
678 (
679  const scalarField& Prms2,
680  const scalarField& f
681 ) const
682 {
683  tmp<scalarField> tspl(10*log10(Prms2/sqr(2e-5)));
684  scalarField& spl = tspl.ref();
685 
686  switch (SPLweighting_)
687  {
688  case weightingType::none:
689  {
690  break;
691  }
692  case weightingType::dBA:
693  {
694  forAll(spl, i)
695  {
696  spl[i] += gainA(f[i]);
697  }
698  break;
699  }
700  case weightingType::dBB:
701  {
702  forAll(spl, i)
703  {
704  spl[i] += gainB(f[i]);
705  }
706  break;
707  }
708  case weightingType::dBC:
709  {
710  forAll(spl, i)
711  {
712  spl[i] += gainC(f[i]);
713  }
714  break;
715  }
716  case weightingType::dBD:
717  {
718  forAll(spl, i)
719  {
720  spl[i] += gainD(f[i]);
721  }
722  break;
723  }
724  default:
725  {
727  << "Unknown weighting " << weightingTypeNames_[SPLweighting_]
728  << abort(FatalError);
729  }
730  }
731 
732  return tspl;
733 }
734 
735 
737 {
738  if (planInfo_.active)
739  {
740  planInfo_.active = false;
741  fftw_destroy_plan(planInfo_.plan);
742  fftw_cleanup();
743  }
744 }
745 
746 
748 {
749  scalar f0 = 10;
750  scalar f1 = 20000;
751 
752  OFstream osA("noiseModel-weight-A");
753  OFstream osB("noiseModel-weight-B");
754  OFstream osC("noiseModel-weight-C");
755  OFstream osD("noiseModel-weight-D");
756 
757  for (label f = f0; f <= f1; ++f)
758  {
759  osA << f << " " << gainA(f) << nl;
760  osB << f << " " << gainB(f) << nl;
761  osC << f << " " << gainC(f) << nl;
762  osD << f << " " << gainD(f) << nl;
763  }
764 }
765 
766 
767 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
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
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
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::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::argList::envGlobalPath
static fileName envGlobalPath()
Global case (directory) from environment variable.
Definition: argList.C:537
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< scalar >
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::noiseModel::RBf
scalar RBf(const scalar f) const
B weighting function.
Definition: noiseModel.C:441
Foam::noiseModel::octaves
tmp< scalarField > octaves(const scalarField &data, const scalarField &f, const labelUList &freqBandIDs) const
Generate octave data.
Definition: noiseModel.C:260
Foam::FatalIOError
IOerror FatalIOError
Foam::noiseModel::findStartTimeIndex
label findStartTimeIndex(const instantList &allTimes, const scalar startTime) const
Find and return start time index.
Definition: noiseModel.C:205
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::windowModel::New
static autoPtr< windowModel > New(const dictionary &dict, const label nSamples)
Return a reference to the selected window model.
Definition: windowModelNew.C:33
Foam::HashSet< label, Hash< label > >
Foam::noiseModel::RCf
scalar RCf(const scalar f) const
C weighting function.
Definition: noiseModel.C:463
noiseModel.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::noiseModel::gainB
scalar gainB(const scalar f) const
B weighting as gain in dB.
Definition: noiseModel.C:457
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::noiseModel::minPressure_
scalar minPressure_
Min pressure value.
Definition: noiseModel.H:254
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Foam::Field< scalar >
Foam::noiseModel::planInfo::active
bool active
Definition: noiseModel.H:188
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::log10
dimensionedScalar log10(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
argList.H
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::functionObject::outputPrefix
static word outputPrefix
Directory prefix.
Definition: functionObject.H:352
Foam::noiseModel::writeWeightings
void writeWeightings() const
Helper function to check weightings.
Definition: noiseModel.C:747
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
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::uniformFrequencies
tmp< scalarField > uniformFrequencies(const scalar deltaT) const
Definition: noiseModel.C:239
Foam::constant::atomic::re
const dimensionedScalar re
Classical electron radius: default SI units: [m].
Foam::FatalError
error FatalError
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::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::constant::physicoChemical::c2
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
Foam::BitOps::none
bool none(const UList< bool > &bools)
True if no entries are 'true'.
Definition: BitOps.H:98
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:608
Foam::noiseModel::RDf
scalar RDf(const scalar f) const
D weighting function.
Definition: noiseModel.C:480
Foam::noiseModel::weightingTypeNames_
static const Enum< weightingType > weightingTypeNames_
Definition: noiseModel.H:182
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
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::RMSmeanPf
tmp< scalarField > RMSmeanPf(const scalarField &p) const
Definition: noiseModel.C:364
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::noiseModel::cleanFFTW
void cleanFFTW()
Clean up the FFTW.
Definition: noiseModel.C:736
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::noiseModel::SPL
tmp< scalarField > SPL(const scalarField &Prms2, const scalar f) const
SPL [dB].
Definition: noiseModel.C:631
f
labelList f(nPoints)
Foam::List< label >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::noiseModel::baseFileDir
fileName baseFileDir(const label dataseti) const
Return the base output directory.
Definition: noiseModel.C:224
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Instant::value
scalar value() const
The value (const access)
Definition: Instant.H:99
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
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
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::UList::size
void size(const label n) noexcept
Override size to be inconsistent with allocated storage.
Definition: UListI.H:360
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:401
startTime
Foam::label startTime
Definition: checkTimeOptions.H:1
Foam::instant
An instant of time. Contains the time value and name.
Definition: instant.H:52
N
const Vector< label > N(dict.get< Vector< label >>("N"))
fft.H
functionObject.H
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::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
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::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:417
Foam::noiseModel::Pf
tmp< scalarField > Pf(const scalarField &p) const
Return the fft of the given pressure data.
Definition: noiseModel.C:300
Foam::fft::realTransform1D
static tmp< complexField > realTransform1D(const scalarField &field)
Transform real-value data.
Definition: fft.C:120
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:328
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