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