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-2021 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  DynamicList<scalar> missedBins;
82 
83  // Convert to lower band limit
84  fTest /= fRatioL2C;
85  while (fTest < fLower)
86  {
87  fTest *= fRatio;
88  }
89 
90  forAll(f, i)
91  {
92  if (f[i] >= fTest)
93  {
94  // Advance band if appropriate
95  label stepi = 0;
96  while (f[i] > fTest)
97  {
98  if (stepi) missedBins.append(fTest/fRatio*fRatioL2C);
99  fTest *= fRatio;
100  ++stepi;
101  }
102  fTest /= fRatio;
103 
104  if (bandIDs.insert(i))
105  {
106  // Also store (next) centre frequency
107  fc.append(fTest*fRatioL2C);
108  }
109  fTest *= fRatio;
110 
111  if (fTest > fUpper)
112  {
113  break;
114  }
115  }
116  }
117 
118  fBandIDs = bandIDs.sortedToc();
119 
120  if (missedBins.size())
121  {
122  label nMiss = missedBins.size();
123  label nTotal = nMiss + fc.size() - 1;
125  << "Empty bands found: " << nMiss << " of " << nTotal
126  << " with centre frequencies " << flatOutput(missedBins)
127  << endl;
128  }
129 
130  if (fc.size())
131  {
132  // Remove the last centre frequency (beyond upper frequency limit)
133  fc.remove();
134 
135  fCentre.transfer(fc);
136  }
137 }
138 
139 
140 namespace Foam
141 {
143  {
144  auto tresult = tmp<scalarField>::New(fld.size(), -GREAT);
145  auto& result = tresult.ref();
146 
147  forAll(result, i)
148  {
149  if (fld[i] > 0)
150  {
151  result[i] = log10(fld[i]);
152  }
153  }
154 
155  return tresult;
156  }
157 }
158 
159 
160 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
161 
163 (
164  const dictionary& dict,
165  const word& lookup,
166  bool& option
167 ) const
168 {
169  dict.readIfPresent(lookup, option);
170 
171  // Only writing on the master process
172  option = option && Pstream::master();
173 
174  if (option)
175  {
176  Info<< " " << lookup << ": " << "yes" << endl;
177  }
178  else
179  {
180  Info<< " " << lookup << ": " << "no" << endl;
181  }
182 }
183 
184 
186 (
187  const scalarList& times
188 ) const
189 {
190  scalar deltaT = -1.0;
191 
192  if (times.size() > 1)
193  {
194  // Uniform time step
195  deltaT = (times.last() - times.first())/scalar(times.size() - 1);
196 
197  for (label i = 1; i < times.size(); ++i)
198  {
199  scalar dT = times[i] - times[i-1];
200 
201  if (mag(dT/deltaT - 1) > 1e-8)
202  {
204  << "Unable to process data with a variable time step"
205  << exit(FatalError);
206  }
207  }
208  }
209  else
210  {
212  << "Unable to create FFT with a single value"
213  << exit(FatalError);
214  }
215 
216  return deltaT;
217 }
218 
219 
221 {
222  forAll(p, i)
223  {
224  if ((p[i] < minPressure_) || (p[i] > maxPressure_))
225  {
227  << "Pressure data at position " << i
228  << " is outside of permitted bounds:" << nl
229  << " pressure: " << p[i] << nl
230  << " minimum pressure: " << minPressure_ << nl
231  << " maximum pressure: " << maxPressure_ << nl
232  << endl;
233 
234  return false;
235  }
236  }
237 
238  return true;
239 }
240 
241 
243 (
244  const instantList& allTimes,
245  const scalar startTime
246 ) const
247 {
248  forAll(allTimes, timeI)
249  {
250  const instant& i = allTimes[timeI];
251 
252  if (i.value() >= startTime)
253  {
254  return timeI;
255  }
256  }
257 
258  return 0;
259 }
260 
261 
263 {
264  return
265  (
268  / "noise"
269  / outputPrefix_
270  / type()
271  / ("input" + Foam::name(dataseti))
272  );
273 }
274 
275 
277 (
278  const scalar deltaT,
279  const bool check
280 ) const
281 {
282  const auto& window = windowModelPtr_();
283  const label N = window.nSamples();
284 
285  auto tf = tmp<scalarField>::New(N/2 + 1, Zero);
286  auto& f = tf.ref();
287 
288  const scalar deltaf = 1.0/(N*deltaT);
289 
290  label nFreq = 0;
291  forAll(f, i)
292  {
293  f[i] = i*deltaf;
294 
295  if (f[i] > fLower_ && f[i] < fUpper_)
296  {
297  ++nFreq;
298  }
299  }
300 
301  if (check && nFreq == 0)
302  {
304  << "No frequencies found in range "
305  << fLower_ << " to " << fUpper_
306  << endl;
307  }
308 
309  return tf;
310 }
311 
312 
314 (
315  const scalarField& data,
316  const scalarField& f,
317  const labelUList& freqBandIDs
318 ) const
319 {
320  if (freqBandIDs.size() < 2)
321  {
323  << "Octave frequency bands are not defined "
324  << "- skipping octaves calculation"
325  << endl;
326 
327  return tmp<scalarField>::New();
328  }
329 
330  auto toctData = tmp<scalarField>::New(freqBandIDs.size() - 1, Zero);
331  auto& octData = toctData.ref();
332 
333  bitSet bandUsed(freqBandIDs.size() - 1);
334  for (label bandI = 0; bandI < freqBandIDs.size() - 1; ++bandI)
335  {
336  label fb0 = freqBandIDs[bandI];
337  label fb1 = freqBandIDs[bandI+1];
338 
339  if (fb0 == fb1) continue;
340 
341  for (label freqI = fb0; freqI < fb1; ++freqI)
342  {
343  label f0 = f[freqI];
344  label f1 = f[freqI + 1];
345  scalar dataAve = 0.5*(data[freqI] + data[freqI + 1]);
346  octData[bandI] += dataAve*(f1 - f0);
347 
348  bandUsed.set(bandI);
349  }
350  }
351 
352  bandUsed.flip();
353  labelList bandUnused = bandUsed.sortedToc();
354  if (bandUnused.size())
355  {
357  << "Empty bands found: " << bandUnused.size() << " of "
358  << bandUsed.size() << endl;
359  }
360 
361  return toctData;
362 }
363 
364 
366 {
367  if (planInfo_.active)
368  {
369  if (p.size() != planInfo_.windowSize)
370  {
372  << "Expected pressure data to have " << planInfo_.windowSize
373  << " values, but received " << p.size() << " values"
374  << abort(FatalError);
375  }
376 
377  List<double>& in = planInfo_.in;
378  const List<double>& out = planInfo_.out;
379  forAll(in, i)
380  {
381  in[i] = p[i];
382  }
383 
384  ::fftw_execute(planInfo_.plan);
385 
386  const label n = planInfo_.windowSize;
387  const label nBy2 = n/2;
388  auto tresult = tmp<scalarField>::New(nBy2 + 1);
389  auto& result = tresult.ref();
390 
391  // 0 th value = DC component
392  // nBy2 th value = real only if n is even
393  result[0] = out[0];
394  for (label i = 1; i <= nBy2; ++i)
395  {
396  const auto re = out[i];
397  const auto im = out[n - i];
398  result[i] = sqrt(re*re + im*im);
399  }
400 
401  return tresult;
402  }
403 
404  return mag(fft::realTransform1D(p));
405 }
406 
407 
409 {
410  const auto& window = windowModelPtr_();
411  const label N = window.nSamples();
412  const label nWindow = window.nWindow();
413 
414  auto tmeanPf = tmp<scalarField>::New(N/2 + 1, Zero);
415  auto& meanPf = tmeanPf.ref();
416 
417  for (label windowI = 0; windowI < nWindow; ++windowI)
418  {
419  meanPf += Pf(window.apply<scalar>(p, windowI));
420  }
421 
422  meanPf /= scalar(nWindow);
423 
424  return tmeanPf;
425 }
426 
427 
429 (
430  const scalarField& p
431 ) const
432 {
433  const auto& window = windowModelPtr_();
434  const label N = window.nSamples();
435  const label nWindow = window.nWindow();
436 
437  scalarField RMSMeanPf(N/2 + 1, Zero);
438  for (label windowI = 0; windowI < nWindow; ++windowI)
439  {
440  RMSMeanPf += sqr(Pf(window.apply<scalar>(p, windowI)));
441  }
442 
443  return sqrt(RMSMeanPf/scalar(nWindow))/scalar(N);
444 }
445 
446 
448 (
449  const scalarField& p,
450  const scalar deltaT
451 ) const
452 {
453  const auto& window = windowModelPtr_();
454  const label N = window.nSamples();
455  const label nWindow = window.nWindow();
456 
457  auto tpsd = tmp<scalarField>::New(N/2 + 1, Zero);
458  auto& psd = tpsd.ref();
459 
460  for (label windowI = 0; windowI < nWindow; ++windowI)
461  {
462  psd += sqr(Pf(window.apply<scalar>(p, windowI)));
463  }
464 
465  scalar fs = 1.0/deltaT;
466  psd /= scalar(nWindow)*fs*N;
467 
468  // Scaling due to use of 1-sided FFT
469  // Note: do not scale DC component
470  psd *= 2;
471  psd.first() /= 2;
472  psd.last() /= 2;
473 
474  if (debug)
475  {
476  Pout<< "Average PSD: " << average(psd) << endl;
477  }
478 
479  return tpsd;
480 }
481 
482 
483 Foam::scalar Foam::noiseModel::RAf(const scalar f) const
484 {
485  const scalar c1 = sqr(12194.0);
486  const scalar c2 = sqr(20.6);
487  const scalar c3 = sqr(107.7);
488  const scalar c4 = sqr(739.9);
489 
490  const scalar f2 = f*f;
491 
492  return
493  c1*sqr(f2)
494  /(
495  (f2 + c2)*sqrt((f2 + c3)*(f2 + c4))*(f2 + c1)
496  );
497 }
498 
499 
500 Foam::scalar Foam::noiseModel::gainA(const scalar f) const
501 {
502  if (f < SMALL)
503  {
504  return 0;
505  }
506 
507  return 20*log10(RAf(f)) - 20*log10(RAf(1000));
508 }
509 
510 
511 Foam::scalar Foam::noiseModel::RBf(const scalar f) const
512 {
513  const scalar c1 = sqr(12194.0);
514  const scalar c2 = sqr(20.6);
515  const scalar c3 = sqr(158.5);
516 
517  const scalar f2 = f*f;
518 
519  return
520  c1*f2*f
521  /(
522  (f2 + c2)*sqrt(f2 + c3)*(f2 + c1)
523  );
524 }
525 
526 
527 Foam::scalar Foam::noiseModel::gainB(const scalar f) const
528 {
529  if (f < SMALL)
530  {
531  return 0;
532  }
533 
534  return 20*log10(RBf(f)) - 20*log10(RBf(1000));
535 }
536 
537 
538 Foam::scalar Foam::noiseModel::RCf(const scalar f) const
539 {
540  const scalar c1 = sqr(12194.0);
541  const scalar c2 = sqr(20.6);
542 
543  const scalar f2 = f*f;
544 
545  return c1*f2/((f2 + c2)*(f2 + c1));
546 }
547 
548 
549 Foam::scalar Foam::noiseModel::gainC(const scalar f) const
550 {
551  if (f < SMALL)
552  {
553  return 0;
554  }
555 
556  return 20*log10(RCf(f)) - 20*log10(RCf(1000));
557 }
558 
559 
560 Foam::scalar Foam::noiseModel::RDf(const scalar f) const
561 {
562  const scalar f2 = f*f;
563 
564  const scalar hf =
565  (sqr(1037918.48 - f2) + 1080768.16*f2)
566  /(sqr(9837328 - f2) + 11723776*f2);
567 
568  return
569  f/6.8966888496476e-5*Foam::sqrt(hf/((f2 + 79919.29)*(f2 + 1345600)));
570 }
571 
572 
573 Foam::scalar Foam::noiseModel::gainD(const scalar f) const
574 {
575  if (f < SMALL)
576  {
577  return 0;
578  }
579 
580  return 20*log10(RDf(f));
581 }
582 
583 
584 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
585 
586 Foam::noiseModel::noiseModel(const dictionary& dict, const bool readFields)
587 :
588  dict_(dict),
589  rhoRef_(1),
590  nSamples_(65536),
591  fLower_(25),
592  fUpper_(10000),
593  startTime_(0),
594  windowModelPtr_(),
595  graphFormat_("raw"),
596  SPLweighting_(weightingType::none),
597  dBRef_(2e-5),
598  minPressure_(-0.5*VGREAT),
599  maxPressure_(0.5*VGREAT),
600  outputPrefix_(),
601  writePrmsf_(true),
602  writeSPL_(true),
603  writePSD_(true),
604  writePSDf_(true),
605  writeOctaves_(true),
606  planInfo_()
607 {
608  planInfo_.active = false;
609 
610  if (readFields)
611  {
612  read(dict);
613  }
614 
615  if (debug)
616  {
617  writeWeightings();
618  }
619 }
620 
621 
622 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
623 
625 {
626  dict.readIfPresent("rhoRef", rhoRef_);
627  dict.readIfPresent("N", nSamples_);
628  dict.readIfPresent("fl", fLower_);
629  dict.readIfPresent("fu", fUpper_);
630  dict.readIfPresent("startTime", startTime_);
631  dict.readIfPresent("graphFormat", graphFormat_);
632  dict.readIfPresent("minPressure", minPressure_);
633  dict.readIfPresent("maxPressure", maxPressure_);
634  dict.readIfPresent("outputPrefix", outputPrefix_);
635 
636  if (fLower_ < 0)
637  {
639  << "fl: lower frequency bound must be greater than zero"
640  << exit(FatalIOError);
641  }
642 
643  if (fUpper_ < 0)
644  {
646  << "fu: upper frequency bound must be greater than zero"
647  << exit(FatalIOError);
648  }
649 
650  if (fUpper_ < fLower_)
651  {
653  << "fu: upper frequency bound must be greater than lower "
654  << "frequency bound (fl)"
655  << exit(FatalIOError);
656 
657  }
658 
659  Info<< " Frequency bounds:" << nl
660  << " lower: " << fLower_ << nl
661  << " upper: " << fUpper_ << endl;
662 
663  weightingTypeNames_.readIfPresent("SPLweighting", dict, SPLweighting_);
664 
665  Info<< " Weighting: " << weightingTypeNames_[SPLweighting_] << endl;
666 
667  if (dict.readIfPresent("dBRef", dBRef_))
668  {
669  Info<< " Reference for dB calculation: " << dBRef_ << endl;
670  }
671 
672  Info<< " Write options:" << endl;
673  dictionary optDict(dict.subOrEmptyDict("writeOptions"));
674  readWriteOption(optDict, "writePrmsf", writePrmsf_);
675  readWriteOption(optDict, "writeSPL", writeSPL_);
676  readWriteOption(optDict, "writePSD", writePSD_);
677  readWriteOption(optDict, "writePSDf", writePSDf_);
678  readWriteOption(optDict, "writeOctaves", writeOctaves_);
679 
680 
681  windowModelPtr_ = windowModel::New(dict, nSamples_);
682 
683  cleanFFTW();
684 
685  const label windowSize = windowModelPtr_->nSamples();
686 
687  if (windowSize > 1)
688  {
689  planInfo_.active = true;
690  planInfo_.windowSize = windowSize;
691  planInfo_.in.setSize(windowSize);
692  planInfo_.out.setSize(windowSize);
693 
694  // Using real to half-complex fftw 'kind'
695  planInfo_.plan =
696  fftw_plan_r2r_1d
697  (
698  windowSize,
699  planInfo_.in.data(),
700  planInfo_.out.data(),
701  FFTW_R2HC,
702  windowSize <= 8192 ? FFTW_MEASURE : FFTW_ESTIMATE
703  );
704  }
705 
706  Info<< nl << endl;
707 
708  return true;
709 }
710 
711 
713 (
714  const scalarField& PSDf
715 ) const
716 {
717  return 10*safeLog10(PSDf/sqr(dBRef_));
718 }
719 
720 
722 (
723  const scalarField& Prms2,
724  const scalar f
725 ) const
726 {
727  tmp<scalarField> tspl(10*safeLog10(Prms2/sqr(dBRef_)));
728  scalarField& spl = tspl.ref();
729 
730  switch (SPLweighting_)
731  {
732  case weightingType::none:
733  {
734  break;
735  }
736  case weightingType::dBA:
737  {
738  spl += gainA(f);
739  break;
740  }
741  case weightingType::dBB:
742  {
743  spl += gainB(f);
744  break;
745  }
746  case weightingType::dBC:
747  {
748  spl += gainC(f);
749  break;
750  }
751  case weightingType::dBD:
752  {
753  spl += gainD(f);
754  break;
755  }
756  default:
757  {
759  << "Unknown weighting " << weightingTypeNames_[SPLweighting_]
760  << abort(FatalError);
761  }
762  }
763 
764  return tspl;
765 }
766 
767 
769 (
770  const scalarField& Prms2,
771  const scalarField& f
772 ) const
773 {
774  tmp<scalarField> tspl(10*safeLog10(Prms2/sqr(dBRef_)));
775  scalarField& spl = tspl.ref();
776 
777  switch (SPLweighting_)
778  {
779  case weightingType::none:
780  {
781  break;
782  }
783  case weightingType::dBA:
784  {
785  forAll(spl, i)
786  {
787  spl[i] += gainA(f[i]);
788  }
789  break;
790  }
791  case weightingType::dBB:
792  {
793  forAll(spl, i)
794  {
795  spl[i] += gainB(f[i]);
796  }
797  break;
798  }
799  case weightingType::dBC:
800  {
801  forAll(spl, i)
802  {
803  spl[i] += gainC(f[i]);
804  }
805  break;
806  }
807  case weightingType::dBD:
808  {
809  forAll(spl, i)
810  {
811  spl[i] += gainD(f[i]);
812  }
813  break;
814  }
815  default:
816  {
818  << "Unknown weighting " << weightingTypeNames_[SPLweighting_]
819  << abort(FatalError);
820  }
821  }
822 
823  return tspl;
824 }
825 
826 
828 {
829  if (planInfo_.active)
830  {
831  planInfo_.active = false;
832  fftw_destroy_plan(planInfo_.plan);
833  fftw_cleanup();
834  }
835 }
836 
837 
839 {
840  scalar f0 = 10;
841  scalar f1 = 20000;
842 
843  OFstream osA("noiseModel-weight-A");
844  OFstream osB("noiseModel-weight-B");
845  OFstream osC("noiseModel-weight-C");
846  OFstream osD("noiseModel-weight-D");
847 
848  for (label f = f0; f <= f1; ++f)
849  {
850  osA << f << " " << gainA(f) << nl;
851  osB << f << " " << gainB(f) << nl;
852  osC << f << " " << gainC(f) << nl;
853  osD << f << " " << gainD(f) << nl;
854  }
855 }
856 
857 
858 // ************************************************************************* //
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:186
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:624
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:163
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::argList::envGlobalPath
static fileName envGlobalPath()
Global case (directory) from environment variable.
Definition: argList.C:549
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:511
Foam::noiseModel::octaves
tmp< scalarField > octaves(const scalarField &data, const scalarField &f, const labelUList &freqBandIDs) const
Generate octave data.
Definition: noiseModel.C:314
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:243
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::safeLog10
tmp< scalarField > safeLog10(const scalarField &fld)
Definition: noiseModel.C:142
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:538
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:527
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::noiseModel::minPressure_
scalar minPressure_
Min pressure value.
Definition: noiseModel.H:263
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Foam::check
static void check(const int retVal, const char *what)
Definition: ptscotchDecomp.C:80
Foam::Field< scalar >
Foam::noiseModel::planInfo::active
bool active
Definition: noiseModel.H:194
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::log10
dimensionedScalar log10(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
argList.H
Foam::noiseModel::weightingType
weightingType
Definition: noiseModel.H:179
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:376
Foam::noiseModel::writeWeightings
void writeWeightings() const
Helper function to check weightings.
Definition: noiseModel.C:838
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
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::noiseModel::gainC
scalar gainC(const scalar f) const
C weighting as gain in dB.
Definition: noiseModel.C:549
dict
dictionary dict
Definition: searchingEngine.H:14
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:123
Foam::noiseModel::validateBounds
bool validateBounds(const scalarList &p) const
Return true if all pressure data is within min/max bounds.
Definition: noiseModel.C:220
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::noiseModel::maxPressure_
scalar maxPressure_
Min pressure value.
Definition: noiseModel.H:266
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:216
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:540
Foam::noiseModel::RDf
scalar RDf(const scalar f) const
D weighting function.
Definition: noiseModel.C:560
Foam::noiseModel::weightingTypeNames_
static const Enum< weightingType > weightingTypeNames_
Definition: noiseModel.H:188
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:429
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::noiseModel::cleanFFTW
void cleanFFTW()
Clean up the FFTW.
Definition: noiseModel.C:827
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::noiseModel::SPL
tmp< scalarField > SPL(const scalarField &Prms2, const scalar f) const
SPL [dB].
Definition: noiseModel.C:722
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:262
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:573
Foam::noiseModel::PSDf
tmp< scalarField > PSDf(const scalarField &p, const scalar deltaT) const
Definition: noiseModel.C:448
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::noiseModel::uniformFrequencies
tmp< scalarField > uniformFrequencies(const scalar deltaT, const bool check) const
Definition: noiseModel.C:277
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
startTime
Foam::label startTime
Definition: checkTimeOptions.H:1
Foam::instant
An instant of time. Contains the time value and name.
Definition: instant.H:52
Foam::UList::size
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
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:500
Foam::noiseModel::PSD
tmp< Foam::scalarField > PSD(const scalarField &PSDf) const
PSD [dB/Hz].
Definition: noiseModel.C:713
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::noiseModel::RAf
scalar RAf(const scalar f) const
A weighting function.
Definition: noiseModel.C:483
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:293
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::noiseModel::Pf
tmp< scalarField > Pf(const scalarField &p) const
Return the fft of the given pressure data.
Definition: noiseModel.C:365
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:408