noise.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  noise
26 
27 Description
28  Utility to perform noise analysis of pressure data using the noiseFFT
29  library.
30 
31  Control settings are read from the $FOAM_CASE/system/noiseDict dictionary,
32  or user-specified dictionary using the -dict option. Pressure data is
33  read using a Table Function1:
34 
35 Usage
36  \verbatim
37  pRef 101325;
38  N 65536;
39  nw 100;
40  f1 25;
41  fU 10000;
42  graphFormat raw;
43 
44  pressureData
45  {
46  file "pressureData";
47  nHeaderLine 1; // number of header lines
48  columns (0 1); // column indices
49  separator " "; // optional (defaults to ",")
50  mergeSeparators no; // merge multiple separators
51  outOfBounds clamp; // optional out-of-bounds handling
52  interpolationScheme linear; // optional interpolation scheme
53  }
54  \endverbatim
55 
56  where
57  \table
58  Property | Description | Required | Default value
59  pRef | Reference pressure | no | 0
60  N | Number of samples in sampling window | no | 65536
61  nw | Number of sampling windows | no | 100
62  fl | Lower frequency band | no | 25
63  fU | Upper frequency band | no | 10000
64  graphFormat | Output graph format | no | raw
65  \endtable
66 
67  Current graph outputs include:
68  - FFT of the pressure data
69  - narrow-band PFL (pressure-fluctuation level) spectrum
70  - one-third-octave-band PFL spectrum
71  - one-third-octave-band pressure spectrum
72 
73 See also
74  Table.H
75  noiseFFT.H
76 
77 \*---------------------------------------------------------------------------*/
78 
79 #include "noiseFFT.H"
80 #include "argList.H"
81 #include "Time.H"
82 #include "Table.H"
83 #include "systemDict.H"
84 #include "setWriter.H"
85 #include "writeFile.H"
86 
87 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
88 
89 using namespace Foam;
90 
91 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
92 
93 Foam::scalar checkUniformTimeStep(const scalarField& t)
94 {
95  // check that a uniform time step has been applied
96  scalar deltaT = -1.0;
97  if (t.size() > 1)
98  {
99  for (label i = 1; i < t.size(); i++)
100  {
101  const scalar dT = t[i] - t[i-1];
102  if (deltaT < 0)
103  {
104  deltaT = dT;
105  }
106 
107  if (mag(deltaT - dT) > rootSmall)
108  {
110  << "Unable to process data with a variable time step"
111  << exit(FatalError);
112  }
113  }
114  }
115  else
116  {
118  << "Unable to create FFT with a single value"
119  << exit(FatalError);
120  }
121 
122  return deltaT;
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
127 
128 int main(int argc, char *argv[])
129 {
131  #include "addDictOption.H"
132  #include "setRootCase.H"
133  #include "createTime.H"
134 
135  const word dictName("noiseDict");
136 
137  IOdictionary dict(systemDict(dictName, args, runTime));
138 
139  // Reference pressure
140  const scalar pRef = dict.lookupOrDefault("pRef", 0.0);
141 
142  // Number of samples in sampling window
143  const label N = dict.lookupOrDefault("N", 65536);
144 
145  // Number of sampling windows
146  const label nw = dict.lookupOrDefault("nw", 100);
147 
148  // Lower frequency of frequency band
149  const scalar f1 = dict.lookupOrDefault("f1", 25.0);
150 
151  // Upper frequency of frequency band
152  const scalar fU = dict.lookupOrDefault("fU", 10000.0);
153 
154  // Graph format
155  const word graphFormat
156  (
157  dict.lookupOrDefault<word>("graphFormat", "raw")
158  );
159 
160  Info<< "Reading data file" << endl;
162  (
163  "pressure",
164  dimTime,
165  dimPressure,
166  dict.subDict("pressureData")
167  );
168 
169  // time history data
170  const scalarField t(pData.x());
171 
172  // pressure data
173  const scalarField p(pData.y());
174 
175  if (t.size() < N)
176  {
178  << "Block size N = " << N
179  << " is larger than number of data = " << t.size()
180  << exit(FatalError);
181  }
182 
183  Info<< " read " << t.size() << " values" << nl << endl;
184 
185 
186  Info<< "Creating noise FFT" << endl;
187  noiseFFT nfft(checkUniformTimeStep(t), p);
188 
189  nfft -= pRef;
190 
191  const fileName pDateFileName(dict.subDict("pressureData").lookup("file"));
192  const fileName baseFileName(pDateFileName.lessExt());
193  const fileName outputPath
194  (
195  runTime.path()
197  /baseFileName
198  );
199 
200  autoPtr<setWriter> writer(setWriter::New(graphFormat));
201 
202  const Pair<scalarField> Pf(nfft.RMSmeanPf(N, min(nfft.size()/N, nw)));
203  Info<< " Creating graph for P(f)" << endl;
204  writer->write
205  (
206  outputPath,
207  "Pf",
208  coordSet(true, "f [Hz]", Pf.first()),
209  "P(f) [Pa]",
210  Pf.second()
211  );
212 
213  const Pair<scalarField> Lf(nfft.Lf(Pf));
214  Info<< " Creating graph for L(f)" << endl;
215  writer->write
216  (
217  outputPath,
218  "Lf",
219  coordSet(true, "f [Hz]", Lf.first()),
220  "L(f) [dB]",
221  Lf.second()
222  );
223 
224  const Pair<scalarField> Ldelta(nfft.Ldelta(Lf, f1, fU));
225  Info<< " Creating graph for Ldelta" << endl;
226  writer->write
227  (
228  outputPath,
229  "Ldelta",
230  coordSet(true, "fm [Hz]", Ldelta.first()),
231  "Ldelta(f) [dB]",
232  Ldelta.second()
233  );
234 
235  const Pair<scalarField> Pdelta(nfft.Pdelta(Pf, f1, fU));
236  Info<< " Creating graph for Pdelta" << endl;
237  writer->write
238  (
239  outputPath,
240  "Pdelta",
241  coordSet(true, "fm [Hz]", Pdelta.first()),
242  "P(f) [dB]",
243  Pdelta.second()
244  );
245 
246  Info<< nl << "End\n" << endl;
247 
248  return 0;
249 }
250 
251 
252 // ************************************************************************* //
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual Ostream & write(const char)=0
Write character.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:65
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Holds list of sampling positions.
Definition: coordSet.H:51
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
A class for handling file names.
Definition: fileName.H:82
static const word outputPrefix
Directory prefix.
Definition: writeFile.H:72
FFT of the pressure field.
Definition: noiseFFT.H:51
static autoPtr< setWriter > New(const word &writeType, const IOstream::streamFormat writeFormat=IOstream::ASCII, const IOstream::compressionType writeCompression=IOstream::UNCOMPRESSED)
Select given write options.
Definition: setWriter.C:286
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
int main(int argc, char *argv[])
Definition: financialFoam.C:44
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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:59
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
const dimensionSet dimPressure
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
static const char nl
Definition: Ostream.H:266
dictionary dict
Foam::argList args(argc, argv)
volScalarField & p