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-2019 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 "noiseFFT.H"
32 #include "argList.H"
33 #include "graph.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace noiseModels
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(surfaceNoise, 0);
46 addToRunTimeSelectionTable(noiseModel, surfaceNoise, dictionary);
47 
48 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
49 
51 {
52  Info<< "Reading data file " << fName << endl;
53 
54  label nAvailableTimes = 0;
55 
56  // All reading performed on the master processor only
57  if (Pstream::master())
58  {
59  // Create the surface reader
61 
62  // Find the index of the pressure data
63  const List<word> fieldNames(readerPtr_->fieldNames(0));
64  pIndex_ = fieldNames.find(pName_);
65  if (pIndex_ == -1)
66  {
68  << "Unable to find pressure field name " << pName_
69  << " in list of available fields: " << fieldNames
70  << exit(FatalError);
71  }
72 
73  // Create the surface writer
74  // - Could be done later, but since this utility can process a lot of
75  // data we can ensure that the user-input is correct prior to doing
76  // the heavy lifting
77 
78  // Set the time range
79  const instantList allTimes = readerPtr_->times();
81 
82  // Determine the windowing
83  nAvailableTimes = allTimes.size() - startTimeIndex_;
84  }
85 
88  Pstream::scatter(nAvailableTimes);
89 
90 
91  // Note: all processors should call the windowing validate function
92  label nRequiredTimes = windowModelPtr_->validate(nAvailableTimes);
93 
94  if (Pstream::master())
95  {
96  // Restrict times
97  const instantList allTimes = readerPtr_->times();
98 
99  times_.setSize(nRequiredTimes);
100  forAll(times_, timeI)
101  {
102  times_[timeI] = allTimes[timeI + startTimeIndex_].value();
103  }
105 
106  // Read the surface geometry
107  // Note: hard-coded to read mesh from first time index
108  const meshedSurface& surf = readerPtr_->geometry(0);
109  nFace_ = surf.size();
110  }
111 
115 }
116 
117 
119 (
120  const labelList& procFaceOffset,
121  List<scalarField>& pData
122 )
123 {
124  // Data is stored as pressure on surface at a given time. Now we need to
125  // pivot the data so that we have the complete pressure time history per
126  // surface face. In serial mode, this results in all pressure data being
127  // loaded into memory (!)
128 
129  Info << "Reading pressure data" << endl;
130 
131  if (Pstream::parRun())
132  {
134 
135  // Procedure:
136  // 1. Master processor reads pressure data for all faces for all times
137  // 2. Master sends each processor a subset of faces
138  // 3. Local processor reconstructs the full time history for the subset
139  // of faces
140  // Note: reading all data on master to avoid potential NFS problems...
141 
142  const label myProcI = Pstream::myProcNo();
143  const label nLocalFace =
144  procFaceOffset[myProcI + 1] - procFaceOffset[myProcI];
145 
146  // Complete pressure time history data for subset of faces
147  pData.setSize(nLocalFace);
148  const label nTimes = times_.size();
149  for (scalarField& pf : pData)
150  {
151  pf.setSize(nTimes);
152  }
153 
154  // Read and send data
155  forAll(times_, i)
156  {
157  pBufs.clear();
158 
159  if (Pstream::master())
160  {
161  label timeI = i + startTimeIndex_;
162 
163  Info<< " time: " << times_[i] << endl;
164 
165  // Read pressure at all faces for time timeI
166  scalarField p(readerPtr_->field(timeI, pIndex_, scalar(0)));
167 
168  // Apply conversions
169  p *= rhoRef_;
170 
171  // Send subset of faces to each processor
172  for (label procI = 0; procI < Pstream::nProcs(); ++procI)
173  {
174  label face0 = procFaceOffset[procI];
175  label nLocalFace =
176  procFaceOffset[procI + 1] - procFaceOffset[procI];
177 
178  UOPstream toProc(procI, pBufs);
179  toProc << SubList<scalar>(p, nLocalFace, face0);
180  }
181  }
182 
183  pBufs.finishedSends();
184 
185  // Receive data from the master
186  UIPstream fromProc(0, pBufs);
187 
188  scalarList pSlice(fromProc);
189 
190  forAll(pSlice, faceI)
191  {
192  pData[faceI][i] = pSlice[faceI];
193  }
194  }
195  }
196  else
197  {
198  const label nLocalFace = procFaceOffset[0];
199 
200  pData.setSize(nLocalFace);
201  for (scalarField& pf : pData)
202  {
203  pf.setSize(times_.size());
204  }
205 
206  forAll(times_, i)
207  {
208  label timeI = i + startTimeIndex_;
209 
210  Info<< " time: " << times_[i] << endl;
211  const scalarField p(readerPtr_->field(timeI, pIndex_, scalar(0)));
212  forAll(p, faceI)
213  {
214  pData[faceI][i] = p[faceI]*rhoRef_;
215  }
216  }
217  }
218 
219  Info<< "Read "
220  << returnReduce(pData.size(), sumOp<label>())
221  << " pressure traces with "
222  << times_.size()
223  << " time values" << nl << endl;
224 }
225 
226 
228 (
229  const fileName& outDirBase,
230  const word& fName,
231  const word& title,
232  const scalar freq,
233  const scalarField& data,
234  const labelList& procFaceOffset,
235  const bool writeSurface
236 ) const
237 {
238  Info<< " processing " << title << " for frequency " << freq << endl;
239 
240  const instant freqInst(freq, Foam::name(freq));
241 
242  if (Pstream::parRun())
243  {
244  // Collect the surface data so that we can output the surfaces
245 
247 
248  if (!Pstream::master())
249  {
250  UOPstream toProc(0, pBufs);
251  toProc << data;
252  }
253 
254  pBufs.finishedSends();
255 
256  scalar areaAverage = 0;
257  if (Pstream::master())
258  {
259  // Note: hard-coded to read mesh from first time index
260  const meshedSurface& surf = readerPtr_->geometry(0);
261 
262  scalarField allData(surf.size());
263 
264  forAll(data, faceI)
265  {
266  // Master procFaceOffset is zero...
267  allData[faceI] = data[faceI];
268  }
269 
270  for (label procI = 1; procI < Pstream::nProcs(); ++procI)
271  {
272  UIPstream fromProc(procI, pBufs);
273  scalarList dataSlice(fromProc);
274  forAll(dataSlice, i)
275  {
276  label faceI = procFaceOffset[procI] + i;
277  allData[faceI] = dataSlice[i];
278  }
279  }
280 
281  if (writeSurface)
282  {
283  // Time-aware, with time spliced into the output path
284  writerPtr_->beginTime(freqInst);
285 
286  writerPtr_->open
287  (
288  surf.points(),
289  surf.surfFaces(),
290  (outDirBase / fName),
291  false // serial - already merged
292  );
293 
294  writerPtr_->write(title, allData);
295 
296  writerPtr_->endTime();
297  writerPtr_->clear();
298  }
299 
300  // TO BE VERIFIED: area-averaged values
301  // areaAverage = sum(allData*surf.magSf())/sum(surf.magSf());
302  areaAverage = sum(allData)/(allData.size() + ROOTVSMALL);
303  }
304  Pstream::scatter(areaAverage);
305 
306  return areaAverage;
307  }
308  else
309  {
310  // Note: hard-coded to read mesh from first time index
311  const meshedSurface& surf = readerPtr_->geometry(0);
312 
313  if (writeSurface)
314  {
315  // Time-aware, with time spliced into the output path
316  writerPtr_->beginTime(freqInst);
317 
318  writerPtr_->open
319  (
320  surf.points(),
321  surf.surfFaces(),
322  (outDirBase / fName),
323  false // serial - already merged
324  );
325 
326  writerPtr_->write(title, data);
327 
328  writerPtr_->endTime();
329  writerPtr_->clear();
330  }
331 
332  // TO BE VERIFIED: area-averaged values
333  // return sum(data*surf.magSf())/sum(surf.magSf());
334  return sum(data)/(data.size() + ROOTVSMALL);
335  }
336 }
337 
338 
339 Foam::scalar surfaceNoise::surfaceAverage
340 (
341  const scalarField& data,
342  const labelList& procFaceOffset
343 ) const
344 {
345  if (Pstream::parRun())
346  {
347  // Collect the surface data so that we can output the surfaces
348 
350 
351  if (!Pstream::master())
352  {
353  UOPstream toProc(0, pBufs);
354  toProc << data;
355  }
356 
357  pBufs.finishedSends();
358 
359  scalar areaAverage = 0;
360  if (Pstream::master())
361  {
362  // Note: hard-coded to read mesh from first time index
363  const meshedSurface& surf = readerPtr_->geometry(0);
364 
365  scalarField allData(surf.size());
366 
367  forAll(data, faceI)
368  {
369  // Master procFaceOffset is zero...
370  allData[faceI] = data[faceI];
371  }
372 
373  for (label procI = 1; procI < Pstream::nProcs(); procI++)
374  {
375  UIPstream fromProc(procI, pBufs);
376  scalarList dataSlice(fromProc);
377  forAll(dataSlice, i)
378  {
379  label faceI = procFaceOffset[procI] + i;
380  allData[faceI] = dataSlice[i];
381  }
382  }
383 
384  // TO BE VERIFIED: area-averaged values
385  // areaAverage = sum(allData*surf.magSf())/sum(surf.magSf());
386  areaAverage = sum(allData)/allData.size();
387  }
388  Pstream::scatter(areaAverage);
389 
390  return areaAverage;
391  }
392  else
393  {
394 
395  // TO BE VERIFIED: area-averaged values
396  // const meshedSurface& surf = readerPtr_->geometry();
397  // return sum(data*surf.magSf())/sum(surf.magSf());
398  return sum(data)/data.size();
399  }
400 }
401 
402 
403 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
404 
406 :
407  noiseModel(dict, false),
408  inputFileNames_(),
409  pName_("p"),
410  pIndex_(0),
411  times_(),
412  deltaT_(0),
413  startTimeIndex_(0),
414  nFace_(0),
415  fftWriteInterval_(1),
416  readerType_(word::null),
417  readerPtr_(nullptr),
418  writerPtr_(nullptr)
419 {
420  if (readFields)
421  {
422  read(dict);
423  }
424 }
425 
426 
427 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
428 
430 {}
431 
432 
433 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
434 
436 {
437  if (noiseModel::read(dict))
438  {
439  if (!dict.readIfPresent("files", inputFileNames_))
440  {
442  dict.readEntry("file", inputFileNames_.first());
443  }
444 
445  dict.readIfPresent("fftWriteInterval", fftWriteInterval_);
446  dict.readIfPresent("p", pName_);
447 
448  readerType_ = dict.get<word>("reader");
449 
450  const word writerType(dict.get<word>("writer"));
451 
453  (
454  writerType,
455  dict.subOrEmptyDict("writeOptions").subOrEmptyDict(writerType)
456  );
457 
458  // Use outputDir/TIME/surface-name
459  writerPtr_->useTimeDir() = true;
460 
461  return true;
462  }
463 
464  return false;
465 }
466 
467 
469 {
470  forAll(inputFileNames_, filei)
471  {
472  fileName fName = inputFileNames_[filei];
473  fName.expand();
474 
475  if (!fName.isAbsolute())
476  {
477  fName = argList::envGlobalPath()/fName;
478  }
479 
480  initialise(fName);
481 
482  // Container for pressure time history data per face
483  List<scalarField> pData;
484 
485  // Processor procFaceOffsets
486  labelList procFaceOffset;
487  if (Pstream::parRun())
488  {
489  const label nProcs = Pstream::nProcs();
490  const label nFacePerProc = floor(nFace_/nProcs) + 1;
491 
492  procFaceOffset.setSize(nProcs + 1, 0);
493  for (label i = 1; i < procFaceOffset.size(); ++i)
494  {
495  procFaceOffset[i] = min(i*nFacePerProc, nFace_);
496  }
497  }
498  else
499  {
500  procFaceOffset.setSize(1, nFace_);
501  }
502 
503  // Read pressure data from file
504  readSurfaceData(procFaceOffset, pData);
505 
506  // Process the pressure data, and store results as surface values per
507  // frequency so that it can be output using the surface writer
508 
509  Info<< "Creating noise FFTs" << endl;
510 
512 
513  // Reset desired frequency range if outside actual frequency range
514  fLower_ = min(fLower_, max(freq1));
515  fUpper_ = min(fUpper_, max(freq1));
516 
517  // Storage for FFT data
518  const label nLocalFace = pData.size();
519  const label nFFT = ceil(freq1.size()/scalar(fftWriteInterval_));
520 
521  List<scalarField> surfPrmsf(nFFT);
522  List<scalarField> surfPSDf(nFFT);
523  forAll(surfPrmsf, freqI)
524  {
525  surfPrmsf[freqI].setSize(nLocalFace);
526  surfPSDf[freqI].setSize(nLocalFace);
527  }
528 
529  // Storage for 1/3 octave data
530  labelList octave13BandIDs;
531  scalarField octave13FreqCentre;
533  (
534  freq1,
535  fLower_,
536  fUpper_,
537  3,
538  octave13BandIDs,
539  octave13FreqCentre
540  );
541 
542  label bandSize = 0;
543  if (octave13BandIDs.empty())
544  {
546  << "Octave band calculation failed (zero sized). "
547  << "Please check your input data"
548  << endl;
549  }
550  else
551  {
552  bandSize = octave13BandIDs.size() - 1;
553  }
554 
555  List<scalarField> surfPrms13f(bandSize);
556  forAll(surfPrms13f, freqI)
557  {
558  surfPrms13f[freqI].setSize(nLocalFace);
559  }
560 
561  const windowModel& win = windowModelPtr_();
562 
563  {
564  noiseFFT nfft(deltaT_, win.nSamples());
565 
566  forAll(pData, faceI)
567  {
568  // Note: passes storage from pData to noiseFFT
569  nfft.setData(pData[faceI]);
570 
571  // Generate the FFT-based data
572  graph Prmsf(nfft.RMSmeanPf(win));
573  graph PSDf(nfft.PSDf(win));
574 
575  // Store the frequency results in slot for face of surface
576  forAll(surfPrmsf, i)
577  {
578  label freqI = i*fftWriteInterval_;
579  surfPrmsf[i][faceI] = Prmsf.y()[freqI];
580  surfPSDf[i][faceI] = PSDf.y()[freqI];
581  }
582 
583  // Integrated PSD = P(rms)^2 [Pa^2]
584  graph Prms13f(nfft.octaves(PSDf, octave13BandIDs));
585 
586  // Store the 1/3 octave results in slot for face of surface
587  forAll(surfPrms13f, freqI)
588  {
589  surfPrms13f[freqI][faceI] = Prms13f.y()[freqI];
590  }
591  }
592  }
593 
594  const word fNameBase = fName.nameLessExt();
595 
596  // Output directory for graphs
597  fileName outDirBase(baseFileDir(filei)/fNameBase);
598 
599  const scalar deltaf = 1.0/(deltaT_*win.nSamples());
600  Info<< "Writing fft surface data";
601  if (fftWriteInterval_ == 1)
602  {
603  Info<< endl;
604  }
605  else
606  {
607  Info<< " at every " << fftWriteInterval_ << " frequency points"
608  << endl;
609  }
610 
611  {
612  fileName outDir(outDirBase/"fft");
613 
614  // Determine frequency range of interest
615  // Note: frequencies have fixed interval, and are in the range
616  // 0 to fftWriteInterval_*(n-1)*deltaf
617  label f0 = ceil(fLower_/deltaf/scalar(fftWriteInterval_));
618  label f1 = floor(fUpper_/deltaf/scalar(fftWriteInterval_));
619  label nFreq = f1 - f0;
620 
621  scalarField PrmsfAve(nFreq, Zero);
622  scalarField PSDfAve(nFreq, Zero);
623  scalarField fOut(nFreq, Zero);
624 
625  if (nFreq == 0)
626  {
628  << "No surface data available using a fftWriteInterval of "
629  << fftWriteInterval_ << endl;
630  }
631  else
632  {
633  forAll(fOut, i)
634  {
635  label freqI = (i + f0)*fftWriteInterval_;
636  fOut[i] = freq1[freqI];
637 
638 
639  PrmsfAve[i] = writeSurfaceData
640  (
641  outDir,
642  fNameBase,
643  "Prmsf",
644  freq1[freqI],
645  surfPrmsf[i + f0],
646  procFaceOffset,
648  );
649 
650  PSDfAve[i] = writeSurfaceData
651  (
652  outDir,
653  fNameBase,
654  "PSDf",
655  freq1[freqI],
656  surfPSDf[i + f0],
657  procFaceOffset,
658  writePSDf_
659  );
661  (
662  outDir,
663  fNameBase,
664  "PSD",
665  freq1[freqI],
666  noiseFFT::PSD(surfPSDf[i + f0]),
667  procFaceOffset,
668  writePSD_
669  );
671  (
672  outDir,
673  fNameBase,
674  "SPL",
675  freq1[freqI],
676  noiseFFT::SPL(surfPSDf[i + f0]*deltaf),
677  procFaceOffset,
678  writeSPL_
679  );
680  }
681  }
682 
683  if (Pstream::master())
684  {
685  graph Prmsfg
686  (
687  "Average Prms(f)",
688  "f [Hz]",
689  "P(f) [Pa]",
690  fOut,
691  PrmsfAve
692  );
693  Prmsfg.write
694  (
695  outDir,
696  graph::wordify(Prmsfg.title()),
698  );
699 
700  graph PSDfg
701  (
702  "Average PSD_f(f)",
703  "f [Hz]",
704  "PSD(f) [PaPa_Hz]",
705  fOut,
706  PSDfAve
707  );
708  PSDfg.write
709  (
710  outDir,
711  graph::wordify(PSDfg.title()),
713  );
714 
715  graph PSDg
716  (
717  "Average PSD_dB_Hz(f)",
718  "f [Hz]",
719  "PSD(f) [dB_Hz]",
720  fOut,
721  noiseFFT::PSD(PSDfAve)
722  );
723  PSDg.write
724  (
725  outDir,
726  graph::wordify(PSDg.title()),
728  );
729 
730  graph SPLg
731  (
732  "Average SPL_dB(f)",
733  "f [Hz]",
734  "SPL(f) [dB]",
735  fOut,
736  noiseFFT::SPL(PSDfAve*deltaf)
737  );
738  SPLg.write
739  (
740  outDir,
741  graph::wordify(SPLg.title()),
743  );
744  }
745  }
746 
747 
748  Info<< "Writing one-third octave surface data" << endl;
749  {
750  fileName outDir(outDirBase/"oneThirdOctave");
751 
752  scalarField PSDfAve(surfPrms13f.size(), Zero);
753  scalarField Prms13fAve(surfPrms13f.size(), Zero);
754 
755  forAll(surfPrms13f, i)
756  {
758  (
759  outDir,
760  fNameBase,
761  "SPL13",
762  octave13FreqCentre[i],
763  noiseFFT::SPL(surfPrms13f[i]),
764  procFaceOffset,
766  );
767 
768  Prms13fAve[i] =
769  surfaceAverage(surfPrms13f[i], procFaceOffset);
770  }
771 
772  if (Pstream::master())
773  {
774  graph SPL13g
775  (
776  "Average SPL13_dB(fm)",
777  "fm [Hz]",
778  "SPL(fm) [dB]",
779  octave13FreqCentre,
780  noiseFFT::SPL(Prms13fAve)
781  );
782  SPL13g.write
783  (
784  outDir,
785  graph::wordify(SPL13g.title()),
787  );
788  }
789  }
790  }
791 }
792 
793 
794 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
795 
796 } // End namespace noiseModels
797 } // End namespace Foam
798 
799 // ************************************************************************* //
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:66
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:189
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::noiseFFT::octaveBandInfo
static void octaveBandInfo(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.
Definition: noiseFFT.C:75
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::noiseModel
Base class for noise models.
Definition: noiseModel.H:158
Foam::noiseFFT
Performs FFT of pressure field to generate noise data.
Definition: noiseFFT.H:69
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:523
Foam::graph::wordify
static word wordify(const string &sname)
Helper function to convert string name into appropriate word.
Definition: graph.C:48
Foam::noiseModels::surfaceNoise::surfaceAverage
scalar surfaceAverage(const scalarField &data, const labelList &procFaceOffset) const
Calculate the area average value.
Definition: surfaceNoise.C:340
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:177
Foam::noiseModels::surfaceNoise::times_
scalarList times_
Sample times.
Definition: surfaceNoise.H:145
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:427
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:142
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:415
noiseFFT.H
Foam::MeshedSurface::surfFaces
const List< Face > & surfFaces() const
Return const access to the faces.
Definition: MeshedSurface.H:411
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:468
Foam::noiseModel::findStartTimeIndex
label findStartTimeIndex(const instantList &allTimes, const scalar startTime) const
Find and return start time index.
Definition: noiseModel.C:125
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:435
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:174
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:228
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:405
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::noiseFFT::PSD
static tmp< scalarField > PSD(const scalarField &PSDf)
Return the PSD [dB/Hz].
Definition: noiseFFT.C:62
Foam::noiseModels::surfaceNoise::~surfaceNoise
virtual ~surfaceNoise()
Destructor.
Definition: surfaceNoise.C:429
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:186
Foam::noiseModel::writeOctaves_
bool writeOctaves_
Write writeOctaves; default = yes.
Definition: noiseModel.H:219
argList.H
Foam::noiseModel::startTime_
scalar startTime_
Start time, default = 0s.
Definition: noiseModel.H:183
Foam::noiseModel::writePSD_
bool writePSD_
Write PSD; default = yes.
Definition: noiseModel.H:213
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::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
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::noiseFFT::frequencies
static tmp< scalarField > frequencies(const label N, const scalar deltaT)
Return the FFT frequencies.
Definition: noiseFFT.C:44
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::noiseModel::writePSDf_
bool writePSDf_
Write PSDf; default = yes.
Definition: noiseModel.H:216
Foam::noiseModel::writePrmsf_
bool writePrmsf_
Write Prmsf; default = yes.
Definition: noiseModel.H:207
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::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:445
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:439
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:372
Foam::nl
constexpr char nl
Definition: Ostream.H:385
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:144
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::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:119
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:189
Foam::MeshedSurface< face >
Foam::noiseModels::surfaceNoise::initialise
void initialise(const fileName &fName)
Initialise.
Definition: surfaceNoise.C:50
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:298
Foam::noiseModel::writeSPL_
bool writeSPL_
Write SPL; default = yes.
Definition: noiseModel.H:210
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:54
Foam::fileName::isAbsolute
static bool isAbsolute(const std::string &str)
Return true if string starts with a '/'.
Definition: fileNameI.H:136
Foam::noiseModel::nSamples_
label nSamples_
Number of samples in sampling window, default = 2^16.
Definition: noiseModel.H:171
Foam::noiseFFT::SPL
static tmp< scalarField > SPL(const scalarField &Prms2)
Return the SPL [dB].
Definition: noiseFFT.C:68
Foam::noiseModels::surfaceNoise::writerPtr_
autoPtr< surfaceWriter > writerPtr_
Pointer to the surface writer.
Definition: surfaceNoise.H:168