Loading...
Searching...
No Matches
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-2023 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
36namespace Foam
37{
40}
41
42
45({
46 {weightingType::none, "dB"},
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
140namespace 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
160// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
161
162namespace Foam
163{
164
165// Get bool option (eg, read/write) and provide Info feedback
166static 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;
177
178} // End namespace Foam
179
180
181// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
182
184(
185 const scalarList& times
186) const
187{
188 scalar deltaT = 1;
189
190 if (times.size() > 1)
191 {
192 // Uniform time step
193 deltaT = (times.back() - times.front())/scalar(times.size() - 1);
194
195 bool nonUniform = false;
196 for (label i = 1; i < times.size() && !nonUniform; ++i)
197 {
198 const scalar dT = times[i] - times[i-1];
199
200 nonUniform = (mag(dT/deltaT - 1) > 1e-8);
201 }
202
203 if (nonUniform)
204 {
206 << "Detected non-uniform time step:"
207 << " resulting FFT may be incorrect"
208 << " (or the deltaT is just very small): " << deltaT
209 << endl;
210 }
211 }
212 else
213 {
215 << "Unable to create FFT with 0 or 1 values"
217 }
218
219 return deltaT;
220}
221
222
224{
225 forAll(p, i)
226 {
227 if ((p[i] < minPressure_) || (p[i] > maxPressure_))
228 {
230 << "Pressure data at position " << i
231 << " is outside of permitted bounds:" << nl
232 << " pressure: " << p[i] << nl
233 << " minimum pressure: " << minPressure_ << nl
234 << " maximum pressure: " << maxPressure_ << nl
235 << endl;
236
237 return false;
239 }
240
241 return true;
242}
243
244
245Foam::fileName Foam::noiseModel::baseFileDir(const label dataseti) const
246{
247 return
248 (
250 / type()
251 / ("input" + Foam::name(dataseti))
252 );
253}
254
255
257(
258 Ostream& os,
259 const string& x,
260 const string& y,
261 const UList<Tuple2<string, token>>& headerValues
262) const
263{
264 writeHeader(os, x + " vs " + y);
265 writeHeaderValue(os, "Lower frequency", fLower_);
266 writeHeaderValue(os, "Upper frequency", fUpper_);
267 writeHeaderValue(os, "Window model", windowModelPtr_->type());
268 writeHeaderValue(os, "Window number", windowModelPtr_->nWindow());
269 writeHeaderValue(os, "Window samples", windowModelPtr_->nSamples());
270 writeHeaderValue(os, "Window overlap %", windowModelPtr_->overlapPercent());
271 writeHeaderValue(os, "dBRef", dBRef_);
272
273 for (const auto& hv : headerValues)
274 {
275 writeHeaderValue(os, hv.first(), hv.second());
276 }
278 writeCommented(os, x.substr(0, x.find(' ')));
279 writeTabbed(os, y.substr(0, y.find(' ')));
280 os << nl;
281}
282
283
285(
286 Ostream& os,
287 const scalarField& f,
288 const scalarField& fx
289) const
290{
291 forAll(f, i)
292 {
293 if (f[i] >= fLower_ && f[i] <= fUpper_)
295 os << f[i] << tab << fx[i] << nl;
296 }
297 }
298}
299
300
301Foam::tmp<Foam::scalarField> Foam::noiseModel::uniformFrequencies
302(
303 const scalar deltaT,
304 const bool check
305) const
306{
307 const auto& window = windowModelPtr_();
308 const label N = window.nSamples();
309
310 auto tf = tmp<scalarField>::New(N/2 + 1, Zero);
311 auto& f = tf.ref();
312
313 const scalar deltaf = 1.0/(N*deltaT);
314
315 label nFreq = 0;
316 forAll(f, i)
317 {
318 f[i] = i*deltaf;
319
320 if (f[i] > fLower_ && f[i] < fUpper_)
321 {
322 ++nFreq;
323 }
324 }
325
326 if (check && nFreq == 0)
327 {
329 << "No frequencies found in range "
330 << fLower_ << " to " << fUpper_
332 }
333
334 return tf;
335}
336
337
339(
340 const scalarField& data,
341 const scalarField& f,
342 const labelUList& freqBandIDs
343) const
344{
345 if (freqBandIDs.size() < 2)
346 {
348 << "Octave frequency bands are not defined "
349 << "- skipping octaves calculation"
350 << endl;
351
352 return tmp<scalarField>::New();
353 }
354
355 auto toctData = tmp<scalarField>::New(freqBandIDs.size() - 1, Zero);
356 auto& octData = toctData.ref();
357
358 bitSet bandUsed(freqBandIDs.size() - 1);
359 for (label bandI = 0; bandI < freqBandIDs.size() - 1; ++bandI)
360 {
361 label fb0 = freqBandIDs[bandI];
362 label fb1 = freqBandIDs[bandI+1];
363
364 if (fb0 == fb1) continue;
365
366 for (label freqI = fb0; freqI < fb1; ++freqI)
367 {
368 label f0 = f[freqI];
369 label f1 = f[freqI + 1];
370 scalar dataAve = 0.5*(data[freqI] + data[freqI + 1]);
371 octData[bandI] += dataAve*(f1 - f0);
372
373 bandUsed.set(bandI);
374 }
375 }
376
377 bandUsed.flip();
378 labelList bandUnused = bandUsed.sortedToc();
379 if (bandUnused.size())
380 {
382 << "Empty bands found: " << bandUnused.size() << " of "
383 << bandUsed.size() << endl;
384 }
385
386 return toctData;
387}
388
389
391{
392 if (planInfo_.active)
393 {
394 if (p.size() != planInfo_.windowSize)
395 {
397 << "Expected pressure data to have " << planInfo_.windowSize
398 << " values, but received " << p.size() << " values"
399 << abort(FatalError);
400 }
401
402 List<double>& in = planInfo_.in;
403 const List<double>& out = planInfo_.out;
404 forAll(in, i)
405 {
406 in[i] = p[i];
407 }
408
409 ::fftw_execute(planInfo_.plan);
410
411 const label n = planInfo_.windowSize;
412 const label nBy2 = n/2;
413 auto tresult = tmp<scalarField>::New(nBy2 + 1);
414 auto& result = tresult.ref();
415
416 // 0 th value = DC component
417 // nBy2 th value = real only if n is even
418 result[0] = out[0];
419 for (label i = 1; i <= nBy2; ++i)
420 {
421 const auto re = out[i];
422 const auto im = out[n - i];
423 result[i] = sqrt(re*re + im*im);
424 }
425
426 return tresult;
427 }
428
429 return mag(fft::realTransform1D(p));
430}
431
432
434{
435 const auto& window = windowModelPtr_();
436 const label N = window.nSamples();
437 const label nWindow = window.nWindow();
438
439 auto tmeanPf = tmp<scalarField>::New(N/2 + 1, Zero);
440 auto& meanPf = tmeanPf.ref();
441
442 for (label windowI = 0; windowI < nWindow; ++windowI)
443 {
444 meanPf += Pf(window.apply<scalar>(p, windowI));
445 }
447 meanPf /= scalar(nWindow);
448
449 return tmeanPf;
450}
451
452
454(
455 const scalarField& p
456) const
457{
458 const auto& window = windowModelPtr_();
459 const label N = window.nSamples();
460 const label nWindow = window.nWindow();
461
462 scalarField RMSMeanPf(N/2 + 1, Zero);
463 for (label windowI = 0; windowI < nWindow; ++windowI)
464 {
465 RMSMeanPf += sqr(Pf(window.apply<scalar>(p, windowI)));
466 }
467
468 return sqrt(RMSMeanPf/scalar(nWindow))/scalar(N);
469}
470
471
473(
474 const scalarField& p,
475 const scalar deltaT
476) const
477{
478 const auto& window = windowModelPtr_();
479 const label N = window.nSamples();
480 const label nWindow = window.nWindow();
481
482 auto tpsd = tmp<scalarField>::New(N/2 + 1, Zero);
483 auto& psd = tpsd.ref();
484
485 for (label windowI = 0; windowI < nWindow; ++windowI)
486 {
487 psd += sqr(Pf(window.apply<scalar>(p, windowI)));
488 }
489
490 scalar fs = 1.0/deltaT;
491 psd /= scalar(nWindow)*fs*N;
492
493 // Scaling due to use of 1-sided FFT
494 // Note: do not scale DC component
495 psd *= 2;
496 psd.first() /= 2;
497 psd.last() /= 2;
498
499 if (debug)
500 {
501 Pout<< "Average PSD: " << average(psd) << endl;
502 }
503
504 return tpsd;
505}
506
507
508Foam::scalar Foam::noiseModel::RAf(const scalar f) const
509{
510 const scalar c1 = sqr(12194.0);
511 const scalar c2 = sqr(20.6);
512 const scalar c3 = sqr(107.7);
513 const scalar c4 = sqr(739.9);
514
515 const scalar f2 = f*f;
516
517 return
518 c1*sqr(f2)
519 /(
520 (f2 + c2)*sqrt((f2 + c3)*(f2 + c4))*(f2 + c1)
521 );
522}
523
524
525Foam::scalar Foam::noiseModel::gainA(const scalar f) const
526{
527 if (f < SMALL)
528 {
529 return 0;
530 }
531
532 return 20*log10(RAf(f)) - 20*log10(RAf(1000));
533}
534
535
536Foam::scalar Foam::noiseModel::RBf(const scalar f) const
537{
538 const scalar c1 = sqr(12194.0);
539 const scalar c2 = sqr(20.6);
540 const scalar c3 = sqr(158.5);
541
542 const scalar f2 = f*f;
543
544 return
545 c1*f2*f
546 /(
547 (f2 + c2)*sqrt(f2 + c3)*(f2 + c1)
548 );
549}
550
551
552Foam::scalar Foam::noiseModel::gainB(const scalar f) const
553{
554 if (f < SMALL)
555 {
556 return 0;
557 }
558
559 return 20*log10(RBf(f)) - 20*log10(RBf(1000));
560}
561
562
563Foam::scalar Foam::noiseModel::RCf(const scalar f) const
564{
565 const scalar c1 = sqr(12194.0);
566 const scalar c2 = sqr(20.6);
568 const scalar f2 = f*f;
569
570 return c1*f2/((f2 + c2)*(f2 + c1));
571}
572
573
574Foam::scalar Foam::noiseModel::gainC(const scalar f) const
575{
576 if (f < SMALL)
577 {
578 return 0;
579 }
580
581 return 20*log10(RCf(f)) - 20*log10(RCf(1000));
582}
583
584
585Foam::scalar Foam::noiseModel::RDf(const scalar f) const
586{
587 const scalar f2 = f*f;
588
589 const scalar hf =
590 (sqr(1037918.48 - f2) + 1080768.16*f2)
591 /(sqr(9837328 - f2) + 11723776*f2);
592
593 return
594 f/6.8966888496476e-5*Foam::sqrt(hf/((f2 + 79919.29)*(f2 + 1345600)));
595}
596
597
598Foam::scalar Foam::noiseModel::gainD(const scalar f) const
599{
600 if (f < SMALL)
601 {
602 return 0;
603 }
605 return 20*log10(RDf(f));
606}
607
608
609// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
610
612(
613 const dictionary& dict,
614 const objectRegistry& obr,
615 const word& name,
616 const bool readFields
617)
618:
619 functionObjects::writeFile(obr, "noise"),
620 dict_(dict),
621 rhoRef_(1),
622 nSamples_(65536),
623 fLower_(25),
624 fUpper_(10000),
625 sampleFreq_(0),
626 startTime_(0),
627 windowModelPtr_(),
628 SPLweighting_(weightingType::none),
629 dBRef_(2e-5),
630 minPressure_(-0.5*VGREAT),
631 maxPressure_(0.5*VGREAT),
632 outputPrefix_(),
633 writePrmsf_(true),
634 writeSPL_(true),
635 writePSD_(true),
636 writePSDf_(true),
637 writeOctaves_(true),
638 planInfo_()
639{
640 planInfo_.active = false;
641
642 if (readFields)
643 {
644 read(dict);
645 }
646
647 if (debug)
648 {
650 }
651}
652
653
654// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
655
657{
659 {
660 return false;
661 }
662
663 // Re-assign defaults (to avoid stickiness)
664 fLower_ = 25;
665 fUpper_ = 10000;
666 sampleFreq_ = 0;
667
668 dict.readIfPresent("rhoRef", rhoRef_);
669 dict.readIfPresent("N", nSamples_);
670 dict.readIfPresentCompat("minFreq", {{"fl", 2312}}, fLower_);
671 dict.readIfPresentCompat("maxFreq", {{"fu", 2312}}, fUpper_);
672 dict.readIfPresent("sampleFreq", sampleFreq_);
673 dict.readIfPresent("startTime", startTime_);
674 dict.readIfPresent("minPressure", minPressure_);
675 dict.readIfPresent("maxPressure", maxPressure_);
676 dict.readIfPresent("outputPrefix", outputPrefix_);
677
678 if (fLower_ < 0)
679 {
681 << "Lower frequency bound must be greater than zero"
682 << exit(FatalIOError);
683 }
684
685 if (fUpper_ < 0)
686 {
688 << "Upper frequency bound must be greater than zero"
689 << exit(FatalIOError);
690 }
691
692 if (fUpper_ < fLower_)
693 {
695 << "Upper frequency bound (" << fUpper_
696 << ") must be greater than lower frequency bound ("
697 << fLower_ << ")"
698 << exit(FatalIOError);
699 }
700
701 Info<< " Frequency bounds:" << nl
702 << " lower: " << fLower_ << nl
703 << " upper: " << fUpper_ << nl
704 << " sample: ";
705
706 if (sampleFreq_ > 0)
707 {
708 Info<< sampleFreq_ << nl;
709 }
710 else
711 {
712 Info<< "auto" << nl;
713 }
714
715 weightingTypeNames_.readIfPresent("SPLweighting", dict, SPLweighting_);
716
717 Info<< " Weighting: " << weightingTypeNames_[SPLweighting_] << endl;
718
719 if (dict.readIfPresent("dBRef", dBRef_))
720 {
721 Info<< " Reference for dB calculation: " << dBRef_ << endl;
722 }
723
724 Info<< " Write options:" << endl;
725 dictionary optDict(dict.subOrEmptyDict("writeOptions"));
726 readWriteOption(optDict, "writePrmsf", writePrmsf_);
727 readWriteOption(optDict, "writeSPL", writeSPL_);
728 readWriteOption(optDict, "writePSD", writePSD_);
729 readWriteOption(optDict, "writePSDf", writePSDf_);
730 readWriteOption(optDict, "writeOctaves", writeOctaves_);
731
732
733 windowModelPtr_ = windowModel::New(dict, nSamples_);
734
735 cleanFFTW();
736
737 const label windowSize = windowModelPtr_->nSamples();
738
739 if (windowSize > 1)
740 {
741 planInfo_.active = true;
742 planInfo_.windowSize = windowSize;
743 planInfo_.in.setSize(windowSize);
744 planInfo_.out.setSize(windowSize);
745
746 // Using real to half-complex fftw 'kind'
747 planInfo_.plan =
748 fftw_plan_r2r_1d
749 (
750 windowSize,
751 planInfo_.in.data(),
752 planInfo_.out.data(),
753 FFTW_R2HC,
754 windowSize <= 8192 ? FFTW_MEASURE : FFTW_ESTIMATE
755 );
756 }
758 Info<< endl;
759
760 return true;
761}
762
763
765(
767) const
768{
769 return 10*safeLog10(PSDf/sqr(dBRef_));
770}
771
772
774(
775 const scalarField& Prms2,
776 const scalar f
777) const
778{
779 tmp<scalarField> tspl(10*safeLog10(Prms2/sqr(dBRef_)));
780 scalarField& spl = tspl.ref();
781
782 switch (SPLweighting_)
783 {
784 case weightingType::none:
785 {
786 break;
787 }
788 case weightingType::dBA:
789 {
790 spl += gainA(f);
791 break;
792 }
793 case weightingType::dBB:
794 {
795 spl += gainB(f);
796 break;
797 }
798 case weightingType::dBC:
799 {
800 spl += gainC(f);
801 break;
802 }
803 case weightingType::dBD:
804 {
805 spl += gainD(f);
806 break;
807 }
808 default:
809 {
811 << "Unknown weighting " << weightingTypeNames_[SPLweighting_]
812 << abort(FatalError);
814 }
815
816 return tspl;
817}
818
819
820Foam::tmp<Foam::scalarField> Foam::noiseModel::SPL
821(
822 const scalarField& Prms2,
823 const scalarField& f
824) const
825{
826 tmp<scalarField> tspl(10*safeLog10(Prms2/sqr(dBRef_)));
827 scalarField& spl = tspl.ref();
828
829 switch (SPLweighting_)
830 {
831 case weightingType::none:
832 {
833 break;
834 }
835 case weightingType::dBA:
836 {
837 forAll(spl, i)
838 {
839 spl[i] += gainA(f[i]);
840 }
841 break;
842 }
843 case weightingType::dBB:
844 {
845 forAll(spl, i)
846 {
847 spl[i] += gainB(f[i]);
848 }
849 break;
850 }
851 case weightingType::dBC:
852 {
853 forAll(spl, i)
854 {
855 spl[i] += gainC(f[i]);
856 }
857 break;
858 }
859 case weightingType::dBD:
860 {
861 forAll(spl, i)
862 {
863 spl[i] += gainD(f[i]);
864 }
865 break;
866 }
867 default:
868 {
870 << "Unknown weighting " << weightingTypeNames_[SPLweighting_]
871 << abort(FatalError);
873 }
874
875 return tspl;
876}
877
878
880{
881 if (planInfo_.active)
882 {
883 planInfo_.active = false;
884 fftw_destroy_plan(planInfo_.plan);
885 fftw_cleanup();
886 }
887}
888
889
891{
892 scalar f0 = 10;
893 scalar f1 = 20000;
894
895 OFstream osA("noiseModel-weight-A");
896 OFstream osB("noiseModel-weight-B");
897 OFstream osC("noiseModel-weight-C");
898 OFstream osD("noiseModel-weight-D");
899
900 for (label f = f0; f <= f1; ++f)
901 {
902 osA << f << " " << gainA(f) << nl;
903 osB << f << " " << gainB(f) << nl;
904 osC << f << " " << gainC(f) << nl;
905 osD << f << " " << gainD(f) << nl;
906 }
907}
908
909
910// ************************************************************************* //
scalar y
label n
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void pop_back(label n=1)
Reduce size by 1 or more elements. Can be called on an empty list.
void append(const T &val)
Copy append an element to the end of this list.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition HashTable.C:156
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
bool set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition List.H:469
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
label size() const noexcept
Number of entries.
Definition PackedList.H:392
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition Tuple2.H:51
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
T & back()
Access last element of the list, position [size()-1].
Definition UListI.H:253
T & front()
Access first element of the list, position [0].
Definition UListI.H:239
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
void flip()
Invert all bits in the addressable region.
Definition bitSetI.H:546
labelList sortedToc() const
The indices of the on bits as a sorted labelList.
Definition bitSetI.H:441
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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:521
bool readIfPresentCompat(const word &keyword, std::initializer_list< std::pair< const char *, int > > compat, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val using any compatibility names if needed....
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...
static tmp< complexField > realTransform1D(const scalarField &field)
Transform real-value data.
Definition fft.C:113
A class for handling file names.
Definition fileName.H:75
fileName baseFileDir() const
Return the base directory for output.
Definition writeFile.C:43
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
virtual bool read(const dictionary &dict)
Read.
Definition writeFile.C:240
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
Base class for noise models.
Definition noiseModel.H:178
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:295
scalar RAf(const scalar f) const
A weighting function.
Definition noiseModel.C:501
scalar rhoRef_
Reference density (to convert from kinematic to static pressure).
Definition noiseModel.H:233
label nSamples_
Number of samples in sampling window, default = 2^16.
Definition noiseModel.H:238
tmp< scalarField > octaves(const scalarField &data, const scalarField &f, const labelUList &freqBandIDs) const
Generate octave data.
Definition noiseModel.C:332
tmp< scalarField > RMSmeanPf(const scalarField &p) const
Return the multi-window RMS mean fft of the complete pressure data [Pa].
Definition noiseModel.C:447
scalar fUpper_
Upper frequency limit, default = 10kHz.
Definition noiseModel.H:248
scalar gainC(const scalar f) const
C weighting as gain in dB.
Definition noiseModel.C:567
tmp< Foam::scalarField > PSD(const scalarField &PSDf) const
PSD [dB/Hz].
Definition noiseModel.C:758
tmp< scalarField > meanPf(const scalarField &p) const
Return the multi-window mean fft of the complete pressure data [Pa].
Definition noiseModel.C:426
scalar maxPressure_
Min pressure value.
Definition noiseModel.H:286
bool writePSD_
Write PSD; default = yes.
Definition noiseModel.H:309
scalar gainA(const scalar f) const
A weighting as gain in dB.
Definition noiseModel.C:518
scalar minPressure_
Min pressure value.
Definition noiseModel.H:281
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
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition noiseModel.C:649
bool writeOctaves_
Write writeOctaves; default = yes.
Definition noiseModel.H:319
fileName outputPrefix_
Output file prefix, default = ''.
Definition noiseModel.H:294
scalar RDf(const scalar f) const
D weighting function.
Definition noiseModel.C:578
tmp< scalarField > Pf(const scalarField &p) const
Return the fft of the given pressure data.
Definition noiseModel.C:383
scalar RCf(const scalar f) const
C weighting function.
Definition noiseModel.C:556
noiseModel(const noiseModel &)=delete
No copy construct.
tmp< scalarField > SPL(const scalarField &Prms2, const scalar f) const
SPL [dB].
Definition noiseModel.C:767
planInfo planInfo_
Plan information for FFTW.
Definition noiseModel.H:327
autoPtr< windowModel > windowModelPtr_
Window model.
Definition noiseModel.H:263
static const Enum< weightingType > weightingTypeNames_
Definition noiseModel.H:190
bool writeSPL_
Write SPL; default = yes.
Definition noiseModel.H:304
const dictionary dict_
Copy of dictionary used for construction.
Definition noiseModel.H:228
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:250
bool writePrmsf_
Write Prmsf; default = yes.
Definition noiseModel.H:299
weightingType SPLweighting_
Weighting.
Definition noiseModel.H:268
scalar checkUniformTimeStep(const scalarList &times) const
Check and return uniform time step.
Definition noiseModel.C:177
bool writePSDf_
Write PSDf; default = yes.
Definition noiseModel.H:314
scalar gainD(const scalar f) const
D weighting as gain in dB.
Definition noiseModel.C:591
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:466
void writeWeightings() const
Helper function to check weightings.
Definition noiseModel.C:883
scalar startTime_
Start time, default = 0s.
Definition noiseModel.H:258
scalar sampleFreq_
Prescribed sample frequency.
Definition noiseModel.H:253
void cleanFFTW()
Clean up the FFTW.
Definition noiseModel.C:872
scalar dBRef_
Reference for dB calculation, default = 2e-5.
Definition noiseModel.H:273
void writeFreqDataToFile(Ostream &os, const scalarField &f, const scalarField &fx) const
Definition noiseModel.C:278
bool validateBounds(const scalarList &p) const
Return true if all pressure data is within min/max bounds.
Definition noiseModel.C:216
scalar fLower_
Lower frequency limit, default = 25Hz.
Definition noiseModel.H:243
scalar gainB(const scalar f) const
B weighting as gain in dB.
Definition noiseModel.C:545
scalar RBf(const scalar f) const
B weighting function.
Definition noiseModel.C:529
Registry of regIOobjects.
Lookup type of boundary radiation properties.
Definition lookup.H:60
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
static autoPtr< windowModel > New(const dictionary &dict, const label nSamples)
Return a reference to the selected window model.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar re
Classical electron radius: default SI units: [m].
Namespace for handling debugging switches.
Definition debug.C:45
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
tmp< scalarField > safeLog10(const scalarField &fld)
Definition noiseModel.C:135
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
static void writeHeader(Ostream &os, const word &fieldName)
messageStream Info
Information stream (stdout output on master, null elsewhere).
static void check(const int retVal, const char *what)
dimensionedScalar log10(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
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:519
static void readWriteOption(const dictionary &dict, const word &lookup, bool &option)
Definition noiseModel.C:160
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
errorManip< error > abort(error &err)
Definition errorManip.H:139
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
labelList f(nPoints)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const Vector< label > N(dict.get< Vector< label > >("N"))