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-2023 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  refColumn 0; // reference column index
49  componentColumns (1); // component column indices
50  separator " "; // optional (defaults to ",")
51  mergeSeparators no; // merge multiple separators
52  outOfBounds clamp; // optional out-of-bounds handling
53  interpolationScheme linear; // optional interpolation scheme
54  }
55  \endverbatim
56 
57  where
58  \table
59  Property | Description | Required | Default value
60  pRef | Reference pressure | no | 0
61  N | Number of samples in sampling window | no | 65536
62  nw | Number of sampling windows | no | 100
63  fl | Lower frequency band | no | 25
64  fU | Upper frequency band | no | 10000
65  graphFormat | Output graph format | no | raw
66  \endtable
67 
68  Current graph outputs include:
69  - FFT of the pressure data
70  - narrow-band PFL (pressure-fluctuation level) spectrum
71  - one-third-octave-band PFL spectrum
72  - one-third-octave-band pressure spectrum
73 
74 See also
75  Table.H
76  noiseFFT.H
77 
78 \*---------------------------------------------------------------------------*/
79 
80 #include "noiseFFT.H"
81 #include "argList.H"
82 #include "Time.H"
83 #include "Table.H"
84 #include "systemDict.H"
85 #include "setWriter.H"
86 #include "writeFile.H"
87 
88 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
89 
90 using namespace Foam;
91 
92 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
93 
94 Foam::scalar checkUniformTimeStep(const scalarField& t)
95 {
96  // check that a uniform time step has been applied
97  scalar deltaT = -1.0;
98  if (t.size() > 1)
99  {
100  for (label i = 1; i < t.size(); i++)
101  {
102  const scalar dT = t[i] - t[i-1];
103  if (deltaT < 0)
104  {
105  deltaT = dT;
106  }
107 
108  if (mag(deltaT - dT) > rootSmall)
109  {
111  << "Unable to process data with a variable time step"
112  << exit(FatalError);
113  }
114  }
115  }
116  else
117  {
119  << "Unable to create FFT with a single value"
120  << exit(FatalError);
121  }
122 
123  return deltaT;
124 }
125 
126 
127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128 
129 int main(int argc, char *argv[])
130 {
132  #include "addDictOption.H"
133  #include "setRootCase.H"
134  #include "createTime.H"
135 
136  const word dictName("noiseDict");
137 
138  IOdictionary dict(systemDict(dictName, args, runTime));
139 
140  // Reference pressure
141  const scalar pRef = dict.lookupOrDefault("pRef", 0.0);
142 
143  // Number of samples in sampling window
144  const label N = dict.lookupOrDefault("N", 65536);
145 
146  // Number of sampling windows
147  const label nw = dict.lookupOrDefault("nw", 100);
148 
149  // Lower frequency of frequency band
150  const scalar f1 = dict.lookupOrDefault("f1", 25.0);
151 
152  // Upper frequency of frequency band
153  const scalar fU = dict.lookupOrDefault("fU", 10000.0);
154 
155  // Graph format
156  const word graphFormat
157  (
158  dict.lookupOrDefault<word>("graphFormat", "raw")
159  );
160 
161  Info<< "Reading data file" << endl;
163  (
164  "pressure",
165  dict.subDict("pressureData")
166  );
167 
168  // time history data
169  const scalarField t(pData.x());
170 
171  // pressure data
172  const scalarField p(pData.y());
173 
174  if (t.size() < N)
175  {
177  << "Block size N = " << N
178  << " is larger than number of data = " << t.size()
179  << exit(FatalError);
180  }
181 
182  Info<< " read " << t.size() << " values" << nl << endl;
183 
184 
185  Info<< "Creating noise FFT" << endl;
186  noiseFFT nfft(checkUniformTimeStep(t), p);
187 
188  nfft -= pRef;
189 
190  const fileName pDateFileName(dict.subDict("pressureData").lookup("file"));
191  const fileName baseFileName(pDateFileName.lessExt());
192  const fileName outputPath
193  (
194  runTime.path()
196  /baseFileName
197  );
198 
199  autoPtr<setWriter> writer(setWriter::New(graphFormat));
200 
201  const Pair<scalarField> Pf(nfft.RMSmeanPf(N, min(nfft.size()/N, nw)));
202  Info<< " Creating graph for P(f)" << endl;
203  writer->write
204  (
205  outputPath,
206  "Pf",
207  coordSet(true, "f [Hz]", Pf.first()),
208  "P(f) [Pa]",
209  Pf.second()
210  );
211 
212  const Pair<scalarField> Lf(nfft.Lf(Pf));
213  Info<< " Creating graph for L(f)" << endl;
214  writer->write
215  (
216  outputPath,
217  "Lf",
218  coordSet(true, "f [Hz]", Lf.first()),
219  "L(f) [dB]",
220  Lf.second()
221  );
222 
223  const Pair<scalarField> Ldelta(nfft.Ldelta(Lf, f1, fU));
224  Info<< " Creating graph for Ldelta" << endl;
225  writer->write
226  (
227  outputPath,
228  "Ldelta",
229  coordSet(true, "fm [Hz]", Ldelta.first()),
230  "Ldelta(f) [dB]",
231  Ldelta.second()
232  );
233 
234  const Pair<scalarField> Pdelta(nfft.Pdelta(Pf, f1, fU));
235  Info<< " Creating graph for Pdelta" << endl;
236  writer->write
237  (
238  outputPath,
239  "Pdelta",
240  coordSet(true, "fm [Hz]", Pdelta.first()),
241  "P(f) [dB]",
242  Pdelta.second()
243  );
244 
245  Info<< nl << "End\n" << endl;
246 
247  return 0;
248 }
249 
250 
251 // ************************************************************************* //
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,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
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:306
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
static const char nl
Definition: Ostream.H:260
dictionary dict
Foam::argList args(argc, argv)
volScalarField & p