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-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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 "argList.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace noiseModels
37{
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
43
44// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
45
47(
48 const scalarField& t0,
49 const scalarField& p0,
50 scalarField& t,
52) const
53{
56
57 forAll(t0, timeI)
58 {
59 if (t0[timeI] >= startTime_)
60 {
61 tf.append(t0[timeI]);
62 pf.append(p0[timeI]);
63 }
64 }
65
66 t.transfer(tf);
67 p.transfer(pf);
68}
69
70
72(
73 const label dataseti,
75)
76{
77 Info<< "Reading data file " << data.fName() << endl;
78
79 const word fNameBase(data.fName().nameLessExt());
80
81 // Time and pressure history data
82 scalarField t, p;
83 filterTimeData(data.x(), data.y(), t, p);
84
85 Info<< " read " << t.size() << " values" << nl << endl;
86
87 if (!validateBounds(p))
88 {
89 Info<< "No noise data generated" << endl;
90 return;
91 }
92
93 Info<< "Creating noise FFT" << endl;
94
95 const scalar deltaT = checkUniformTimeStep(t);
96
97 // Apply conversions
98 p *= rhoRef_;
99 p -= average(p);
100
101 // Determine the windowing
102 windowModelPtr_->validate(t.size());
103 const windowModel& win = windowModelPtr_();
104 const scalar deltaf = 1.0/(deltaT*win.nSamples());
105 const fileName outDir(baseFileDir(dataseti)/fNameBase);
106
107
108 // Narrow band data
109 // ----------------
110
111 scalarField f(uniformFrequencies(deltaT, true));
112
113 // RMS pressure [Pa]
114 if (writePrmsf_)
115 {
116 graph g
117 (
118 "Prms(f)",
119 "f [Hz]",
120 "Prms(f) [Pa]",
121 f,
122 RMSmeanPf(p)
123 );
124
125 g.setXRange(fLower_, fUpper_);
126
127 Info<< " Creating graph for " << g.title() << endl;
128 g.write(outDir, graph::wordify(g.title()), graphFormat_);
129 }
130
131 // PSD [Pa^2/Hz]
132 const scalarField PSDf(this->PSDf(p, deltaT));
133
134 if (writePSDf_)
135 {
136 graph g
137 (
138 "PSD(f)",
139 "f [Hz]",
140 "PSD(f) [PaPa_Hz]",
141 f,
142 PSDf
143 );
144
145 g.setXRange(fLower_, fUpper_);
146
147 Info<< " Creating graph for " << g.title() << endl;
148 g.write(outDir, graph::wordify(g.title()), graphFormat_);
149 }
150
151 // PSD [dB/Hz]
152 if (writePSD_)
153 {
154 graph g
155 (
156 "PSD_dB_Hz(f)",
157 "f [Hz]",
158 "PSD(f) [dB_Hz]",
159 f,
160 PSD(PSDf)
161 );
162
163 g.setXRange(fLower_, fUpper_);
164
165 Info<< " Creating graph for " << g.title() << endl;
166 g.write(outDir, graph::wordify(g.title()), graphFormat_);
167 }
168
169 // SPL [dB]
170 if (writeSPL_)
171 {
172 graph g
173 (
174 "SPL_dB(f)",
175 "f [Hz]",
176 "SPL(f) [" + weightingTypeNames_[SPLweighting_] + "]",
177 f,
178 SPL(PSDf*deltaf, f)
179 );
180
181 g.setXRange(fLower_, fUpper_);
182
183 Info<< " Creating graph for " << g.title() << endl;
184 g.write(outDir, graph::wordify(g.title()), graphFormat_);
185 }
186
187 if (writeOctaves_)
188 {
189 labelList octave13BandIDs;
190 scalarField octave13FreqCentre;
192 (
193 f,
194 fLower_,
195 fUpper_,
196 3,
197 octave13BandIDs,
198 octave13FreqCentre
199 );
200
201
202 // 1/3 octave data
203 // ---------------
204
205 // Integrated PSD = P(rms)^2 [Pa^2]
206 scalarField Prms13f(octaves(PSDf, f, octave13BandIDs));
207
208 graph g
209 (
210 "SPL13_dB(fm)",
211 "fm [Hz]",
212 "SPL(fm) [" + weightingTypeNames_[SPLweighting_] + "]",
213 octave13FreqCentre,
214 SPL(Prms13f, octave13FreqCentre)
215 );
216
217 Info<< " Creating graph for " << g.title() << endl;
218 g.write(outDir, graph::wordify(g.title()), graphFormat_);
219 }
220}
221
222
223// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
224
226:
227 noiseModel(dict, false)
228{
229 if (readFields)
230 {
231 read(dict);
232 }
233}
234
235
236// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
237
239{
240 // Point data only handled by master
241 if (!Pstream::master())
242 {
243 return;
244 }
245
246 forAll(inputFileNames_, filei)
247 {
248 fileName fName = inputFileNames_[filei];
249 fName.expand();
250
251 if (!fName.isAbsolute())
252 {
253 fName = argList::envGlobalPath()/fName;
254 }
255
256 Function1Types::CSV<scalar> data("pressure", dict_, nullptr, fName);
257 processData(filei, data);
258 }
259}
260
261
263{
265 {
266 if (!dict.readIfPresent("files", inputFileNames_))
267 {
269
270 // Note: lookup uses same keyword as used by the CSV constructor
272 }
273
274 return true;
275 }
276
277 return false;
278}
279
280
281// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282
283} // End namespace noiseModels
284} // End namespace Foam
285
286// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const uniformDimensionedVectorField & g
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
Templated CSV function.
Definition: CSV.H:78
void transfer(List< T > &list)
Definition: List.C:447
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
virtual bool read()
Re-read model coefficients if they have changed.
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static fileName envGlobalPath()
Global case (directory) from environment variable.
Definition: argList.C:549
Database for solution data, solver performance and other reduced data.
Definition: data.H:58
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
A class for handling file names.
Definition: fileName.H:76
static bool isAbsolute(const std::string &str)
Definition: fileNameI.H:136
Class to create, store and output qgraph files.
Definition: graph.H:62
static word wordify(const string &sname)
Helper function to convert string name into appropriate word.
Definition: graph.C:47
Base class for noise models.
Definition: noiseModel.H:176
tmp< scalarField > uniformFrequencies(const scalar deltaT, const bool check) const
Definition: noiseModel.C:277
scalar rhoRef_
Reference density (to convert from kinematic to static pressure)
Definition: noiseModel.H:233
tmp< scalarField > octaves(const scalarField &data, const scalarField &f, const labelUList &freqBandIDs) const
Generate octave data.
Definition: noiseModel.C:314
tmp< scalarField > RMSmeanPf(const scalarField &p) const
Definition: noiseModel.C:429
scalar fUpper_
Upper frequency limit, default = 10kHz.
Definition: noiseModel.H:242
tmp< Foam::scalarField > PSD(const scalarField &PSDf) const
PSD [dB/Hz].
Definition: noiseModel.C:713
bool writePSD_
Write PSD; default = yes.
Definition: noiseModel.H:281
word graphFormat_
Graph format.
Definition: noiseModel.H:251
static void setOctaveBands(const scalarField &f, const scalar fLower, const scalar fUpper, const scalar octave, labelList &fBandIDs, scalarField &fCentre)
Definition: noiseModel.C:55
fileName baseFileDir(const label dataseti) const
Return the base output directory.
Definition: noiseModel.C:262
bool writeOctaves_
Write writeOctaves; default = yes.
Definition: noiseModel.H:287
tmp< scalarField > SPL(const scalarField &Prms2, const scalar f) const
SPL [dB].
Definition: noiseModel.C:722
autoPtr< windowModel > windowModelPtr_
Window model.
Definition: noiseModel.H:248
static const Enum< weightingType > weightingTypeNames_
Definition: noiseModel.H:188
bool writeSPL_
Write SPL; default = yes.
Definition: noiseModel.H:278
const dictionary dict_
Copy of dictionary used for construction.
Definition: noiseModel.H:230
bool writePrmsf_
Write Prmsf; default = yes.
Definition: noiseModel.H:275
weightingType SPLweighting_
Weighting.
Definition: noiseModel.H:254
scalar checkUniformTimeStep(const scalarList &times) const
Check and return uniform time step.
Definition: noiseModel.C:186
bool writePSDf_
Write PSDf; default = yes.
Definition: noiseModel.H:284
tmp< scalarField > PSDf(const scalarField &p, const scalar deltaT) const
Definition: noiseModel.C:448
scalar startTime_
Start time, default = 0s.
Definition: noiseModel.H:245
bool validateBounds(const scalarList &p) const
Return true if all pressure data is within min/max bounds.
Definition: noiseModel.C:220
scalar fLower_
Lower frequency limit, default = 25Hz.
Definition: noiseModel.H:239
Perform noise analysis on point-based pressure data.
Definition: pointNoise.H:97
void filterTimeData(const scalarField &t0, const scalarField &p0, scalarField &t, scalarField &p) const
Definition: pointNoise.C:47
void processData(const label dataseti, const Function1Types::CSV< scalar > &data)
Process the CSV data.
Definition: pointNoise.C:72
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: pointNoise.C:262
List< fileName > inputFileNames_
Input file names - optional.
Definition: pointNoise.H:104
virtual void calculate()
Calculate.
Definition: pointNoise.C:238
virtual bool write(const bool valid=true) const
Write using setting from DB.
splitCell * master() const
Definition: splitCell.H:113
string & expand(const bool allowEmpty=false)
Definition: string.C:173
Base class for windowing models.
Definition: windowModel.H:55
label nSamples() const
Return the number of samples in the window.
Definition: windowModel.C:56
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const volScalarField & p0
Definition: EEqn.H:36
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333