Public Member Functions | Private Attributes | List of all members
pd_monitor::PDWaveform Class Reference
Inheritance diagram for pd_monitor::PDWaveform:
art::EDAnalyzer art::detail::Analyzer art::detail::LegacyModule art::Observer art::ModuleBase

Public Member Functions

 PDWaveform (fhicl::ParameterSet const &pset)
 
virtual ~PDWaveform ()
 
void setRootObjects ()
 
void beginJob ()
 
void endJob ()
 
void beginRun (const art::Run &run)
 
void reconfigure (fhicl::ParameterSet const &pset)
 
void analyze (const art::Event &evt)
 
- Public Member Functions inherited from art::EDAnalyzer
 EDAnalyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDAnalyzer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Analyzer
virtual ~Analyzer () noexcept
 
 Analyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 Analyzer (Table< Config > const &config)
 
void doBeginJob (SharedResources const &resources)
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Observer
 ~Observer () noexcept
 
 Observer (Observer const &)=delete
 
 Observer (Observer &&)=delete
 
Observeroperator= (Observer const &)=delete
 
Observeroperator= (Observer &&)=delete
 
void registerProducts (ProductDescriptions &, ModuleDescription const &)
 
void fillDescriptions (ModuleDescription const &)
 
fhicl::ParameterSetID selectorConfig () const
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Private Attributes

std::string fSSPInput
 
std::string fSSPInstance
 
int fSSP_m1 =10
 
int fSSP_m2 =10
 
int fSSP_i1 =40
 
int fSSP_i2 =1200
 
int fSSP_readout_pretrigger =50
 
int fSSP_disc_width =20
 
int fSSP_win =20
 
unsigned int fSSP_wfm_verbose =0
 
unsigned int fPDwaveform_fft =0
 
unsigned int fSSP_corrchan1 =205
 
unsigned int fSSP_corrchan2 =134
 
unsigned int fSSP_smoothing =0
 
float fSSP_rarenum =3.0
 
float fSSP_peak_sense =0.5
 
unsigned int fSSP_nump =10
 
unsigned int fSSP_peak =1
 
bool fIsSSP
 
art::ServiceHandle< geo::GeometryfGeom
 
TH2F * PDchanMean
 
TH2F * PDchanRMS
 
TH2F * PDchanRMSwide
 
TH2F * PDchanFFT
 
TH2F * PDchanPED
 
TH2F * PDchanMax
 
TH2F * PDchanMaxPed
 
TH2F * PDCalibInt
 
TH2F * PDchanCorr
 
TH2F * PDchanCorPerTrace [288]
 
TH2F * PDchanRawPerTrace [288]
 
TH2F * PDchanPEDRough [288]
 
TH1F * PDchanWaveInt [288]
 
TH1I * PDtrigs
 
TH1F * PDPEDhist
 
TH1F * PDchanThres
 
TH2F * PDchanCurPeak [288]
 

Additional Inherited Members

- Public Types inherited from art::EDAnalyzer
using WorkerType = WorkerT< EDAnalyzer >
 
using ModuleType = EDAnalyzer
 
- Protected Member Functions inherited from art::Observer
std::string const & processName () const
 
bool wantAllEvents () const noexcept
 
bool wantEvent (ScheduleID id, Event const &e) const
 
Handle< TriggerResultsgetTriggerResults (Event const &e) const
 
 Observer (fhicl::ParameterSet const &config)
 
 Observer (std::vector< std::string > const &select_paths, std::vector< std::string > const &reject_paths, fhicl::ParameterSet const &config)
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

Definition at line 74 of file PDWaveform_module.cc.

Constructor & Destructor Documentation

pd_monitor::PDWaveform::PDWaveform ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 132 of file PDWaveform_module.cc.

133  : EDAnalyzer(parameterSet){
134  reconfigure(parameterSet);
135  }
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
void reconfigure(fhicl::ParameterSet const &pset)
pd_monitor::PDWaveform::~PDWaveform ( )
virtual

Definition at line 174 of file PDWaveform_module.cc.

174 {}

Member Function Documentation

void pd_monitor::PDWaveform::analyze ( const art::Event evt)

