Base class for noise models. More...
#include <noiseModel.H>


Classes | |
| struct | planInfo |
| FFTW planner information. More... | |
| struct | octaveBandInfo |
| Octave band information. More... | |
Public Types | |
| enum class | weightingType { none , dBA , dBB , dBC , dBD } |
Public Member Functions | |
| 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. | |
| virtual bool | read (const dictionary &dict) |
| Read from dictionary. | |
| virtual void | calculate ()=0 |
| Abstract call to calculate. | |
| 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. | |
Static Public Member Functions | |
| static autoPtr< noiseModel > | New (const dictionary &dict, const objectRegistry &obr) |
| Selector. | |
Static Public Attributes | |
| static const Enum< weightingType > | weightingTypeNames_ |
| Static Public Attributes inherited from writeFile | |
| static label | addChars = 8 |
| Additional characters for writing. | |
Protected Member Functions | |
| 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. | |
Static Protected Member Functions | |
| 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. | |
Protected Attributes | |
| 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. | |
Base class for noise models.
Data is read from a dictionary, e.g.
rhoRef 1;
N 4096;
minFreq 25;
maxFreq 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) |
minFreq | Lower frequency bounds (fl) | no | 25 |
maxFreq | Upper frequency bounds (fu) | no | 10000 |
sampleFreq | Sample frequency | no | 0 |
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 |
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 |
If provided, the sampleFreq is used to define the deltaT between samples. Otherwise the deltaT is determined from the time samples themselves.
Definition at line 175 of file noiseModel.H.
|
strong |
| Enumerator | |
|---|---|
| none | |
| dBA | |
| dBB | |
| dBC | |
| dBD | |
Definition at line 181 of file noiseModel.H.
|
delete |
No copy construct.
References noiseModel().
Referenced by declareRunTimeSelectionTable(), noiseModel(), operator=(), pointNoise::pointNoise(), and surfaceNoise::surfaceNoise().


| noiseModel | ( | const dictionary & | dict, |
| const objectRegistry & | obr, | ||
| const word & | name, | ||
| const bool | readFields = true ) |
Constructor.
Definition at line 604 of file noiseModel.C.
References dBRef_, dict, dict_, e, fLower_, fUpper_, maxPressure_, minPressure_, Foam::name(), nSamples_, outputPrefix_, planInfo_, read(), Foam::readFields(), rhoRef_, sampleFreq_, SPLweighting_, startTime_, windowModelPtr_, writeFile::writeFile(), writeOctaves_, writePrmsf_, writePSD_, writePSDf_, writeSPL_, and writeWeightings().

|
virtualdefault |
Destructor.
References dict.
|
protected |
Check and return uniform time step.
Definition at line 176 of file noiseModel.C.
References UList< T >::back(), e, Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, UList< T >::front(), Foam::mag(), UList< T >::size(), and WarningInFunction.
Referenced by surfaceNoise::initialise(), and pointNoise::processData().


|
protected |
Return true if all pressure data is within min/max bounds.
Definition at line 216 of file noiseModel.C.
References Foam::endl(), forAll, maxPressure_, minPressure_, Foam::nl, p, and WarningInFunction.
Referenced by pointNoise::processData().


|
inlinestaticprotected |
Find and return start time index.
Definition at line 349 of file noiseModel.H.
References instant::findStart(), findStartTimeIndex(), and startTime.
Referenced by findStartTimeIndex().


|
protected |
Return the base output directory.
Definition at line 238 of file noiseModel.C.
References Foam::name(), outputPrefix_, and Foam::type().

|
protected |
Write output file header.
Definition at line 249 of file noiseModel.C.
References dBRef_, fLower_, fUpper_, Foam::nl, os(), windowModelPtr_, writeFile::writeCommented(), Foam::writeHeader(), writeFile::writeHeaderValue(), writeFile::writeTabbed(), x, and y.
Referenced by surfaceNoise::calculate(), and pointNoise::processData().


|
protected |
Definition at line 277 of file noiseModel.C.
References f(), fLower_, forAll, fUpper_, Foam::nl, os(), and Foam::tab.
Referenced by surfaceNoise::calculate(), and pointNoise::processData().


|
protected |
Create a field of equally spaced frequencies for the current set of data - assumes a constant time step.
Definition at line 294 of file noiseModel.C.
References Foam::check(), Foam::endl(), f(), fLower_, forAll, fUpper_, N(), tmp< T >::New(), WarningInFunction, windowModelPtr_, and Foam::Zero.
Referenced by surfaceNoise::calculate(), and pointNoise::processData().


|
staticprotected |
Return a list of the frequency indices wrt f field that correspond to the bands limits for a given octave.
Definition at line 47 of file noiseModel.C.
References DynamicList< T, SizeMin >::append(), Foam::endl(), f(), Foam::flatOutput(), forAll, HashSet< Key, Hash >::insert(), DynamicList< T, SizeMin >::pop_back(), Foam::pow(), UList< T >::size(), HashTable< T, Key, Hash >::sortedToc(), List< T >::transfer(), and WarningInFunction.
Referenced by surfaceNoise::calculate(), and pointNoise::processData().


