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-2020 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_->write(title, allData);
308 
309  writerPtr_->endTime();
310  writerPtr_->clear();
311  }
312 
313  // TO BE VERIFIED: area-averaged values
314  // areaAverage = sum(allData*surf.magSf())/sum(surf.magSf());
315  areaAverage = sum(allData)/(allData.size() + ROOTVSMALL);
316  }
317  Pstream::scatter(areaAverage);
318 
319  return areaAverage;
320  }
321  else
322  {
323  // Note: hard-coded to read mesh from first time index
324  const meshedSurface& surf = readerPtr_->geometry(0);
325 
326  if (writeSurface)
327  {
328  // Time-aware, with time spliced into the output path
329  writerPtr_->beginTime(freqInst);
330 
331  writerPtr_->open
332  (
333  surf.points(),
334  surf.surfFaces(),
335  (outDirBase / fName),
336  false // serial - already merged
337  );
338 
339  writerPtr_->write(title, data);
340 
341  writerPtr_->endTime();
342  writerPtr_->clear();
343  }
344 
345  // TO BE VERIFIED: area-averaged values
346  // return sum(data*surf.magSf())/sum(surf.magSf());
347  return sum(data)/(data.size() + ROOTVSMALL);
348  }
349 }
350 
351 
353 (
354  const scalarField& data,
355  const labelList& procFaceOffset
356 ) const
357 {
358  if (Pstream::parRun())
359  {
360  // Collect the surface data so that we can output the surfaces
361 
363 
364  if (!Pstream::master())
365  {
366  UOPstream toProc(0, pBufs);
367  toProc << data;
368  }
369 
370  pBufs.finishedSends();
371 
372  scalar areaAverage = 0;
373  if (Pstream::master())
374  {
375  // Note: hard-coded to read mesh from first time index
376  const meshedSurface& surf = readerPtr_->geometry(0);
377 
378  scalarField allData(surf.size());
379 
380  forAll(data, faceI)
381  {
382  // Master procFaceOffset is zero...
383  allData[faceI] = data[faceI];
384  }
385 
386  for (const int procI : Pstream::subProcs())
387  {
388  UIPstream fromProc(procI, pBufs);
389  scalarList dataSlice(fromProc);
390  forAll(dataSlice, i)
391  {
392  label faceI = procFaceOffset[procI] + i;
393  allData[faceI] = dataSlice[i];
394  }
395  }
396 
397  // TO BE VERIFIED: area-averaged values
398  // areaAverage = sum(allData*surf.magSf())/sum(surf.magSf());
399  areaAverage = sum(allData)/allData.size();
400  }
401  Pstream::scatter(areaAverage);
402 
403  return areaAverage;
404  }
405  else
406  {
407 
408  // TO BE VERIFIED: area-averaged values
409  // const meshedSurface& surf = readerPtr_->geometry();
410  // return sum(data*surf.magSf())/sum(surf.magSf());
411  return sum(data)/data.size();
412  }
413 }
414 
415 
416 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
417 
419 :
420  noiseModel(dict, false),
421  inputFileNames_(),
422  pName_("p"),
423  pIndex_(0),
424  times_(),
425  deltaT_(0),
426  startTimeIndex_(0),
427  nFace_(0),
428  fftWriteInterval_(1),
429  readerType_(word::null),
430  readerPtr_(nullptr),
431  writerPtr_(nullptr)
432 {
433  if (readFields)
434  {
435  read(dict);
436  }
437 }
438 
439 
440 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
441 
443 {
444  if (noiseModel::read(dict))
445  {
446  if (!dict.readIfPresent("files", inputFileNames_))
447  {
449  dict.readEntry("file", inputFileNames_.first());
450  }
451 
452  dict.readIfPresent("fftWriteInterval", fftWriteInterval_);
453  dict.readIfPresent("p", pName_);
454 
455  readerType_ = dict.get<word>("reader");
456 
457  const word writerType(dict.get<word>("writer"));
458 
460  (
461  writerType,
462  dict.subOrEmptyDict("writeOptions").subOrEmptyDict(writerType)
463  );
464 
465  // Use outputDir/TIME/surface-name
466  writerPtr_->useTimeDir(true);
467 
468  return true;
469  }
470 
471  return false;
472 }
473 
474 
476 {
477  forAll(inputFileNames_, filei)
478  {
479  fileName fName = inputFileNames_[filei];
480  fName.expand();
481 
482  if (!fName.isAbsolute())
483  {
484  fName = argList::envGlobalPath()/fName;
485  }
486 
487  initialise(fName);
488 
489  // Container for pressure time history data per face
490  List<scalarField> pData;
491 
492  // Processor procFaceOffsets
493  labelList procFaceOffset;
494  if (Pstream::parRun())
495  {
496  const label nProcs = Pstream::nProcs();
497  const label nFacePerProc = floor(nFace_/nProcs) + 1;
498 
499  procFaceOffset.setSize(nProcs + 1, 0);
500  for (label i = 1; i < procFaceOffset.size(); ++i)
501  {
502  procFaceOffset[i] = min(i*nFacePerProc, nFace_);
503  }
504  }
505  else
506  {
507  procFaceOffset.setSize(1, nFace_);
508  }
509 
510  // Read pressure data from file
511  readSurfaceData(procFaceOffset, pData);
512 
513  // Process the pressure data, and store results as surface values per
514  // frequency so that it can be output using the surface writer
515 
516  Info<< "Creating noise FFTs" << endl;
517 
518  const scalarField freq1(uniformFrequencies(deltaT_));
519 
520  // Reset desired frequency range if outside actual frequency range
521  fLower_ = min(fLower_, max(freq1));
522  fUpper_ = min(fUpper_, max(freq1));
523 
524  // Storage for FFT data
525  const label nLocalFace = pData.size();
526  const label nFFT = ceil(freq1.size()/scalar(fftWriteInterval_));
527 
528  List<scalarField> surfPrmsf(nFFT);
529  List<scalarField> surfPSDf(nFFT);
530  forAll(surfPrmsf, freqI)
531  {
532  surfPrmsf[freqI].setSize(nLocalFace);
533  surfPSDf[freqI].setSize(nLocalFace);
534  }
535 
536  // Storage for 1/3 octave data
537  labelList octave13BandIDs;
538  scalarField octave13FreqCentre;
540  (
541  freq1,
542  fLower_,
543  fUpper_,
544  3,
545  octave13BandIDs,
546  octave13FreqCentre
547  );
548 
549  label bandSize = 0;
550  if (octave13BandIDs.empty())
551  {
553  << "Octave band calculation failed (zero sized). "
554  << "Please check your input data"
555  << endl;
556  }
557  else
558  {
559  bandSize = octave13BandIDs.size() - 1;
560  }
561 
562  List<scalarField> surfPrms13f(bandSize);
563  forAll(surfPrms13f, freqI)
564  {
565  surfPrms13f[freqI].setSize(nLocalFace);
566  }
567 
568  const windowModel& win = windowModelPtr_();
569 
570  {
571  forAll(pData, faceI)
572  {
573  const scalarField& p = pData[faceI];
574 
575  // Generate the FFT-based data
576  const scalarField Prmsf(RMSmeanPf(p));
577  const scalarField PSDf(this->PSDf(p, deltaT_));
578 
579  // Store the frequency results in slot for face of surface
580  forAll(surfPrmsf, i)
581  {
582  label freqI = i*fftWriteInterval_;
583  surfPrmsf[i][faceI] = Prmsf[freqI];
584  surfPSDf[i][faceI] = PSDf[freqI];
585  }
586 
587  // Integrated PSD = P(rms)^2 [Pa^2]
588  const scalarField Prms13f
589  (
590  octaves
591  (
592  PSDf,
593  freq1,
594  octave13BandIDs
595  )
596  );
597 
598  // Store the 1/3 octave results in slot for face of surface
599  forAll(surfPrms13f, freqI)
600  {
601  surfPrms13f[freqI][faceI] = Prms13f[freqI];
602  }
603  }
604  }
605 
606  const word fNameBase = fName.nameLessExt();
607 
608  // Output directory for graphs
609  fileName outDirBase(baseFileDir(filei)/fNameBase);
610 
611  const scalar deltaf = 1.0/(deltaT_*win.nSamples());
612  Info<< "Writing fft surface data";
613  if (fftWriteInterval_ == 1)
614  {
615  Info<< endl;
616  }
617  else
618  {
619  Info<< " at every " << fftWriteInterval_ << " frequency points"
620  << endl;
621  }
622 
623  {
624  fileName outDir(outDirBase/"fft");
625 
626  // Determine frequency range of interest
627  // Note: frequencies have fixed interval, and are in the range
628  // 0 to fftWriteInterval_*(n-1)*deltaf
629  label f0 = ceil(fLower_/deltaf/scalar(fftWriteInterval_));
630  label f1 = floor(fUpper_/deltaf/scalar(fftWriteInterval_));
631  label nFreq = f1 - f0;
632 
633  scalarField PrmsfAve(nFreq, Zero);
634  scalarField PSDfAve(nFreq, Zero);
635  scalarField fOut(nFreq, Zero);
636 
637  if (nFreq == 0)
638  {
640  << "No surface data available using a fftWriteInterval of "
641  << fftWriteInterval_ << endl;
642  }
643  else
644  {
645  forAll(fOut, i)
646  {
647  label freqI = (i + f0)*fftWriteInterval_;
648  fOut[i] = freq1[freqI];
649 
650 
651  PrmsfAve[i] = writeSurfaceData
652  (
653  outDir,
654  fNameBase,
655  "Prmsf",
656  freq1[freqI],
657  surfPrmsf[i + f0],
658  procFaceOffset,
660  );
661 
662  PSDfAve[i] = writeSurfaceData
663  (
664  outDir,
665  fNameBase,
666  "PSDf",
667  freq1[freqI],
668  surfPSDf[i + f0],
669  procFaceOffset,
670  writePSDf_
671  );
673  (
674  outDir,
675  fNameBase,
676  "PSD",
677  freq1[freqI],
678  PSD(surfPSDf[i + f0]),
679  procFaceOffset,
680  writePSD_
681  );
683  (
684  outDir,
685  fNameBase,
686  "SPL",
687  freq1[freqI],
688  SPL(surfPSDf[i + f0]*deltaf, freq1[freqI]),
689  procFaceOffset,
690  writeSPL_
691  );
692  }
693  }
694 
695  if (Pstream::master())
696  {
697  graph Prmsfg
698  (
699  "Average Prms(f)",
700  "f [Hz]",
701  "P(f) [Pa]",
702  fOut,
703  PrmsfAve
704  );
705  Prmsfg.write
706  (
707  outDir,
708  graph::wordify(Prmsfg.title()),
710  );
711 
712  graph PSDfg
713  (
714  "Average PSD_f(f)",
715  "f [Hz]",
716  "PSD(f) [PaPa_Hz]",
717  fOut,
718  PSDfAve
719  );
720  PSDfg.write
721  (
722  outDir,
723  graph::wordify(PSDfg.title()),
725  );
726 
727  graph PSDg
728  (
729  "Average PSD_dB_Hz(f)",
730  "f [Hz]",
731  "PSD(f) [dB_Hz]",
732  fOut,
733  PSD(PSDfAve)
734  );
735  PSDg.write
736  (
737  outDir,
738  graph::wordify(PSDg.title()),
740  );
741 
742  graph SPLg
743  (
744  "Average SPL_dB(f)",
745  "f [Hz]",
746  "SPL(f) [dB]",
747  fOut,
748  SPL(PSDfAve*deltaf, fOut)
749  );
750  SPLg.write
751  (
752  outDir,
753  graph::wordify(SPLg.title()),
755  );
756  }
757  }
758 
759 
760  Info<< "Writing one-third octave surface data" << endl;
761  {
762  fileName outDir(outDirBase/"oneThirdOctave");
763 
764  scalarField PSDfAve(surfPrms13f.size(), Zero);
765  scalarField Prms13fAve(surfPrms13f.size(), Zero);
766 
767  forAll(surfPrms13f, i)
768  {
770  (
771  outDir,
772  fNameBase,
773  "SPL13",
774  octave13FreqCentre[i],
775  SPL(surfPrms13f[i], octave13FreqCentre[i]),
776  procFaceOffset,
778  );
779 
780  Prms13fAve[i] =
781  surfaceAverage(surfPrms13f[i], procFaceOffset);
782  }
783 
784  if (Pstream::master())
785  {
786  graph SPL13g
787  (
788  "Average SPL13_dB(fm)",
789  "fm [Hz]",
790  "SPL(fm) [dB]",
791  octave13FreqCentre,
792  SPL(Prms13fAve, octave13FreqCentre)
793  );
794  SPL13g.write
795  (
796  outDir,
797  graph::wordify(SPL13g.title()),
799  );
800  }
801  }
802  }
803 }
804 
805 
806 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
807 
808 } // End namespace noiseModels
809 } // End namespace Foam
810 
811 // ************************************************************************* //
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:148
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::UPstream::commsTypes::nonBlocking
Foam::noiseModel::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: noiseModel.C:538
Foam::noiseModels::surfaceNoise::nFace_
label nFace_
Number of surface faces.
Definition: surfaceNoise.H:154
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:402
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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:69
Foam::noiseModel
Base class for noise models.
Definition: noiseModel.H:169
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::argList::envGlobalPath
static fileName envGlobalPath()
Global case (directory) from environment variable.
Definition: argList.C:537
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:236
Foam::noiseModels::surfaceNoise::times_
scalarList times_
Sample times.
Definition: surfaceNoise.H:145
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:87
Foam::noiseModels::surfaceNoise::pIndex_
label pIndex_
Index of pressure field in reader field list.
Definition: surfaceNoise.H:142
Foam::UPstream::parRun
static bool & parRun()
Test if this a parallel run, or allow modify access.
Definition: UPstream.H:434
Foam::MeshedSurface::surfFaces
const List< Face > & surfFaces() const
Return const access to the faces.
Definition: MeshedSurface.H:411
Foam::noiseModel::octaves
tmp< scalarField > octaves(const scalarField &data, const scalarField &f, const labelUList &freqBandIDs) const
Generate octave data.
Definition: noiseModel.C:260
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:458
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:475
Foam::noiseModel::findStartTimeIndex
label findStartTimeIndex(const instantList &allTimes, const scalar startTime) const
Find and return start time index.
Definition: noiseModel.C:205
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::noiseModels::surfaceNoise::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: surfaceNoise.C:442
surfaceWriter.H
Foam::PstreamBuffers::clear
void clear()
Clear storage and reset.
Definition: PstreamBuffers.C:130
surfaceReader.H
Foam::noiseModel::fLower_
scalar fLower_
Lower frequency limit, default = 25Hz.
Definition: noiseModel.H:233
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:139
Foam::noiseModels::surfaceNoise::surfaceNoise
surfaceNoise(const dictionary &dict, const bool readFields=true)
Constructor.
Definition: surfaceNoise.C:418
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:516
Foam::noiseModels::surfaceNoise::deltaT_
scalar deltaT_
Time step (constant)
Definition: surfaceNoise.H:148
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::noiseModel::windowModelPtr_
autoPtr< windowModel > windowModelPtr_
Window model.
Definition: noiseModel.H:242
Foam::noiseModel::writeOctaves_
bool writeOctaves_
Write writeOctaves; default = yes.
Definition: noiseModel.H:278
argList.H
Foam::noiseModel::startTime_
scalar startTime_
Start time, default = 0s.
Definition: noiseModel.H:239
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::noiseModel::writePSD_
bool writePSD_
Write PSD; default = yes.
Definition: noiseModel.H:272
Foam::noiseModels::surfaceNoise::startTimeIndex_
label startTimeIndex_
Start time index.
Definition: surfaceNoise.H:151
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::List::resize
void resize(const label newSize)
Adjust allocated size of list.
Definition: ListI.H:139
Foam::PstreamBuffers::finishedSends
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
Definition: PstreamBuffers.C:80
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::noiseModel::uniformFrequencies
tmp< scalarField > uniformFrequencies(const scalar deltaT) const
Definition: noiseModel.C:239
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:121
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
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:275
Foam::noiseModel::writePrmsf_
bool writePrmsf_
Write Prmsf; default = yes.
Definition: noiseModel.H:266
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:364
Foam::noiseModels::surfaceNoise::fftWriteInterval_
label fftWriteInterval_
Frequency data output interval, default = 1.
Definition: surfaceNoise.H:159
Foam::noiseModels::surfaceNoise::readerPtr_
autoPtr< surfaceReader > readerPtr_
Pointer to the surface reader.
Definition: surfaceNoise.H:165
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::UPstream::allProcs
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition: UPstream.H:509
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:464
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::noiseModel::SPL
tmp< scalarField > SPL(const scalarField &Prms2, const scalar f) const
SPL [dB].
Definition: noiseModel.C:631
Foam::noiseModels::surfaceNoise::readerType_
word readerType_
Reader type.
Definition: surfaceNoise.H:162
Foam::surfaceWriter::New
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:64
Foam::List< word >
Foam::noiseModel::baseFileDir
fileName baseFileDir(const label dataseti) const
Return the base output directory.
Definition: noiseModel.C:224
Foam::string::expand
string & expand(const bool allowEmpty=false)
Definition: string.C:186
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:115
Foam::noiseModel::PSDf
tmp< scalarField > PSDf(const scalarField &p, const scalar deltaT) const
Definition: noiseModel.C:383
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
graph.H
Foam::windowModel::nSamples
label nSamples() const
Return the number of samples in the window.
Definition: windowModel.C:56
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::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::noiseModel::graphFormat_
word graphFormat_
Graph format.
Definition: noiseModel.H:245
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:622
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:405
Foam::noiseModels::surfaceNoise::inputFileNames_
List< fileName > inputFileNames_
Input file names.
Definition: surfaceNoise.H:136
Foam::noiseModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(noiseModel, pointNoise, dictionary)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
Foam::noiseModel::writeSPL_
bool writeSPL_
Write SPL; default = yes.
Definition: noiseModel.H:269
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)
Return true if string starts with a '/'.
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:353
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:446
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:168