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-2022 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
80  DynamicList<scalar> fc;
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.pop_back();
134 
135  fCentre.transfer(fc);
136  }
137 }
138 
139 
140 namespace Foam
141 {
142  tmp<scalarField> safeLog10(const scalarField& fld)
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 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
161 
162 namespace Foam
163 {
164 
165 // Get bool option (eg, read/write) and provide Info feedback
166 static void readWriteOption
167 (
168  const dictionary& dict,
169  const word& lookup,
170  bool& option
171 )
172 {
173  dict.readIfPresent(lookup, option);
174 
175  Info<< " " << lookup << ": " << (option ? "yes" : "no") << endl;
176 }
178 } // End namespace Foam
179 
180 
181 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
182 
184 (
185  const scalarList& times
186 ) const
187 {
188  scalar deltaT = -1.0;
189 
190  if (times.size() > 1)
191  {
192  // Uniform time step
193  deltaT = (times.last() - times.first())/scalar(times.size() - 1);
194 
195  for (label i = 1; i < times.size(); ++i)
196  {
197  scalar dT = times[i] - times[i-1];
198 
199  if (mag(dT/deltaT - 1) > 1e-8)
200  {
202  << "Unable to process data with a variable time step"
203  << exit(FatalError);
204  }
205  }
206  }
207  else
208  {
210  << "Unable to create FFT with a single value"
212  }
213 
214  return deltaT;
215 }
216 
217 
219 {
220  forAll(p, i)
221  {
222  if ((p[i] < minPressure_) || (p[i] > maxPressure_))
223  {
225  << "Pressure data at position " << i
226  << " is outside of permitted bounds:" << nl
227  << " pressure: " << p[i] << nl
228  << " minimum pressure: " << minPressure_ << nl
229  << " maximum pressure: " << maxPressure_ << nl
230  << endl;
231 
232  return false;
233  }
234  }
235 
236  return true;
237 }
238 
239 
240 Foam::fileName Foam::noiseModel::baseFileDir(const label dataseti) const
241 {
242  return
243  (
244  outputPrefix_
245  / type()
246  / ("input" + Foam::name(dataseti))
247  );
248 }
249 
250 
252 (
253  Ostream& os,
254  const string& x,
255  const string& y,
256  const UList<Tuple2<string, token>>& headerValues
257 ) const
258 {
259  writeHeader(os, x + " vs " + y);
260  writeHeaderValue(os, "Lower frequency", fLower_);
261  writeHeaderValue(os, "Upper frequency", fUpper_);
262  writeHeaderValue(os, "Window model", windowModelPtr_->type());
263  writeHeaderValue(os, "Window number", windowModelPtr_->nWindow());
264  writeHeaderValue(os, "Window samples", windowModelPtr_->nSamples());
265  writeHeaderValue(os, "Window overlap %", windowModelPtr_->overlapPercent());
266  writeHeaderValue(os, "dBRef", dBRef_);
267 
268  for (const auto& hv : headerValues)
269  {
270  writeHeaderValue(os, hv.first(), hv.second());
271  }
272 
273  writeCommented(os, x.substr(0, x.find(' ')));
274  writeTabbed(os, y.substr(0, y.find(' ')));
275  os << nl;
276 }
277 
278 
280 (
281  Ostream& os,
282  const scalarField& f,
283  const scalarField& fx
284 ) const
285 {
286  forAll(f, i)
287  {
288  if (f[i] >= fLower_ && f[i] <= fUpper_)
289  {
290  os << f[i] << tab << fx[i] << nl;
291  }
292  }
293 }
294 
295 
297 (
298  const scalar deltaT,
299  const bool check
300 ) const
301 {
302  const auto& window = windowModelPtr_();
303  const label N = window.nSamples();
304 
305  auto tf = tmp<scalarField>::New(N/2 + 1, Zero);
306  auto& f = tf.ref();
307 
308  const scalar deltaf = 1.0/(N*deltaT);
309 
310  label nFreq = 0;
311  forAll(f, i)
312  {
313  f[i] = i*deltaf;
314 
315  if (f[i] > fLower_ && f[i] < fUpper_)
316  {
317  ++nFreq;
318  }
319  }
320 
321  if (check && nFreq == 0)
322  {
324  << "No frequencies found in range "
325  << fLower_ << " to " << fUpper_
326  << endl;
327  }
328 
329  return tf;
330 }
331 
332 
334 (
335  const scalarField& data,
336  const scalarField& f,
337  const labelUList& freqBandIDs
338 ) const
339 {
340  if (freqBandIDs.size() < 2)
341  {
343  << "Octave frequency bands are not defined "
344  << "- skipping octaves calculation"
345  << endl;
346 
347  return tmp<scalarField>::New();
348  }
349 
350  auto toctData = tmp<scalarField>::New(freqBandIDs.size() - 1, Zero);
351  auto& octData = toctData.ref();
352 
353  bitSet bandUsed(freqBandIDs.size() - 1);
354  for (label bandI = 0; bandI < freqBandIDs.size() - 1; ++bandI)
355  {
356  label fb0 = freqBandIDs[bandI];
357  label fb1 = freqBandIDs[bandI+1];
358 
359  if (fb0 == fb1) continue;
360 
361  for (label freqI = fb0; freqI < fb1; ++freqI)
362  {
363  label f0 = f[freqI];
364  label f1 = f[freqI + 1];
365  scalar dataAve = 0.5*(data[freqI] + data[freqI + 1]);
366  octData[bandI] += dataAve*(f1 - f0);
367 
368  bandUsed.set(bandI);
369  }
370  }
371 
372  bandUsed.flip();
373  labelList bandUnused = bandUsed.sortedToc();
374  if (bandUnused.size())
375  {
377  << "Empty bands found: " << bandUnused.size() << " of "
378  << bandUsed.size() << endl;
379  }
380 
381  return toctData;
382 }
383 
384 
386 {
387  if (planInfo_.active)
388  {
389  if (p.size() != planInfo_.windowSize)
390  {
392  << "Expected pressure data to have " << planInfo_.windowSize
393  << " values, but received " << p.size() << " values"
394  << abort(FatalError);
395  }
396 
397  List<double>& in = planInfo_.in;
398  const List<double>& out = planInfo_.out;
399  forAll(in, i)
400  {
401  in[i] = p[i];
402  }
403 
404  ::fftw_execute(planInfo_.plan);
405 
406  const label n = planInfo_.windowSize;
407  const label nBy2 = n/2;
408  auto tresult = tmp<scalarField>::New(nBy2 + 1);
409  auto& result = tresult.ref();
410 
411  // 0 th value = DC component
412  // nBy2 th value = real only if n is even
413  result[0] = out[0];
414  for (label i = 1; i <= nBy2; ++i)
415  {
416  const auto re = out[i];
417  const auto im = out[n - i];
418  result[i] = sqrt(re*re + im*im);
419  }
420 
421  return tresult;
422  }
423 
424  return mag(fft::realTransform1D(p));
425 }
426 
427 
429 {
430  const auto& window = windowModelPtr_();
431  const label N = window.nSamples();
432  const label nWindow = window.nWindow();
433 
434  auto tmeanPf = tmp<scalarField>::New(N/2 + 1, Zero);
435  auto& meanPf = tmeanPf.ref();
436 
437  for (label windowI = 0; windowI < nWindow; ++windowI)
438  {
439  meanPf += Pf(window.apply<scalar>(p, windowI));
440  }
441 
442  meanPf /= scalar(nWindow);
443 
444  return tmeanPf;
445 }
446 
447 
449 (
450  const scalarField& p
451 ) const
452 {
453  const auto& window = windowModelPtr_();
454  const label N = window.nSamples();
455  const label nWindow = window.nWindow();
456 
457  scalarField RMSMeanPf(N/2 + 1, Zero);
458  for (label windowI = 0; windowI < nWindow; ++windowI)
459  {
460  RMSMeanPf += sqr(Pf(window.apply<scalar>(p, windowI)));
461  }
462 
463  return sqrt(RMSMeanPf/scalar(nWindow))/scalar(N);
464 }
465 
466 
468 (
469  const scalarField& p,
470  const scalar deltaT
471 ) const
472 {
473  const auto& window = windowModelPtr_();
474  const label N = window.nSamples();
475  const label nWindow = window.nWindow();
476 
477  auto tpsd = tmp<scalarField>::New(N/2 + 1, Zero);
478  auto& psd = tpsd.ref();
479 
480  for (label windowI = 0; windowI < nWindow; ++windowI)
481  {
482  psd += sqr(Pf(window.apply<scalar>(p, windowI)));
483  }
484 
485  scalar fs = 1.0/deltaT;
486  psd /= scalar(nWindow)*fs*N;
487 
488  // Scaling due to use of 1-sided FFT
489  // Note: do not scale DC component
490  psd *= 2;
491  psd.first() /= 2;
492  psd.last() /= 2;
493 
494  if (debug)
495  {
496  Pout<< "Average PSD: " << average(psd) << endl;
497  }
498 
499  return tpsd;
500 }
501 
502 
503 Foam::scalar Foam::noiseModel::RAf(const scalar f) const
504 {
505  const scalar c1 = sqr(12194.0);
506  const scalar c2 = sqr(20.6);
507  const scalar c3 = sqr(107.7);
508  const scalar c4 = sqr(739.9);
509 
510  const scalar f2 = f*f;
511 
512  return
513  c1*sqr(f2)
514  /(
515  (f2 + c2)*sqrt((f2 + c3)*(f2 + c4))*(f2 + c1)
516  );
517 }
518 
519 
520 Foam::scalar Foam::noiseModel::gainA(const scalar f) const
521 {
522  if (f < SMALL)
523  {
524  return 0;
525  }
526 
527  return 20*log10(RAf(f)) - 20*log10(RAf(1000));
528 }
529 
530 
531 Foam::scalar Foam::noiseModel::RBf(const scalar f) const
532 {
533  const scalar c1 = sqr(12194.0);
534  const scalar c2 = sqr(20.6);
535  const scalar c3 = sqr(158.5);
536 
537  const scalar f2 = f*f;
538 
539  return
540  c1*f2*f
541  /(
542  (f2 + c2)*sqrt(f2 + c3)*(f2 + c1)
543  );
544 }
545 
546 
547 Foam::scalar Foam::noiseModel::gainB(const scalar f) const
548 {
549  if (f < SMALL)
550  {
551  return 0;
552  }
553 
554  return 20*log10(RBf(f)) - 20*log10(RBf(1000));
555 }
556 
557 
558 Foam::scalar Foam::noiseModel::RCf(const scalar f) const
559 {
560  const scalar c1 = sqr(12194.0);
561  const scalar c2 = sqr(20.6);
563  const scalar f2 = f*f;
564 
565  return c1*f2/((f2 + c2)*(f2 + c1));
566 }
567 
568 
569 Foam::scalar Foam::noiseModel::gainC(const scalar f) const
570 {
571  if (f < SMALL)
572  {
573  return 0;
574  }
575 
576  return 20*log10(RCf(f)) - 20*log10(RCf(1000));
577 }
578 
579 
580 Foam::scalar Foam::noiseModel::RDf(const scalar f) const
581 {
582  const scalar f2 = f*f;
583 
584  const scalar hf =
585  (sqr(1037918.48 - f2) + 1080768.16*f2)
586  /(sqr(9837328 - f2) + 11723776*f2);
587 
588  return
589  f/6.8966888496476e-5*Foam::sqrt(hf/((f2 + 79919.29)*(f2 + 1345600)));
590 }
591 
592 
593 Foam::scalar Foam::noiseModel::gainD(const scalar f) const
594 {
595  if (f < SMALL)
596  {
597  return 0;
598  }
599 
600  return 20*log10(RDf(f));
601 }
602 
603 
604 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
605 
607 (
608  const dictionary& dict,
609  const objectRegistry& obr,
610  const word& name,
611  const bool readFields
612 )
613 :
614  functionObjects::writeFile(obr, "noise"),
615  dict_(dict),
616  rhoRef_(1),
617  nSamples_(65536),
618  fLower_(25),
619  fUpper_(10000),
620  startTime_(0),
621  windowModelPtr_(),
622  SPLweighting_(weightingType::none),
623  dBRef_(2e-5),
624  minPressure_(-0.5*VGREAT),
625  maxPressure_(0.5*VGREAT),
626  outputPrefix_(),
627  writePrmsf_(true),
628  writeSPL_(true),
629  writePSD_(true),
630  writePSDf_(true),
631  writeOctaves_(true),
632  planInfo_()
633 {
634  planInfo_.active = false;
635 
636  if (readFields)
637  {
638  read(dict);
639  }
640 
641  if (debug)
642  {
644  }
645 }
646 
647 
648 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
649 
650 bool Foam::noiseModel::read(const dictionary& dict)
651 {
653  {
654  return false;
655  }
656 
657  dict.readIfPresent("rhoRef", rhoRef_);
658  dict.readIfPresent("N", nSamples_);
659  dict.readIfPresent("fl", fLower_);
660  dict.readIfPresent("fu", fUpper_);
661  dict.readIfPresent("startTime", startTime_);
662  dict.readIfPresent("minPressure", minPressure_);
663  dict.readIfPresent("maxPressure", maxPressure_);
664  dict.readIfPresent("outputPrefix", outputPrefix_);
665 
666  if (fLower_ < 0)
667  {
669  << "fl: lower frequency bound must be greater than zero"
670  << exit(FatalIOError);
671  }
672 
673  if (fUpper_ < 0)
674  {
676  << "fu: upper frequency bound must be greater than zero"
677  << exit(FatalIOError);
678  }
679 
680  if (fUpper_ < fLower_)
681  {
683  << "fu: upper frequency bound must be greater than lower "
684  << "frequency bound (fl)"
685  << exit(FatalIOError);
686 
687  }
688 
689  Info<< " Frequency bounds:" << nl
690  << " lower: " << fLower_ << nl
691  << " upper: " << fUpper_ << endl;
692 
693  weightingTypeNames_.readIfPresent("SPLweighting", dict, SPLweighting_);
694 
695  Info<< " Weighting: " << weightingTypeNames_[SPLweighting_] << endl;
696 
697  if (dict.readIfPresent("dBRef", dBRef_))
698  {
699  Info<< " Reference for dB calculation: " << dBRef_ << endl;
700  }
701 
702  Info<< " Write options:" << endl;
703  dictionary optDict(dict.subOrEmptyDict("writeOptions"));
704  readWriteOption(optDict, "writePrmsf", writePrmsf_);
705  readWriteOption(optDict, "writeSPL", writeSPL_);
706  readWriteOption(optDict, "writePSD", writePSD_);
707  readWriteOption(optDict, "writePSDf", writePSDf_);
708  readWriteOption(optDict, "writeOctaves", writeOctaves_);
709 
710 
711  windowModelPtr_ = windowModel::New(dict, nSamples_);
712 
713  cleanFFTW();
714 
715  const label windowSize = windowModelPtr_->nSamples();
716 
717  if (windowSize > 1)
718  {
719  planInfo_.active = true;
720  planInfo_.windowSize = windowSize;
721  planInfo_.in.setSize(windowSize);
722  planInfo_.out.setSize(windowSize);
723 
724  // Using real to half-complex fftw 'kind'
725  planInfo_.plan =
726  fftw_plan_r2r_1d
727  (
728  windowSize,
729  planInfo_.in.data(),
730  planInfo_.out.data(),
731  FFTW_R2HC,
732  windowSize <= 8192 ? FFTW_MEASURE : FFTW_ESTIMATE
733  );
734  }
735 
737 
738  return true;
739 }
740 
741 
743 (
744  const scalarField& PSDf
745 ) const
746 {
747  return 10*safeLog10(PSDf/sqr(dBRef_));
748 }
749 
750 
752 (
753  const scalarField& Prms2,
754  const scalar f
755 ) const
756 {
757  tmp<scalarField> tspl(10*safeLog10(Prms2/sqr(dBRef_)));
758  scalarField& spl = tspl.ref();
759 
760  switch (SPLweighting_)
761  {
762  case weightingType::none:
763  {
764  break;
765  }
766  case weightingType::dBA:
767  {
768  spl += gainA(f);
769  break;
770  }
771  case weightingType::dBB:
772  {
773  spl += gainB(f);
774  break;
775  }
776  case weightingType::dBC:
777  {
778  spl += gainC(f);
779  break;
780  }
781  case weightingType::dBD:
782  {
783  spl += gainD(f);
784  break;
785  }
786  default:
787  {
789  << "Unknown weighting " << weightingTypeNames_[SPLweighting_]
790  << abort(FatalError);
791  }
792  }
793 
794  return tspl;
795 }
796 
797 
799 (
800  const scalarField& Prms2,
801  const scalarField& f
802 ) const
803 {
804  tmp<scalarField> tspl(10*safeLog10(Prms2/sqr(dBRef_)));
805  scalarField& spl = tspl.ref();
806 
807  switch (SPLweighting_)
808  {
809  case weightingType::none:
810  {
811  break;
812  }
813  case weightingType::dBA:
814  {
815  forAll(spl, i)
816  {
817  spl[i] += gainA(f[i]);
818  }
819  break;
820  }
821  case weightingType::dBB:
822  {
823  forAll(spl, i)
824  {
825  spl[i] += gainB(f[i]);
826  }
827  break;
828  }
829  case weightingType::dBC:
830  {
831  forAll(spl, i)
832  {
833  spl[i] += gainC(f[i]);
834  }
835  break;
836  }
837  case weightingType::dBD:
838  {
839  forAll(spl, i)
840  {
841  spl[i] += gainD(f[i]);
842  }
843  break;
844  }
845  default:
846  {
848  << "Unknown weighting " << weightingTypeNames_[SPLweighting_]
849  << abort(FatalError);
850  }
851  }
852 
853  return tspl;
854 }
855 
856 
858 {
859  if (planInfo_.active)
860  {
861  planInfo_.active = false;
862  fftw_destroy_plan(planInfo_.plan);
863  fftw_cleanup();
864  }
865 }
866 
867 
869 {
870  scalar f0 = 10;
871  scalar f1 = 20000;
872 
873  OFstream osA("noiseModel-weight-A");
874  OFstream osB("noiseModel-weight-B");
875  OFstream osC("noiseModel-weight-C");
876  OFstream osD("noiseModel-weight-D");
877 
878  for (label f = f0; f <= f1; ++f)
879  {
880  osA << f << " " << gainA(f) << nl;
881  osB << f << " " << gainB(f) << nl;
882  osC << f << " " << gainC(f) << nl;
883  osD << f << " " << gainD(f) << nl;
884  }
885 }
886 
887 
888 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
void cleanFFTW()
Clean up the FFTW.
Definition: noiseModel.C:850
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1)
tmp< scalarField > meanPf(const scalarField &p) const
Return the multi-window mean fft of the complete pressure data [Pa].
Definition: noiseModel.C:421
scalar gainD(const scalar f) const
D weighting as gain in dB.
Definition: noiseModel.C:586
A class for handling file names.
Definition: fileName.H:71
static void writeHeader(Ostream &os, const word &fieldName)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
bool none(const UList< bool > &bools)
True if no entries are &#39;true&#39;.
Definition: BitOps.H:103
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:54
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:472
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionedScalar re
Classical electron radius: default SI units: [m].
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
T & first()
Access first element of the list, position [0].
Definition: UList.H:814
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:48
scalar gainC(const scalar f) const
C weighting as gain in dB.
Definition: noiseModel.C:562
tmp< scalarField > Pf(const scalarField &p) const
Return the fft of the given pressure data.
Definition: noiseModel.C:378
Lookup type of boundary radiation properties.
Definition: lookup.H:57
scalar RBf(const scalar f) const
B weighting function.
Definition: noiseModel.C:524
bool validateBounds(const scalarList &p) const
Return true if all pressure data is within min/max bounds.
Definition: noiseModel.C:211
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
scalar checkUniformTimeStep(const scalarList &times) const
Check and return uniform time step.
Definition: noiseModel.C:177
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
scalar y
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: noiseModel.C:643
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
tmp< scalarField > octaves(const scalarField &data, const scalarField &f, const labelUList &freqBandIDs) const
Generate octave data.
Definition: noiseModel.C:327
noiseModel(const noiseModel &)=delete
No copy construct.
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:245
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
scalar RDf(const scalar f) const
D weighting function.
Definition: noiseModel.C:573
static const Enum< weightingType > weightingTypeNames_
Definition: noiseModel.H:181
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:461
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:212
scalar RAf(const scalar f) const
A weighting function.
Definition: noiseModel.C:496
void writeFreqDataToFile(Ostream &os, const scalarField &f, const scalarField &fx) const
Definition: noiseModel.C:273
errorManip< error > abort(error &err)
Definition: errorManip.H:139
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:99
tmp< Foam::scalarField > PSD(const scalarField &PSDf) const
PSD [dB/Hz].
Definition: noiseModel.C:736
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
int debug
Static debugging option.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
planInfo planInfo_
Plan information for FFTW.
Definition: noiseModel.H:313
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
static void check(const int retVal, const char *what)
labelList f(nPoints)
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Find and return a sub-dictionary as a copy, otherwise return an empty dictionary. ...
Definition: dictionary.C:533
T & last()
Access last element of the list, position [size()-1].
Definition: UList.H:828
const Vector< label > N(dict.get< Vector< label >>("N"))
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))
static tmp< complexField > realTransform1D(const scalarField &field)
Transform real-value data.
Definition: fft.C:113
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
scalar gainB(const scalar f) const
B weighting as gain in dB.
Definition: noiseModel.C:540
static void readWriteOption(const dictionary &dict, const word &lookup, bool &option)
Definition: noiseModel.C:160
fileName baseFileDir() const
Return the base directory for output.
Definition: writeFile.C:43
#define WarningInFunction
Report a warning using Foam::Warning.
static autoPtr< windowModel > New(const dictionary &dict, const label nSamples)
Return a reference to the selected window model.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:607
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
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
scalar gainA(const scalar f) const
A weighting as gain in dB.
Definition: noiseModel.C:513
tmp< scalarField > SPL(const scalarField &Prms2, const scalar f) const
SPL [dB].
Definition: noiseModel.C:745
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
void writeWeightings() const
Helper function to check weightings.
Definition: noiseModel.C:861
tmp< scalarField > safeLog10(const scalarField &fld)
Definition: noiseModel.C:135
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Registry of regIOobjects.
dimensionedScalar log10(const dimensionedScalar &ds)
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:290
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
tmp< scalarField > RMSmeanPf(const scalarField &p) const
Return the multi-window RMS mean fft of the complete pressure data [Pa].
Definition: noiseModel.C:442
Namespace for OpenFOAM.
scalar RCf(const scalar f) const
C weighting function.
Definition: noiseModel.C:551
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133