All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
icarus_tool::SBNNoise Class Reference
Inheritance diagram for icarus_tool::SBNNoise:
icarus_tool::IGenNoise

Public Member Functions

 SBNNoise (const fhicl::ParameterSet &pset)
 
 ~SBNNoise ()
 
void configure (const fhicl::ParameterSet &pset) override
 
void nextEvent () override
 
void generateNoise (CLHEP::HepRandomEngine &noise_engine, CLHEP::HepRandomEngine &cornoise_engine, icarusutil::TimeVec &noise, detinfo::DetectorPropertiesData const &, double noise_factor, const geo::PlaneID &, unsigned int board) override
 

Private Types

using WaveformTools = icarus_signal_processing::WaveformTools< icarusutil::SigProcPrecision >
 

Private Member Functions

void GenerateCorrelatedNoise (CLHEP::HepRandomEngine &, icarusutil::TimeVec &, double, unsigned int)
 
void GenerateUncorrelatedNoise (CLHEP::HepRandomEngine &, icarusutil::TimeVec &, double)
 
void GenNoise (std::function< void(double[])> &, const icarusutil::TimeVec &, icarusutil::TimeVec &, float)
 
void ComputeRMSs ()
 
void makeHistograms ()
 
- Private Member Functions inherited from icarus_tool::IGenNoise
virtual ~IGenNoise () noexcept=default
 

Private Attributes

size_t fPlane
 
int fMedianNumBins
 
float fNoiseRand
 
long fCorrelatedSeed
 
long fUncorrelatedSeed
 
float fIncoherentNoiseFrac
 
bool fStoreHistograms
 
std::string fInputNoiseHistFileName
 
std::string fHistogramName
 
std::string fCorrelatedHistogramName
 
std::string fUncorrelatedHistogramName
 
std::string fCorrelatedRMSHistoName
 
std::string fUncorrelatedRMSHistoName
 
std::string fTotalRMSHistoName
 
WaveformTools fWaveformTool
 
icarusutil::TimeVec fCoherentNoiseVec
 
icarusutil::TimeVec fIncoherentNoiseVec
 
std::vector< double > fCorrAmpDistVec
 
double fIncoherentNoiseRMS
 
double fCoherentNoiseRMS
 
icarusutil::FrequencyVec fNoiseFrequencyVec
 
bool fNeedFirstSeed =true
 
TProfile * fInputNoiseHist
 
TProfile * fMedianNoiseHist
 
TProfile * fPeakNoiseHist
 
TProfile * fCorAmpDistHist
 
TH1D * corrRMSHistPtr
 
TH1D * uncorrRMSHistPtr
 
TH1D * totalRMSHistPtr
 
float totalRMS
 
Eigen::FFT< double > fEigenFFT
 

Detailed Description

Definition at line 40 of file SBNNoise_tool.cc.

Member Typedef Documentation

using icarus_tool::SBNNoise::WaveformTools = icarus_signal_processing::WaveformTools<icarusutil::SigProcPrecision>
private

Definition at line 82 of file SBNNoise_tool.cc.

Constructor & Destructor Documentation

icarus_tool::SBNNoise::SBNNoise ( const fhicl::ParameterSet &  pset)
explicit

Definition at line 121 of file SBNNoise_tool.cc.

122 {
123  // Recover the configuration of the tool from the input fhicl file and set up
124  configure(pset);
125 ComputeRMSs();
126 
127  // Output some histograms to catalogue what's been done
128  makeHistograms();
129 }
void configure(const fhicl::ParameterSet &pset) override
icarus_tool::SBNNoise::~SBNNoise ( )

Definition at line 131 of file SBNNoise_tool.cc.

132 {
133 }

Member Function Documentation

void icarus_tool::SBNNoise::ComputeRMSs ( )
private

Definition at line 342 of file SBNNoise_tool.cc.

