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-2018 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 CSV reader:
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  CSV.H
76  noiseFFT.H
77 
78 \*---------------------------------------------------------------------------*/
79 
80 
81 #include "noiseFFT.H"
82 #include "argList.H"
83 #include "Time.H"
84 #include "CSV.H"
85 #include "IOdictionary.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  scalar dT = t[i] - t[i-1];
102  if (deltaT < 0)
103  {
104  deltaT = dT;
105  }
106 
107  if (mag(deltaT - dT) > small)
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  #include "createFields.H"
135 
136  Info<< "Reading data file" << endl;
138  (
139  "pressure",
140  dict.subDict("pressureData")
141  );
142 
143  // time history data
144  const scalarField t(pData.x());
145 
146  // pressure data
147  const scalarField p(pData.y());
148 
149  if (t.size() < N)
150  {
152  << "Block size N = " << N
153  << " is larger than number of data = " << t.size()
154  << exit(FatalError);
155  }
156 
157  Info<< " read " << t.size() << " values" << nl << endl;
158 
159 
160  Info<< "Creating noise FFT" << endl;
161  noiseFFT nfft(checkUniformTimeStep(t), p);
162 
163  nfft -= pRef;
164 
165  fileName baseFileName(pData.fName().lessExt());
166 
167  graph Pf(nfft.RMSmeanPf(N, min(nfft.size()/N, nw)));
168  Info<< " Creating graph for " << Pf.title() << endl;
169  Pf.write(baseFileName + graph::wordify(Pf.title()), graphFormat);
170 
171  graph Lf(nfft.Lf(Pf));
172  Info<< " Creating graph for " << Lf.title() << endl;
173  Lf.write(baseFileName + graph::wordify(Lf.title()), graphFormat);
174 
175  graph Ldelta(nfft.Ldelta(Lf, f1, fU));
176  Info<< " Creating graph for " << Ldelta.title() << endl;
177  Ldelta.write(baseFileName + graph::wordify(Ldelta.title()), graphFormat);
178 
179  graph Pdelta(nfft.Pdelta(Pf, f1, fU));
180  Info<< " Creating graph for " << Pdelta.title() << endl;
181  Pdelta.write(baseFileName + graph::wordify(Pdelta.title()), graphFormat);
182 
183  Info<< nl << "End\n" << endl;
184 
185  return 0;
186 }
187 
188 
189 // ************************************************************************* //
dictionary dict
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
A class for handling file names.
Definition: fileName.H:69
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
word graphFormat
Definition: createFields.H:34
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
const string & title() const
Return the title of this error type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
static void noParallel()
Remove the parallel options.
Definition: argList.C:167
scalar f1
Definition: createFields.H:28
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
Class to create, store and output qgraph files.
Definition: graph.H:58
static word wordify(const string &sname)
Helper function to convert string name into appropriate word.
Definition: graph.C:44
static const char nl
Definition: Ostream.H:265
scalar pRef
Definition: createFields.H:19
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
messageStream Info
label N
Definition: createFields.H:22
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar fU
Definition: createFields.H:31
virtual Ostream & write(const token &)=0
Write next token to stream.
volScalarField & p
FFT of the pressure field.
Definition: noiseFFT.H:49
label nw
Definition: createFields.H:25
Namespace for OpenFOAM.