Loading...
Searching...
No Matches
surfaceNoise.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-2025 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 "surfaceNoise.H"
29#include "surfaceReader.H"
30#include "surfaceWriter.H"
31#include "globalIndex.H"
32#include "argList.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
39namespace noiseModels
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
46
47// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
48
49void surfaceNoise::initialise(const fileName& fName)
50{
51 Info<< "Reading data file: "
52 << fileObr_.time().relativePath(fName) << endl;
53
54 instantList allTimes;
55 label nAvailableTimes = 0;
56
57 // All reading performed on the master processor only
58 if (Pstream::master())
59 {
60 // Create the surface reader
62
63 // Find the index of the pressure data
64 const wordList fieldNames(readerPtr_->fieldNames(0));
65 pIndex_ = fieldNames.find(pName_);
66 if (pIndex_ == -1)
67 {
69 << "Unable to find pressure field name " << pName_
70 << " in list of available fields: " << fieldNames
71 << exit(FatalError);
72 }
73
74 // Create the surface writer
75 // - Could be done later, but since this utility can process a lot of
76 // data we can ensure that the user-input is correct prior to doing
77 // the heavy lifting
78
79 // Set the time range
80 allTimes = readerPtr_->times();
82
83 // Determine the windowing
84 nAvailableTimes = allTimes.size() - startTimeIndex_;
85 }
86
88 (
90 pIndex_,
92 nAvailableTimes
93 );
94
95
96 // Note: all processors should call the windowing validate function
97 label nRequiredTimes = windowModelPtr_->validate(nAvailableTimes);
98
99 if (UPstream::master())
100 {
101 // Restrict times
102 times_.resize_nocopy(nRequiredTimes);
103 forAll(times_, timei)
104 {
105 times_[timei] = allTimes[timei + startTimeIndex_].value();
106 }
107
108 deltaT_ =
109 (
110 sampleFreq_ > 0
111 ? (1.0/sampleFreq_)
113 );
114
115 // Read the surface geometry
116 // Note: hard-coded to read mesh from first time index
117 const meshedSurface& surf = readerPtr_->geometry(0);
118
119 nFaces_ = surf.nFaces();
120 }
121
123 (
126 deltaT_,
127 nFaces_
128 );
129}
130
131
133(
134 const globalIndex& procFaceAddr,
135 List<scalarField>& pData
136)
137{
138 // Data is stored as pressure on surface at a given time. Now we need to
139 // pivot the data so that we have the complete pressure time history per
140 // surface face. In serial mode, this results in all pressure data being
141 // loaded into memory (!)
142
143 const label nLocalFace = procFaceAddr.localSize();
144
145 // Complete pressure time history data for subset of faces
146 pData.resize_nocopy(nLocalFace);
147 const label nTimes = times_.size();
148 for (scalarField& pf : pData)
149 {
150 pf.resize_nocopy(nTimes);
151 }
152
153 Info<< "Reading pressure data" << endl;
154
155 // Master only
156 scalarField allData;
157
158 if (Pstream::parRun())
159 {
160 // Procedure:
161 // 1. Master processor reads pressure data for all faces for all times
162 // 2. Master sends each processor a subset of faces
163 // 3. Local processor reconstructs the full time history for the subset
164 // of faces
165 // Note: reading all data on master to avoid potential NFS problems...
166
167 scalarField scratch;
168
169 if (!useBroadcast_)
170 {
171 scratch.resize(nLocalFace);
172 }
173
174 // Read data and send to sub-ranks
175 forAll(times_, timei)
176 {
177 const label fileTimeIndex = timei + startTimeIndex_;
178
179 if (Pstream::master())
180 {
181 Info<< " time: " << times_[timei] << endl;
182
183 // Read pressure at all faces for time timeI
184 allData = readerPtr_->field(fileTimeIndex, pIndex_, scalar(0));
185 }
186
187 if (useBroadcast_)
188 {
189 Pstream::broadcast(allData);
190 }
191 else
192 {
193 procFaceAddr.scatter
194 (
195 allData,
196 scratch,
198 commType_,
200 );
201 }
202
203 scalarField::subField procData =
204 (
206 ? allData.slice(procFaceAddr.range())
207 : scratch.slice(0, nLocalFace)
208 );
209
210 // Apply conversions
211 procData *= rhoRef_;
212
213 // Transcribe this time snapshot (transpose)
214 forAll(procData, facei)
215 {
216 pData[facei][timei] = procData[facei];
217 }
218 }
219 }
220 else
221 {
222 // Read data - no sub-ranks
223 forAll(times_, timei)
224 {
225 const label fileTimeIndex = timei + startTimeIndex_;
226
227 Info<< " time: " << times_[timei] << endl;
228
229 allData = readerPtr_->field(fileTimeIndex, pIndex_, scalar(0));
230
231 auto& procData = allData;
232
233 // Apply conversions
234 procData *= rhoRef_;
235
236 // Transcribe this time snapshot (transpose)
237 forAll(procData, facei)
238 {
239 pData[facei][timei] = procData[facei];
240 }
241 }
242 }
243
244 forAll(pData, facei)
245 {
246 pData[facei] -= average(pData[facei]);
247 }
248
249
250 Info<< "Read "
251 << returnReduce(pData.size(), sumOp<label>())
252 << " pressure traces with "
253 << times_.size()
254 << " time values" << nl << endl;
255}
256
257
259(
260 const scalarField& data,
261 const globalIndex& procFaceAddr
262) const
263{
264 if (!nFaces_)
265 {
266 // Already reduced, can use as sanity check
267 return 0;
268 }
269
270 scalar areaAverage = 0;
271
272 if (areaAverage_)
273 {
274 if (Pstream::parRun())
275 {
276 // Collect the surface data so that we can output the surfaces
277 scalarField allData = procFaceAddr.gather
278 (
279 data,
281 commType_,
283 );
284
285 if (Pstream::master())
286 {
287 // Note: hard-coded to read mesh from first time index
288 const meshedSurface& surf = readerPtr_->geometry(0);
289
290 areaAverage = sum(allData*surf.magSf())/sum(surf.magSf());
291 }
292 }
293 else
294 {
295 // Note: hard-coded to read mesh from first time index
296 const meshedSurface& surf = readerPtr_->geometry(0);
297
298 areaAverage = sum(data*surf.magSf())/sum(surf.magSf());
299 }
300
301 Pstream::broadcast(areaAverage);
302 }
303 else
304 {
305 // Ensemble averaged
306 // - same as gAverage, but already know number of faces
307
308 areaAverage = sum(data);
309 reduce(areaAverage, sumOp<scalar>());
310
311 areaAverage /= (scalar(nFaces_) + ROOTVSMALL);
312 }
313
314 return areaAverage;
315}
316
317
319(
320 const fileName& outDirBase,
321 const word& fName,
322 const word& title,
323 const scalar freq,
324 const scalarField& data,
325 const globalIndex& procFaceAddr,
326 const bool writeSurface
327) const
328{
329 Info<< " processing " << title << " for frequency " << freq << endl;
330
331 const instant freqInst(freq, Foam::name(freq));
332
333 if (!writeSurface)
334 {
335 return surfaceAverage(data, procFaceAddr);
336 }
337
338 scalar areaAverage = 0;
339
340 if (Pstream::parRun())
341 {
342 // Collect the surface data so that we can output the surfaces
343 scalarField allData = procFaceAddr.gather
344 (
345 data,
347 commType_,
349 );
350
351 if (Pstream::master())
352 {
353 // Note: hard-coded to read mesh from first time index
354 const meshedSurface& surf = readerPtr_->geometry(0);
355
356 if (areaAverage_)
357 {
358 areaAverage = sum(allData*surf.magSf())/sum(surf.magSf());
359 }
360 else
361 {
362 areaAverage = sum(allData)/(allData.size() + ROOTVSMALL);
363 }
364
365 // (writeSurface == true)
366 {
367 // Time-aware, with time spliced into the output path
368 writerPtr_->beginTime(freqInst);
369
370 writerPtr_->open
371 (
372 surf.points(),
373 surf.surfFaces(),
374 (outDirBase / fName),
375 false // serial - already merged
376 );
377
378 writerPtr_->nFields(1); // Legacy VTK
379 writerPtr_->write(title, allData);
380
381 writerPtr_->endTime();
382 writerPtr_->clear();
383 }
384 }
385 }
386 else
387 {
388 // Note: hard-coded to read mesh from first time index
389 const meshedSurface& surf = readerPtr_->geometry(0);
390
391 if (areaAverage_)
392 {
393 areaAverage = sum(data*surf.magSf())/sum(surf.magSf());
394 }
395 else
396 {
397 areaAverage = sum(data)/(data.size() + ROOTVSMALL);
398 }
399
400 // (writeSurface == true)
401 {
402 // Time-aware, with time spliced into the output path
403 writerPtr_->beginTime(freqInst);
404
405 writerPtr_->open
406 (
407 surf.points(),
408 surf.surfFaces(),
409 (outDirBase / fName),
410 false // serial - already merged
411 );
412
413 writerPtr_->nFields(1); // Legacy VTK
414 writerPtr_->write(title, data);
415
416 writerPtr_->endTime();
417 writerPtr_->clear();
418 }
419 }
420
421 Pstream::broadcast(areaAverage);
422 return areaAverage;
423}
424
425
426// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
427
429(
430 const dictionary& dict,
431 const objectRegistry& obr,
432 const word& name,
433 const bool readFields
434)
435:
436 noiseModel(dict, obr, name, false),
437 inputFileNames_(),
438 pName_("p"),
439 pIndex_(0),
440 times_(),
441 deltaT_(0),
442 startTimeIndex_(0),
443 nFaces_(0),
444 fftWriteInterval_(1),
445 areaAverage_(false),
446 useBroadcast_(false),
447 commType_(UPstream::commsTypes::scheduled),
448 readerType_(),
449 readerPtr_(nullptr),
450 writerPtr_(nullptr)
451{
452 if (readFields)
453 {
455 }
456}
457
458
459// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
460
462{
464 {
465 if (!dict.readIfPresent("files", inputFileNames_))
466 {
468 dict.readEntry("file", inputFileNames_.first());
469 }
470
471 dict.readIfPresent("p", pName_);
472 dict.readIfPresent("fftWriteInterval", fftWriteInterval_);
473
474 Info<< this->type() << nl
475 << " Pressure field name: " << pName_ << nl
476 << " FFT write interval: " << fftWriteInterval_ << nl;
477
478 dict.readIfPresent("areaAverage", areaAverage_);
479
480 if (areaAverage_)
481 {
482 Info<< " Averaging: area weighted" << endl;
483 }
484 else
485 {
486 Info<< " Averaging: ensemble" << endl;
487 }
488
489 useBroadcast_ = false;
490 dict.readIfPresent("broadcast", useBroadcast_);
492
493 if (Pstream::parRun())
494 {
495 Info<< " Distribute fields: "
497
498 if (useBroadcast_)
499 {
500 Info<< " (broadcast all)";
501 }
502 Info<< endl;
503 }
504
505 readerType_ = dict.get<word>("reader");
506
507 // Surface writer (keywords: writer, writeOptions)
508
509 const word writerType = dict.get<word>("writer");
510
512 (
513 writerType,
514 surfaceWriter::formatOptions(dict, writerType, "writeOptions")
515 );
516
517 // Use outputDir/TIME/surface-name
518 writerPtr_->useTimeDir(true);
519
520 Info<< endl;
521
522 return true;
523 }
524
525 return false;
526}
527
528
530{
531 forAll(inputFileNames_, filei)
532 {
533 fileName fName = inputFileNames_[filei];
534 fName.expand();
535
536 if (!fName.isAbsolute())
537 {
538 fName = argList::envGlobalPath()/fName;
539 }
540
541 initialise(fName);
542
543 // Processor face addressing
544 globalIndex procFaceAddr;
545
546 if (Pstream::parRun())
547 {
548 // Calculate face/proc offsets manually
549 labelList procFaceOffsets(Pstream::nProcs() + 1);
550 const label nFacePerProc = floor(nFaces_/Pstream::nProcs()) + 1;
551
552 procFaceOffsets[0] = 0;
553 for (label proci = 1; proci < procFaceOffsets.size(); ++proci)
554 {
555 procFaceOffsets[proci] = min(proci*nFacePerProc, nFaces_);
556 }
557
558 procFaceAddr.offsets() = std::move(procFaceOffsets);
559
560 // Don't need to broadcast. Already identical on all ranks
561 }
562 else
563 {
564 // Local data only
565 procFaceAddr.reset(globalIndex::gatherNone{}, nFaces_);
566 }
567
568 // Pressure time history data per face
569 List<scalarField> pData;
570
571 // Read pressure data from file
572 readSurfaceData(procFaceAddr, pData);
573
574 // Process the pressure data, and store results as surface values per
575 // frequency so that it can be output using the surface writer
576
577 Info<< "Creating noise FFTs" << endl;
578
579 const scalarField freq1(uniformFrequencies(deltaT_, true));
580
581 const scalar maxFreq1 = max(freq1);
582
583 // Reset desired frequency range if outside actual frequency range
584 fLower_ = min(fLower_, maxFreq1);
585 fUpper_ = min(fUpper_, maxFreq1);
586
587 // Storage for FFT data
588 const label nLocalFace = pData.size();
589 const label nFFT = ceil(freq1.size()/scalar(fftWriteInterval_));
590
591 List<scalarField> surfPrmsf(nFFT);
592 List<scalarField> surfPSDf(nFFT);
593 forAll(surfPrmsf, freqI)
594 {
595 surfPrmsf[freqI].setSize(nLocalFace);
596 surfPSDf[freqI].setSize(nLocalFace);
597 }
598
599 // Storage for 1/3 octave data
600 labelList octave13BandIDs;
601 scalarField octave13FreqCentre;
603 (
604 freq1,
605 fLower_,
606 fUpper_,
607 3,
608 octave13BandIDs,
609 octave13FreqCentre
610 );
611
612 label bandSize = 0;
613 if (octave13BandIDs.empty())
614 {
616 << "Octave band calculation failed (zero sized). "
617 << "Please check your input data"
618 << endl;
619 }
620 else
621 {
622 bandSize = octave13BandIDs.size() - 1;
623 }
624
625 List<scalarField> surfPrms13f(bandSize);
626 forAll(surfPrms13f, freqI)
627 {
628 surfPrms13f[freqI].setSize(nLocalFace);
629 }
630
631 const windowModel& win = windowModelPtr_();
632
633 {
634 forAll(pData, faceI)
635 {
636 const scalarField& p = pData[faceI];
637
638 // Generate the FFT-based data
639 const scalarField Prmsf(RMSmeanPf(p));
640 const scalarField PSDf(this->PSDf(p, deltaT_));
641
642 // Store the frequency results in slot for face of surface
643 forAll(surfPrmsf, i)
644 {
645 label freqI = i*fftWriteInterval_;
646 surfPrmsf[i][faceI] = Prmsf[freqI];
647 surfPSDf[i][faceI] = PSDf[freqI];
648 }
649
650 // Integrated PSD = P(rms)^2 [Pa^2]
651 const scalarField Prms13f
652 (
653 octaves
654 (
655 PSDf,
656 freq1,
657 octave13BandIDs
658 )
659 );
660
661 // Store the 1/3 octave results in slot for face of surface
662 forAll(surfPrms13f, freqI)
663 {
664 surfPrms13f[freqI][faceI] = Prms13f[freqI];
665 }
666 }
667 }
668
669 const word fNameBase = fName.stem();
670
671 // Output directory
672 fileName outDirBase(baseFileDir(filei)/fNameBase);
673
674 const scalar deltaf = 1.0/(deltaT_*win.nSamples());
675 Info<< "Writing fft surface data";
676 if (fftWriteInterval_ == 1)
677 {
678 Info<< endl;
679 }
680 else
681 {
682 Info<< " at every " << fftWriteInterval_ << " frequency points"
683 << endl;
684 }
685
686 // Common output information
687 // Note: hard-coded to read mesh from first time index
688 scalar surfArea = 0;
689 label surfSize = 0;
690 if (Pstream::master())
691 {
692 const meshedSurface& surf = readerPtr_->geometry(0);
693 surfArea = sum(surf.magSf());
694 surfSize = surf.size();
695 }
697 (
699 surfArea,
700 surfSize
701 );
702
703 List<Tuple2<string, token>> commonInfo
704 ({
705 {"Area average", token(word(Switch::name(areaAverage_)))},
706 {"Area sum", token(surfArea)},
707 {"Number of faces", token(surfSize)}
708 });
709
710 {
711 fileName outDir(outDirBase/"fft");
712 fileName outSurfDir(filePath(outDir));
713
714 // Determine frequency range of interest
715 // Note: frequencies have fixed interval, and are in the range
716 // 0 to fftWriteInterval_*(n-1)*deltaf
717 label f0 = ceil(fLower_/deltaf/scalar(fftWriteInterval_));
718 label f1 = floor(fUpper_/deltaf/scalar(fftWriteInterval_));
719 label nFreq = f1 - f0;
720
721 scalarField PrmsfAve(nFreq, Zero);
722 scalarField PSDfAve(nFreq, Zero);
723 scalarField fOut(nFreq, Zero);
724
725 if (nFreq == 0)
726 {
728 << "No surface data available using a fftWriteInterval of "
730 }
731 else
732 {
733 forAll(fOut, i)
734 {
735 label freqI = (i + f0)*fftWriteInterval_;
736 fOut[i] = freq1[freqI];
737
738
739 PrmsfAve[i] = writeSurfaceData
740 (
741 outSurfDir,
742 fNameBase,
743 "Prmsf",
744 freq1[freqI],
745 surfPrmsf[i + f0],
746 procFaceAddr,
748 );
749
750 PSDfAve[i] = writeSurfaceData
751 (
752 outSurfDir,
753 fNameBase,
754 "PSDf",
755 freq1[freqI],
756 surfPSDf[i + f0],
757 procFaceAddr,
759 );
761 (
762 outSurfDir,
763 fNameBase,
764 "PSD",
765 freq1[freqI],
766 PSD(surfPSDf[i + f0]),
767 procFaceAddr,
769 );
771 (
772 outSurfDir,
773 fNameBase,
774 "SPL",
775 freq1[freqI],
776 SPL(surfPSDf[i + f0]*deltaf, freq1[freqI]),
777 procFaceAddr,
779 );
780 }
781 }
782
783 if (Pstream::master())
784 {
785 {
786 auto filePtr = newFile(outDir/"Average_Prms_f");
787 auto& os = filePtr();
788
789 Info<< " Writing " << os.relativeName() << endl;
790
791 writeFileHeader(os, "f [Hz]", "P(f) [Pa]", commonInfo);
792 writeFreqDataToFile(os, fOut, PrmsfAve);
793 }
794 {
795 auto filePtr = newFile(outDir/"Average_PSD_f_f");
796 auto& os = filePtr();
797
798 Info<< " Writing " << os.relativeName() << endl;
799
801 (
802 os,
803 "f [Hz]",
804 "PSD(f) [PaPa_Hz]",
805 commonInfo
806 );
807 writeFreqDataToFile(os, fOut, PSDfAve);
808 }
809 {
810 auto filePtr = newFile(outDir/"Average_PSD_dB_Hz_f");
811 auto& os = filePtr();
812
813 Info<< " Writing " << os.relativeName() << endl;
814
816 (
817 os,
818 "f [Hz]",
819 "PSD(f) [dB_Hz]",
820 commonInfo
821 );
822 writeFreqDataToFile(os, fOut, PSD(PSDfAve));
823 }
824 {
825 auto filePtr = newFile(outDir/"Average_SPL_dB_f");
826 auto& os = filePtr();
827
828 Info<< " Writing " << os.relativeName() << endl;
829
831 (
832 os,
833 "f [Hz]",
834 "SPL(f) [dB]",
835 commonInfo
836 );
837 writeFreqDataToFile(os, fOut, SPL(PSDfAve*deltaf, fOut));
838 }
839 }
840 }
841
842
843 Info<< "Writing one-third octave surface data" << endl;
844 {
845 fileName outDir(outDirBase/"oneThirdOctave");
846 fileName outSurfDir(filePath(outDir));
847
848 scalarField PSDfAve(surfPrms13f.size(), Zero);
849 scalarField Prms13fAve(surfPrms13f.size(), Zero);
850
851 forAll(surfPrms13f, i)
852 {
854 (
855 outSurfDir,
856 fNameBase,
857 "SPL13",
858 octave13FreqCentre[i],
859 SPL(surfPrms13f[i], octave13FreqCentre[i]),
860 procFaceAddr,
862 );
863
864 Prms13fAve[i] =
865 surfaceAverage(surfPrms13f[i], procFaceAddr);
866 }
867
868 if (Pstream::master())
869 {
870 auto filePtr = newFile(outDir/"Average_SPL13_dB_fm");
871 auto& os = filePtr();
872
873 Info<< " Writing " << os.relativeName() << endl;
874
876 (
877 os,
878 "fm [Hz]",
879 "SPL(fm) [dB]",
880 commonInfo
881 );
883 (
884 os,
885 octave13FreqCentre,
886 SPL(Prms13fAve, octave13FreqCentre)
887 );
888 }
889 }
890 }
891}
892
893
894// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
895
896} // End namespace noiseModels
897} // End namespace Foam
898
899// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
bool readIfPresent(const word &key, const dictionary &dict, EnumType &val, const bool warnOnly=false) const
Find an entry if present, and assign to T val.
Definition EnumI.H:111
SubField< scalar > subField
Definition Field.H:183
SubField< Type > slice(const label pos, label len=-1)
Return SubField slice (non-const access) - no range checking.
Definition SubField.H:230
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 resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
void setSize(label n)
Alias for resize().
Definition List.H:536
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
const scalarField & magSf() const
Face area magnitudes.
label size() const
The surface size is the number of faces.
const List< Face > & surfFaces() const
Return const access to the faces.
label nFaces() const noexcept
Number of faces in the patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
static void broadcasts(const int communicator, Type &value, Args &&... values)
Broadcast multiple items to all communicator ranks. Does nothing in non-parallel.
static const char * name(bool b) noexcept
A string representation of bool as "false" / "true".
T & first()
Access first element of the list, position [0].
Definition UList.H:957
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
Inter-processor communications stream.
Definition UPstream.H:69
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static const Enum< commsTypes > commsTypeNames
Enumerated names for the communication types.
Definition UPstream.H:92
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition UPstream.H:1069
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
static fileName envGlobalPath()
Global case (directory) from environment variable.
Definition argList.C:668
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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...
A class for handling file names.
Definition fileName.H:75
static std::string stem(const std::string &str)
Return the basename, without extension.
Definition fileName.C:391
static bool isAbsolute(const std::string &str)
Return true if filename starts with a '/' or '\' or (windows-only) with a filesystem-root.
Definition fileNameI.H:129
fileName baseFileDir() const
Return the base directory for output.
Definition writeFile.C:43
virtual autoPtr< OFstream > newFile(const fileName &fName) const
Return autoPtr to a new file using file name.
Definition writeFile.C:82
fileName filePath(const fileName &fName) const
Return the full path for the supplied file name.
Definition writeFile.C:73
const objectRegistry & fileObr_
Reference to the region objectRegistry.
Definition writeFile.H:121
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
static void gather(const labelUList &offsets, const label comm, const ProcIDsContainer &procIDs, const UList< Type > &fld, UList< Type > &allFld, const int tag=UPstream::msgType(), UPstream::commsTypes commsType=UPstream::commsTypes::nonBlocking)
Collect data in processor order on master (== procIDs[0]).
labelRange range(label proci) const noexcept
Return start/size range of proci data.
static void scatter(const labelUList &offsets, const label comm, const ProcIDsContainer &procIDs, const UList< Type > &allFld, UList< Type > &fld, const int tag=UPstream::msgType(), UPstream::commsTypes commsType=UPstream::commsTypes::nonBlocking)
Distribute data in processor order.
void reset(label localSize, const label comm=UPstream::worldComm, const bool parallel=UPstream::parRun())
Reset from local size, using gather/broadcast with default/specified communicator if parallel.
label localSize(const label proci) const
Size of proci data.
const labelList & offsets() const noexcept
Const-access to the offsets.
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name.
Definition instant.H:56
static label findStart(const UList< instant > &times, const scalar timeVal)
Find and return index of given start time (linear search).
Definition instant.C:37
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 rhoRef_
Reference density (to convert from kinematic to static pressure).
Definition noiseModel.H:233
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
tmp< Foam::scalarField > PSD(const scalarField &PSDf) const
PSD [dB/Hz].
Definition noiseModel.C:758
bool writePSD_
Write PSD; default = yes.
Definition noiseModel.H:309
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
noiseModel(const noiseModel &)=delete
No copy construct.
tmp< scalarField > SPL(const scalarField &Prms2, const scalar f) const
SPL [dB].
Definition noiseModel.C:767
autoPtr< windowModel > windowModelPtr_
Window model.
Definition noiseModel.H:263
bool writeSPL_
Write SPL; default = yes.
Definition noiseModel.H:304
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
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
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
scalar startTime_
Start time, default = 0s.
Definition noiseModel.H:258
scalar sampleFreq_
Prescribed sample frequency.
Definition noiseModel.H:253
void writeFreqDataToFile(Ostream &os, const scalarField &f, const scalarField &fx) const
Definition noiseModel.C:278
scalar fLower_
Lower frequency limit, default = 25Hz.
Definition noiseModel.H:243
Perform noise analysis on surface-based pressure data.
scalarList times_
Sample times.
word pName_
Name of pressure field.
void initialise(const fileName &fName)
Initialise.
label nFaces_
Global number of surface faces.
surfaceNoise(const dictionary &dict, const objectRegistry &obr, const word &name=typeName, const bool readFields=true)
Constructor.
scalar surfaceAverage(const scalarField &data, const globalIndex &procFaceAddr) const
Calculate the area average value.
scalar deltaT_
Time step (constant).
virtual bool read(const dictionary &dict)
Read from dictionary.
List< fileName > inputFileNames_
Input file names.
UPstream::commsTypes commType_
Communication type (for sending/receiving fields).
autoPtr< surfaceWriter > writerPtr_
Pointer to the surface writer.
bool useBroadcast_
Use broadcast to send entire field to sub-ranks.
label startTimeIndex_
Start time index.
label fftWriteInterval_
Frequency data output interval, default = 1.
void readSurfaceData(const globalIndex &procFaceAddr, List< scalarField > &pData)
Read surface data.
autoPtr< surfaceReader > readerPtr_
Pointer to the surface reader.
bool areaAverage_
Apply area average; default = no (ensemble average) for backwards compatibility.
label pIndex_
Index of pressure field in reader field list.
scalar writeSurfaceData(const fileName &outDirBase, const word &fName, const word &title, const scalar freq, const scalarField &data, const globalIndex &procFaceAddr, const bool writeSurface) const
Write surface data to file.
virtual void calculate()
Calculate.
Registry of regIOobjects.
string & expand(const bool allowEmpty=false)
Inplace expand initial tags, tildes, and all occurrences of environment variables as per stringOps::e...
Definition string.C:166
static autoPtr< surfaceReader > New(const word &readType, const fileName &fName, const dictionary &options=dictionary())
Return a reference to the selected surfaceReader.
static autoPtr< surfaceWriter > New(const word &writeType)
Select construct a surfaceWriter.
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
A token holds an item read from Istream.
Definition token.H:70
Base class for windowing models.
Definition windowModel.H:51
label nSamples() const
Return the number of samples in the window.
Definition windowModel.C:45
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 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.
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
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.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< instant > instantList
List of instants.
Definition instantList.H:41
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.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
MeshedSurface< face > meshedSurface
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Dispatch tag: Construct with a single (local size) entry, no communication.
autoPtr< OFstream > filePtr