diff --git a/src/OpenFOAM/graph/graph.C b/src/OpenFOAM/graph/graph.C index a7698a20f450fb4190188e5c7fda59475fb0a96d..5eb1599d6e4658a6c76a4068dd6fa7632ca9b86f 100644 --- a/src/OpenFOAM/graph/graph.C +++ b/src/OpenFOAM/graph/graph.C @@ -28,6 +28,7 @@ License #include "IOmanip.H" #include "Pair.H" #include "OSspecific.H" +#include "SubField.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // @@ -160,6 +161,45 @@ Foam::scalarField& Foam::graph::y() } +void Foam::graph::setXRange(const scalar x0, const scalar x1) +{ + if (x1 < x0) + { + FatalErrorInFunction + << "When setting limits, x1 must be greater than x0" << nl + << " x0: " << x0 << nl + << " x1: " << x1 << nl + << abort(FatalError); + } + + label i0 = 0; + label i1 = 0; + + forAll(x_, i) + { + if (x_[i] < x0) + { + i0 = i + 1; + } + if (x_[i] < x1) + { + i1 = i; + } + } + + label nX = i1 - i0 + 1; + scalarField xNew(SubField<scalar>(x_, nX, i0)); + x_.transfer(xNew); + + forAllIter(HashPtrTable<curve>, *this, iter) + { + curve* c = iter(); + scalarField cNew(SubField<scalar>(*c, nX, i0)); + c->transfer(cNew); + } +} + + Foam::autoPtr<Foam::graph::writer> Foam::graph::writer::New ( const word& graphFormat diff --git a/src/OpenFOAM/graph/graph.H b/src/OpenFOAM/graph/graph.H index b3ff00c37d5938a826bc3cc889d88622fc78ba8f..d77d468d34cec2c52d0aa1b3e61e0f454c55e2c1 100644 --- a/src/OpenFOAM/graph/graph.H +++ b/src/OpenFOAM/graph/graph.H @@ -178,6 +178,9 @@ public: scalarField& y(); + // Limit the data for all axes based on x-limits + void setXRange(const scalar x0, const scalar x1); + // Write diff --git a/src/randomProcesses/noise/noiseModels/noiseModel/noiseModel.C b/src/randomProcesses/noise/noiseModels/noiseModel/noiseModel.C index ae80887fe89f8d15fff01d275f050dd86d837e3b..a0a088c2c59e2b74cef0a2f2271835854675229d 100644 --- a/src/randomProcesses/noise/noiseModels/noiseModel/noiseModel.C +++ b/src/randomProcesses/noise/noiseModels/noiseModel/noiseModel.C @@ -35,6 +35,26 @@ namespace Foam // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * // +void Foam::noiseModel::readWriteOption +( + const dictionary& dict, + const word& lookup, + bool& option +) const +{ + dict.readIfPresent(lookup, option); + + if (option) + { + Info<< " " << lookup << ": " << "yes" << endl; + } + else + { + Info<< " " << lookup << ": " << "no" << endl; + } +} + + Foam::scalar Foam::noiseModel::checkUniformTimeStep ( const scalarList& times @@ -92,6 +112,14 @@ Foam::label Foam::noiseModel::findStartTimeIndex } +Foam::fileName Foam::noiseModel::baseFileDir() const +{ + fileName baseDir("$FOAM_CASE"); + baseDir = baseDir.expand()/"postProcessing"/"noise"; + return baseDir; +} + + // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::noiseModel::noiseModel(const dictionary& dict, const bool readFields) @@ -101,9 +129,15 @@ Foam::noiseModel::noiseModel(const dictionary& dict, const bool readFields) nSamples_(65536), fLower_(25), fUpper_(10000), + customBounds_(false), startTime_(0), windowModelPtr_(), - graphFormat_("raw") + graphFormat_("raw"), + writePrmsf_(true), + writeSPL_(true), + writePSD_(true), + writePSDf_(true), + writeOctaves_(true) { if (readFields) { @@ -124,8 +158,11 @@ bool Foam::noiseModel::read(const dictionary& dict) { dict.readIfPresent("rhoRef", rhoRef_); dict.readIfPresent("N", nSamples_); - dict.readIfPresent("fl", fLower_); - dict.readIfPresent("fu", fUpper_); + customBounds_ = false; + if (dict.readIfPresent("fl", fLower_) || dict.readIfPresent("fu", fUpper_)) + { + customBounds_ = true; + } dict.readIfPresent("startTime", startTime_); dict.readIfPresent("graphFormat", graphFormat_); @@ -164,6 +201,14 @@ bool Foam::noiseModel::read(const dictionary& dict) } + Info<< " Write options:" << endl; + dictionary optDict(dict.subOrEmptyDict("writeOptions")); + readWriteOption(optDict, "writePrmsf", writePrmsf_); + readWriteOption(optDict, "writeSPL", writeSPL_); + readWriteOption(optDict, "writePSD", writePSD_); + readWriteOption(optDict, "writeOctaves", writeOctaves_); + + windowModelPtr_ = windowModel::New(dict, nSamples_); return true; diff --git a/src/randomProcesses/noise/noiseModels/noiseModel/noiseModel.H b/src/randomProcesses/noise/noiseModels/noiseModel/noiseModel.H index 89c41edb84c927554f505235f19fe83668d1ced5..1db4efb6135d88cb5023fc96ab2c00c941394515 100644 --- a/src/randomProcesses/noise/noiseModels/noiseModel/noiseModel.H +++ b/src/randomProcesses/noise/noiseModels/noiseModel/noiseModel.H @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2015-2016 OpenCFD Ltd. + \\ / A nd | Copyright (C) 2015-2017 OpenCFD Ltd. \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -35,6 +35,16 @@ Description fl 25; fu 25; startTime 0; + + // Optional write options dictionary + writeOptions + { + writePrmsf no; + writeSPL yes; + writePSD yes; + writePSDf no; + writeOctaves yes; + } \endverbatim where @@ -46,6 +56,11 @@ Description fu | Upper frequency bounds | no | 10000 startTime | Start time | no | 0 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 \endtable Note @@ -108,6 +123,9 @@ protected: //- Upper frequency limit, default = 10kHz scalar fUpper_; + //- Flagto indicate that custom frequenct bounds are being used + bool customBounds_; + //- Start time, default = 0s scalar startTime_; @@ -118,8 +136,34 @@ protected: word graphFormat_; + // Write options + + //- Write Prmsf; default = yes + bool writePrmsf_; + + //- Write SPL; default = yes + bool writeSPL_; + + //- Write PSD; default = yes + bool writePSD_; + + //- Write PSDf; default = yes + bool writePSDf_; + + //- Write writeOctaves; default = yes + bool writeOctaves_; + + // Protected Member Functions + //- Helper function to read write options and provide info feedback + void readWriteOption + ( + const dictionary& dict, + const word& lookup, + bool& option + ) const; + //- Check and return uniform time step scalar checkUniformTimeStep ( @@ -133,6 +177,9 @@ protected: const scalar startTime ) const; + //- Return the base output directory + fileName baseFileDir() const; + public: diff --git a/src/randomProcesses/noise/noiseModels/pointNoise/pointNoise.C b/src/randomProcesses/noise/noiseModels/pointNoise/pointNoise.C index b711a55bac82a9c5ccc60e3b361aa141ecf0260e..7af0c1c9e06317a449f466548c9c82e16fd7bed8 100644 --- a/src/randomProcesses/noise/noiseModels/pointNoise/pointNoise.C +++ b/src/randomProcesses/noise/noiseModels/pointNoise/pointNoise.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2015-2016 OpenCFD Ltd. + \\ / A nd | Copyright (C) 2015-2017 OpenCFD Ltd. \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -86,7 +86,7 @@ void pointNoise::processData(const Function1Types::CSV<scalar>& data) windowModelPtr_->validate(t.size()); const windowModel& win = windowModelPtr_(); const scalar deltaf = 1.0/(deltaT*win.nSamples()); - fileName outDir(fileName("postProcessing")/"noise"/typeName/fNameBase); + fileName outDir(baseFileDir()/typeName/fNameBase); // Create the fft noiseFFT nfft(deltaT, p); @@ -97,13 +97,27 @@ void pointNoise::processData(const Function1Types::CSV<scalar>& data) // RMS pressure [Pa] graph Prmsf(nfft.RMSmeanPf(win)); - Info<< " Creating graph for " << Prmsf.title() << endl; - Prmsf.write(outDir, graph::wordify(Prmsf.title()), graphFormat_); + if (customBounds_) + { + Prmsf.setXRange(fLower_, fUpper_); + } + if (writePrmsf_) + { + Info<< " Creating graph for " << Prmsf.title() << endl; + Prmsf.write(outDir, graph::wordify(Prmsf.title()), graphFormat_); + } // PSD [Pa^2/Hz] graph PSDf(nfft.PSDf(win)); - Info<< " Creating graph for " << PSDf.title() << endl; - PSDf.write(outDir, graph::wordify(PSDf.title()), graphFormat_); + if (customBounds_) + { + PSDf.setXRange(fLower_, fUpper_); + } + if (writePSDf_) + { + Info<< " Creating graph for " << PSDf.title() << endl; + PSDf.write(outDir, graph::wordify(PSDf.title()), graphFormat_); + } // PSD [dB/Hz] graph PSDg @@ -114,8 +128,12 @@ void pointNoise::processData(const Function1Types::CSV<scalar>& data) Prmsf.x(), noiseFFT::PSD(PSDf.y()) ); - Info<< " Creating graph for " << PSDg.title() << endl; - PSDg.write(outDir, graph::wordify(PSDg.title()), graphFormat_); + + if (writePSD_) + { + Info<< " Creating graph for " << PSDg.title() << endl; + PSDg.write(outDir, graph::wordify(PSDg.title()), graphFormat_); + } // SPL [dB] graph SPLg @@ -126,52 +144,59 @@ void pointNoise::processData(const Function1Types::CSV<scalar>& data) Prmsf.x(), noiseFFT::SPL(PSDf.y()*deltaf) ); - Info<< " Creating graph for " << SPLg.title() << endl; - SPLg.write(outDir, graph::wordify(SPLg.title()), graphFormat_); - - labelList octave13BandIDs; - scalarField octave13FreqCentre; - noiseFFT::octaveBandInfo - ( - Prmsf.x(), - fLower_, - fUpper_, - 3, - octave13BandIDs, - octave13FreqCentre - ); - - - // 1/3 octave data - // --------------- - // PSD [Pa^2/Hz] - graph PSD13f(nfft.octaves(PSDf, octave13BandIDs, false)); - - // Integrated PSD = P(rms)^2 [Pa^2] - graph Prms13f2(nfft.octaves(PSDf, octave13BandIDs, true)); - - graph PSD13g - ( - "PSD13_dB_Hz(fm)", - "fm [Hz]", - "PSD(fm) [dB_Hz]", - octave13FreqCentre, - noiseFFT::PSD(PSD13f.y()) - ); - Info<< " Creating graph for " << PSD13g.title() << endl; - PSD13g.write(outDir, graph::wordify(PSD13g.title()), graphFormat_); + if (writeSPL_) + { + Info<< " Creating graph for " << SPLg.title() << endl; + SPLg.write(outDir, graph::wordify(SPLg.title()), graphFormat_); + } - graph SPL13g - ( - "SPL13_dB(fm)", - "fm [Hz]", - "SPL(fm) [dB]", - octave13FreqCentre, - noiseFFT::SPL(Prms13f2.y()) - ); - Info<< " Creating graph for " << SPL13g.title() << endl; - SPL13g.write(outDir, graph::wordify(SPL13g.title()), graphFormat_); + if (writeOctaves_) + { + labelList octave13BandIDs; + scalarField octave13FreqCentre; + noiseFFT::octaveBandInfo + ( + Prmsf.x(), + fLower_, + fUpper_, + 3, + octave13BandIDs, + octave13FreqCentre + ); + + + // 1/3 octave data + // --------------- + + // PSD [Pa^2/Hz] + graph PSD13f(nfft.octaves(PSDf, octave13BandIDs, false)); + + // Integrated PSD = P(rms)^2 [Pa^2] + graph Prms13f2(nfft.octaves(PSDf, octave13BandIDs, true)); + + graph PSD13g + ( + "PSD13_dB_Hz(fm)", + "fm [Hz]", + "PSD(fm) [dB_Hz]", + octave13FreqCentre, + noiseFFT::PSD(PSD13f.y()) + ); + Info<< " Creating graph for " << PSD13g.title() << endl; + PSD13g.write(outDir, graph::wordify(PSD13g.title()), graphFormat_); + + graph SPL13g + ( + "SPL13_dB(fm)", + "fm [Hz]", + "SPL(fm) [dB]", + octave13FreqCentre, + noiseFFT::SPL(Prms13f2.y()) + ); + Info<< " Creating graph for " << SPL13g.title() << endl; + SPL13g.write(outDir, graph::wordify(SPL13g.title()), graphFormat_); + } } @@ -186,18 +211,15 @@ void pointNoise::calculate() } - if (inputFileNames_.size()) + forAll(inputFileNames_, i) { - forAll(inputFileNames_, i) + fileName fName = inputFileNames_[i]; + if (!fName.isAbsolute()) { - const fileName fName = inputFileNames_[i].expand(); - Function1Types::CSV<scalar> data("pressure", dict_, "Data", fName); - processData(data); + fName = "$FOAM_CASE"/fName; } - } - else - { - Function1Types::CSV<scalar> data("pressure", dict_, "Data"); + fName.expand(); + Function1Types::CSV<scalar> data("pressure", dict_, "Data", fName); processData(data); } } @@ -226,7 +248,15 @@ bool pointNoise::read(const dictionary& dict) { if (noiseModel::read(dict)) { - dict.readIfPresent("inputFiles", inputFileNames_); + if (!dict.readIfPresent("inputFiles", inputFileNames_)) + { + inputFileNames_.setSize(1); + + // Note: unsafe lookup of file name from pressureData sub-dict! + dict.subDict("pressureData").lookup("fileName") + >> inputFileNames_[0]; + } + return true; } diff --git a/src/randomProcesses/noise/noiseModels/surfaceNoise/surfaceNoise.C b/src/randomProcesses/noise/noiseModels/surfaceNoise/surfaceNoise.C index faed14eea6d9cb10253b93b666237274a2446954..23a661df2ab67ddebe660e636659155aa44d016f 100644 --- a/src/randomProcesses/noise/noiseModels/surfaceNoise/surfaceNoise.C +++ b/src/randomProcesses/noise/noiseModels/surfaceNoise/surfaceNoise.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2015-2016 OpenCFD Ltd. + \\ / A nd | Copyright (C) 2015-2017 OpenCFD Ltd. \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -231,15 +231,13 @@ Foam::scalar surfaceNoise::writeSurfaceData const word& title, const scalar freq, const scalarField& data, - const labelList& procFaceOffset + const labelList& procFaceOffset, + const bool writeSurface ) const { Info<< " processing " << title << " for frequency " << freq << endl; - fileName outDir - ( - fileName("postProcessing")/"noise"/groupName/Foam::name(freq) - ); + fileName outDir(baseFileDir()/groupName/Foam::name(freq)); if (Pstream::parRun()) { @@ -279,20 +277,23 @@ Foam::scalar surfaceNoise::writeSurfaceData } } - // could also have meshedSurface implement meshedSurf - fileName outFileName = writerPtr_->write - ( - outDir, - fName, - meshedSurfRef + // Could also have meshedSurface implement meshedSurf + if (writeSurface) + { + fileName outFileName = writerPtr_->write ( - surf.points(), - surf.surfFaces() - ), - title, - allData, - false - ); + outDir, + fName, + meshedSurfRef + ( + surf.points(), + surf.surfFaces() + ), + title, + allData, + false + ); + } // TO BE VERIFIED: area-averaged values // areaAverage = sum(allData*surf.magSf())/sum(surf.magSf()); @@ -306,20 +307,23 @@ Foam::scalar surfaceNoise::writeSurfaceData { const meshedSurface& surf = readerPtr_->geometry(); - // could also have meshedSurface implement meshedSurf - writerPtr_->write - ( - outDir, - fName, - meshedSurfRef + // Could also have meshedSurface implement meshedSurf + if (writeSurface) + { + writerPtr_->write ( - surf.points(), - surf.surfFaces() - ), - title, - data, - false - ); + outDir, + fName, + meshedSurfRef + ( + surf.points(), + surf.surfFaces() + ), + title, + data, + false + ); + } // TO BE VERIFIED: area-averaged values // return sum(data*surf.magSf())/sum(surf.magSf()); @@ -463,6 +467,11 @@ void surfaceNoise::calculate() { fileName fName = inputFileNames_[i]; + if (!fName.isAbsolute()) + { + fName = "$FOAM_CASE"/fName; + } + initialise(fName.expand()); // Container for pressure time history data per face @@ -575,19 +584,23 @@ void surfaceNoise::calculate() const word& fNameBase = fName.name(true); // Output directory for graphs - fileName outDir - ( - fileName("postProcessing")/"noise"/typeName/fNameBase - ); + fileName outDir(baseFileDir()/typeName/fNameBase); const scalar deltaf = 1.0/(deltaT_*win.nSamples()); Info<< "Writing fft surface data" << endl; { - scalarField PrmsfAve(surfPrmsf.size(), 0); - scalarField PSDfAve(surfPrmsf.size(), 0); - scalarField fOut(surfPrmsf.size(), 0); - - forAll(surfPrmsf, i) + // Determine frequency range of interest + // Note: freqencies have fixed interval, and are in the range + // 0 to (n-1)*deltaf + label f0 = ceil(fLower_/deltaf); + label f1 = floor(fUpper_/deltaf); + label nFreq = f1 - f0 + 1; + + scalarField PrmsfAve(nFreq, 0); + scalarField PSDfAve(nFreq, 0); + scalarField fOut(nFreq, 0); + + for (label i = f0; i <= f1; ++i) { label freqI = (i + 1)*fftWriteInterval_ - 1; fOut[i] = freq1[freqI]; @@ -599,7 +612,8 @@ void surfaceNoise::calculate() "Prmsf", freq1[freqI], surfPrmsf[i], - procFaceOffset + procFaceOffset, + writePrmsf_ ); PSDfAve[i] = writeSurfaceData @@ -609,7 +623,8 @@ void surfaceNoise::calculate() "PSDf", freq1[freqI], surfPSDf[i], - procFaceOffset + procFaceOffset, + writePSDf_ ); writeSurfaceData ( @@ -618,7 +633,8 @@ void surfaceNoise::calculate() "PSD", freq1[freqI], noiseFFT::PSD(surfPSDf[i]), - procFaceOffset + procFaceOffset, + writePSD_ ); writeSurfaceData ( @@ -627,7 +643,8 @@ void surfaceNoise::calculate() "SPL", freq1[freqI], noiseFFT::SPL(surfPSDf[i]*deltaf), - procFaceOffset + procFaceOffset, + writeSPL_ ); } @@ -688,7 +705,8 @@ void surfaceNoise::calculate() "PSD13f", octave13FreqCentre[i], surfPSD13f[i], - procFaceOffset + procFaceOffset, + writeOctaves_ ); writeSurfaceData ( @@ -697,7 +715,8 @@ void surfaceNoise::calculate() "PSD13", octave13FreqCentre[i], noiseFFT::PSD(surfPSD13f[i]), - procFaceOffset + procFaceOffset, + writeOctaves_ ); writeSurfaceData ( @@ -706,7 +725,8 @@ void surfaceNoise::calculate() "SPL13", octave13FreqCentre[i], noiseFFT::SPL(surfPrms13f2[i]), - procFaceOffset + procFaceOffset, + writeOctaves_ ); Prms13f2Ave[i] = diff --git a/src/randomProcesses/noise/noiseModels/surfaceNoise/surfaceNoise.H b/src/randomProcesses/noise/noiseModels/surfaceNoise/surfaceNoise.H index 165a1cefd77963b882d10c4d44285eaa79c6258e..b513af162fb77c87cb13edd867ccc1ea2188003b 100644 --- a/src/randomProcesses/noise/noiseModels/surfaceNoise/surfaceNoise.H +++ b/src/randomProcesses/noise/noiseModels/surfaceNoise/surfaceNoise.H @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2015-2016 OpenCFD Ltd. + \\ / A nd | Copyright (C) 2015-2017 OpenCFD Ltd. \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -71,7 +71,23 @@ Description { collateTimes true; } + + // Write Prmsf; default = yes + writePrmsf no; + + // Write SPL; default = yes + writeSPL yes; + + // Write PSD; default = yes + writePSD yes; + + // Write PSDf; default = yes + writePSDf no; + + // Write writeOctaves; default = yes + writeOctaves yes; } + \endverbatim SourceFiles @@ -172,7 +188,8 @@ protected: const word& title, const scalar freq, const scalarField& data, - const labelList& procFaceOffset + const labelList& procFaceOffset, + const bool writeSurface ) const; //- Calculate the area average value