Perform noise analysis on point-based pressure data. More...
#include <pointNoise.H>


Public Member Functions | |
| TypeName ("pointNoise") | |
| Runtime type information. | |
| pointNoise (const dictionary &dict, const objectRegistry &obr, const word &name=typeName, const bool readFields=true) | |
| Constructor. | |
| virtual | ~pointNoise ()=default |
| Destructor. | |
| virtual bool | read (const dictionary &dict) |
| Read from dictionary. | |
| virtual void | calculate () |
| Calculate. | |
| Public Member Functions inherited from noiseModel | |
| TypeName ("noiseModel") | |
| Runtime type information. | |
| declareRunTimeSelectionTable (autoPtr, noiseModel, dictionary,(const dictionary &dict, const objectRegistry &obr),(dict, obr)) | |
| Run time selection table. | |
| noiseModel (const noiseModel &)=delete | |
| No copy construct. | |
| void | operator= (const noiseModel &)=delete |
| No copy assignment. | |
| noiseModel (const dictionary &dict, const objectRegistry &obr, const word &name, const bool readFields=true) | |
| Constructor. | |
| virtual | ~noiseModel ()=default |
| Destructor. | |
| tmp< Foam::scalarField > | PSD (const scalarField &PSDf) const |
| PSD [dB/Hz]. | |
| tmp< scalarField > | SPL (const scalarField &Prms2, const scalar f) const |
| SPL [dB]. | |
| tmp< scalarField > | SPL (const scalarField &Prms2, const scalarField &f) const |
| SPL [dB]. | |
| void | cleanFFTW () |
| Clean up the FFTW. | |
| void | writeWeightings () const |
| Helper function to check weightings. | |
| Public Member Functions inherited from writeFile | |
| writeFile (const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat") | |
| Construct from objectRegistry, prefix, fileName. | |
| writeFile (const objectRegistry &obr, const fileName &prefix, const word &name, const dictionary &dict, const bool writeToFile=true, const string &ext=".dat") | |
| Construct from objectRegistry, prefix, fileName and read options from dictionary. | |
| writeFile (const writeFile &wf) | |
| Construct copy. | |
| virtual | ~writeFile ()=default |
| Destructor. | |
| virtual const string & | setExt (const string &ext) |
| Set extension. | |
| virtual OFstream & | file () |
| Return access to the file (if only 1). | |
| virtual bool | writeToFile () const |
| Flag to allow writing to file. | |
| virtual bool | canWriteToFile () const |
| Flag to allow writing to the file. | |
| virtual bool | canResetFile () const |
| Flag to allow resetting the file. | |
| virtual bool | canWriteHeader () const |
| Flag to allow writing the header. | |
| virtual label | charWidth () const |
| Return width of character stream output. | |
| virtual void | writeCommented (Ostream &os, const string &str) const |
| Write a commented string to stream. | |
| virtual void | writeTabbed (Ostream &os, const string &str) const |
| Write a tabbed string to stream. | |
| virtual void | writeHeader (Ostream &os, const string &str) const |
| Write a commented header to stream. | |
| virtual void | writeCurrentTime (Ostream &os) const |
| Write the current time to stream. | |
| virtual void | writeBreak (Ostream &os) const |
| Write a break marker to the stream. | |
| template<class Type> | |
| void | writeHeaderValue (Ostream &os, const string &property, const Type &value) const |
| Write a (commented) header property and value pair. | |
| template<class Type> | |
| void | writeValue (Ostream &os, const Type &val) const |
| Write a given value to stream with the space delimiter. | |
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. | |
| Protected Member Functions inherited from noiseModel | |
| scalar | checkUniformTimeStep (const scalarList ×) const |
| Check and return uniform time step. | |
| bool | validateBounds (const scalarList &p) const |
| Return true if all pressure data is within min/max bounds. | |
| fileName | baseFileDir (const label dataseti) const |
| Return the base output directory. | |
| void | writeFileHeader (Ostream &os, const string &x, const string &y, const UList< Tuple2< string, token > > &headerValues=UList< Tuple2< string, token > >::null()) const |
| Write output file header. | |
| void | writeFreqDataToFile (Ostream &os, const scalarField &f, const scalarField &fx) const |
| tmp< scalarField > | uniformFrequencies (const scalar deltaT, const bool check) const |
| Create a field of equally spaced frequencies for the current set of data - assumes a constant time step. | |
| tmp< scalarField > | octaves (const scalarField &data, const scalarField &f, const labelUList &freqBandIDs) const |
| Generate octave data. | |
| tmp< scalarField > | Pf (const scalarField &p) const |
| Return the fft of the given pressure data. | |
| tmp< scalarField > | meanPf (const scalarField &p) const |
| Return the multi-window mean fft of the complete pressure data [Pa]. | |
| tmp< scalarField > | RMSmeanPf (const scalarField &p) const |
| Return the multi-window RMS mean fft of the complete pressure data [Pa]. | |
| tmp< scalarField > | PSDf (const scalarField &p, const scalar deltaT) const |
| Return the multi-window Power Spectral Density (PSD) of the complete pressure data [Pa^2/Hz]. | |
| scalar | RAf (const scalar f) const |
| A weighting function. | |
| scalar | gainA (const scalar f) const |
| A weighting as gain in dB. | |
| scalar | RBf (const scalar f) const |
| B weighting function. | |
| scalar | gainB (const scalar f) const |
| B weighting as gain in dB. | |
| scalar | RCf (const scalar f) const |
| C weighting function. | |
| scalar | gainC (const scalar f) const |
| C weighting as gain in dB. | |
| scalar | RDf (const scalar f) const |
| D weighting function. | |
| scalar | gainD (const scalar f) const |
| D weighting as gain in dB. | |
| Protected Member Functions inherited from writeFile | |
| void | initStream (Ostream &os) const |
| Initialise the output stream for writing. | |
| fileName | baseFileDir () const |
| Return the base directory for output. | |
| fileName | baseTimeDir () const |
| Return the base directory for the current time value. | |
| fileName | filePath (const fileName &fName) const |
| Return the full path for the supplied file name. | |
| virtual autoPtr< OFstream > | newFile (const fileName &fName) const |
| Return autoPtr to a new file using file name. | |
| virtual autoPtr< OFstream > | newFileAtTime (const word &name, scalar timeValue) const |
| Return autoPtr to a new file for a given time. | |
| virtual autoPtr< OFstream > | newFileAtStartTime (const word &name) const |
| Return autoPtr to a new file using the simulation start time. | |
| virtual void | resetFile (const word &name) |
| Reset internal file pointer to new file with new name. | |
| Omanip< int > | valueWidth (const label offset=0) const |
| Return the value width when writing to stream with optional offset. | |
| void | operator= (const writeFile &)=delete |
| No copy assignment. | |
| virtual autoPtr< OFstream > | createFile (const word &name, scalar timeValue) const |
| Deprecated(2022-09) Return autoPtr to a new file for a given time. | |
| virtual autoPtr< OFstream > | createFile (const word &name) const |
| Deprecated(2022-09) Return autoPtr to a new file using the simulation start time. | |
Protected Attributes | |
| List< fileName > | inputFileNames_ |
| Input file names - optional. | |
| Protected Attributes inherited from noiseModel | |
| const dictionary | dict_ |
| Copy of dictionary used for construction. | |
| scalar | rhoRef_ |
| Reference density (to convert from kinematic to static pressure). | |
| label | nSamples_ |
| Number of samples in sampling window, default = 2^16. | |
| scalar | fLower_ |
| Lower frequency limit, default = 25Hz. | |
| scalar | fUpper_ |
| Upper frequency limit, default = 10kHz. | |
| scalar | sampleFreq_ |
| Prescribed sample frequency. | |
| scalar | startTime_ |
| Start time, default = 0s. | |
| autoPtr< windowModel > | windowModelPtr_ |
| Window model. | |
| weightingType | SPLweighting_ |
| Weighting. | |
| scalar | dBRef_ |
| Reference for dB calculation, default = 2e-5. | |
| scalar | minPressure_ |
| Min pressure value. | |
| scalar | maxPressure_ |
| Min pressure value. | |
| fileName | outputPrefix_ |
| Output file prefix, default = ''. | |
| bool | writePrmsf_ |
| Write Prmsf; default = yes. | |
| bool | writeSPL_ |
| Write SPL; default = yes. | |
| bool | writePSD_ |
| Write PSD; default = yes. | |
| bool | writePSDf_ |
| Write PSDf; default = yes. | |
| bool | writeOctaves_ |
| Write writeOctaves; default = yes. | |
| planInfo | planInfo_ |
| Plan information for FFTW. | |
| Protected Attributes inherited from writeFile | |
| const objectRegistry & | fileObr_ |
| Reference to the region objectRegistry. | |
| const fileName | prefix_ |
| Prefix. | |
| word | fileName_ |
| Name of file. | |
| autoPtr< OFstream > | filePtr_ |
| File pointer. | |
| label | writePrecision_ |
| Write precision. | |
| bool | writeToFile_ |
| Flag to enable/disable writing to file. | |
| bool | updateHeader_ |
| Flag to update the header, e.g. on mesh changes. Default is true. | |
| bool | writtenHeader_ |
| Flag to identify whether the header has been written. | |
| bool | useUserTime_ |
| Flag to use the specified user time, e.g. CA deg instead of seconds. Default = true. | |
| scalar | startTime_ |
| Start time value. | |
| string | ext_ |
| File extension; default = .dat. | |
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, const objectRegistry &obr) |
| Selector. | |
| Static Public Attributes inherited from noiseModel | |
| static const Enum< weightingType > | weightingTypeNames_ |
| Static Public Attributes inherited from writeFile | |
| static label | addChars = 8 |
| Additional characters for writing. | |
| Static Protected Member Functions inherited from noiseModel | |
| static label | findStartTimeIndex (const instantList &allTimes, const scalar startTime) |
| Find and return start time index. | |
| static void | setOctaveBands (const scalarField &f, const scalar fLower, const scalar fUpper, const scalar octave, labelList &fBandIDs, scalarField &fCentre) |
| Return a list of the frequency indices wrt f field that correspond to the bands limits for a given octave. | |
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;
Definition at line 88 of file pointNoise.H.
| pointNoise | ( | const dictionary & | dict, |
| const objectRegistry & | obr, | ||
| const word & | name = typeName, | ||
| const bool | readFields = true ) |
Constructor.
Definition at line 219 of file pointNoise.C.
References dict, Foam::name(), noiseModel::noiseModel(), read(), and Foam::readFields().

