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-2022 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 "argList.H"
32#include "graph.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace noiseModels
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
46
47// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
48
50{
51 Info<< "Reading data file " << fName << endl;
52
53 label nAvailableTimes = 0;
54
55 // All reading performed on the master processor only
56 if (Pstream::master())
57 {
58 // Create the surface reader
60
61 // Find the index of the pressure data
62 const List<word> fieldNames(readerPtr_->fieldNames(0));
63 pIndex_ = fieldNames.find(pName_);
64 if (pIndex_ == -1)
65 {
67 << "Unable to find pressure field name " << pName_
68 << " in list of available fields: " << fieldNames
69 << exit(FatalError);
70 }
71
72 // Create the surface writer
73 // - Could be done later, but since this utility can process a lot of
74 // data we can ensure that the user-input is correct prior to doing
75 // the heavy lifting
76
77 // Set the time range
78 const instantList allTimes = readerPtr_->times();
80
81 // Determine the windowing
82 nAvailableTimes = allTimes.size() - startTimeIndex_;
83 }
84
86 (
88 pIndex_,
90 nAvailableTimes
91 );
92
93
94 // Note: all processors should call the windowing validate function
95 label nRequiredTimes = windowModelPtr_->validate(nAvailableTimes);
96
97 if (Pstream::master())
98 {
99 // Restrict times
100 const instantList allTimes = readerPtr_->times();
101
102 times_.setSize(nRequiredTimes);
103 forAll(times_, timeI)
104 {
105 times_[timeI] = allTimes[timeI + startTimeIndex_].value();
106 }
108
109 // Read the surface geometry
110 // Note: hard-coded to read mesh from first time index
111 const meshedSurface& surf = readerPtr_->geometry(0);
112 nFace_ = surf.size();
113 }
114
116 (
118 times_,
119 deltaT_,
120 nFace_
121 );
122}
123
124
126(
127 const labelList& procFaceOffset,
128 List<scalarField>& pData
129)
130{
131 // Data is stored as pressure on surface at a given time. Now we need to
132 // pivot the data so that we have the complete pressure time history per
133 // surface face. In serial mode, this results in all pressure data being
134 // loaded into memory (!)
135
136 Info << "Reading pressure data" << endl;
137
138 if (Pstream::parRun())
139 {
141
142 // Procedure:
143 // 1. Master processor reads pressure data for all faces for all times
144 // 2. Master sends each processor a subset of faces
145 // 3. Local processor reconstructs the full time history for the subset
146 // of faces
147 // Note: reading all data on master to avoid potential NFS problems...
148
149 const label myProcI = Pstream::myProcNo();
150 const label nLocalFace =
151 procFaceOffset[myProcI + 1] - procFaceOffset[myProcI];
152
153 // Complete pressure time history data for subset of faces
154 pData.setSize(nLocalFace);
155 const label nTimes = times_.size();
156 for (scalarField& pf : pData)
157 {
158 pf.setSize(nTimes);
159 }
160
161 // Read and send data
162 forAll(times_, i)
163 {
164 pBufs.clear();
165
166 if (Pstream::master())
167 {
168 label timeI = i + startTimeIndex_;
169
170 Info<< " time: " << times_[i] << endl;
171
172 // Read pressure at all faces for time timeI
173 scalarField p(readerPtr_->field(timeI, pIndex_, scalar(0)));
174
175 // Apply conversions
176 p *= rhoRef_;
177
178 // Send subset of faces to each processor
179 for (const int procI : Pstream::allProcs())
180 {
181 label face0 = procFaceOffset[procI];
182 label nLocalFace =
183 procFaceOffset[procI + 1] - procFaceOffset[procI];
184
185 UOPstream toProc(procI, pBufs);
186 toProc << SubList<scalar>(p, nLocalFace, face0);
187 }
188 }
189
190 pBufs.finishedScatters();
191
192 // Receive data from the master
193 UIPstream fromProc(Pstream::masterNo(), pBufs);
194 scalarList pSlice(fromProc);
195
196 forAll(pSlice, faceI)
197 {
198 pData[faceI][i] = pSlice[faceI];
199 }
200 }
201
202 forAll(pData, faceI)
203 {
204 pData[faceI] -= average(pData[faceI]);
205 }
206 }
207 else
208 {
209 const label nLocalFace = procFaceOffset[0];
210
211 pData.setSize(nLocalFace);
212 for (scalarField& pf : pData)
213 {
214 pf.setSize(times_.size());
215 }
216
217 forAll(times_, i)
218 {
219 label timeI = i + startTimeIndex_;
220
221 Info<< " time: " << times_[i] << endl;
222 scalarField p(readerPtr_->field(timeI, pIndex_, scalar(0)));
223
224 // Apply conversions
225 p *= rhoRef_;
226
227 forAll(p, faceI)
228 {
229 pData[faceI][i] = p[faceI];
230 }
231 }
232
233 forAll(pData, faceI)
234 {
235 pData[faceI] -= average(pData[faceI]);
236 }
237 }
238
239 Info<< "Read "
240 << returnReduce(pData.size(), sumOp<label>())
241 << " pressure traces with "
242 << times_.size()
243 << " time values" << nl << endl;
244}
245
246
248(
249 const fileName& outDirBase,
250 const word& fName,
251 const word& title,
252 const scalar freq,
253 const scalarField& data,
254 const labelList& procFaceOffset,
255 const bool writeSurface
256) const
257{
258 Info<< " processing " << title << " for frequency " << freq << endl;
259
260 const instant freqInst(freq, Foam::name(freq));
261
262 if (Pstream::parRun())
263 {
264 // Collect the surface data so that we can output the surfaces
265
267
268 if (!Pstream::master())
269 {
270 UOPstream toProc(Pstream::masterNo(), pBufs);
271 toProc << data;
272 }
273
274 pBufs.finishedGathers();
275
276 scalar areaAverage = 0;
277 if (Pstream::master())
278 {
279 // Note: hard-coded to read mesh from first time index
280 const meshedSurface& surf = readerPtr_->geometry(0);
281
282 scalarField allData(surf.size());
283
284 forAll(data, faceI)
285 {
286 // Master procFaceOffset is zero...
287 allData[faceI] = data[faceI];
288 }
289
290 for (const int procI : Pstream::subProcs())
291 {
292 UIPstream fromProc(procI, pBufs);
293 scalarList dataSlice(fromProc);
294 forAll(dataSlice, i)
295 {
296 label faceI = procFaceOffset[procI] + i;
297 allData[faceI] = dataSlice[i];
298 }
299 }
300
301 if (writeSurface)
302 {
303 // Time-aware, with time spliced into the output path
304 writerPtr_->beginTime(freqInst);
305
306 writerPtr_->open
307 (
308 surf.points(),
309 surf.surfFaces(),
310 (outDirBase / fName),
311 false // serial - already merged
312 );
313
314 writerPtr_->nFields(1); // Legacy VTK
315 writerPtr_->write(title, allData);
316
317 writerPtr_->endTime();
318 writerPtr_->clear();
319 }
320
321 if (areaAverage_)
322 {
323 areaAverage = sum(allData*surf.magSf())/sum(surf.magSf());
324 }
325 else
326 {
327 areaAverage = sum(allData)/(allData.size() + ROOTVSMALL);
328 }
329 }
330 Pstream::broadcast(areaAverage);
331
332 return areaAverage;
333 }
334 else
335 {
336 // Note: hard-coded to read mesh from first time index
337 const meshedSurface& surf = readerPtr_->geometry(0);
338
339 if (writeSurface)
340 {
341 // Time-aware, with time spliced into the output path
342 writerPtr_->beginTime(freqInst);
343
344 writerPtr_->open
345 (
346 surf.points(),
347 surf.surfFaces(),
348 (outDirBase / fName),
349 false // serial - already merged
350 );
351
352 writerPtr_->nFields(1); // Legacy VTK
353 writerPtr_->write(title, data);
354
355 writerPtr_->endTime();
356 writerPtr_->clear();
357 }
358
359 if (areaAverage_)
360 {
361 return sum(data*surf.magSf())/sum(surf.magSf());
362 }
363 else
364 {
365 return sum(data)/(data.size() + ROOTVSMALL);
366 }
367 }
368}
369
370
372(
373 const scalarField& data,
374 const labelList& procFaceOffset
375) const
376{
377 if (Pstream::parRun())
378 {
379 // Collect the surface data so that we can output the surfaces
380
382
383 if (!Pstream::master())
384 {
385 UOPstream toProc(Pstream::masterNo(), pBufs);
386 toProc << data;
387 }
388
389 pBufs.finishedGathers();
390
391 scalar areaAverage = 0;
392 if (Pstream::master())
393 {
394 // Note: hard-coded to read mesh from first time index
395 const meshedSurface& surf = readerPtr_->geometry(0);
396
397 scalarField allData(surf.size());
398
399 forAll(data, faceI)
400 {
401 // Master procFaceOffset is zero...
402 allData[faceI] = data[faceI];
403 }
404
405 for (const int procI : Pstream::subProcs())
406 {
407 UIPstream fromProc(procI, pBufs);
408 scalarList dataSlice(fromProc);
409 forAll(dataSlice, i)
410 {
411 label faceI = procFaceOffset[procI] + i;
412 allData[faceI] = dataSlice[i];
413 }
414 }
415
416 if (areaAverage_)
417 {
418 areaAverage = sum(allData*surf.magSf())/sum(surf.magSf());
419 }
420 else
421 {
422 areaAverage = sum(allData)/allData.size();
423 }
424 }
425 Pstream::broadcast(areaAverage);
426
427 return areaAverage;
428 }
429 else
430 {
431 if (areaAverage_)
432 {
433 // Note: hard-coded to read mesh from first time index
434 const meshedSurface& surf = readerPtr_->geometry(0);
435 return sum(data*surf.magSf())/sum(surf.magSf());
436 }
437
438 return sum(data)/data.size();
439 }
440}
441
442
443// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
444
446:
447 noiseModel(dict, false),
448 inputFileNames_(),
449 pName_("p"),
450 pIndex_(0),
451 times_(),
452 deltaT_(0),
453 startTimeIndex_(0),
454 nFace_(0),
455 fftWriteInterval_(1),
456 areaAverage_(false),
457 readerType_(word::null),
458 readerPtr_(nullptr),
459 writerPtr_(nullptr)
460{
461 if (readFields)
462 {
463 read(dict);
464 }
465}
466
467
468// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
469
471{
473 {
474 if (!dict.readIfPresent("files", inputFileNames_))
475 {
478 }
479
481 dict.readIfPresent("fftWriteInterval", fftWriteInterval_);
482
483 Info<< " Pressure field name: " << pName_ << nl
484 << " FFT write interval: " << fftWriteInterval_ << nl;
485
486 dict.readIfPresent("areaAverage", areaAverage_);
487
488 if (areaAverage_)
489 {
490 Info<< " Averaging: area weighted" << endl;
491 }
492 else
493 {
494 Info<< " Averaging: ensemble" << endl;
495 }
496
497 readerType_ = dict.get<word>("reader");
498
499 const word writerType(dict.get<word>("writer"));
500
502 (
503 writerType,
504 dict.subOrEmptyDict("writeOptions").subOrEmptyDict(writerType)
505 );
506
507 // Use outputDir/TIME/surface-name
508 writerPtr_->useTimeDir(true);
509
510 Info << endl;
511
512 return true;
513 }
514
515 return false;
516}
517
518
520{
521 forAll(inputFileNames_, filei)
522 {
523 fileName fName = inputFileNames_[filei];
524 fName.expand();
525
526 if (!fName.isAbsolute())
527 {
528 fName = argList::envGlobalPath()/fName;
529 }
530
531 initialise(fName);
532
533 // Container for pressure time history data per face
534 List<scalarField> pData;
535
536 // Processor procFaceOffsets
537 labelList procFaceOffset;
538 if (Pstream::parRun())
539 {
540 const label nProcs = Pstream::nProcs();
541 const label nFacePerProc = floor(nFace_/nProcs) + 1;
542
543 procFaceOffset.setSize(nProcs + 1, 0);
544 for (label i = 1; i < procFaceOffset.size(); ++i)
545 {
546 procFaceOffset[i] = min(i*nFacePerProc, nFace_);
547 }
548 }
549 else
550 {
551 procFaceOffset.setSize(1, nFace_);
552 }
553
554 // Read pressure data from file
555 readSurfaceData(procFaceOffset, pData);
556
557 // Process the pressure data, and store results as surface values per
558 // frequency so that it can be output using the surface writer
559
560 Info<< "Creating noise FFTs" << endl;
561
562 const scalarField freq1(uniformFrequencies(deltaT_, true));
563
564 // Reset desired frequency range if outside actual frequency range
565 fLower_ = min(fLower_, max(freq1));
566 fUpper_ = min(fUpper_, max(freq1));
567
568 // Storage for FFT data
569 const label nLocalFace = pData.size();
570 const label nFFT = ceil(freq1.size()/scalar(fftWriteInterval_));
571
572 List<scalarField> surfPrmsf(nFFT);
573 List<scalarField> surfPSDf(nFFT);
574 forAll(surfPrmsf, freqI)
575 {
576 surfPrmsf[freqI].setSize(nLocalFace);
577 surfPSDf[freqI].setSize(nLocalFace);
578 }
579
580 // Storage for 1/3 octave data
581 labelList octave13BandIDs;
582 scalarField octave13FreqCentre;
584 (
585 freq1,
586 fLower_,
587 fUpper_,
588 3,
589 octave13BandIDs,
590 octave13FreqCentre
591 );
592
593 label bandSize = 0;
594 if (octave13BandIDs.empty())
595 {
597 << "Octave band calculation failed (zero sized). "
598 << "Please check your input data"
599 << endl;
600 }
601 else
602 {
603 bandSize = octave13BandIDs.size() - 1;
604 }
605
606 List<scalarField> surfPrms13f(bandSize);
607 forAll(surfPrms13f, freqI)
608 {
609 surfPrms13f[freqI].setSize(nLocalFace);
610 }
611
612 const windowModel& win = windowModelPtr_();
613
614 {
615 forAll(pData, faceI)
616 {
617 const scalarField& p = pData[faceI];
618
619 // Generate the FFT-based data
620 const scalarField Prmsf(RMSmeanPf(p));
621 const scalarField PSDf(this->PSDf(p, deltaT_));
622
623 // Store the frequency results in slot for face of surface
624 forAll(surfPrmsf, i)
625 {
626 label freqI = i*fftWriteInterval_;
627 surfPrmsf[i][faceI] = Prmsf[freqI];
628 surfPSDf[i][faceI] = PSDf[freqI];
629 }
630
631 // Integrated PSD = P(rms)^2 [Pa^2]
632 const scalarField Prms13f
633 (
634 octaves
635 (
636 PSDf,
637 freq1,
638 octave13BandIDs
639 )
640 );
641
642 // Store the 1/3 octave results in slot for face of surface
643 forAll(surfPrms13f, freqI)
644 {
645 surfPrms13f[freqI][faceI] = Prms13f[freqI];
646 }
647 }
648 }
649
650 const word fNameBase = fName.nameLessExt();
651
652 // Output directory for graphs
653 fileName outDirBase(baseFileDir(filei)/fNameBase);
654
655 const scalar deltaf = 1.0/(deltaT_*win.nSamples());
656 Info<< "Writing fft surface data";
657 if (fftWriteInterval_ == 1)
658 {
659 Info<< endl;
660 }
661 else
662 {
663 Info<< " at every " << fftWriteInterval_ << " frequency points"
664 << endl;
665 }
666
667 {
668 fileName outDir(outDirBase/"fft");
669
670 // Determine frequency range of interest
671 // Note: frequencies have fixed interval, and are in the range
672 // 0 to fftWriteInterval_*(n-1)*deltaf
673 label f0 = ceil(fLower_/deltaf/scalar(fftWriteInterval_));
674 label f1 = floor(fUpper_/deltaf/scalar(fftWriteInterval_));
675 label nFreq = f1 - f0;
676
677 scalarField PrmsfAve(nFreq, Zero);
678 scalarField PSDfAve(nFreq, Zero);
679 scalarField fOut(nFreq, Zero);
680
681 if (nFreq == 0)
682 {
684 << "No surface data available using a fftWriteInterval of "
686 }
687 else
688 {
689 forAll(fOut, i)
690 {
691 label freqI = (i + f0)*fftWriteInterval_;
692 fOut[i] = freq1[freqI];
693
694
695 PrmsfAve[i] = writeSurfaceData
696 (
697 outDir,
698 fNameBase,
699 "Prmsf",
700 freq1[freqI],
701 surfPrmsf[i + f0],
702 procFaceOffset,
704 );
705
706 PSDfAve[i] = writeSurfaceData
707 (
708 outDir,
709 fNameBase,
710 "PSDf",
711 freq1[freqI],
712 surfPSDf[i + f0],
713 procFaceOffset,
715 );
717 (
718 outDir,
719 fNameBase,
720 "PSD",
721 freq1[freqI],
722 PSD(surfPSDf[i + f0]),
723 procFaceOffset,
725 );
727 (
728 outDir,
729 fNameBase,
730 "SPL",
731 freq1[freqI],
732 SPL(surfPSDf[i + f0]*deltaf, freq1[freqI]),
733 procFaceOffset,
735 );
736 }
737 }
738
739 if (Pstream::master())
740 {
741 graph Prmsfg
742 (
743 "Average Prms(f)",
744 "f [Hz]",
745 "P(f) [Pa]",
746 fOut,
747 PrmsfAve
748 );
749 Prmsfg.write
750 (
751 outDir,
752 graph::wordify(Prmsfg.title()),
754 );
755
756 graph PSDfg
757 (
758 "Average PSD_f(f)",
759 "f [Hz]",
760 "PSD(f) [PaPa_Hz]",
761 fOut,
762 PSDfAve
763 );
764 PSDfg.write
765 (
766 outDir,
767 graph::wordify(PSDfg.title()),
769 );
770
771 graph PSDg
772 (
773 "Average PSD_dB_Hz(f)",
774 "f [Hz]",
775 "PSD(f) [dB_Hz]",
776 fOut,
777 PSD(PSDfAve)
778 );
779 PSDg.write
780 (
781 outDir,
782 graph::wordify(PSDg.title()),
784 );
785
786 graph SPLg
787 (
788 "Average SPL_dB(f)",
789 "f [Hz]",
790 "SPL(f) [dB]",
791 fOut,
792 SPL(PSDfAve*deltaf, fOut)
793 );
794 SPLg.write
795 (
796 outDir,
797 graph::wordify(SPLg.title()),
799 );
800 }
801 }
802
803
804 Info<< "Writing one-third octave surface data" << endl;
805 {
806 fileName outDir(outDirBase/"oneThirdOctave");
807
808 scalarField PSDfAve(surfPrms13f.size(), Zero);
809 scalarField Prms13fAve(surfPrms13f.size(), Zero);
810
811 forAll(surfPrms13f, i)
812 {
814 (
815 outDir,
816 fNameBase,
817 "SPL13",
818 octave13FreqCentre[i],
819 SPL(surfPrms13f[i], octave13FreqCentre[i]),
820 procFaceOffset,
822 );
823
824 Prms13fAve[i] =
825 surfaceAverage(surfPrms13f[i], procFaceOffset);
826 }
827
828 if (Pstream::master())
829 {
830 graph SPL13g
831 (
832 "Average SPL13_dB(fm)",
833 "fm [Hz]",
834 "SPL(fm) [dB]",
835 octave13FreqCentre,
836 SPL(Prms13fAve, octave13FreqCentre)
837 );
838 SPL13g.write
839 (
840 outDir,
841 graph::wordify(SPL13g.title()),
843 );
844 }
845 }
846 }
847}
848
849
850// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
851
852} // End namespace noiseModels
853} // End namespace Foam
854
855// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
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.
const Field< point_type > & points() const noexcept
Return reference to global points.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
UPstream::rangeType allProcs() const noexcept
Range of ranks indices associated with PstreamBuffers.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
void finishedScatters(const bool wait=true)
Mark all sends to sub-procs as done.
void clear()
Clear individual buffers and reset states.
void finishedGathers(const bool wait=true)
Mark all sends to master as done.
UPstream::rangeType subProcs() const noexcept
Range of sub-processes indices associated with PstreamBuffers.
static void broadcast(Type &value, const label comm=UPstream::worldComm)
static void broadcasts(const label comm, Type &arg1, Args &&... args)
Broadcast multiple items to all processes in communicator.
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
T & first()
Return the first element of the list.
Definition: UListI.H:202
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:212
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
static constexpr int masterNo() noexcept
Process index of the master (always 0)
Definition: UPstream.H:451
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:293
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
static fileName envGlobalPath()
Global case (directory) from environment variable.
Definition: argList.C:549
Database for solution data, solver performance and other reduced data.
Definition: data.H:58
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
A class for handling file names.
Definition: fileName.H:76
static std::string nameLessExt(const std::string &str)
Return basename, without extension.
Definition: fileName.C:396
static bool isAbsolute(const std::string &str)
Definition: fileNameI.H:136
Class to create, store and output qgraph files.
Definition: graph.H:62
const string & title() const
Definition: graph.H:147
void write(Ostream &, const word &format) const
Write graph to stream in given format.
Definition: graph.C:266
static word wordify(const string &sname)
Helper function to convert string name into appropriate word.
Definition: graph.C:47
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name.
Definition: instant.H:56
Base class for noise models.
Definition: noiseModel.H:176
tmp< scalarField > uniformFrequencies(const scalar deltaT, const bool check) const
Definition: noiseModel.C:277
label findStartTimeIndex(const instantList &allTimes, const scalar startTime) const
Find and return start time index.
Definition: noiseModel.C:243
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:314
tmp< scalarField > RMSmeanPf(const scalarField &p) const
Definition: noiseModel.C:429
scalar fUpper_
Upper frequency limit, default = 10kHz.
Definition: noiseModel.H:242
tmp< Foam::scalarField > PSD(const scalarField &PSDf) const
PSD [dB/Hz].
Definition: noiseModel.C:713
bool writePSD_
Write PSD; default = yes.
Definition: noiseModel.H:281
word graphFormat_
Graph format.
Definition: noiseModel.H:251
static void setOctaveBands(const scalarField &f, const scalar fLower, const scalar fUpper, const scalar octave, labelList &fBandIDs, scalarField &fCentre)
Definition: noiseModel.C:55
fileName baseFileDir(const label dataseti) const
Return the base output directory.
Definition: noiseModel.C:262
bool writeOctaves_
Write writeOctaves; default = yes.
Definition: noiseModel.H:287
tmp< scalarField > SPL(const scalarField &Prms2, const scalar f) const
SPL [dB].
Definition: noiseModel.C:722
autoPtr< windowModel > windowModelPtr_
Window model.
Definition: noiseModel.H:248
bool writeSPL_
Write SPL; default = yes.
Definition: noiseModel.H:278
bool writePrmsf_
Write Prmsf; default = yes.
Definition: noiseModel.H:275
scalar checkUniformTimeStep(const scalarList &times) const
Check and return uniform time step.
Definition: noiseModel.C:186
bool writePSDf_
Write PSDf; default = yes.
Definition: noiseModel.H:284
tmp< scalarField > PSDf(const scalarField &p, const scalar deltaT) const
Definition: noiseModel.C:448
scalar startTime_
Start time, default = 0s.
Definition: noiseModel.H:245
scalar fLower_
Lower frequency limit, default = 25Hz.
Definition: noiseModel.H:239
Perform noise analysis on surface-based pressure data.
Definition: surfaceNoise.H:133
scalarList times_
Sample times.
Definition: surfaceNoise.H:149
scalar surfaceAverage(const scalarField &data, const labelList &procFaceOffset) const
Calculate the area average value.
Definition: surfaceNoise.C:372
word pName_
Name of pressure field.
Definition: surfaceNoise.H:143
void initialise(const fileName &fName)
Initialise.
Definition: surfaceNoise.C:49
scalar deltaT_
Time step (constant)
Definition: surfaceNoise.H:152
label nFace_
Number of surface faces.
Definition: surfaceNoise.H:158
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: surfaceNoise.C:470
List< fileName > inputFileNames_
Input file names.
Definition: surfaceNoise.H:140
autoPtr< surfaceWriter > writerPtr_
Pointer to the surface writer.
Definition: surfaceNoise.H:176
label startTimeIndex_
Start time index.
Definition: surfaceNoise.H:155
label fftWriteInterval_
Frequency data output interval, default = 1.
Definition: surfaceNoise.H:163
scalar writeSurfaceData(const fileName &outDirBase, const word &fName, const word &title, const scalar freq, const scalarField &data, const labelList &procFaceOffset, const bool writeSurface) const
Write surface data to file.
Definition: surfaceNoise.C:248
autoPtr< surfaceReader > readerPtr_
Pointer to the surface reader.
Definition: surfaceNoise.H:173
label pIndex_
Index of pressure field in reader field list.
Definition: surfaceNoise.H:146
void readSurfaceData(const labelList &procFaceOffset, List< scalarField > &pData)
Read surface data.
Definition: surfaceNoise.C:126
virtual void calculate()
Calculate.
Definition: surfaceNoise.C:519
int myProcNo() const noexcept
Return processor number.
splitCell * master() const
Definition: splitCell.H:113
string & expand(const bool allowEmpty=false)
Definition: string.C:173
void initialise()
Initialise integral-scale box properties.
Base class for windowing models.
Definition: windowModel.H:55
label nSamples() const
Return the number of samples in the window.
Definition: windowModel.C:56
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333