Perform noise analysis on point-based pressure data. More...
Public Member Functions | |
TypeName ("pointNoise") | |
Runtime type information. More... | |
pointNoise (const dictionary &dict, const bool readFields=true) | |
Constructor. More... | |
virtual | ~pointNoise ()=default |
Destructor. More... | |
virtual bool | read (const dictionary &dict) |
Read from dictionary. More... | |
virtual void | calculate () |
Calculate. More... | |
Public Member Functions inherited from noiseModel | |
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... | |
Protected Member Functions | |
void | filterTimeData (const scalarField &t0, const scalarField &p0, scalarField &t, scalarField &p) const |
void | processData (const label dataseti, const Function1Types::CSV< scalar > &data) |
Process the CSV data. More... | |
Protected Member Functions inherited from noiseModel | |
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... | |
Protected Attributes | |
List< fileName > | inputFileNames_ |
Input file names - optional. More... | |
Protected Attributes inherited from noiseModel | |
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... | |
Additional Inherited Members | |
Public Types inherited from noiseModel | |
enum class | weightingType { none , dBA , dBB , dBC , dBD } |
Static Public Member Functions inherited from noiseModel | |
static autoPtr< noiseModel > | New (const dictionary &dict) |
Selector. More... | |
Static Public Attributes inherited from noiseModel | |
static const Enum< weightingType > | weightingTypeNames_ |
Static Protected Member Functions inherited from noiseModel | |
static void | setOctaveBands (const scalarField &f, const scalar fLower, const scalar fUpper, const scalar octave, labelList &fBandIDs, scalarField &fCentre) |
Perform noise analysis on point-based pressure data.
Input data is read from a dictionary, e.g.
// Pressure reference pRef 0; // Number of samples in sampling window, default = 2^16 (=65536) N 4096; // Lower frequency bounds fl 25; // Upper frequency bounds fu 25; // Start time startTime 0; windowModel <modelType> <modelType>Coeffs { ... } // Pressure data supplied in CSV file format file "pressureData"; //files ("pressureData1" "pressureData2"); nHeaderLine 1; refColumn 0; componentColumns (1); separator " "; mergeSeparators yes; graphFormat raw;
Definition at line 94 of file pointNoise.H.
pointNoise | ( | const dictionary & | dict, |
const bool | readFields = true |
||
) |
Constructor.
Definition at line 225 of file pointNoise.C.
References dict, pointNoise::read(), and Foam::readFields().
|
virtualdefault |
Destructor.
|
protected |
Definition at line 46 of file pointNoise.C.
References DynamicList< T, SizeMin >::append(), forAll, p, p0, UList< T >::size(), noiseModel::startTime_, and List< T >::transfer().
Referenced by pointNoise::processData().
|
protected |
Process the CSV data.
Definition at line 71 of file pointNoise.C.
References Foam::average(), noiseModel::baseFileDir(), noiseModel::checkUniformTimeStep(), Foam::endl(), f(), pointNoise::filterTimeData(), noiseModel::fLower_, noiseModel::fUpper_, g, noiseModel::graphFormat_, Foam::Info, Foam::nl, windowModel::nSamples(), noiseModel::octaves(), p, noiseModel::PSD(), noiseModel::PSDf(), noiseModel::rhoRef_, noiseModel::RMSmeanPf(), noiseModel::setOctaveBands(), UList< T >::size(), noiseModel::SPL(), noiseModel::SPLweighting_, noiseModel::uniformFrequencies(), noiseModel::validateBounds(), noiseModel::weightingTypeNames_, noiseModel::windowModelPtr_, graph::wordify(), regIOobject::write(), noiseModel::writeOctaves_, noiseModel::writePrmsf_, noiseModel::writePSD_, noiseModel::writePSDf_, and noiseModel::writeSPL_.
Referenced by pointNoise::calculate().
TypeName | ( | "pointNoise" | ) |
Runtime type information.
|
virtual |
Read from dictionary.
Reimplemented from noiseModel.
Definition at line 262 of file pointNoise.C.
References dict, UList< T >::first(), pointNoise::inputFileNames_, kEpsilonLopesdaCosta< BasicTurbulenceModel >::read(), dictionary::readEntry(), dictionary::readIfPresent(), and List< T >::resize().
Referenced by pointNoise::pointNoise().
|
virtual |
Calculate.
Implements noiseModel.
Definition at line 238 of file pointNoise.C.
References noiseModel::dict_, argList::envGlobalPath(), string::expand(), forAll, pointNoise::inputFileNames_, fileName::isAbsolute(), splitCell::master(), and pointNoise::processData().
Input file names - optional.
Definition at line 104 of file pointNoise.H.
Referenced by pointNoise::calculate(), and pointNoise::read().