|
virtualdefault |
Destructor.
References dict.
|
protected |
Definition at line 39 of file pointNoise.C.
References DynamicList< T, SizeMin >::append(), forAll, p, p0, UList< T >::size(), noiseModel::startTime_, and List< T >::transfer().
Referenced by processData().


|
protected |
Process the CSV data.
Definition at line 64 of file pointNoise.C.
References Foam::average(), writeFile::baseFileDir(), noiseModel::checkUniformTimeStep(), Foam::endl(), f(), writeFile::fileObr_, filePtr, filterTimeData(), noiseModel::fLower_, CSV< Type >::fName(), noiseModel::fUpper_, Foam::Info, UPstream::master(), writeFile::newFile(), Foam::nl, windowModel::nSamples(), noiseModel::octaves(), os(), p, noiseModel::PSD(), noiseModel::PSDf(), noiseModel::rhoRef_, noiseModel::RMSmeanPf(), noiseModel::sampleFreq_, noiseModel::setOctaveBands(), UList< T >::size(), noiseModel::SPL(), noiseModel::SPLweighting_, fileName::stem(), noiseModel::uniformFrequencies(), noiseModel::validateBounds(), WarningInFunction, noiseModel::weightingTypeNames_, noiseModel::windowModelPtr_, noiseModel::writeFileHeader(), noiseModel::writeFreqDataToFile(), noiseModel::writeOctaves_, noiseModel::writePrmsf_, noiseModel::writePSD_, noiseModel::writePSDf_, noiseModel::writeSPL_, TableBase< Type >::x(), and TableBase< Type >::y().
Referenced by calculate().


| TypeName | ( | "pointNoise" | ) |
Runtime type information.
References dict, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, Foam::name(), and Foam::readFields().

|
virtual |
Read from dictionary.
Reimplemented from noiseModel.
Definition at line 262 of file pointNoise.C.
References dict, inputFileNames_, and noiseModel::read().
Referenced by pointNoise().


|
virtual |
Calculate.
Implements noiseModel.
Definition at line 238 of file pointNoise.C.
References noiseModel::dict_, argList::envGlobalPath(), string::expand(), forAll, inputFileNames_, fileName::isAbsolute(), UPstream::master(), and processData().

Input file names - optional.
Definition at line 100 of file pointNoise.H.
Referenced by calculate(), and read().