noiseModel.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-2018 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 "noiseModel.H"
29 #include "functionObject.H"
30 #include "argList.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(noiseModel, 0);
37  defineRunTimeSelectionTable(noiseModel, dictionary);
38 }
39 
40 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
41 
43 (
44  const dictionary& dict,
45  const word& lookup,
46  bool& option
47 ) const
48 {
49  dict.readIfPresent(lookup, option);
50 
51  // Only writing on the master process
52  option = option && Pstream::master();
53 
54  if (option)
55  {
56  Info<< " " << lookup << ": " << "yes" << endl;
57  }
58  else
59  {
60  Info<< " " << lookup << ": " << "no" << endl;
61  }
62 }
63 
64 
66 (
67  const scalarList& times
68 ) const
69 {
70  scalar deltaT = -1.0;
71 
72  if (times.size() > 1)
73  {
74  for (label i = 1; i < times.size(); i++)
75  {
76  scalar dT = times[i] - times[i-1];
77 
78  if (deltaT < 0)
79  {
80  deltaT = dT;
81  }
82 
83  if (mag(dT/deltaT - 1) > 1e-8)
84  {
86  << "Unable to process data with a variable time step"
87  << exit(FatalError);
88  }
89  }
90  }
91  else
92  {
94  << "Unable to create FFT with a single value"
95  << exit(FatalError);
96  }
97 
98  return deltaT;
99 }
100 
101 
103 {
104  forAll(p, i)
105  {
106  if ((p[i] < minPressure_) || (p[i] > maxPressure_))
107  {
109  << "Pressure data at position " << i
110  << " is outside of permitted bounds:" << nl
111  << " pressure: " << p[i] << nl
112  << " minimum pressure: " << minPressure_ << nl
113  << " maximum pressure: " << maxPressure_ << nl
114  << endl;
115 
116  return false;
117  }
118  }
119 
120  return true;
121 }
122 
123 
125 (
126  const instantList& allTimes,
127  const scalar startTime
128 ) const
129 {
130  forAll(allTimes, timeI)
131  {
132  const instant& i = allTimes[timeI];
133 
134  if (i.value() >= startTime)
135  {
136  return timeI;
137  }
138  }
139 
140  return 0;
141 }
142 
143 
145 {
146  return
147  (
150  / "noise"
151  / outputPrefix_
152  / type()
153  / ("input" + Foam::name(dataseti))
154  );
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
159 
161 :
162  dict_(dict),
163  rhoRef_(1),
164  nSamples_(65536),
165  fLower_(25),
166  fUpper_(10000),
167  customBounds_(false),
168  startTime_(0),
169  windowModelPtr_(),
170  graphFormat_("raw"),
171  minPressure_(-0.5*VGREAT),
172  maxPressure_(0.5*VGREAT),
173  outputPrefix_(),
174  writePrmsf_(true),
175  writeSPL_(true),
176  writePSD_(true),
177  writePSDf_(true),
178  writeOctaves_(true)
179 {
180  if (readFields)
181  {
182  read(dict);
183  }
184 }
185 
186 
187 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
188 
190 {
191  dict.readIfPresent("rhoRef", rhoRef_);
192  dict.readIfPresent("N", nSamples_);
193  customBounds_ = false;
194  if (dict.readIfPresent("fl", fLower_))
195  {
196  customBounds_ = true;
197  }
198  if (dict.readIfPresent("fu", fUpper_))
199  {
200  customBounds_ = true;
201  }
202  dict.readIfPresent("startTime", startTime_);
203  dict.readIfPresent("graphFormat", graphFormat_);
204  dict.readIfPresent("minPressure", minPressure_);
205  dict.readIfPresent("maxPressure", maxPressure_);
206  dict.readIfPresent("outputPrefix", outputPrefix_);
207 
208  if (fLower_ < 0)
209  {
211  << "fl: lower frequency bound must be greater than zero"
212  << exit(FatalIOError);
213 
214  }
215 
216  if (fUpper_ < 0)
217  {
219  << "fu: upper frequency bound must be greater than zero"
220  << exit(FatalIOError);
221 
222  }
223 
224  if (fUpper_ < fLower_)
225  {
227  << "fu: upper frequency bound must be greater than lower "
228  << "frequency bound (fl)"
229  << exit(FatalIOError);
230 
231  }
232 
233  Info<< " Write options:" << endl;
234  dictionary optDict(dict.subOrEmptyDict("writeOptions"));
235  readWriteOption(optDict, "writePrmsf", writePrmsf_);
236  readWriteOption(optDict, "writeSPL", writeSPL_);
237  readWriteOption(optDict, "writePSD", writePSD_);
238  readWriteOption(optDict, "writePSDf", writePSDf_);
239  readWriteOption(optDict, "writeOctaves", writeOctaves_);
240 
241 
242  windowModelPtr_ = windowModel::New(dict, nSamples_);
243 
244  Info<< nl << endl;
245 
246  return true;
247 }
248 
249 
250 // ************************************************************************* //
Foam::noiseModel::checkUniformTimeStep
scalar checkUniformTimeStep(const scalarList &times) const
Check and return uniform time step.
Definition: noiseModel.C:66
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::noiseModel::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: noiseModel.C:189
Foam::noiseModel::readWriteOption
void readWriteOption(const dictionary &dict, const word &lookup, bool &option) const
Helper function to read write options and provide info feedback.
Definition: noiseModel.C:43
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::argList::envGlobalPath
static fileName envGlobalPath()
Global case (directory) from environment variable.
Definition: argList.C:513
Foam::noiseModel::noiseModel
noiseModel(const noiseModel &)=delete
No copy construct.
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::FatalIOError
IOerror FatalIOError
Foam::noiseModel::findStartTimeIndex
label findStartTimeIndex(const instantList &allTimes, const scalar startTime) const
Find and return start time index.
Definition: noiseModel.C:125
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::windowModel::New
static autoPtr< windowModel > New(const dictionary &dict, const label nSamples)
Return a reference to the selected window model.
Definition: windowModelNew.C:33
noiseModel.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::noiseModel::minPressure_
scalar minPressure_
Min pressure value.
Definition: noiseModel.H:195
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:62
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
argList.H
Foam::functionObject::outputPrefix
static word outputPrefix
Directory prefix.
Definition: functionObject.H:259
Foam::radiation::lookup
Lookup type of boundary radiation properties.
Definition: lookup.H:63
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::noiseModel::validateBounds
bool validateBounds(const scalarList &p) const
Return true if all pressure data is within min/max bounds.
Definition: noiseModel.C:102
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::noiseModel::maxPressure_
scalar maxPressure_
Min pressure value.
Definition: noiseModel.H:198
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:603
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::List< scalar >
Foam::noiseModel::baseFileDir
fileName baseFileDir(const label dataseti) const
Return the base output directory.
Definition: noiseModel.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Instant::value
scalar value() const
The value (const access)
Definition: Instant.H:99
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
startTime
Foam::label startTime
Definition: checkTimeOptions.H:1
Foam::instant
An instant of time. Contains the time value and name.
Definition: instant.H:52
functionObject.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:417