mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: noiseModels - enable models to accept lists of file names
This commit is contained in:
@ -204,7 +204,8 @@ Foam::Function1Types::CSV<Type>::CSV
|
|||||||
(
|
(
|
||||||
const word& entryName,
|
const word& entryName,
|
||||||
const dictionary& dict,
|
const dictionary& dict,
|
||||||
const word& ext
|
const word& ext,
|
||||||
|
const fileName& fName
|
||||||
)
|
)
|
||||||
:
|
:
|
||||||
TableBase<Type>(entryName, dict.subDict(entryName + ext)),
|
TableBase<Type>(entryName, dict.subDict(entryName + ext)),
|
||||||
@ -214,7 +215,7 @@ Foam::Function1Types::CSV<Type>::CSV
|
|||||||
componentColumns_(coeffs_.lookup("componentColumns")),
|
componentColumns_(coeffs_.lookup("componentColumns")),
|
||||||
separator_(coeffs_.lookupOrDefault<string>("separator", string(","))[0]),
|
separator_(coeffs_.lookupOrDefault<string>("separator", string(","))[0]),
|
||||||
mergeSeparators_(readBool(coeffs_.lookup("mergeSeparators"))),
|
mergeSeparators_(readBool(coeffs_.lookup("mergeSeparators"))),
|
||||||
fName_(coeffs_.lookup("fileName"))
|
fName_(fName != fileName::null ? fName : coeffs_.lookup("fileName"))
|
||||||
{
|
{
|
||||||
if (componentColumns_.size() != pTraits<Type>::nComponents)
|
if (componentColumns_.size() != pTraits<Type>::nComponents)
|
||||||
{
|
{
|
||||||
|
|||||||
@ -122,7 +122,8 @@ public:
|
|||||||
(
|
(
|
||||||
const word& entryName,
|
const word& entryName,
|
||||||
const dictionary& dict,
|
const dictionary& dict,
|
||||||
const word& ext = "Coeffs"
|
const word& ext = "Coeffs",
|
||||||
|
const fileName& fName = fileName::null
|
||||||
);
|
);
|
||||||
|
|
||||||
//- Copy constructor
|
//- Copy constructor
|
||||||
|
|||||||
@ -94,17 +94,41 @@ Foam::label Foam::noiseModel::findStartTimeIndex
|
|||||||
|
|
||||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
Foam::noiseModel::noiseModel(const dictionary& dict)
|
Foam::noiseModel::noiseModel(const dictionary& dict, const bool readFields)
|
||||||
:
|
:
|
||||||
dict_(dict),
|
dict_(dict),
|
||||||
rhoRef_(dict.lookupOrDefault<scalar>("rhoRef", 1)),
|
rhoRef_(1),
|
||||||
nSamples_(dict.lookupOrDefault<label>("N", 65536)),
|
nSamples_(65536),
|
||||||
fLower_(dict.lookupOrDefault<scalar>("fl", 25)),
|
fLower_(25),
|
||||||
fUpper_(dict.lookupOrDefault<scalar>("fu", 10000)),
|
fUpper_(10000),
|
||||||
startTime_(dict.lookupOrDefault<scalar>("startTime", 0)),
|
startTime_(0),
|
||||||
windowModelPtr_(windowModel::New(dict, nSamples_)),
|
windowModelPtr_(),
|
||||||
graphFormat_(dict.lookupOrDefault<word>("graphFormat", "raw"))
|
graphFormat_("raw")
|
||||||
{
|
{
|
||||||
|
if (readFields)
|
||||||
|
{
|
||||||
|
read(dict);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::noiseModel::~noiseModel()
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
bool Foam::noiseModel::read(const dictionary& dict)
|
||||||
|
{
|
||||||
|
dict.readIfPresent("rhoRef", rhoRef_);
|
||||||
|
dict.readIfPresent("N", nSamples_);
|
||||||
|
dict.readIfPresent("fl", fLower_);
|
||||||
|
dict.readIfPresent("fu", fUpper_);
|
||||||
|
dict.readIfPresent("startTime", startTime_);
|
||||||
|
dict.readIfPresent("graphFormat", graphFormat_);
|
||||||
|
|
||||||
// Check number of samples - must be a power of 2 for our FFT
|
// Check number of samples - must be a power of 2 for our FFT
|
||||||
bool powerOf2 = ((nSamples_ != 0) && !(nSamples_ & (nSamples_ - 1)));
|
bool powerOf2 = ((nSamples_ != 0) && !(nSamples_ & (nSamples_ - 1)));
|
||||||
if (!powerOf2)
|
if (!powerOf2)
|
||||||
@ -139,13 +163,11 @@ Foam::noiseModel::noiseModel(const dictionary& dict)
|
|||||||
<< exit(FatalIOError);
|
<< exit(FatalIOError);
|
||||||
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
windowModelPtr_ = windowModel::New(dict, nSamples_);
|
||||||
|
|
||||||
|
return true;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
|
||||||
|
|
||||||
Foam::noiseModel::~noiseModel()
|
|
||||||
{}
|
|
||||||
|
|
||||||
|
|
||||||
// ************************************************************************* //
|
// ************************************************************************* //
|
||||||
|
|||||||
@ -124,7 +124,7 @@ protected:
|
|||||||
scalar checkUniformTimeStep
|
scalar checkUniformTimeStep
|
||||||
(
|
(
|
||||||
const scalarList& times
|
const scalarList& times
|
||||||
) const;
|
) const;
|
||||||
|
|
||||||
//- Find and return start time index
|
//- Find and return start time index
|
||||||
label findStartTimeIndex
|
label findStartTimeIndex
|
||||||
@ -139,7 +139,7 @@ public:
|
|||||||
//- Runtime type information
|
//- Runtime type information
|
||||||
TypeName("noiseModel");
|
TypeName("noiseModel");
|
||||||
|
|
||||||
//- Run time selection table
|
//- Run time selection table
|
||||||
declareRunTimeSelectionTable
|
declareRunTimeSelectionTable
|
||||||
(
|
(
|
||||||
autoPtr,
|
autoPtr,
|
||||||
@ -155,7 +155,7 @@ public:
|
|||||||
static autoPtr<noiseModel> New(const dictionary& dict);
|
static autoPtr<noiseModel> New(const dictionary& dict);
|
||||||
|
|
||||||
//- Constructor
|
//- Constructor
|
||||||
noiseModel(const dictionary& dict);
|
noiseModel(const dictionary& dict, const bool readFields = true);
|
||||||
|
|
||||||
//- Destructor
|
//- Destructor
|
||||||
virtual ~noiseModel();
|
virtual ~noiseModel();
|
||||||
@ -163,6 +163,9 @@ public:
|
|||||||
|
|
||||||
// Public Member Functions
|
// Public Member Functions
|
||||||
|
|
||||||
|
//- Read from dictionary
|
||||||
|
virtual bool read(const dictionary& dict);
|
||||||
|
|
||||||
//- Abstract call to calculate
|
//- Abstract call to calculate
|
||||||
virtual void calculate() = 0;
|
virtual void calculate() = 0;
|
||||||
};
|
};
|
||||||
|
|||||||
@ -43,14 +43,12 @@ addToRunTimeSelectionTable(noiseModel, pointNoise, dictionary);
|
|||||||
|
|
||||||
void pointNoise::filterTimeData
|
void pointNoise::filterTimeData
|
||||||
(
|
(
|
||||||
const Function1Types::CSV<scalar>& pData,
|
const scalarField& t0,
|
||||||
|
const scalarField& p0,
|
||||||
scalarField& t,
|
scalarField& t,
|
||||||
scalarField& p
|
scalarField& p
|
||||||
)
|
) const
|
||||||
{
|
{
|
||||||
const scalarField t0(pData.x());
|
|
||||||
const scalarField p0(pData.y());
|
|
||||||
|
|
||||||
DynamicList<scalar> tf(t0.size());
|
DynamicList<scalar> tf(t0.size());
|
||||||
DynamicList<scalar> pf(t0.size());
|
DynamicList<scalar> pf(t0.size());
|
||||||
|
|
||||||
@ -68,23 +66,15 @@ void pointNoise::filterTimeData
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
void pointNoise::processData(const Function1Types::CSV<scalar>& data)
|
||||||
|
|
||||||
void pointNoise::calculate()
|
|
||||||
{
|
{
|
||||||
// Point data only handled by master
|
Info<< "Reading data file " << data.fName() << endl;
|
||||||
if (!Pstream::master())
|
|
||||||
{
|
|
||||||
return;
|
|
||||||
}
|
|
||||||
|
|
||||||
Info<< "Reading data file" << endl;
|
const fileName& fNameBase = data.fName()(true);
|
||||||
|
|
||||||
Function1Types::CSV<scalar> pData("pressure", dict_, "Data");
|
|
||||||
|
|
||||||
// Time and pressure history data
|
// Time and pressure history data
|
||||||
scalarField t, p;
|
scalarField t, p;
|
||||||
filterTimeData(pData, t, p);
|
filterTimeData(data.x(), data.y(), t, p);
|
||||||
p *= rhoRef_;
|
p *= rhoRef_;
|
||||||
|
|
||||||
Info<< " read " << t.size() << " values" << nl << endl;
|
Info<< " read " << t.size() << " values" << nl << endl;
|
||||||
@ -96,7 +86,7 @@ void pointNoise::calculate()
|
|||||||
windowModelPtr_->validate(t.size());
|
windowModelPtr_->validate(t.size());
|
||||||
const windowModel& win = windowModelPtr_();
|
const windowModel& win = windowModelPtr_();
|
||||||
const scalar deltaf = 1.0/(deltaT*win.nSamples());
|
const scalar deltaf = 1.0/(deltaT*win.nSamples());
|
||||||
fileName outDir(fileName("postProcessing")/"noise"/typeName);
|
fileName outDir(fileName("postProcessing")/"noise"/typeName/fNameBase);
|
||||||
|
|
||||||
// Create the fft
|
// Create the fft
|
||||||
noiseFFT nfft(deltaT, p);
|
noiseFFT nfft(deltaT, p);
|
||||||
@ -185,12 +175,45 @@ void pointNoise::calculate()
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
void pointNoise::calculate()
|
||||||
|
{
|
||||||
|
// Point data only handled by master
|
||||||
|
if (!Pstream::master())
|
||||||
|
{
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
if (inputFileNames_.size())
|
||||||
|
{
|
||||||
|
forAll(inputFileNames_, i)
|
||||||
|
{
|
||||||
|
const fileName fName = inputFileNames_[i].expand();
|
||||||
|
Function1Types::CSV<scalar> data("pressure", dict_, "Data", fName);
|
||||||
|
processData(data);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
Function1Types::CSV<scalar> data("pressure", dict_, "Data");
|
||||||
|
processData(data);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
pointNoise::pointNoise(const dictionary& dict)
|
pointNoise::pointNoise(const dictionary& dict, const bool readFields)
|
||||||
:
|
:
|
||||||
noiseModel(dict)
|
noiseModel(dict)
|
||||||
{}
|
{
|
||||||
|
if (readFields)
|
||||||
|
{
|
||||||
|
read(dict);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||||
@ -199,6 +222,18 @@ pointNoise::~pointNoise()
|
|||||||
{}
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
bool pointNoise::read(const dictionary& dict)
|
||||||
|
{
|
||||||
|
if (noiseModel::read(dict))
|
||||||
|
{
|
||||||
|
dict.readIfPresent("inputFiles", inputFileNames_);
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
} // End namespace noiseModels
|
} // End namespace noiseModels
|
||||||
|
|||||||
@ -76,8 +76,8 @@ SeeAlso
|
|||||||
|
|
||||||
\*---------------------------------------------------------------------------*/
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
#ifndef pointNoise_H
|
#ifndef noiseModels_pointNoise_H
|
||||||
#define pointNoise_H
|
#define noiseModels_pointNoise_H
|
||||||
|
|
||||||
#include "noiseModel.H"
|
#include "noiseModel.H"
|
||||||
#include "CSV.H"
|
#include "CSV.H"
|
||||||
@ -100,14 +100,24 @@ class pointNoise
|
|||||||
|
|
||||||
protected:
|
protected:
|
||||||
|
|
||||||
|
// Protected data
|
||||||
|
|
||||||
|
//- Input file names - optional
|
||||||
|
List<fileName> inputFileNames_;
|
||||||
|
|
||||||
|
|
||||||
// Protected Member Functions
|
// Protected Member Functions
|
||||||
|
|
||||||
void filterTimeData
|
void filterTimeData
|
||||||
(
|
(
|
||||||
const Function1Types::CSV<scalar>& pData,
|
const scalarField& t0,
|
||||||
|
const scalarField& p0,
|
||||||
scalarField& t,
|
scalarField& t,
|
||||||
scalarField& p
|
scalarField& p
|
||||||
);
|
) const;
|
||||||
|
|
||||||
|
//- Process the CSV data
|
||||||
|
void processData(const Function1Types::CSV<scalar>& data);
|
||||||
|
|
||||||
|
|
||||||
public:
|
public:
|
||||||
@ -116,7 +126,7 @@ public:
|
|||||||
TypeName("pointNoise");
|
TypeName("pointNoise");
|
||||||
|
|
||||||
//- Constructor
|
//- Constructor
|
||||||
pointNoise(const dictionary& dict);
|
pointNoise(const dictionary& dict, const bool readFields = true);
|
||||||
|
|
||||||
//- Destructor
|
//- Destructor
|
||||||
virtual ~pointNoise();
|
virtual ~pointNoise();
|
||||||
@ -124,6 +134,9 @@ public:
|
|||||||
|
|
||||||
// Public Member Functions
|
// Public Member Functions
|
||||||
|
|
||||||
|
//- Read from dictionary
|
||||||
|
virtual bool read(const dictionary& dict);
|
||||||
|
|
||||||
//- Calculate
|
//- Calculate
|
||||||
virtual void calculate();
|
virtual void calculate();
|
||||||
|
|
||||||
|
|||||||
@ -44,12 +44,9 @@ addToRunTimeSelectionTable(noiseModel, surfaceNoise, dictionary);
|
|||||||
|
|
||||||
// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
|
// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
|
||||||
|
|
||||||
void surfaceNoise::initialise(const dictionary& dict)
|
void surfaceNoise::initialise(const fileName& fName)
|
||||||
{
|
{
|
||||||
dict.lookup("inputFile") >> inputFileName_;
|
Info<< "Reading data file " << fName << endl;
|
||||||
inputFileName_.expand();
|
|
||||||
|
|
||||||
dict.readIfPresent("fftWriteInterval", fftWriteInterval_);
|
|
||||||
|
|
||||||
label nAvailableTimes = 0;
|
label nAvailableTimes = 0;
|
||||||
|
|
||||||
@ -57,17 +54,15 @@ void surfaceNoise::initialise(const dictionary& dict)
|
|||||||
if (Pstream::master())
|
if (Pstream::master())
|
||||||
{
|
{
|
||||||
// Create the surface reader
|
// Create the surface reader
|
||||||
const word readerType(dict.lookup("reader"));
|
readerPtr_ = surfaceReader::New(readerType_, fName);
|
||||||
readerPtr_.reset(surfaceReader::New(readerType, inputFileName_).ptr());
|
|
||||||
|
|
||||||
// Find the index of the pressure data
|
// Find the index of the pressure data
|
||||||
const word pName(dict.lookupOrDefault<word>("p", "p"));
|
|
||||||
const List<word> fieldNames(readerPtr_->fieldNames(0));
|
const List<word> fieldNames(readerPtr_->fieldNames(0));
|
||||||
pIndex_ = findIndex(fieldNames, pName);
|
pIndex_ = findIndex(fieldNames, pName_);
|
||||||
if (pIndex_ == -1)
|
if (pIndex_ == -1)
|
||||||
{
|
{
|
||||||
FatalErrorInFunction
|
FatalErrorInFunction
|
||||||
<< "Unable to find pressure field name " << pName
|
<< "Unable to find pressure field name " << pName_
|
||||||
<< " in list of available fields: " << fieldNames
|
<< " in list of available fields: " << fieldNames
|
||||||
<< exit(FatalError);
|
<< exit(FatalError);
|
||||||
}
|
}
|
||||||
@ -76,12 +71,6 @@ void surfaceNoise::initialise(const dictionary& dict)
|
|||||||
// - Could be done later, but since this utility can process a lot of
|
// - Could be done later, but since this utility can process a lot of
|
||||||
// data we can ensure that the user-input is correct prior to doing
|
// data we can ensure that the user-input is correct prior to doing
|
||||||
// the heavy lifting
|
// the heavy lifting
|
||||||
const word writerType(dict.lookup("writer"));
|
|
||||||
dictionary optDict
|
|
||||||
(
|
|
||||||
dict.subOrEmptyDict("writeOptions").subOrEmptyDict(writerType)
|
|
||||||
);
|
|
||||||
writerPtr_.reset(surfaceWriter::New(writerType, optDict).ptr());
|
|
||||||
|
|
||||||
// Set the time range
|
// Set the time range
|
||||||
const instantList allTimes = readerPtr_->times();
|
const instantList allTimes = readerPtr_->times();
|
||||||
@ -404,18 +393,25 @@ Foam::scalar surfaceNoise::surfaceAverage
|
|||||||
|
|
||||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
surfaceNoise::surfaceNoise(const dictionary& dict)
|
surfaceNoise::surfaceNoise(const dictionary& dict, const bool readFields)
|
||||||
:
|
:
|
||||||
noiseModel(dict),
|
noiseModel(dict, false),
|
||||||
inputFileName_("unknown-inputFile"),
|
inputFileNames_(),
|
||||||
|
pName_("p"),
|
||||||
pIndex_(0),
|
pIndex_(0),
|
||||||
times_(),
|
times_(),
|
||||||
deltaT_(0),
|
deltaT_(0),
|
||||||
startTimeIndex_(0),
|
startTimeIndex_(0),
|
||||||
nFace_(0),
|
nFace_(0),
|
||||||
fftWriteInterval_(1)
|
fftWriteInterval_(1),
|
||||||
|
readerType_(word::null),
|
||||||
|
readerPtr_(nullptr),
|
||||||
|
writerPtr_(nullptr)
|
||||||
{
|
{
|
||||||
initialise(dict);
|
if (readFields)
|
||||||
|
{
|
||||||
|
read(dict);
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
@ -427,271 +423,307 @@ surfaceNoise::~surfaceNoise()
|
|||||||
|
|
||||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
bool surfaceNoise::read(const dictionary& dict)
|
||||||
|
{
|
||||||
|
if (noiseModel::read(dict))
|
||||||
|
{
|
||||||
|
dict.lookup("inputFiles") >> inputFileNames_;
|
||||||
|
dict.readIfPresent("fftWriteInterval", fftWriteInterval_);
|
||||||
|
dict.readIfPresent("p", pName_);
|
||||||
|
|
||||||
|
dict.lookup("reader") >> readerType_;
|
||||||
|
|
||||||
|
word writerType(dict.lookup("writer"));
|
||||||
|
dictionary optDict
|
||||||
|
(
|
||||||
|
dict.subOrEmptyDict("writeOptions").subOrEmptyDict(writerType)
|
||||||
|
);
|
||||||
|
|
||||||
|
writerPtr_ = surfaceWriter::New(writerType, optDict);
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
void surfaceNoise::calculate()
|
void surfaceNoise::calculate()
|
||||||
{
|
{
|
||||||
// Container for pressure time history data per face
|
forAll(inputFileNames_, i)
|
||||||
List<scalarField> pData;
|
|
||||||
|
|
||||||
// Processor procFaceOffsets
|
|
||||||
labelList procFaceOffset;
|
|
||||||
if (Pstream::parRun())
|
|
||||||
{
|
{
|
||||||
const label nProcs = Pstream::nProcs();
|
fileName fName = inputFileNames_[i];
|
||||||
const label nFacePerProc = floor(nFace_/nProcs) + 1;
|
|
||||||
|
|
||||||
procFaceOffset.setSize(nProcs + 1, 0);
|
initialise(fName.expand());
|
||||||
for (label i = 1; i < procFaceOffset.size(); i++)
|
|
||||||
|
// Container for pressure time history data per face
|
||||||
|
List<scalarField> pData;
|
||||||
|
|
||||||
|
// Processor procFaceOffsets
|
||||||
|
labelList procFaceOffset;
|
||||||
|
if (Pstream::parRun())
|
||||||
{
|
{
|
||||||
procFaceOffset[i] = min(i*nFacePerProc, nFace_);
|
const label nProcs = Pstream::nProcs();
|
||||||
|
const label nFacePerProc = floor(nFace_/nProcs) + 1;
|
||||||
|
|
||||||
|
procFaceOffset.setSize(nProcs + 1, 0);
|
||||||
|
for (label i = 1; i < procFaceOffset.size(); i++)
|
||||||
|
{
|
||||||
|
procFaceOffset[i] = min(i*nFacePerProc, nFace_);
|
||||||
|
}
|
||||||
}
|
}
|
||||||
}
|
else
|
||||||
else
|
|
||||||
{
|
|
||||||
procFaceOffset.setSize(1, nFace_);
|
|
||||||
}
|
|
||||||
|
|
||||||
// Read pressure data from file
|
|
||||||
readSurfaceData(procFaceOffset, pData);
|
|
||||||
|
|
||||||
// Process the pressure data, and store results as surface values per
|
|
||||||
// frequency so that it can be output using the surface writer
|
|
||||||
|
|
||||||
Info<< "Creating noise FFTs" << endl;
|
|
||||||
|
|
||||||
// Storage for FFT data
|
|
||||||
const label nLocalFace = pData.size();
|
|
||||||
const scalarField freq1(noiseFFT::frequencies(nSamples_, deltaT_));
|
|
||||||
const label nFFT = freq1.size()/fftWriteInterval_;
|
|
||||||
List<scalarField> surfPrmsf(nFFT);
|
|
||||||
List<scalarField> surfPSDf(nFFT);
|
|
||||||
forAll(surfPrmsf, freqI)
|
|
||||||
{
|
|
||||||
surfPrmsf[freqI].setSize(nLocalFace);
|
|
||||||
surfPSDf[freqI].setSize(nLocalFace);
|
|
||||||
}
|
|
||||||
|
|
||||||
// Storage for 1/3 octave data
|
|
||||||
labelList octave13BandIDs;
|
|
||||||
scalarField octave13FreqCentre;
|
|
||||||
noiseFFT::octaveBandInfo
|
|
||||||
(
|
|
||||||
freq1,
|
|
||||||
fLower_,
|
|
||||||
fUpper_,
|
|
||||||
3,
|
|
||||||
octave13BandIDs,
|
|
||||||
octave13FreqCentre
|
|
||||||
);
|
|
||||||
|
|
||||||
label bandSize = 0;
|
|
||||||
if (octave13BandIDs.empty())
|
|
||||||
{
|
|
||||||
WarningInFunction
|
|
||||||
<< "Ocatve band calculation failed (zero sized). "
|
|
||||||
<< "please check your input data"
|
|
||||||
<< endl;
|
|
||||||
}
|
|
||||||
else
|
|
||||||
{
|
|
||||||
bandSize = octave13BandIDs.size() - 1;
|
|
||||||
}
|
|
||||||
|
|
||||||
List<scalarField> surfPSD13f(bandSize);
|
|
||||||
List<scalarField> surfPrms13f2(bandSize);
|
|
||||||
forAll(surfPSD13f, freqI)
|
|
||||||
{
|
|
||||||
surfPSD13f[freqI].setSize(nLocalFace);
|
|
||||||
surfPrms13f2[freqI].setSize(nLocalFace);
|
|
||||||
}
|
|
||||||
|
|
||||||
const windowModel& win = windowModelPtr_();
|
|
||||||
|
|
||||||
forAll(pData, faceI)
|
|
||||||
{
|
|
||||||
const scalarField& p = pData[faceI];
|
|
||||||
|
|
||||||
noiseFFT nfft(deltaT_, p);
|
|
||||||
graph Prmsf(nfft.RMSmeanPf(win));
|
|
||||||
graph PSDf(nfft.PSDf(win));
|
|
||||||
|
|
||||||
// Store the frequency results in slot for face of surface
|
|
||||||
forAll(surfPrmsf, i)
|
|
||||||
{
|
{
|
||||||
label freqI = (i + 1)*fftWriteInterval_ - 1;
|
procFaceOffset.setSize(1, nFace_);
|
||||||
surfPrmsf[i][faceI] = Prmsf.y()[freqI];
|
|
||||||
surfPSDf[i][faceI] = PSDf.y()[freqI];
|
|
||||||
}
|
}
|
||||||
|
|
||||||
// PSD [Pa^2/Hz]
|
// Read pressure data from file
|
||||||
graph PSD13f(nfft.octaves(PSDf, octave13BandIDs, false));
|
readSurfaceData(procFaceOffset, pData);
|
||||||
|
|
||||||
// Integrated PSD = P(rms)^2 [Pa^2]
|
// Process the pressure data, and store results as surface values per
|
||||||
graph Prms13f2(nfft.octaves(PSDf, octave13BandIDs, true));
|
// frequency so that it can be output using the surface writer
|
||||||
|
|
||||||
// Store the 1/3 octave results in slot for face of surface
|
Info<< "Creating noise FFTs" << endl;
|
||||||
|
|
||||||
|
// Storage for FFT data
|
||||||
|
const label nLocalFace = pData.size();
|
||||||
|
const scalarField freq1(noiseFFT::frequencies(nSamples_, deltaT_));
|
||||||
|
const label nFFT = freq1.size()/fftWriteInterval_;
|
||||||
|
List<scalarField> surfPrmsf(nFFT);
|
||||||
|
List<scalarField> surfPSDf(nFFT);
|
||||||
|
forAll(surfPrmsf, freqI)
|
||||||
|
{
|
||||||
|
surfPrmsf[freqI].setSize(nLocalFace);
|
||||||
|
surfPSDf[freqI].setSize(nLocalFace);
|
||||||
|
}
|
||||||
|
|
||||||
|
// Storage for 1/3 octave data
|
||||||
|
labelList octave13BandIDs;
|
||||||
|
scalarField octave13FreqCentre;
|
||||||
|
noiseFFT::octaveBandInfo
|
||||||
|
(
|
||||||
|
freq1,
|
||||||
|
fLower_,
|
||||||
|
fUpper_,
|
||||||
|
3,
|
||||||
|
octave13BandIDs,
|
||||||
|
octave13FreqCentre
|
||||||
|
);
|
||||||
|
|
||||||
|
label bandSize = 0;
|
||||||
|
if (octave13BandIDs.empty())
|
||||||
|
{
|
||||||
|
WarningInFunction
|
||||||
|
<< "Ocatve band calculation failed (zero sized). "
|
||||||
|
<< "please check your input data"
|
||||||
|
<< endl;
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
bandSize = octave13BandIDs.size() - 1;
|
||||||
|
}
|
||||||
|
|
||||||
|
List<scalarField> surfPSD13f(bandSize);
|
||||||
|
List<scalarField> surfPrms13f2(bandSize);
|
||||||
forAll(surfPSD13f, freqI)
|
forAll(surfPSD13f, freqI)
|
||||||
{
|
{
|
||||||
surfPSD13f[freqI][faceI] = PSD13f.y()[freqI];
|
surfPSD13f[freqI].setSize(nLocalFace);
|
||||||
surfPrms13f2[freqI][faceI] = Prms13f2.y()[freqI];
|
surfPrms13f2[freqI].setSize(nLocalFace);
|
||||||
}
|
}
|
||||||
}
|
|
||||||
|
|
||||||
// Output directory for graphs
|
const windowModel& win = windowModelPtr_();
|
||||||
fileName outDir(fileName("postProcessing")/"noise"/typeName);
|
|
||||||
|
|
||||||
const scalar deltaf = 1.0/(deltaT_*win.nSamples());
|
forAll(pData, faceI)
|
||||||
Info<< "Writing fft surface data" << endl;
|
|
||||||
{
|
|
||||||
scalarField PrmsfAve(surfPrmsf.size(), 0);
|
|
||||||
scalarField PSDfAve(surfPrmsf.size(), 0);
|
|
||||||
scalarField fOut(surfPrmsf.size(), 0);
|
|
||||||
|
|
||||||
forAll(surfPrmsf, i)
|
|
||||||
{
|
{
|
||||||
label freqI = i*fftWriteInterval_;
|
const scalarField& p = pData[faceI];
|
||||||
fOut[i] = freq1[freqI];
|
|
||||||
const word& fName = inputFileName_.name(true);
|
|
||||||
const word gName = "fft";
|
|
||||||
PrmsfAve[i] = writeSurfaceData
|
|
||||||
(
|
|
||||||
fName,
|
|
||||||
gName,
|
|
||||||
"Prmsf",
|
|
||||||
freq1[freqI],
|
|
||||||
surfPrmsf[i],
|
|
||||||
procFaceOffset
|
|
||||||
);
|
|
||||||
|
|
||||||
PSDfAve[i] = writeSurfaceData
|
noiseFFT nfft(deltaT_, p);
|
||||||
(
|
graph Prmsf(nfft.RMSmeanPf(win));
|
||||||
fName,
|
graph PSDf(nfft.PSDf(win));
|
||||||
gName,
|
|
||||||
"PSDf",
|
// Store the frequency results in slot for face of surface
|
||||||
freq1[freqI],
|
forAll(surfPrmsf, i)
|
||||||
surfPSDf[i],
|
{
|
||||||
procFaceOffset
|
label freqI = (i + 1)*fftWriteInterval_ - 1;
|
||||||
);
|
surfPrmsf[i][faceI] = Prmsf.y()[freqI];
|
||||||
writeSurfaceData
|
surfPSDf[i][faceI] = PSDf.y()[freqI];
|
||||||
(
|
}
|
||||||
fName,
|
|
||||||
gName,
|
// PSD [Pa^2/Hz]
|
||||||
"PSD",
|
graph PSD13f(nfft.octaves(PSDf, octave13BandIDs, false));
|
||||||
freq1[freqI],
|
|
||||||
noiseFFT::PSD(surfPSDf[i]),
|
// Integrated PSD = P(rms)^2 [Pa^2]
|
||||||
procFaceOffset
|
graph Prms13f2(nfft.octaves(PSDf, octave13BandIDs, true));
|
||||||
);
|
|
||||||
writeSurfaceData
|
// Store the 1/3 octave results in slot for face of surface
|
||||||
(
|
forAll(surfPSD13f, freqI)
|
||||||
fName,
|
{
|
||||||
gName,
|
surfPSD13f[freqI][faceI] = PSD13f.y()[freqI];
|
||||||
"SPL",
|
surfPrms13f2[freqI][faceI] = Prms13f2.y()[freqI];
|
||||||
freq1[freqI],
|
}
|
||||||
noiseFFT::SPL(surfPSDf[i]*deltaf),
|
|
||||||
procFaceOffset
|
|
||||||
);
|
|
||||||
}
|
}
|
||||||
|
|
||||||
graph Prmsfg
|
const word& fNameBase = fName.name(true);
|
||||||
|
|
||||||
|
// Output directory for graphs
|
||||||
|
fileName outDir
|
||||||
(
|
(
|
||||||
"Average Prms(f)",
|
fileName("postProcessing")/"noise"/typeName/fNameBase
|
||||||
"f [Hz]",
|
|
||||||
"P(f) [Pa]",
|
|
||||||
fOut,
|
|
||||||
PrmsfAve
|
|
||||||
);
|
);
|
||||||
Prmsfg.write(outDir, graph::wordify(Prmsfg.title()), graphFormat_);
|
|
||||||
|
|
||||||
graph PSDfg
|
const scalar deltaf = 1.0/(deltaT_*win.nSamples());
|
||||||
(
|
Info<< "Writing fft surface data" << endl;
|
||||||
"Average PSD_f(f)",
|
|
||||||
"f [Hz]",
|
|
||||||
"PSD(f) [PaPa_Hz]",
|
|
||||||
fOut,
|
|
||||||
PSDfAve
|
|
||||||
);
|
|
||||||
PSDfg.write(outDir, graph::wordify(PSDfg.title()), graphFormat_);
|
|
||||||
|
|
||||||
graph PSDg
|
|
||||||
(
|
|
||||||
"Average PSD_dB_Hz(f)",
|
|
||||||
"f [Hz]",
|
|
||||||
"PSD(f) [dB_Hz]",
|
|
||||||
fOut,
|
|
||||||
noiseFFT::PSD(PSDfAve)
|
|
||||||
);
|
|
||||||
PSDg.write(outDir, graph::wordify(PSDg.title()), graphFormat_);
|
|
||||||
|
|
||||||
graph SPLg
|
|
||||||
(
|
|
||||||
"Average SPL_dB(f)",
|
|
||||||
"f [Hz]",
|
|
||||||
"SPL(f) [dB]",
|
|
||||||
fOut,
|
|
||||||
noiseFFT::SPL(PSDfAve*deltaf)
|
|
||||||
);
|
|
||||||
SPLg.write(outDir, graph::wordify(SPLg.title()), graphFormat_);
|
|
||||||
}
|
|
||||||
|
|
||||||
|
|
||||||
Info<< "Writing one-third octave surface data" << endl;
|
|
||||||
{
|
|
||||||
scalarField PSDfAve(surfPSD13f.size(), 0);
|
|
||||||
scalarField Prms13f2Ave(surfPSD13f.size(), 0);
|
|
||||||
|
|
||||||
forAll(surfPSD13f, i)
|
|
||||||
{
|
{
|
||||||
const word& fName = inputFileName_.name(true);
|
scalarField PrmsfAve(surfPrmsf.size(), 0);
|
||||||
const word gName = "oneThirdOctave";
|
scalarField PSDfAve(surfPrmsf.size(), 0);
|
||||||
PSDfAve[i] = writeSurfaceData
|
scalarField fOut(surfPrmsf.size(), 0);
|
||||||
(
|
|
||||||
fName,
|
|
||||||
gName,
|
|
||||||
"PSD13f",
|
|
||||||
octave13FreqCentre[i],
|
|
||||||
surfPSD13f[i],
|
|
||||||
procFaceOffset
|
|
||||||
);
|
|
||||||
writeSurfaceData
|
|
||||||
(
|
|
||||||
fName,
|
|
||||||
gName,
|
|
||||||
"PSD13",
|
|
||||||
octave13FreqCentre[i],
|
|
||||||
noiseFFT::PSD(surfPSD13f[i]),
|
|
||||||
procFaceOffset
|
|
||||||
);
|
|
||||||
writeSurfaceData
|
|
||||||
(
|
|
||||||
fName,
|
|
||||||
gName,
|
|
||||||
"SPL13",
|
|
||||||
octave13FreqCentre[i],
|
|
||||||
noiseFFT::SPL(surfPrms13f2[i]),
|
|
||||||
procFaceOffset
|
|
||||||
);
|
|
||||||
|
|
||||||
Prms13f2Ave[i] = surfaceAverage(surfPrms13f2[i], procFaceOffset);
|
forAll(surfPrmsf, i)
|
||||||
|
{
|
||||||
|
label freqI = i*fftWriteInterval_;
|
||||||
|
fOut[i] = freq1[freqI];
|
||||||
|
const word gName = "fft";
|
||||||
|
PrmsfAve[i] = writeSurfaceData
|
||||||
|
(
|
||||||
|
fNameBase,
|
||||||
|
gName,
|
||||||
|
"Prmsf",
|
||||||
|
freq1[freqI],
|
||||||
|
surfPrmsf[i],
|
||||||
|
procFaceOffset
|
||||||
|
);
|
||||||
|
|
||||||
|
PSDfAve[i] = writeSurfaceData
|
||||||
|
(
|
||||||
|
fNameBase,
|
||||||
|
gName,
|
||||||
|
"PSDf",
|
||||||
|
freq1[freqI],
|
||||||
|
surfPSDf[i],
|
||||||
|
procFaceOffset
|
||||||
|
);
|
||||||
|
writeSurfaceData
|
||||||
|
(
|
||||||
|
fNameBase,
|
||||||
|
gName,
|
||||||
|
"PSD",
|
||||||
|
freq1[freqI],
|
||||||
|
noiseFFT::PSD(surfPSDf[i]),
|
||||||
|
procFaceOffset
|
||||||
|
);
|
||||||
|
writeSurfaceData
|
||||||
|
(
|
||||||
|
fNameBase,
|
||||||
|
gName,
|
||||||
|
"SPL",
|
||||||
|
freq1[freqI],
|
||||||
|
noiseFFT::SPL(surfPSDf[i]*deltaf),
|
||||||
|
procFaceOffset
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
graph Prmsfg
|
||||||
|
(
|
||||||
|
"Average Prms(f)",
|
||||||
|
"f [Hz]",
|
||||||
|
"P(f) [Pa]",
|
||||||
|
fOut,
|
||||||
|
PrmsfAve
|
||||||
|
);
|
||||||
|
Prmsfg.write(outDir, graph::wordify(Prmsfg.title()), graphFormat_);
|
||||||
|
|
||||||
|
graph PSDfg
|
||||||
|
(
|
||||||
|
"Average PSD_f(f)",
|
||||||
|
"f [Hz]",
|
||||||
|
"PSD(f) [PaPa_Hz]",
|
||||||
|
fOut,
|
||||||
|
PSDfAve
|
||||||
|
);
|
||||||
|
PSDfg.write(outDir, graph::wordify(PSDfg.title()), graphFormat_);
|
||||||
|
|
||||||
|
graph PSDg
|
||||||
|
(
|
||||||
|
"Average PSD_dB_Hz(f)",
|
||||||
|
"f [Hz]",
|
||||||
|
"PSD(f) [dB_Hz]",
|
||||||
|
fOut,
|
||||||
|
noiseFFT::PSD(PSDfAve)
|
||||||
|
);
|
||||||
|
PSDg.write(outDir, graph::wordify(PSDg.title()), graphFormat_);
|
||||||
|
|
||||||
|
graph SPLg
|
||||||
|
(
|
||||||
|
"Average SPL_dB(f)",
|
||||||
|
"f [Hz]",
|
||||||
|
"SPL(f) [dB]",
|
||||||
|
fOut,
|
||||||
|
noiseFFT::SPL(PSDfAve*deltaf)
|
||||||
|
);
|
||||||
|
SPLg.write(outDir, graph::wordify(SPLg.title()), graphFormat_);
|
||||||
}
|
}
|
||||||
|
|
||||||
graph PSD13g
|
|
||||||
(
|
|
||||||
"Average PSD13_dB_Hz(fm)",
|
|
||||||
"fm [Hz]",
|
|
||||||
"PSD(fm) [dB_Hz]",
|
|
||||||
octave13FreqCentre,
|
|
||||||
noiseFFT::PSD(PSDfAve)
|
|
||||||
);
|
|
||||||
PSD13g.write(outDir, graph::wordify(PSD13g.title()), graphFormat_);
|
|
||||||
|
|
||||||
graph SPL13g
|
Info<< "Writing one-third octave surface data" << endl;
|
||||||
(
|
{
|
||||||
"Average SPL13_dB(fm)",
|
scalarField PSDfAve(surfPSD13f.size(), 0);
|
||||||
"fm [Hz]",
|
scalarField Prms13f2Ave(surfPSD13f.size(), 0);
|
||||||
"SPL(fm) [dB]",
|
|
||||||
octave13FreqCentre,
|
forAll(surfPSD13f, i)
|
||||||
noiseFFT::SPL(Prms13f2Ave)
|
{
|
||||||
);
|
const word gName = "oneThirdOctave";
|
||||||
SPL13g.write(outDir, graph::wordify(SPL13g.title()), graphFormat_);
|
PSDfAve[i] = writeSurfaceData
|
||||||
|
(
|
||||||
|
fNameBase,
|
||||||
|
gName,
|
||||||
|
"PSD13f",
|
||||||
|
octave13FreqCentre[i],
|
||||||
|
surfPSD13f[i],
|
||||||
|
procFaceOffset
|
||||||
|
);
|
||||||
|
writeSurfaceData
|
||||||
|
(
|
||||||
|
fNameBase,
|
||||||
|
gName,
|
||||||
|
"PSD13",
|
||||||
|
octave13FreqCentre[i],
|
||||||
|
noiseFFT::PSD(surfPSD13f[i]),
|
||||||
|
procFaceOffset
|
||||||
|
);
|
||||||
|
writeSurfaceData
|
||||||
|
(
|
||||||
|
fNameBase,
|
||||||
|
gName,
|
||||||
|
"SPL13",
|
||||||
|
octave13FreqCentre[i],
|
||||||
|
noiseFFT::SPL(surfPrms13f2[i]),
|
||||||
|
procFaceOffset
|
||||||
|
);
|
||||||
|
|
||||||
|
Prms13f2Ave[i] =
|
||||||
|
surfaceAverage(surfPrms13f2[i], procFaceOffset);
|
||||||
|
}
|
||||||
|
|
||||||
|
graph PSD13g
|
||||||
|
(
|
||||||
|
"Average PSD13_dB_Hz(fm)",
|
||||||
|
"fm [Hz]",
|
||||||
|
"PSD(fm) [dB_Hz]",
|
||||||
|
octave13FreqCentre,
|
||||||
|
noiseFFT::PSD(PSDfAve)
|
||||||
|
);
|
||||||
|
PSD13g.write(outDir, graph::wordify(PSD13g.title()), graphFormat_);
|
||||||
|
|
||||||
|
graph SPL13g
|
||||||
|
(
|
||||||
|
"Average SPL13_dB(fm)",
|
||||||
|
"fm [Hz]",
|
||||||
|
"SPL(fm) [dB]",
|
||||||
|
octave13FreqCentre,
|
||||||
|
noiseFFT::SPL(Prms13f2Ave)
|
||||||
|
);
|
||||||
|
SPL13g.write(outDir, graph::wordify(SPL13g.title()), graphFormat_);
|
||||||
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|||||||
@ -53,7 +53,7 @@ Description
|
|||||||
}
|
}
|
||||||
|
|
||||||
// Input file
|
// Input file
|
||||||
inputFile "postProcessing/faceSource1/surface/patch/patch.case";
|
inputFiles ("postProcessing/faceSource1/surface/patch/patch.case");
|
||||||
|
|
||||||
// Write interval for FFT data, default = 1
|
// Write interval for FFT data, default = 1
|
||||||
fftWriteInterval 100;
|
fftWriteInterval 100;
|
||||||
@ -83,8 +83,8 @@ SeeAlso
|
|||||||
|
|
||||||
\*---------------------------------------------------------------------------*/
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
#ifndef surfaceNoise_H
|
#ifndef noiseModels_surfaceNoise_H
|
||||||
#define surfaceNoise_H
|
#define noiseModels_surfaceNoise_H
|
||||||
|
|
||||||
#include "noiseModel.H"
|
#include "noiseModel.H"
|
||||||
#include "labelList.H"
|
#include "labelList.H"
|
||||||
@ -115,8 +115,11 @@ protected:
|
|||||||
|
|
||||||
// Protected Data
|
// Protected Data
|
||||||
|
|
||||||
//- Input file name
|
//- Input file names
|
||||||
fileName inputFileName_;
|
List<fileName> inputFileNames_;
|
||||||
|
|
||||||
|
//- Name of pressure field
|
||||||
|
word pName_;
|
||||||
|
|
||||||
//- Index of pressure field in reader field list
|
//- Index of pressure field in reader field list
|
||||||
label pIndex_;
|
label pIndex_;
|
||||||
@ -138,17 +141,20 @@ protected:
|
|||||||
// result in a very large number of output files (1 per frequency)
|
// result in a very large number of output files (1 per frequency)
|
||||||
label fftWriteInterval_;
|
label fftWriteInterval_;
|
||||||
|
|
||||||
|
//- Reader type
|
||||||
|
word readerType_;
|
||||||
|
|
||||||
//- Pointer to the surface reader
|
//- Pointer to the surface reader
|
||||||
mutable autoPtr<surfaceReader> readerPtr_;
|
mutable autoPtr<surfaceReader> readerPtr_;
|
||||||
|
|
||||||
//- Pointer to the surface writer
|
//- Pointer to the surface writer
|
||||||
autoPtr<surfaceWriter> writerPtr_;
|
mutable autoPtr<surfaceWriter> writerPtr_;
|
||||||
|
|
||||||
|
|
||||||
// Protected Member Functions
|
// Protected Member Functions
|
||||||
|
|
||||||
//- Initialise
|
//- Initialise
|
||||||
void initialise(const dictionary& dict);
|
void initialise(const fileName& fName);
|
||||||
|
|
||||||
//- Read surface data
|
//- Read surface data
|
||||||
void readSurfaceData
|
void readSurfaceData
|
||||||
@ -183,7 +189,7 @@ public:
|
|||||||
TypeName("surfaceNoise");
|
TypeName("surfaceNoise");
|
||||||
|
|
||||||
//- Constructor
|
//- Constructor
|
||||||
surfaceNoise(const dictionary& dict);
|
surfaceNoise(const dictionary& dict, const bool readFields = true);
|
||||||
|
|
||||||
//- Destructor
|
//- Destructor
|
||||||
virtual ~surfaceNoise();
|
virtual ~surfaceNoise();
|
||||||
@ -191,6 +197,9 @@ public:
|
|||||||
|
|
||||||
// Public Member Functions
|
// Public Member Functions
|
||||||
|
|
||||||
|
//- Read from dictionary
|
||||||
|
virtual bool read(const dictionary& dict);
|
||||||
|
|
||||||
//- Calculate
|
//- Calculate
|
||||||
virtual void calculate();
|
virtual void calculate();
|
||||||
};
|
};
|
||||||
|
|||||||
Reference in New Issue
Block a user