noise.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 | Copyright (C) 2011-2015 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  \heading Usage
36 
37  \verbatim
38  pRef 101325;
39  N 65536;
40  nw 100;
41  f1 25;
42  fU 10000;
43  graphFormat raw;
44 
45  pressureData
46  {
47  fileName "pressureData"
48  nHeaderLine 1; // number of header lines
49  refColumn 0; // reference column index
50  componentColumns (1); // component column indices
51  separator " "; // optional (defaults to ",")
52  mergeSeparators no; // merge multiple separators
53  outOfBounds clamp; // optional out-of-bounds handling
54  interpolationScheme linear; // optional interpolation scheme
55  }
56  \endverbatim
57 
58  where
59  \table
60  Property | Description | Required | Default value
61  pRef | Reference pressure | no | 0
62  N | Number of samples in sampling window | no | 65536
63  nw | Number of sampling windows | no | 100
64  fl | Lower frequency band | no | 25
65  fU | Upper frequency band | no | 10000
66  graphFormat | Output graph format | no | raw
67  \endtable
68 
69  Current graph outputs include:
70  - FFT of the pressure data
71  - narrow-band PFL (pressure-fluctuation level) spectrum
72  - one-third-octave-band PFL spectrum
73  - one-third-octave-band pressure spectrum
74 
75 SeeAlso
76  CSV.H
77  noiseFFT.H
78 
79 \*---------------------------------------------------------------------------*/
80 
81 
82 #include "noiseFFT.H"
83 #include "argList.H"
84 #include "Time.H"
85 #include "functionObjectFile.H"
86 #include "CSV.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  scalar dT = t[i] - t[i-1];
103  if (deltaT < 0)
104  {
105  deltaT = dT;
106  }
107 
108  if (mag(deltaT - dT) > SMALL)
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  #include "createFields.H"
136 
137  Info<< "Reading data file" << endl;
138  CSV<scalar> pData("pressure", dict, "Data");
139 
140  // time history data
141  const scalarField t(pData.x());
142 
143  // pressure data
144  const scalarField p(pData.y());
145 
146  if (t.size() < N)
147  {
149  << "Block size N = " << N
150  << " is larger than number of data = " << t.size()
151  << exit(FatalError);
152  }
153 
154  Info<< " read " << t.size() << " values" << nl << endl;
155 
156 
157  Info<< "Creating noise FFT" << endl;
158  noiseFFT nfft(checkUniformTimeStep(t), p);
159 
160  nfft -= pRef;
161 
162  fileName baseFileName(pData.fName().lessExt());
163 
164  graph Pf(nfft.RMSmeanPf(N, min(nfft.size()/N, nw)));
165  Info<< " Creating graph for " << Pf.title() << endl;
166  Pf.write(baseFileName + graph::wordify(Pf.title()), graphFormat);
167 
168  graph Lf(nfft.Lf(Pf));
169  Info<< " Creating graph for " << Lf.title() << endl;
170  Lf.write(baseFileName + graph::wordify(Lf.title()), graphFormat);
171 
172  graph Ldelta(nfft.Ldelta(Lf, f1, fU));
173  Info<< " Creating graph for " << Ldelta.title() << endl;
174  Ldelta.write(baseFileName + graph::wordify(Ldelta.title()), graphFormat);
175 
176  graph Pdelta(nfft.Pdelta(Pf, f1, fU));
177  Info<< " Creating graph for " << Pdelta.title() << endl;
178  Pdelta.write(baseFileName + graph::wordify(Pdelta.title()), graphFormat);
179 
180  Info<< nl << "End\n" << endl;
181 
182  return 0;
183 }
184 
185 
186 // ************************************************************************* //
Foam::graph
Class to create, store and output qgraph files.
Definition: graph.H:58
p
p
Definition: pEqn.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::noiseFFT
FFT of the pressure field.
Definition: noiseFFT.H:49
Foam::graph::wordify
static word wordify(const string &sname)
Helper function to convert string name into appropriate word.
Definition: graph.C:44
noiseFFT.H
nw
label nw
Definition: createFields.H:25
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
addDictOption.H
Foam::CSV
Templated CSV container data entry. Reference column is always a scalar, e.g. time.
Definition: CSV.H:65
fU
scalar fU
Definition: createFields.H:31
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
argList.H
f1
scalar f1
Definition: createFields.H:28
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
pRef
scalar pRef
Definition: createFields.H:19
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
setRootCase.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
functionObjectFile.H
CSV.H
createTime.H
Foam::messageStream::title
const string & title() const
Return the title of this error type.
Definition: messageStream.H:118
Foam::argList::noParallel
static void noParallel()
Remove the parallel options.
Definition: argList.C:161
graphFormat
word graphFormat
Definition: createFields.H:34
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
N
label N
Definition: createFields.H:22