|
protected |
Generate octave data.
Definition at line 331 of file noiseModel.C.
References Foam::endl(), f(), bitSet::flip(), tmp< T >::New(), bitSet::set(), PackedList< Width >::size(), UList< T >::size(), bitSet::sortedToc(), WarningInFunction, and Foam::Zero.
Referenced by surfaceNoise::calculate(), and pointNoise::processData().


|
protected |
Return the fft of the given pressure data.
Definition at line 383 of file noiseModel.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::mag(), n, tmp< T >::New(), p, planInfo_, fft::realTransform1D(), and Foam::sqrt().
Referenced by meanPf(), PSDf(), and RMSmeanPf().


|
protected |
Return the multi-window mean fft of the complete pressure data [Pa].
Definition at line 426 of file noiseModel.C.
References meanPf(), N(), tmp< T >::New(), p, Pf(), windowModelPtr_, and Foam::Zero.
Referenced by meanPf().


|
protected |
Return the multi-window RMS mean fft of the complete pressure data [Pa].
Definition at line 446 of file noiseModel.C.
References N(), p, Pf(), Foam::sqr(), Foam::sqrt(), windowModelPtr_, and Foam::Zero.
Referenced by surfaceNoise::calculate(), and pointNoise::processData().


|
protected |
Return the multi-window Power Spectral Density (PSD) of the complete pressure data [Pa^2/Hz].
Definition at line 465 of file noiseModel.C.
References Foam::average(), Foam::endl(), N(), tmp< T >::New(), p, Pf(), Foam::Pout, Foam::sqr(), windowModelPtr_, and Foam::Zero.
Referenced by calculate(), surfaceNoise::calculate(), pointNoise::processData(), and PSD().


|
protected |
A weighting function.
Definition at line 501 of file noiseModel.C.
References f(), Foam::sqr(), and Foam::sqrt().
Referenced by gainA().


|
protected |
A weighting as gain in dB.
Definition at line 518 of file noiseModel.C.
References f(), Foam::log10(), and RAf().
Referenced by SPL(), SPL(), and writeWeightings().


|
protected |
B weighting function.
Definition at line 529 of file noiseModel.C.
References f(), Foam::sqr(), and Foam::sqrt().
Referenced by gainB().


|
protected |
B weighting as gain in dB.
Definition at line 545 of file noiseModel.C.
References f(), Foam::log10(), and RBf().
Referenced by SPL(), SPL(), and writeWeightings().


|
protected |
C weighting function.
Definition at line 556 of file noiseModel.C.
References f(), and Foam::sqr().
Referenced by gainC().


|
protected |
C weighting as gain in dB.
Definition at line 567 of file noiseModel.C.
References f(), Foam::log10(), and RCf().
Referenced by SPL(), SPL(), and writeWeightings().


|
protected |
D weighting function.
Definition at line 578 of file noiseModel.C.
References f(), Foam::sqr(), and Foam::sqrt().
Referenced by gainD().


|
protected |
D weighting as gain in dB.
Definition at line 591 of file noiseModel.C.
References f(), Foam::log10(), and RDf().
Referenced by SPL(), SPL(), and writeWeightings().


| TypeName | ( | "noiseModel" | ) |
Runtime type information.
| declareRunTimeSelectionTable | ( | autoPtr | , |
| noiseModel | , | ||
| dictionary | , | ||
| (const dictionary &dict, const objectRegistry &obr) | , | ||
| (dict, obr) | ) |
Run time selection table.
References dict, and noiseModel().

|
delete |
No copy assignment.
References dict, Foam::name(), noiseModel(), and Foam::readFields().

|
static |
Selector.
Definition at line 25 of file noiseModelNew.C.
References dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, and Foam::Info.

|
virtual |
Read from dictionary.
Reimplemented from writeFile.
Reimplemented in pointNoise, and surfaceNoise.
Definition at line 649 of file noiseModel.C.
References cleanFFTW(), dBRef_, dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, fLower_, fUpper_, Foam::Info, maxPressure_, minPressure_, windowModel::New(), Foam::nl, nSamples_, outputPrefix_, planInfo_, writeFile::read(), Foam::readWriteOption(), rhoRef_, sampleFreq_, SPLweighting_, startTime_, weightingTypeNames_, windowModelPtr_, writeOctaves_, writePrmsf_, writePSD_, writePSDf_, and writeSPL_.
Referenced by noiseModel(), pointNoise::read(), and surfaceNoise::read().


|
pure virtual |
Abstract call to calculate.
Implemented in pointNoise, and surfaceNoise.

| Foam::tmp< Foam::scalarField > PSD | ( | const scalarField & | PSDf | ) | const |
PSD [dB/Hz].
Definition at line 757 of file noiseModel.C.
References dBRef_, PSDf(), Foam::safeLog10(), and Foam::sqr().
Referenced by surfaceNoise::calculate(), and pointNoise::processData().