343 {
344  // Let's get the rms we expect from the incoherent noise contribution to the waveform
345  // A couple of ways to do this, let's basically invert the frequency spectrum to
346  // produce a waveform and then get the rms from that
347  std::function<void (double[])> randGenFunc = [](double randArray[]){randArray[0]=0.5; randArray[1]=0.5;};
348 
349  icarusutil::TimeVec waveNoise(fIncoherentNoiseVec.size());
350  float scaleFactor = 1.;
351 
352  GenNoise(randGenFunc, fIncoherentNoiseVec, waveNoise, scaleFactor);
353 
354  // Now get the details...
355  double nSig(3.);
356  double mean,rmsTrunc;
357  int nTrunc;
358  int range;
359 
360  // Use the waveform tool to recover the full rms
361  fWaveformTool.getTruncatedMeanRMS(waveNoise, nSig, mean, fIncoherentNoiseRMS, rmsTrunc, nTrunc, range);
362 
363  // Do the same for the coherent term
364  GenNoise(randGenFunc, fCoherentNoiseVec, waveNoise, scaleFactor);
365 
366  fWaveformTool.getTruncatedMeanRMS(waveNoise, nSig, mean, fCoherentNoiseRMS, rmsTrunc, nTrunc, range);
367 
368  if (fStoreHistograms)
369  {
370  for(size_t idx = 0; idx < fCorrAmpDistVec.size(); idx++)
371  fCorAmpDistHist->Fill(idx, fCorrAmpDistVec[idx], 1.);
372  }
373 
374  totalRMS=totalRMSHistPtr->GetMean();
375 double rmsUnc=uncorrRMSHistPtr->GetMean();
376 //double rmsCorr=corrRMSHistPtr->GetMean();
378  return;
379 }
WaveformTools fWaveformTool
icarusutil::TimeVec fIncoherentNoiseVec
void GenNoise(std::function< void(double[])> &, const icarusutil::TimeVec &, icarusutil::TimeVec &, float)
std::vector< SigProcPrecision > TimeVec
std::vector< double > fCorrAmpDistVec
icarusutil::TimeVec fCoherentNoiseVec
double mean(const std::vector< short > &wf, size_t start, size_t nsample)
Definition: UtilFunc.cxx:13
void icarus_tool::SBNNoise::configure ( const fhicl::ParameterSet &  pset)
overridevirtual

Implements icarus_tool::IGenNoise.

Definition at line 135 of file SBNNoise_tool.cc.