Definition at line 220 of file PDWaveform_module.cc.

220  {
221 
222  MF_LOG_INFO("PDWaveform")
223  << "-------------------- Photodetector waveforms -------------------";
224 
225  // Get the data with the correct label and instance from the root file
226  //art::InputTag itag1("ssprawdecoder", "external");
228  auto RawSSP = event.getHandle< std::vector<raw::OpDetWaveform> >(itag1);
229 
230  // Make sure data is collected
231  try { RawSSP->size(); }
232  catch(std::exception e) { fIsSSP = false; }
233  fIsSSP=true;
234 
235  // Put data into a vector of pointers
236  std::vector< art::Ptr<raw::OpDetWaveform> > RawPulses;
237  if (fIsSSP) { art::fill_ptr_vector(RawPulses, RawSSP); }
238 
239  //channel correlation variables
240  double tschan1=0.0, tschan2=0.0, ampchan1=0.0,ampchan2=0.0;
241 
242  // Loop over waveforms
243  for(auto const & RawPulse : RawPulses) {
244 
245  // Get waveform from pointer
246  const raw::OpDetWaveform & PDdigit = *RawPulse;
247  unsigned int CurChannel = PDdigit.ChannelNumber();
248  PDtrigs->Fill(CurChannel);
249 
250  // Loop through individual waveform and print ADC's at each position
251  long int ADCval,sum=0,sum2=0,N=0,ADCMax=0;
252  double sumthres=0,sumpedsub=0;
253  int nBins = PDdigit.size();
254  int nfound = 0;
255  TH1F WfmHist("Waveform","Waveform",nBins,0,nBins),
256  WfmFFT("WfmFFT","WfmFFT",nBins,0,nBins);
257  TH1F *htemp = new TH1F("htemp","htemp",nBins,0,nBins);
258  /*long int presum=0,postsum=0;
259  long int runavg;
260  unsigned int forwin, backwin;
261  std::vector< long int > smoothwave;
262 
263  if(fSSP_smoothing){
264  if(fSSP_smoothing%2==0){
265  forwin = fSSP_smoothing/2;
266  backwin = fSSP_smoothing/2;
267  }
268  else{
269  forwin = (fSSP_smoothing-1)/2;
270  backwin = (fSSP_smoothing+1)/2;
271  }
272  for(size_t i=0;i<PDdigit.size();i++){
273  runavg = 0;
274  if(i<=backwin) {
275  presum = 0;
276  for(unsigned int j=0;j<i+forwin;j++) {
277  presum += PDdigit[j];
278  }
279  smoothwave.push_back(presum/(i+forwin));
280  }
281  else if(i>PDdigit.size()-forwin) {
282  postsum = 0;
283  for(unsigned int j=PDdigit.size()-i-backwin;j<PDdigit.size();j++) {
284  postsum += PDdigit[j];
285  }
286  smoothwave.push_back(postsum/(i+backwin));
287  }
288  else{
289  for(unsigned int k=i-backwin;k<i+forwin;k++) runavg += PDdigit[k];
290  smoothwave.push_back(runavg/fSSP_smoothing);
291  }
292  }
293  }*/
294 
295  if(fSSP_smoothing){
296  for (size_t i = 0; i < PDdigit.size(); i++) htemp->SetBinContent(i+1,PDdigit.at(i));
297  htemp->Smooth(fSSP_smoothing);
298  }
299  for (size_t i = 0; i < PDdigit.size(); i++) {
300  if(fSSP_smoothing) ADCval = htemp->GetBinContent(i+1);
301  else {
302  ADCval = PDdigit.at(i);
303  if(fSSP_peak) htemp->SetBinContent(i+1,ADCval);
304  }
305 
306  PDchanRawPerTrace[CurChannel]->Fill(i+1,ADCval);
307  ADCMax=std::max(ADCMax,ADCval);
308  N=N+1;
309  WfmHist.SetBinContent(i+1,ADCval);
310  sum+=ADCval;
311  sum2+=ADCval*ADCval;
312  PDchanPED->Fill(CurChannel,ADCval);
313 
314  if(i < static_cast<unsigned int>(fSSP_i1)) {
315  PDchanPEDRough[CurChannel]->Fill(i+1,ADCval);
316  sumthres += ADCval;
317  }
318 
319 
320  if(i > static_cast<unsigned int>(fSSP_i1)) PDchanCorPerTrace[CurChannel]->Fill(i+1,(ADCval-(sumthres/static_cast<float>(fSSP_i1))));
321 
322  }
323 
324  if(fSSP_peak){
325  TSpectrum *peakloc = new TSpectrum(fSSP_nump);
326  nfound = peakloc->Search(htemp,fSSP_rarenum,"",fSSP_peak_sense);
327  double *peaks = peakloc->GetPositionX();
328 
329 
330  for(int j=0;j<nfound;j++){
331  bool goodpeak = true;
332  for(int k=0;k<nfound;k++){
333  int peakpos1 = htemp->GetXaxis()->FindBin(peaks[j]);
334  int peakpos2 = htemp->GetXaxis()->FindBin(peaks[k]);
335  if((TMath::Abs(peakpos1-peakpos2)<fSSP_win) && (k!=j)){
336  if((htemp->GetBinContent(htemp->GetXaxis()->FindBin(peaks[j])) < htemp->GetBinContent(htemp->GetXaxis()->FindBin(peaks[k]))) && (fSSP_peak==1)) goodpeak=false;
337  if(fSSP_peak==2) goodpeak=false;
338  }
339  if((peakpos1 < (fSSP_i1+fSSP_disc_width))|| (peakpos1 > (nBins-fSSP_win))) goodpeak = false;
340  if(!goodpeak) break;
341  }
342 
343  if(goodpeak){
344  sumpedsub=0;
345  int intbin = htemp->GetXaxis()->FindBin(peaks[j]);
346  std::vector <int> base;
347  for(int l=intbin-static_cast<int>(fSSP_disc_width)-static_cast<int>(fSSP_i1);l<intbin-static_cast<int>(fSSP_disc_width);l++) base.push_back(htemp->GetBinContent(l));
348  double basemean = TMath::Mean(base.begin(),base.end());
349  double basestddev = TMath::StdDev(base.begin(),base.end());
350  base.clear();
351  if(((basemean+(fSSP_rarenum*basestddev)) < (htemp->GetBinContent(intbin))) && (basestddev < 4.0)) for(int m=intbin-static_cast<int>(fSSP_disc_width);m<intbin+fSSP_win;m++) sumpedsub += (htemp->GetBinContent(m))-basemean;
352  if(sumpedsub>1) {
353  PDchanWaveInt[CurChannel]->Fill(sumpedsub);
354  PDCalibInt->Fill(CurChannel,sumpedsub);
355  PDchanCurPeak[CurChannel]->Fill(htemp->GetBinContent(intbin)-basemean,sumpedsub);
356  }
357  }
358  }
359  htemp->Delete();
360  }
361 
362  float mean = (float)sum/(float)N;
363  float rms = sqrt((float)sum2/(float)N - mean*mean );
364  //std::cout <<"ADC Sum: " << sum << std::endl;
365 
366  PDchanMean->Fill(CurChannel,mean);
367  PDchanMax->Fill(CurChannel,ADCMax);
368  float thres = sumthres/static_cast<float>(fSSP_i1);
369  PDchanMaxPed->Fill(CurChannel,ADCMax-thres);
370 
371  //PDchanMin->Fill(CurChannel,ADCMin);
372  PDchanRMS->Fill(CurChannel,rms);
373  PDchanRMSwide->Fill(CurChannel,rms);
374 
375  if (event.event()%1==0){ // only do FFT on 100% of events (or average them)
376  WfmHist.FFT(&WfmFFT,"MAG");
377  WfmFFT.Scale(1.0/(float)N);
378  for (int freqBin = 2 ; freqBin < nBins/2 ; freqBin++){
379  // Just set the value
380  //PDchanFFT->SetBinContent(CurChannel,freqBin,WfmFFT.GetBinContent(freqBin));
381  PDchanFFT->Fill(CurChannel,((float)freqBin+0.5)*75.0/1000.0,WfmFFT.GetBinContent(freqBin));
382  }
383  }
384 
385  if(ADCMax-thres > 100.0){
386  if(CurChannel == fSSP_corrchan1){
387  tschan1=PDdigit.TimeStamp();
388  ampchan1=(ADCMax-thres);
389  }
390  if(CurChannel == fSSP_corrchan2){
391  tschan2=PDdigit.TimeStamp();
392  ampchan2=(ADCMax-thres);
393  }
394  if(tschan1 != 0.0 && tschan2 != 0.0 && (abs(tschan1-tschan2) < 100.0E-9)){
395  PDchanCorr->Fill(ampchan1,ampchan2);
396  tschan1=0.0;
397  tschan2=0.0;
398  ampchan1=0.0;
399  ampchan2=0.0;
400  }
401  }
402  }
403 
404  return;
405 
406  }
Channel_t ChannelNumber() const
Definition: OpDetWaveform.h:65
double rms(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:40
TimeStamp_t TimeStamp() const
Definition: OpDetWaveform.h:66
static QStrList * l
Definition: config.cpp:1044
T abs(T value)
const double e
static int max(int a, int b)
#define MF_LOG_INFO(category)
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
double mean(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:16
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
Event finding and building.
void pd_monitor::PDWaveform::beginJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 176 of file PDWaveform_module.cc.

176  {
177 
179 
180  PDchanMean = tFileService->make<TH2F>("Mean vs. Channel","Mean vs. Channel",288,0.,288.,1000,1000.,2000.);
181  PDchanMax = tFileService->make<TH2F>("Max vs. Channel","Max vs. Channel",288.,0.,288.,2500,1500.,4000.);
182  PDchanMaxPed = tFileService->make<TH2F>("Max - Pedestal vs. Channel","Max - Pedestal vs. Channel",288.,0.,288.,2500,0,2500.);
183  PDchanPED = tFileService->make<TH2F>("PedVals vs. Channel","PedVals vs. Channel",288,0.,288.,1000,1000.,2000.);
184  PDchanRMS = tFileService->make<TH2F>("RMS vs. Channel","RMS vs. Channel",288,0.,288.,100,0.,10.);
185  PDchanRMSwide = tFileService->make<TH2F>("Coarse RMS vs. Channel","Coarse RMS vs. Channel",288,0.,288.,100,0.,100.);
186  PDchanFFT = tFileService->make<TH2F>("FFTFreq vs. Channel","FFTFreq vs. Channel",288,0.,288.,1000,0.,75.);
187  PDCalibInt = tFileService->make<TH2F>("Integral_cal_int","Integral Calibration by Channel",288,0,288.,100000,0.0,1000000.0);
188  PDtrigs = tFileService->make<TH1I>("Triggers vs. Channel","Triggers vs. Channel",288.,0.,288.);
189  PDPEDhist = tFileService->make<TH1F>("Pedestal vs. Channel","Pedestal vs. Channel",288.,0.,288.);
190  PDchanThres = tFileService->make<TH1F>("Threshold vs. Channel","Threshold vs. Channel",288,0.,288.);
191  PDchanCorr = tFileService->make<TH2F>(Form("ADCMax_chan%d_chan%d",fSSP_corrchan1,fSSP_corrchan2),
192  Form("ADCMax Channel %d vs. ADCMax Channel %d",fSSP_corrchan1,fSSP_corrchan2),4000,0.0,4000.0,4000,0.0,4000.0);
193 
194  for(int i=0;i<288;i++){
195  PDchanRawPerTrace[i] = tFileService->make<TH2F>(Form("raw_per_trace_chan_%d",i),
196  Form("Raw Persistence Traces Channel %d",i),2000.,0.,2000.,1000.,1350.,4000.);
197  PDchanCorPerTrace[i] = tFileService->make<TH2F>(Form("cor_per_trace_chan_%d",i),
198  Form("Pedestal Substracted Persistence Traces Channel %d",i),2000.,0.,2000.,1000.,0.,2000.);
199  PDchanPEDRough[i] = tFileService->make<TH2F>(Form("ped_calc_trace_chan_%d",i),
200  Form("Wave Form Fraction for Pedestal %d",i),40.,0.,40.,1000.,1500.,2500.);
201  PDchanWaveInt[i] = tFileService->make<TH1F>(Form("wave_integrals_pedsub_chan_%d",i),
202  Form("Pedestal Subtracted Wave Integrals Channel %d",i),100000,0.0,1000000.0);
203  PDchanCurPeak[i] = tFileService->make<TH2F>(Form("current_peak_%d",i),
204  Form("current_peak_%d",i),1500,0.0,1500.0,1000,0.0,10000.0);
205  }
206  }
void pd_monitor::PDWaveform::beginRun ( const art::Run run)

Definition at line 218 of file PDWaveform_module.cc.

218 { }
void pd_monitor::PDWaveform::endJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 208 of file PDWaveform_module.cc.

208  {
209  std::cout << "Finalizing Histograms." << std::endl;
210  for(int i=0;i<288;i++) {
211  PDPEDhist->SetBinContent(i+1,PDchanPEDRough[i]->ProjectionY()->GetMean());
212  PDchanThres->SetBinContent(i,PDchanMaxPed->ProjectionY("",i,i)->GetXaxis()->GetBinLowEdge(PDchanMaxPed->ProjectionY("",i,i)->GetMaximumBin()));
213  }
214  }
QTextStream & endl(QTextStream &s)
void pd_monitor::PDWaveform::reconfigure ( fhicl::ParameterSet const &  pset)

Definition at line 137 of file PDWaveform_module.cc.

137  {
138  fSSPInput = p.get< std::string >("SSPInputModule");
139  fSSPInstance = p.get< std::string >("SSPInstanceName");
140 
141  //Try this, reconfigure fhicl settings
142  fSSP_m1=p.get<int>("SSP_m1");
143  fSSP_m2=p.get<int>("SSP_m2");
144  fSSP_i1=p.get<int>("SSP_i1");
145  fSSP_i2=p.get<int>("SSP_i2");
146  fSSP_disc_width=p.get<int>("SSP_disc_width");
147  fSSP_readout_pretrigger=p.get<int>("SSP_readout_pretrigger");
148  fSSP_wfm_verbose=p.get<unsigned int>("SSP_wfm_verbose");
149  fPDwaveform_fft=p.get<unsigned int>("PDwaveform_fft");
150  fSSP_corrchan1=p.get<unsigned int>("SSP_Corr_Chan_1");
151  fSSP_corrchan2=p.get<unsigned int>("SSP_Corr_Chan_2");
152  fSSP_win=p.get<int>("SSP_calib_win");
153 
154  fSSP_smoothing=p.get<unsigned int>("SSP_smoothing");
155  fSSP_rarenum=p.get<float>("SSP_rare_num");
156  fSSP_nump=p.get<unsigned int>("SSP_nump");
157  fSSP_peak_sense=p.get<float>("SSP_peak_sense");
158  fSSP_peak=p.get<unsigned int>("SSP_peak");
159 
160  if( fSSP_wfm_verbose ){
161  std::cout << " fSSP_m1=" << fSSP_m1 <<std::endl;
162  std::cout << " fSSP_m2=" << fSSP_m2 <<std::endl;
163  std::cout << " fSSP_i1=" << fSSP_i1 <<std::endl;
164  std::cout << " fSSP_i2=" << fSSP_i2 <<std::endl;
165  std::cout << " fSSP_disc_width=" << fSSP_disc_width <<std::endl;
166  std::cout << " fSSP_readout_pretrigger=" << fSSP_readout_pretrigger <<std::endl;
167  std::cout << " fSSP_corrchan1="<< fSSP_corrchan1 << std::endl;
168  std::cout << " fSSP_corrchan2="<< fSSP_corrchan2 << std::endl;
169  std::cout << " fSSP_win="<< fSSP_win << std::endl;
170  std::cout << " fSSP_smoothing="<< fSSP_smoothing << std::endl;
171  }
172  }
std::string string
Definition: nybbler.cc:12
p
Definition: test.py:223
QTextStream & endl(QTextStream &s)
void pd_monitor::PDWaveform::setRootObjects ( )

Definition at line 216 of file PDWaveform_module.cc.

216 { }

Member Data Documentation

art::ServiceHandle<geo::Geometry> pd_monitor::PDWaveform::fGeom
private

Definition at line 109 of file PDWaveform_module.cc.

bool pd_monitor::PDWaveform::fIsSSP
private

Definition at line 107 of file PDWaveform_module.cc.

unsigned int pd_monitor::PDWaveform::fPDwaveform_fft =0
private

Definition at line 98 of file PDWaveform_module.cc.

unsigned int pd_monitor::PDWaveform::fSSP_corrchan1 =205
private

Definition at line 99 of file PDWaveform_module.cc.

unsigned int pd_monitor::PDWaveform::fSSP_corrchan2 =134
private

Definition at line 99 of file PDWaveform_module.cc.

int pd_monitor::PDWaveform::fSSP_disc_width =20
private

Definition at line 97 of file PDWaveform_module.cc.

int pd_monitor::PDWaveform::fSSP_i1 =40
private

Definition at line 97 of file PDWaveform_module.cc.

int pd_monitor::PDWaveform::fSSP_i2 =1200
private

Definition at line 97 of file PDWaveform_module.cc.

int pd_monitor::PDWaveform::fSSP_m1 =10
private

Definition at line 97 of file PDWaveform_module.cc.

int pd_monitor::PDWaveform::fSSP_m2 =10
private

Definition at line 97 of file PDWaveform_module.cc.

unsigned int pd_monitor::PDWaveform::fSSP_nump =10
private

Definition at line 104 of file PDWaveform_module.cc.

unsigned int pd_monitor::PDWaveform::fSSP_peak =1
private

Definition at line 105 of file PDWaveform_module.cc.

float pd_monitor::PDWaveform::fSSP_peak_sense =0.5
private

Definition at line 102 of file PDWaveform_module.cc.

float pd_monitor::PDWaveform::fSSP_rarenum =3.0
private

Definition at line 101 of file PDWaveform_module.cc.

int pd_monitor::PDWaveform::fSSP_readout_pretrigger =50
private

Definition at line 97 of file PDWaveform_module.cc.

unsigned int pd_monitor::PDWaveform::fSSP_smoothing =0
private

Definition at line 100 of file PDWaveform_module.cc.

unsigned int pd_monitor::PDWaveform::fSSP_wfm_verbose =0
private

Definition at line 98 of file PDWaveform_module.cc.

int pd_monitor::PDWaveform::fSSP_win =20
private

Definition at line 97 of file PDWaveform_module.cc.

std::string pd_monitor::PDWaveform::fSSPInput
private

Definition at line 95 of file PDWaveform_module.cc.

std::string pd_monitor::PDWaveform::fSSPInstance
private

Definition at line 96 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDCalibInt
private

Definition at line 118 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanCorPerTrace[288]
private

Definition at line 120 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanCorr
private

Definition at line 119 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanCurPeak[288]
private

Definition at line 127 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanFFT
private

Definition at line 114 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanMax
private

Definition at line 116 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanMaxPed
private

Definition at line 117 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanMean
private

Definition at line 111 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanPED
private

Definition at line 115 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanPEDRough[288]
private

Definition at line 122 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanRawPerTrace[288]
private

Definition at line 121 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanRMS
private

Definition at line 112 of file PDWaveform_module.cc.

TH2F* pd_monitor::PDWaveform::PDchanRMSwide
private

Definition at line 113 of file PDWaveform_module.cc.

TH1F* pd_monitor::PDWaveform::PDchanThres
private

Definition at line 126 of file PDWaveform_module.cc.

TH1F* pd_monitor::PDWaveform::PDchanWaveInt[288]
private

Definition at line 123 of file PDWaveform_module.cc.

TH1F* pd_monitor::PDWaveform::PDPEDhist
private

Definition at line 125 of file PDWaveform_module.cc.

TH1I* pd_monitor::PDWaveform::PDtrigs
private

Definition at line 124 of file PDWaveform_module.cc.


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