| Foam::tmp< Foam::scalarField > SPL | ( | const scalarField & | Prms2, |
| const scalar | f ) const |
SPL [dB].
Definition at line 766 of file noiseModel.C.
References Foam::abort(), dBA, dBB, dBC, dBD, dBRef_, f(), Foam::FatalError, FatalErrorInFunction, gainA(), gainB(), gainC(), gainD(), none, tmp< T >::ref(), Foam::safeLog10(), SPLweighting_, Foam::sqr(), and weightingTypeNames_.
Referenced by surfaceNoise::calculate(), and pointNoise::processData().


| Foam::tmp< Foam::scalarField > SPL | ( | const scalarField & | Prms2, |
| const scalarField & | f ) const |
SPL [dB].
Definition at line 813 of file noiseModel.C.
References Foam::abort(), dBA, dBB, dBC, dBD, dBRef_, f(), Foam::FatalError, FatalErrorInFunction, forAll, gainA(), gainB(), gainC(), gainD(), none, tmp< T >::ref(), Foam::safeLog10(), SPLweighting_, Foam::sqr(), and weightingTypeNames_.

| void cleanFFTW | ( | ) |
Clean up the FFTW.
Definition at line 872 of file noiseModel.C.
References planInfo_.
Referenced by read().

| void writeWeightings | ( | ) | const |
Helper function to check weightings.
Definition at line 883 of file noiseModel.C.
References f(), gainA(), gainB(), gainC(), gainD(), and Foam::nl.
Referenced by noiseModel().


|
static |
Definition at line 190 of file noiseModel.H.
Referenced by pointNoise::processData(), read(), SPL(), and SPL().
|
protected |
Copy of dictionary used for construction.
Definition at line 228 of file noiseModel.H.
Referenced by pointNoise::calculate(), and noiseModel().
|
protected |
Reference density (to convert from kinematic to static pressure).
Definition at line 233 of file noiseModel.H.
Referenced by noiseModel(), pointNoise::processData(), read(), and surfaceNoise::readSurfaceData().
|
protected |
Number of samples in sampling window, default = 2^16.
Definition at line 238 of file noiseModel.H.
Referenced by noiseModel(), and read().
|
protected |
Lower frequency limit, default = 25Hz.
Definition at line 243 of file noiseModel.H.
Referenced by surfaceNoise::calculate(), noiseModel(), pointNoise::processData(), read(), uniformFrequencies(), writeFileHeader(), and writeFreqDataToFile().
|
protected |
Upper frequency limit, default = 10kHz.
Definition at line 248 of file noiseModel.H.
Referenced by surfaceNoise::calculate(), noiseModel(), pointNoise::processData(), read(), uniformFrequencies(), writeFileHeader(), and writeFreqDataToFile().
|
protected |
Prescribed sample frequency.
Definition at line 253 of file noiseModel.H.
Referenced by surfaceNoise::initialise(), noiseModel(), pointNoise::processData(), and read().
|
protected |
Start time, default = 0s.
Definition at line 258 of file noiseModel.H.
Referenced by pointNoise::filterTimeData(), surfaceNoise::initialise(), noiseModel(), and read().
|
protected |
Window model.
Definition at line 263 of file noiseModel.H.
Referenced by surfaceNoise::calculate(), surfaceNoise::initialise(), meanPf(), noiseModel(), pointNoise::processData(), PSDf(), read(), RMSmeanPf(), uniformFrequencies(), and writeFileHeader().
|
protected |
Weighting.
Definition at line 268 of file noiseModel.H.
Referenced by noiseModel(), pointNoise::processData(), read(), SPL(), and SPL().
|
protected |
Reference for dB calculation, default = 2e-5.
Definition at line 273 of file noiseModel.H.
Referenced by noiseModel(), PSD(), read(), SPL(), SPL(), and writeFileHeader().
|
protected |
Min pressure value.
Definition at line 281 of file noiseModel.H.
Referenced by noiseModel(), read(), and validateBounds().
|
protected |
Min pressure value.
Definition at line 286 of file noiseModel.H.
Referenced by noiseModel(), read(), and validateBounds().
|
protected |
Output file prefix, default = ''.
Definition at line 294 of file noiseModel.H.
Referenced by baseFileDir(), noiseModel(), and read().
|
protected |
Write Prmsf; default = yes.
Definition at line 299 of file noiseModel.H.
Referenced by surfaceNoise::calculate(), noiseModel(), pointNoise::processData(), and read().
|
protected |
Write SPL; default = yes.
Definition at line 304 of file noiseModel.H.
Referenced by surfaceNoise::calculate(), noiseModel(), pointNoise::processData(), and read().
|
protected |
Write PSD; default = yes.
Definition at line 309 of file noiseModel.H.
Referenced by surfaceNoise::calculate(), noiseModel(), pointNoise::processData(), and read().
|
protected |
Write PSDf; default = yes.
Definition at line 314 of file noiseModel.H.
Referenced by surfaceNoise::calculate(), noiseModel(), pointNoise::processData(), and read().
|
protected |
Write writeOctaves; default = yes.
Definition at line 319 of file noiseModel.H.
Referenced by surfaceNoise::calculate(), noiseModel(), pointNoise::processData(), and read().
|
mutableprotected |
Plan information for FFTW.
Definition at line 327 of file noiseModel.H.
Referenced by cleanFFTW(), noiseModel(), Pf(), and read().