136 {
137  // Recover the histogram used for noise generation
138  fPlane = pset.get< size_t >("Plane");
139  fMedianNumBins = pset.get< int >("MedianNumBins");
140  fNoiseRand = pset.get< float >("NoiseRand");
141  fCorrelatedSeed = pset.get< long >("CorrelatedSeed",1000);
142  fUncorrelatedSeed = pset.get< long >("UncorrelatedSeed",5000);
143  fIncoherentNoiseFrac = pset.get< float >("IncoherentNoiseFraction",0.5);
144  fStoreHistograms = pset.get< bool >("StoreHistograms");
145  fInputNoiseHistFileName = pset.get< std::string >("NoiseHistFileName");
146  fCorrelatedHistogramName = pset.get< std::string >("CorrelatedHistogramName");
147  fUncorrelatedHistogramName = pset.get< std::string >("UncorrelatedHistogramName");
148  fCorrelatedRMSHistoName = pset.get< std::string >("CorrelatedRMSHistoName");
149  fUncorrelatedRMSHistoName = pset.get< std::string >("UncorrelatedRMSHistoName");
150  fTotalRMSHistoName = pset.get< std::string >("TotalRMSHistoName");
151  // Initialize the work vector
152  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataForJob();
153  auto const detProp = art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataForJob(clockData);
154  fNoiseFrequencyVec.resize(detProp.NumberTimeSamples(),std::complex<float>(0.,0.));
155 
156  // Set up to input the histogram with the overall noise spectrum
157  std::string fullFileName;
158  std::string corrAmpFileName;
159 
160  cet::search_path searchPath("FW_SEARCH_PATH");
161  searchPath.find_file(fInputNoiseHistFileName, fullFileName);
162 
163  TFile inputFile(fullFileName.c_str(), "READ");
164 
165  if (!inputFile.IsOpen())
166  throw cet::exception("NoiseFromHist::configure") << "Unable to open input file: " << fInputNoiseHistFileName << std::endl;
167 
168  TH1D* corrHistPtr = (TH1D*)inputFile.Get(fCorrelatedHistogramName.c_str());
169  if (!corrHistPtr)
170  throw cet::exception("NoiseFromHist::configure") << "Unable to recover desired histogram: " << fCorrelatedHistogramName << std::endl;
171 
172  TH1D* uncorrHistPtr = (TH1D*)inputFile.Get(fUncorrelatedHistogramName.c_str());
173  if (!uncorrHistPtr)
174  throw cet::exception("NoiseFromHist::configure") << "Unable to recover desired histogram: " << fUncorrelatedHistogramName << std::endl;
175 
176  corrRMSHistPtr = (TH1D*)inputFile.Get(fCorrelatedRMSHistoName.c_str());
177  if (!corrRMSHistPtr)
178  throw cet::exception("NoiseFromHist::configure") << "Unable to recover desired histogram: " << fCorrelatedRMSHistoName << std::endl;
179 
180  uncorrRMSHistPtr = (TH1D*)inputFile.Get(fUncorrelatedRMSHistoName.c_str());
181  if (!uncorrRMSHistPtr)
182  throw cet::exception("NoiseFromHist::configure") << "Unable to recover desired histogram: " << fUncorrelatedRMSHistoName << std::endl;
183 
184  totalRMSHistPtr = (TH1D*)inputFile.Get(fTotalRMSHistoName.c_str());
185  if (!totalRMSHistPtr)
186  throw cet::exception("NoiseFromHist::configure") << "Unable to recover desired histogram: " << fTotalRMSHistoName << std::endl;
187 
188  // Close the input file
189  inputFile.Close();
190 std::cout << " corr nbins " << corrHistPtr->GetNbinsX() << std::endl;
191 std::cout << " uncorr nbins " << uncorrHistPtr->GetNbinsX() << std::endl;
192  fCoherentNoiseVec.resize(corrHistPtr->GetNbinsX(), 0.);
193 
194  // Recover the bin contents into local vectors
195  for(size_t histIdx = 0; histIdx < size_t(corrHistPtr->GetNbinsX()); histIdx++)
196  fCoherentNoiseVec[histIdx] = corrHistPtr->GetBinContent(histIdx+1);
197 fIncoherentNoiseVec.resize(uncorrHistPtr->GetNbinsX(), 0.);
198 
199  // Recover the bin contents into local vectors
200  for(size_t histIdx = 0; histIdx < size_t(uncorrHistPtr->GetNbinsX()); histIdx++)
201  fIncoherentNoiseVec[histIdx] = uncorrHistPtr->GetBinContent(histIdx+1);
202  // Should we store hists?
203 
204 std::cout << " after filling vectors " << std::endl;
205  if (fStoreHistograms)
206  {
207  // Define histograms
208  art::ServiceHandle<art::TFileService> tfs;
209 
210  art::TFileDirectory* histDirectory = tfs.get();
211 
212  // Make a directory for these histograms
213  art::TFileDirectory dir = histDirectory->mkdir(Form("CorNoisePlane%1zu",fPlane));
214 
215  float sampleRate = sampling_rate(clockData);
216  float readOutSize = detProp.ReadOutWindowSize();
217  float maxFreq = 1.e6 / (2. * sampleRate);
218  float minFreq = 1.e6 / (2. * sampleRate * readOutSize);
219  int numSamples = readOutSize / 2;
220 
221  fInputNoiseHist = dir.make<TProfile>("InNoise", ";freq(kHz)", numSamples, minFreq, maxFreq);
222  fMedianNoiseHist = dir.make<TProfile>("MedNoise", ";freq(kHz)", numSamples, minFreq, maxFreq);;
223  fPeakNoiseHist = dir.make<TProfile>("PeakNoise", ";freq(kHz)", numSamples, minFreq, maxFreq);;
224 
225  fCorAmpDistHist = dir.make<TProfile>("CorAmp", ";Motherboard", fCorrAmpDistVec.size(),0.,fCorrAmpDistVec.size());
226  }
227 
228  return;
229 }
std::string fCorrelatedRMSHistoName
icarusutil::TimeVec fIncoherentNoiseVec
std::string fUncorrelatedRMSHistoName
std::vector< double > fCorrAmpDistVec
std::string fUncorrelatedHistogramName
std::string fCorrelatedHistogramName
icarusutil::TimeVec fCoherentNoiseVec
tuple dir
Definition: dropbox.py:28
std::string fTotalRMSHistoName
std::string fInputNoiseHistFileName
art::ServiceHandle< art::TFileService > tfs
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
icarusutil::FrequencyVec fNoiseFrequencyVec
BEGIN_PROLOG could also be cout
auto const detProp
void icarus_tool::SBNNoise::GenerateCorrelatedNoise ( CLHEP::HepRandomEngine &  engine,
icarusutil::TimeVec noise,
double  noise_factor,
unsigned int  board 
)
private

