Base class for noise models. More...
Classes | |
struct | octaveBandInfo |
Octave band information. More... | |
struct | planInfo |
FFTW planner information. More... | |
Public Types | |
enum | weightingType { none, dBA, dBB, dBC, dBD } |
Public Member Functions | |
TypeName ("noiseModel") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, noiseModel, dictionary,(const dictionary &dict),(dict)) | |
Run time selection table. More... | |
noiseModel (const dictionary &dict, const bool readFields=true) | |
Constructor. More... | |
virtual | ~noiseModel ()=default |
Destructor. More... | |
virtual bool | read (const dictionary &dict) |
Read from dictionary. More... | |
virtual void | calculate ()=0 |
Abstract call to calculate. More... | |
tmp< Foam::scalarField > | PSD (const scalarField &PSDf) const |
PSD [dB/Hz]. More... | |
tmp< scalarField > | SPL (const scalarField &Prms2, const scalar f) const |
SPL [dB]. More... | |
tmp< scalarField > | SPL (const scalarField &Prms2, const scalarField &f) const |
SPL [dB]. More... | |
void | cleanFFTW () |
Clean up the FFTW. More... | |
void | writeWeightings () const |
Helper function to check weightings. More... | |
Static Public Member Functions | |
static autoPtr< noiseModel > | New (const dictionary &dict) |
Selector. More... | |
Static Public Attributes | |
static const Enum< weightingType > | weightingTypeNames_ |
Protected Member Functions | |
void | readWriteOption (const dictionary &dict, const word &lookup, bool &option) const |
Helper function to read write options and provide info feedback. More... | |
scalar | checkUniformTimeStep (const scalarList ×) const |
Check and return uniform time step. More... | |
bool | validateBounds (const scalarList &p) const |
Return true if all pressure data is within min/max bounds. More... | |
label | findStartTimeIndex (const instantList &allTimes, const scalar startTime) const |
Find and return start time index. More... | |
fileName | baseFileDir (const label dataseti) const |
Return the base output directory. More... | |
tmp< scalarField > | uniformFrequencies (const scalar deltaT, const bool check) const |
tmp< scalarField > | octaves (const scalarField &data, const scalarField &f, const labelUList &freqBandIDs) const |
Generate octave data. More... | |
tmp< scalarField > | Pf (const scalarField &p) const |
Return the fft of the given pressure data. More... | |
tmp< scalarField > | meanPf (const scalarField &p) const |
Return the multi-window mean fft of the complete pressure data [Pa]. More... | |
tmp< scalarField > | RMSmeanPf (const scalarField &p) const |
tmp< scalarField > | PSDf (const scalarField &p, const scalar deltaT) const |
scalar | RAf (const scalar f) const |
A weighting function. More... | |
scalar | gainA (const scalar f) const |
A weighting as gain in dB. More... | |
scalar | RBf (const scalar f) const |
B weighting function. More... | |
scalar | gainB (const scalar f) const |
B weighting as gain in dB. More... | |
scalar | RCf (const scalar f) const |
C weighting function. More... | |
scalar | gainC (const scalar f) const |
C weighting as gain in dB. More... | |
scalar | RDf (const scalar f) const |
D weighting function. More... | |
scalar | gainD (const scalar f) const |
D weighting as gain in dB. More... | |
Static Protected Member Functions | |
static void | setOctaveBands (const scalarField &f, const scalar fLower, const scalar fUpper, const scalar octave, labelList &fBandIDs, scalarField &fCentre) |
Protected Attributes | |
const dictionary | dict_ |
Copy of dictionary used for construction. More... | |
scalar | rhoRef_ |
Reference density (to convert from kinematic to static pressure) More... | |
label | nSamples_ |
Number of samples in sampling window, default = 2^16. More... | |
scalar | fLower_ |
Lower frequency limit, default = 25Hz. More... | |
scalar | fUpper_ |
Upper frequency limit, default = 10kHz. More... | |
scalar | startTime_ |
Start time, default = 0s. More... | |
autoPtr< windowModel > | windowModelPtr_ |
Window model. More... | |
word | graphFormat_ |
Graph format. More... | |
weightingType | SPLweighting_ |
Weighting. More... | |
scalar | dBRef_ |
Reference for dB calculation, default = 2e-5. More... | |
scalar | minPressure_ |
Min pressure value. More... | |
scalar | maxPressure_ |
Min pressure value. More... | |
fileName | outputPrefix_ |
Output file prefix, default = ''. More... | |
bool | writePrmsf_ |
Write Prmsf; default = yes. More... | |
bool | writeSPL_ |
Write SPL; default = yes. More... | |
bool | writePSD_ |
Write PSD; default = yes. More... | |
bool | writePSDf_ |
Write PSDf; default = yes. More... | |
bool | writeOctaves_ |
Write writeOctaves; default = yes. More... | |
planInfo | planInfo_ |
Plan information for FFTW. More... | |
Base class for noise models.
Data is read from a dictionary, e.g.
rhoRef 1; N 4096; fl 25; fu 10000; startTime 0; outputPrefix "test1"; SPLweighting dBA; // Optional write options dictionary writeOptions { writePrmsf no; writeSPL yes; writePSD yes; writePSDf no; writeOctaves yes; }
where
Property | Description | Required | Default value |
---|---|---|---|
rhoRef | Reference density | no | 1 |
N | Number of samples in sampling window | no | 65536 (2^16) |
fl | Lower frequency bounds | no | 25 |
fu | Upper frequency bounds | no | 10000 |
startTime | Start time | no | 0 |
outputPrefix | Prefix applied to output files | no | '' |
SPLweighting | Weighting: dBA, dBB, dBC, DBD | no | none |
dBRef | Reference for dB calculation | no | 2e-5 |
graphFormat | Graph format | no | raw |
writePrmsf | Write Prmsf data | no | yes |
writeSPL | Write SPL data | no | yes |
writePSD | Write PSD data | no | yes |
writePSDf | Write PSDf data | no | yes |
writeOctaves | Write octaves data | no | yes |
Definition at line 175 of file noiseModel.H.
|
strong |
Enumerator | |
---|---|
none | |
dBA | |
dBB | |
dBC | |
dBD |
Definition at line 179 of file noiseModel.H.
noiseModel | ( | const dictionary & | dict, |
const bool | readFields = true |
||
) |
Constructor.
Definition at line 586 of file noiseModel.C.
References noiseModel::planInfo::active, Foam::expressions::patchExpr::debug, dict, noiseModel::planInfo_, noiseModel::read(), Foam::readFields(), and noiseModel::writeWeightings().
|
virtualdefault |
Destructor.
|
protected |
Helper function to read write options and provide info feedback.
Definition at line 163 of file noiseModel.C.
References dict, Foam::endl(), and Foam::Info.
|
protected |
Check and return uniform time step.
Definition at line 186 of file noiseModel.C.
References Foam::constant::electromagnetic::e, Foam::exit(), Foam::FatalError, FatalErrorInFunction, and Foam::mag().
Referenced by surfaceNoise::initialise().
|
protected |
Return true if all pressure data is within min/max bounds.
Definition at line 220 of file noiseModel.C.
References Foam::endl(), forAll, noiseModel::maxPressure_, noiseModel::minPressure_, Foam::nl, p, and WarningInFunction.
|
protected |
Find and return start time index.
Definition at line 243 of file noiseModel.C.
References forAll, startTime, and Instant< T >::value().
Referenced by surfaceNoise::initialise().
|
protected |
Return the base output directory.
Definition at line 262 of file noiseModel.C.
References argList::envGlobalPath(), Foam::name(), functionObject::outputPrefix, and Foam::type().
Referenced by surfaceNoise::calculate().
|
protected |
Create a field of equally spaced frequencies for the current set of data - assumes a constant time step
Definition at line 277 of file noiseModel.C.
References Foam::check(), Foam::endl(), f(), forAll, N(), tmp< T >::New(), WarningInFunction, and Foam::Zero.
Referenced by surfaceNoise::calculate().
|
staticprotected |
Return a list of the frequency indices wrt f field that correspond to the bands limits for a given octave
Definition at line 55 of file noiseModel.C.
References DynamicList< T, SizeMin >::append(), Foam::endl(), f(), Foam::flatOutput(), forAll, Foam::pow(), and WarningInFunction.
Referenced by surfaceNoise::calculate().
|
protected |
Generate octave data.
Definition at line 314 of file noiseModel.C.
References Foam::endl(), f(), tmp< T >::New(), UList< T >::size(), WarningInFunction, and Foam::Zero.
Referenced by surfaceNoise::calculate().
|
protected |
Return the fft of the given pressure data.
Definition at line 365 of file noiseModel.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::mag(), n, tmp< T >::New(), p, Foam::constant::atomic::re, fft::realTransform1D(), and Foam::sqrt().
|
protected |
Return the multi-window mean fft of the complete pressure data [Pa].
Definition at line 408 of file noiseModel.C.
References N(), tmp< T >::New(), p, and Foam::Zero.
|
protected |
Return the multi-window RMS mean fft of the complete pressure data [Pa]
Definition at line 429 of file noiseModel.C.
References N(), p, Foam::sqr(), Foam::sqrt(), and Foam::Zero.
Referenced by surfaceNoise::calculate().
|
protected |
Return the multi-window Power Spectral Density (PSD) of the complete pressure data [Pa^2/Hz]
Definition at line 448 of file noiseModel.C.
References Foam::average(), Foam::expressions::patchExpr::debug, Foam::endl(), N(), tmp< T >::New(), p, Foam::Pout, Foam::sqr(), and Foam::Zero.
Referenced by surfaceNoise::calculate().
|
protected |
A weighting function.
Definition at line 483 of file noiseModel.C.
References Foam::constant::physicoChemical::c1, Foam::constant::physicoChemical::c2, f(), Foam::sqr(), and Foam::sqrt().
|
protected |
A weighting as gain in dB.
Definition at line 500 of file noiseModel.C.
References f(), and Foam::log10().
|
protected |
B weighting function.
Definition at line 511 of file noiseModel.C.
References Foam::constant::physicoChemical::c1, Foam::constant::physicoChemical::c2, f(), Foam::sqr(), and Foam::sqrt().
|
protected |
B weighting as gain in dB.
Definition at line 527 of file noiseModel.C.
References f(), and Foam::log10().
|
protected |
C weighting function.
Definition at line 538 of file noiseModel.C.
References Foam::constant::physicoChemical::c1, Foam::constant::physicoChemical::c2, f(), and Foam::sqr().
|
protected |
C weighting as gain in dB.
Definition at line 549 of file noiseModel.C.
References f(), and Foam::log10().
|
protected |
D weighting function.
Definition at line 560 of file noiseModel.C.
References f(), Foam::sqr(), and Foam::sqrt().
|
protected |
D weighting as gain in dB.
Definition at line 573 of file noiseModel.C.
References f(), and Foam::log10().
TypeName | ( | "noiseModel" | ) |
Runtime type information.
declareRunTimeSelectionTable | ( | autoPtr | , |
noiseModel | , | ||
dictionary | , | ||
(const dictionary &dict) | , | ||
(dict) | |||
) |
Run time selection table.
|
static |
Selector.
Definition at line 32 of file noiseModelNew.C.
References dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, and Foam::Info.
|
virtual |
Read from dictionary.
Reimplemented in surfaceNoise, and pointNoise.
Definition at line 624 of file noiseModel.C.
References dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, Foam::Info, windowModel::New(), Foam::nl, dictionary::readIfPresent(), and dictionary::subOrEmptyDict().
Referenced by noiseModel::noiseModel(), pointNoise::read(), and surfaceNoise::read().
|
pure virtual |
Abstract call to calculate.
Implemented in surfaceNoise, and pointNoise.
Foam::tmp< Foam::scalarField > PSD | ( | const scalarField & | PSDf | ) | const |
PSD [dB/Hz].
Definition at line 713 of file noiseModel.C.
References Foam::safeLog10(), and Foam::sqr().
Referenced by surfaceNoise::calculate().
Foam::tmp< Foam::scalarField > SPL | ( | const scalarField & | Prms2, |
const scalar | f | ||
) | const |
SPL [dB].
Definition at line 722 of file noiseModel.C.
References Foam::abort(), f(), Foam::FatalError, FatalErrorInFunction, Foam::BitOps::none(), tmp< T >::ref(), Foam::safeLog10(), and Foam::sqr().
Referenced by surfaceNoise::calculate().
Foam::tmp< Foam::scalarField > SPL | ( | const scalarField & | Prms2, |
const scalarField & | f | ||
) | const |
SPL [dB].
Definition at line 769 of file noiseModel.C.
References Foam::abort(), f(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::BitOps::none(), tmp< T >::ref(), Foam::safeLog10(), and Foam::sqr().
void cleanFFTW | ( | ) |
Clean up the FFTW.
Definition at line 827 of file noiseModel.C.
void writeWeightings | ( | ) | const |
Helper function to check weightings.
Definition at line 838 of file noiseModel.C.
Referenced by noiseModel::noiseModel().
|
static |
Definition at line 188 of file noiseModel.H.
|
protected |
Copy of dictionary used for construction.
Definition at line 230 of file noiseModel.H.
Referenced by pointNoise::calculate().
|
protected |
Reference density (to convert from kinematic to static pressure)
Definition at line 233 of file noiseModel.H.
|
protected |
Number of samples in sampling window, default = 2^16.
Definition at line 236 of file noiseModel.H.
|
protected |
Lower frequency limit, default = 25Hz.
Definition at line 239 of file noiseModel.H.
Referenced by surfaceNoise::calculate().
|
protected |
Upper frequency limit, default = 10kHz.
Definition at line 242 of file noiseModel.H.
Referenced by surfaceNoise::calculate().
|
protected |
Start time, default = 0s.
Definition at line 245 of file noiseModel.H.
Referenced by surfaceNoise::initialise().
|
protected |
Window model.
Definition at line 248 of file noiseModel.H.
Referenced by surfaceNoise::calculate(), and surfaceNoise::initialise().
|
protected |
|
protected |
Weighting.
Definition at line 254 of file noiseModel.H.
|
protected |
Reference for dB calculation, default = 2e-5.
Definition at line 257 of file noiseModel.H.
|
protected |
Min pressure value.
Definition at line 263 of file noiseModel.H.
Referenced by noiseModel::validateBounds().
|
protected |
Min pressure value.
Definition at line 266 of file noiseModel.H.
Referenced by noiseModel::validateBounds().
|
protected |
Output file prefix, default = ''.
Definition at line 272 of file noiseModel.H.
|
protected |
Write Prmsf; default = yes.
Definition at line 275 of file noiseModel.H.
Referenced by surfaceNoise::calculate().
|
protected |
Write SPL; default = yes.
Definition at line 278 of file noiseModel.H.
Referenced by surfaceNoise::calculate().
|
protected |
Write PSD; default = yes.
Definition at line 281 of file noiseModel.H.
Referenced by surfaceNoise::calculate().
|
protected |
Write PSDf; default = yes.
Definition at line 284 of file noiseModel.H.
Referenced by surfaceNoise::calculate().
|
protected |
Write writeOctaves; default = yes.
Definition at line 287 of file noiseModel.H.
Referenced by surfaceNoise::calculate().
|
mutableprotected |
Plan information for FFTW.
Definition at line 293 of file noiseModel.H.
Referenced by noiseModel::noiseModel().