pointNoise.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-2018 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 "pointNoise.H"
29 #include "noiseFFT.H"
30 #include "argList.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace noiseModels
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
44 
45 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
46 
48 (
49  const scalarField& t0,
50  const scalarField& p0,
51  scalarField& t,
52  scalarField& p
53 ) const
54 {
55  DynamicList<scalar> tf(t0.size());
56  DynamicList<scalar> pf(t0.size());
57 
58  forAll(t0, timeI)
59  {
60  if (t0[timeI] >= startTime_)
61  {
62  tf.append(t0[timeI]);
63  pf.append(p0[timeI]);
64  }
65  }
66 
67  t.transfer(tf);
68  p.transfer(pf);
69 }
70 
71 
73 (
74  const label dataseti,
76 )
77 {
78  Info<< "Reading data file " << data.fName() << endl;
79 
80  const word fNameBase = data.fName().nameLessExt();
81 
82  // Time and pressure history data
83  scalarField t, p;
84  filterTimeData(data.x(), data.y(), t, p);
85  p *= rhoRef_;
86 
87  Info<< " read " << t.size() << " values" << nl << endl;
88 
89  Info<< "Creating noise FFT" << endl;
90 
91  const scalar deltaT = checkUniformTimeStep(t);
92 
93  if (!validateBounds(p))
94  {
95  Info<< "No noise data generated" << endl;
96  return;
97  }
98 
99 
100  // Determine the windowing
101  windowModelPtr_->validate(t.size());
102  const windowModel& win = windowModelPtr_();
103  const scalar deltaf = 1.0/(deltaT*win.nSamples());
104  fileName outDir(baseFileDir(dataseti)/fNameBase);
105 
106  // Create the fft
107  noiseFFT nfft(deltaT, win.nSamples());
108  nfft.setData(p);
109 
110 
111  // Narrow band data
112  // ----------------
113 
114  // RMS pressure [Pa]
115  graph Prmsf(nfft.RMSmeanPf(win));
116  if (customBounds_)
117  {
118  Prmsf.setXRange(fLower_, fUpper_);
119  }
120  if (writePrmsf_)
121  {
122  Info<< " Creating graph for " << Prmsf.title() << endl;
123  Prmsf.write(outDir, graph::wordify(Prmsf.title()), graphFormat_);
124  }
125 
126  // PSD [Pa^2/Hz]
127  graph PSDf(nfft.PSDf(win));
128  if (customBounds_)
129  {
130  PSDf.setXRange(fLower_, fUpper_);
131  }
132  if (writePSDf_)
133  {
134  Info<< " Creating graph for " << PSDf.title() << endl;
135  PSDf.write(outDir, graph::wordify(PSDf.title()), graphFormat_);
136  }
137 
138  // PSD [dB/Hz]
139  graph PSDg
140  (
141  "PSD_dB_Hz(f)",
142  "f [Hz]",
143  "PSD(f) [dB_Hz]",
144  Prmsf.x(),
145  noiseFFT::PSD(PSDf.y())
146  );
147 
148  if (writePSD_)
149  {
150  Info<< " Creating graph for " << PSDg.title() << endl;
151  PSDg.write(outDir, graph::wordify(PSDg.title()), graphFormat_);
152  }
153 
154  // SPL [dB]
155  graph SPLg
156  (
157  "SPL_dB(f)",
158  "f [Hz]",
159  "SPL(f) [dB]",
160  Prmsf.x(),
161  noiseFFT::SPL(PSDf.y()*deltaf)
162  );
163 
164  if (writeSPL_)
165  {
166  Info<< " Creating graph for " << SPLg.title() << endl;
167  SPLg.write(outDir, graph::wordify(SPLg.title()), graphFormat_);
168  }
169 
170  if (writeOctaves_)
171  {
172  labelList octave13BandIDs;
173  scalarField octave13FreqCentre;
175  (
176  Prmsf.x(),
177  fLower_,
178  fUpper_,
179  3,
180  octave13BandIDs,
181  octave13FreqCentre
182  );
183 
184 
185  // 1/3 octave data
186  // ---------------
187 
188  // Integrated PSD = P(rms)^2 [Pa^2]
189  graph Prms13f(nfft.octaves(PSDf, octave13BandIDs));
190 
191  graph SPL13g
192  (
193  "SPL13_dB(fm)",
194  "fm [Hz]",
195  "SPL(fm) [dB]",
196  octave13FreqCentre,
197  noiseFFT::SPL(Prms13f.y())
198  );
199  Info<< " Creating graph for " << SPL13g.title() << endl;
200  SPL13g.write(outDir, graph::wordify(SPL13g.title()), graphFormat_);
201  }
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
206 
208 :
209  noiseModel(dict, false)
210 {
211  if (readFields)
212  {
213  read(dict);
214  }
215 }
216 
217 
218 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
219 
221 {
222  // Point data only handled by master
223  if (!Pstream::master())
224  {
225  return;
226  }
227 
228  forAll(inputFileNames_, filei)
229  {
230  fileName fName = inputFileNames_[filei];
231  fName.expand();
232 
233  if (!fName.isAbsolute())
234  {
235  fName = argList::envGlobalPath()/fName;
236  }
237 
238  Function1Types::CSV<scalar> data("pressure", dict_, fName);
239  processData(filei, data);
240  }
241 }
242 
243 
245 {
246  if (noiseModel::read(dict))
247  {
248  if (!dict.readIfPresent("files", inputFileNames_))
249  {
251 
252  // Note: lookup uses same keyword as used by the CSV constructor
253  dict.readEntry("file", inputFileNames_.first());
254  }
255 
256  return true;
257  }
258 
259  return false;
260 }
261 
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 } // End namespace noiseModels
266 } // End namespace Foam
267 
268 // ************************************************************************* //
Foam::noiseFFT::setData
void setData(scalarList &data)
Definition: noiseFFT.C:189
Foam::graph::title
const string & title() const
Definition: graph.H:152
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::noiseModel::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: noiseModel.C:189
Foam::graph
Class to create, store and output qgraph files.
Definition: graph.H:61
Foam::noiseModels::pointNoise::pointNoise
pointNoise(const dictionary &dict, const bool readFields=true)
Constructor.
Definition: pointNoise.C:207
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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::argList::envGlobalPath
static fileName envGlobalPath()
Global case (directory) from environment variable.
Definition: argList.C:513
Foam::graph::wordify
static word wordify(const string &sname)
Helper function to convert string name into appropriate word.
Definition: graph.C:48
Foam::DynamicList< scalar >
Foam::graph::setXRange
void setXRange(const scalar x0, const scalar x1)
Definition: graph.C:171
Foam::Function1Types::CSV
Templated CSV function.
Definition: CSV.H:76
Foam::noiseModels::pointNoise::inputFileNames_
List< fileName > inputFileNames_
Input file names - optional.
Definition: pointNoise.H:104
noiseFFT.H
Foam::noiseModels::defineTypeNameAndDebug
defineTypeNameAndDebug(pointNoise, 0)
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
pointNoise.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::noiseFFT::PSD
static tmp< scalarField > PSD(const scalarField &PSDf)
Return the PSD [dB/Hz].
Definition: noiseFFT.C:62
Foam::noiseModel::dict_
const dictionary dict_
Copy of dictionary used for construction.
Definition: noiseModel.H:165
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
argList.H
Foam::noiseModels::pointNoise::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: pointNoise.C:244
Foam::List::resize
void resize(const label newSize)
Adjust allocated size of list.
Definition: ListI.H:139
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::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::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
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::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::noiseModels::pointNoise
Perform noise analysis on point-based pressure data.
Definition: pointNoise.H:94
Foam::List< label >
Foam::string::expand
string & expand(const bool allowEmpty=false)
Definition: string.C:159
Foam::messageStream::title
const string & title() const
Return the title of this error type.
Definition: messageStream.H:129
Foam::windowModel::nSamples
label nSamples() const
Return the number of samples in the window.
Definition: windowModel.C:56
Foam::noiseModels::pointNoise::calculate
virtual void calculate()
Calculate.
Definition: pointNoise.C:220
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::noiseModels::pointNoise::processData
void processData(const label dataseti, const Function1Types::CSV< scalar > &data)
Process the CSV data.
Definition: pointNoise.C:73
Foam::windowModel
Base class for windowing models.
Definition: windowModel.H:52
Foam::noiseModels::pointNoise::filterTimeData
void filterTimeData(const scalarField &t0, const scalarField &p0, scalarField &t, scalarField &p) const
Definition: pointNoise.C:48
Foam::noiseModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(noiseModel, pointNoise, dictionary)
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::noiseFFT::SPL
static tmp< scalarField > SPL(const scalarField &Prms2)
Return the SPL [dB].
Definition: noiseFFT.C:68