Definition at line 293 of file SBNNoise_tool.cc.

294 {
295 
296  // Set the engine seed to the board being considered
297  engine.setSeed(fCorrelatedSeed+board,0);
298 
299  CLHEP::RandFlat noiseGen(engine,0,1);
300 
301  std::function<void (double[])> randGenFunc = [&noiseGen](double randArray[]){noiseGen.fireArray(2,randArray);};
302 
303  float fraction = std::sqrt(1. - fIncoherentNoiseFrac * fIncoherentNoiseFrac);
304  float scaleFactor = fraction * noise_factor / fCoherentNoiseRMS;
305 
306  GenNoise(randGenFunc, fCoherentNoiseVec, noise, scaleFactor);
307 
308 
309  return;
310 }
void GenNoise(std::function< void(double[])> &, const icarusutil::TimeVec &, icarusutil::TimeVec &, float)
icarusutil::TimeVec fCoherentNoiseVec
void icarus_tool::SBNNoise::generateNoise ( CLHEP::HepRandomEngine &  noise_engine,
CLHEP::HepRandomEngine &  cornoise_engine,
icarusutil::TimeVec noise,
detinfo::DetectorPropertiesData const &  ,
double  noise_factor,
const geo::PlaneID planeID,
unsigned int  board 
)
overridevirtual

Implements icarus_tool::IGenNoise.

Definition at line 239 of file SBNNoise_tool.cc.

246 {
247 noise_factor=totalRMS;
248 //std::cout << " generating noise " << std::endl;
249  // Define a couple of vectors to hold intermediate work
250  icarusutil::TimeVec noise_unc(noise.size(),0.);
251  icarusutil::TimeVec noise_corr(noise.size(),0.);
252 
253  // Make sure the work vector is size right with the output
254  if (fNoiseFrequencyVec.size() != noise.size()) fNoiseFrequencyVec.resize(noise.size(),std::complex<float>(0.,0.));
255  //std::cout << " generating uncorrelated noise " << std::endl;
256  // If applying incoherent noise call the generator
257  if (fIncoherentNoiseFrac > 0.) GenerateUncorrelatedNoise(engine_unc,noise_unc,noise_factor);
258 
259 // int board=channel/64;
260 //std::cout << " generating correlated noise " << std::endl;
261  // If applying coherent noise call the generator
262  if (fIncoherentNoiseFrac < 1.) GenerateCorrelatedNoise(engine_corr, noise_corr, noise_factor, board);
263  // std::cout << " summing noise " << std::endl;
264  // Take the noise as the simple sum of the two contributions
265  std::transform(noise_unc.begin(),noise_unc.end(),noise_corr.begin(),noise.begin(),std::plus<float>());
266 
267  return;
268 }
static constexpr Sample_t transform(Sample_t sample)
void GenerateCorrelatedNoise(CLHEP::HepRandomEngine &, icarusutil::TimeVec &, double, unsigned int)
std::vector< SigProcPrecision > TimeVec
void GenerateUncorrelatedNoise(CLHEP::HepRandomEngine &, icarusutil::TimeVec &, double)
icarusutil::FrequencyVec fNoiseFrequencyVec
void icarus_tool::SBNNoise::GenerateUncorrelatedNoise ( CLHEP::HepRandomEngine &  engine,
icarusutil::TimeVec noise,
double  noise_factor 
)
private

Definition at line 270 of file SBNNoise_tool.cc.

271 {
272  // Here we aim to produce a waveform consisting of incoherent noise
273  // Note that this is expected to be the dominate noise contribution
274  // Check for seed initialization
275  if (fNeedFirstSeed)
276  {
277  engine.setSeed(fUncorrelatedSeed,0);
278  fNeedFirstSeed = false;
279  }
280 
281  // Get the generator
282  CLHEP::RandFlat noiseGen(engine,0,1);
283 
284  std::function<void (double[])> randGenFunc = [&noiseGen](double randArray[]){noiseGen.fireArray(2,randArray);};
285 
286  float scaleFactor = fIncoherentNoiseFrac * noise_factor / fIncoherentNoiseRMS;
287 
288  GenNoise(randGenFunc, fIncoherentNoiseVec, noise, scaleFactor);
289 
290  return;
291 }
icarusutil::TimeVec fIncoherentNoiseVec
void GenNoise(std::function< void(double[])> &, const icarusutil::TimeVec &, icarusutil::TimeVec &, float)
void icarus_tool::SBNNoise::GenNoise ( std::function< void(double[])> &  gen,
const icarusutil::TimeVec freqDist,
icarusutil::TimeVec noise,
float  scaleFactor 
)
private

Definition at line 312 of file SBNNoise_tool.cc.

313 {
314  double rnd_corr[2] = {0.,0.};
315 
316  // Build out the frequency vector
317  for(size_t i=0; i< noise.size()/2; ++i)
318  {
319  // exponential noise spectrum
320  gen(rnd_corr);
321 
322  float pval = freqDist[i] * ((1-fNoiseRand) + 2 * fNoiseRand*rnd_corr[0]) * scaleFactor;
323  float phase = rnd_corr[1] * 2. * M_PI;
324 
325  std::complex<float> tc(pval*cos(phase),pval*sin(phase));
326 
327  fNoiseFrequencyVec[i] = tc;
328  }
329 
330  // inverse FFT MCSignal
331  fEigenFFT.inv(noise, fNoiseFrequencyVec);
332 
333  return;
334 }
icarusutil::FrequencyVec fNoiseFrequencyVec
Eigen::FFT< double > fEigenFFT
void icarus_tool::SBNNoise::makeHistograms ( )
private

Definition at line 336 of file SBNNoise_tool.cc.

337 {
338 
339  return;
340 }
void icarus_tool::SBNNoise::nextEvent ( )
overridevirtual

Implements icarus_tool::IGenNoise.

Definition at line 231 of file SBNNoise_tool.cc.

232 {
233  // We update the correlated seed because we want to see different noise event-by-event
234  fCorrelatedSeed = (333 * fCorrelatedSeed) % 900000000;
235 
236  return;
237 }

Member Data Documentation

TH1D* icarus_tool::SBNNoise::corrRMSHistPtr
private

Definition at line 107 of file SBNNoise_tool.cc.

double icarus_tool::SBNNoise::fCoherentNoiseRMS
private

Definition at line 94 of file SBNNoise_tool.cc.

icarusutil::TimeVec icarus_tool::SBNNoise::fCoherentNoiseVec
private

Definition at line 88 of file SBNNoise_tool.cc.

TProfile* icarus_tool::SBNNoise::fCorAmpDistHist
private

Definition at line 106 of file SBNNoise_tool.cc.

std::vector<double> icarus_tool::SBNNoise::fCorrAmpDistVec
private

Definition at line 91 of file SBNNoise_tool.cc.

std::string icarus_tool::SBNNoise::fCorrelatedHistogramName
private

Definition at line 76 of file SBNNoise_tool.cc.

std::string icarus_tool::SBNNoise::fCorrelatedRMSHistoName
private

Definition at line 78 of file SBNNoise_tool.cc.

long icarus_tool::SBNNoise::fCorrelatedSeed
private

Definition at line 70 of file SBNNoise_tool.cc.

Eigen::FFT<double> icarus_tool::SBNNoise::fEigenFFT
private

Definition at line 115 of file SBNNoise_tool.cc.

std::string icarus_tool::SBNNoise::fHistogramName
private

Definition at line 75 of file SBNNoise_tool.cc.

float icarus_tool::SBNNoise::fIncoherentNoiseFrac
private

Definition at line 72 of file SBNNoise_tool.cc.

double icarus_tool::SBNNoise::fIncoherentNoiseRMS
private

Definition at line 93 of file SBNNoise_tool.cc.

icarusutil::TimeVec icarus_tool::SBNNoise::fIncoherentNoiseVec
private

Definition at line 89 of file SBNNoise_tool.cc.

TProfile* icarus_tool::SBNNoise::fInputNoiseHist
private

Definition at line 103 of file SBNNoise_tool.cc.

std::string icarus_tool::SBNNoise::fInputNoiseHistFileName
private

Definition at line 74 of file SBNNoise_tool.cc.

TProfile* icarus_tool::SBNNoise::fMedianNoiseHist
private

Definition at line 104 of file SBNNoise_tool.cc.

int icarus_tool::SBNNoise::fMedianNumBins
private

Definition at line 68 of file SBNNoise_tool.cc.

bool icarus_tool::SBNNoise::fNeedFirstSeed =true
private

Definition at line 100 of file SBNNoise_tool.cc.

icarusutil::FrequencyVec icarus_tool::SBNNoise::fNoiseFrequencyVec
private

Definition at line 97 of file SBNNoise_tool.cc.

float icarus_tool::SBNNoise::fNoiseRand
private

Definition at line 69 of file SBNNoise_tool.cc.

TProfile* icarus_tool::SBNNoise::fPeakNoiseHist
private

Definition at line 105 of file SBNNoise_tool.cc.

size_t icarus_tool::SBNNoise::fPlane
private

Definition at line 67 of file SBNNoise_tool.cc.

bool icarus_tool::SBNNoise::fStoreHistograms
private

Definition at line 73 of file SBNNoise_tool.cc.

std::string icarus_tool::SBNNoise::fTotalRMSHistoName
private

Definition at line 80 of file SBNNoise_tool.cc.

std::string icarus_tool::SBNNoise::fUncorrelatedHistogramName
private

Definition at line 77 of file SBNNoise_tool.cc.

std::string icarus_tool::SBNNoise::fUncorrelatedRMSHistoName
private

Definition at line 79 of file SBNNoise_tool.cc.

long icarus_tool::SBNNoise::fUncorrelatedSeed
private

Definition at line 71 of file SBNNoise_tool.cc.

WaveformTools icarus_tool::SBNNoise::fWaveformTool
private

Definition at line 84 of file SBNNoise_tool.cc.

float icarus_tool::SBNNoise::totalRMS
private

Definition at line 111 of file SBNNoise_tool.cc.

TH1D* icarus_tool::SBNNoise::totalRMSHistPtr
private

Definition at line 109 of file SBNNoise_tool.cc.

TH1D* icarus_tool::SBNNoise::uncorrRMSHistPtr
private

Definition at line 108 of file SBNNoise_tool.cc.


The documentation for this class was generated from the following file: