Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Private Attributes | List of all members
NBPulse Class Reference

Analyzes the signal of a nanobeacon in a PMT. More...

#include <NBPulse.hh>

Public Member Functions

 NBPulse ()
 Default constructor. More...
 
 NBPulse (TH2D *h)
 Constructor. More...
 
 ~NBPulse ()
 Destructor. More...
 
double Kolmogorov_Smirnov_uniform (TH1D *h)
 Performs a Kolmogorov-Smirnov compatibility test between the hit time distribution and a uniform distribution. More...
 
void initialize ()
 Leaves the pulse ready for analysis by projecting the time vs ToT distribution to the time axis and moving it out of the Root file directory. More...
 
void analyzeFast ()
 Performs a fast analysis of the hit time distribution. More...
 
void computeBaseline (TH1 *h)
 Computes baseline of the hit time distribution measured by the PMT. More...
 
void fit ()
 Fits the hit nanobeacon signal to a model composed by a Landau and a Gauss functions. More...
 
bool IsGood ()
 Check if the peak is good. More...
 
bool IsFitted ()
 Check if the peak is fitted. More...
 
bool IsSaturatedHit ()
 Check if the hit time distribution is saturated. More...
 
bool IsWeak ()
 Check if the hit time distribution is weak. More...
 
TH1D * getHtime_full ()
 Get the hit time distribution. More...
 
TH1D * gettime_peak ()
 Get the hit time distribution. More...
 
TH1D * gettot_peak ()
 Get the hit time distribution. More...
 
TH2D * gettime_vs_tot ()
 Get the hit time vs ToT distribution. More...
 
RooWorkspace getWorkspace ()
 Get all the information from the fit. More...
 
double getPeakTime ()
 Get the arrival time of the nanobeacon pulse. More...
 
double getMeanToT ()
 Get the mean of the ToT distribution of the hits produced by the nanobeacon. More...
 
double getSigmaToT ()
 Get the standard deviation of the ToT distribution of the hits produced by the nanobeacon. More...
 
string getName ()
 Get name of the hit time vs ToT histogram. More...
 

Private Attributes

TH2D * time_vs_tot
 
TH1D * htime_full
 
TH1D * htime_peak
 
TH1D * tot_peak
 
RooWorkspace w
 
pair< double, double > baseline
 
double tlow
 
double thigh
 
double tpeak
 
double nhits
 
double meanToT
 
double sigmaToT
 
double max_ToT
 
int n_bins_right
 
int n_bins_left
 
int n_bins_before_peak
 
int n_bins_afterpulse
 
int max_n_bins_baseline
 
int min_n_bins_baseline
 
int n_sigma_peak
 
int n_sigma_before_peak
 
double fit_region_min
 
double fit_region_max
 
double fit_range_min
 
double fit_range_max
 
double landau_m_min
 
double landau_m_max
 
double landau_s_min
 
double landau_s_max
 
double gaus_m_min
 
double gaus_m_max
 
double gaus_s_min
 
double gaus_s_max
 
double l_g_ratio_min
 
double l_g_ratio_max
 
double KS_threshold
 
bool analyzed
 
bool fitted
 
bool is_good
 
bool saturated_hit
 
bool saturated_tot
 
bool weak
 

Detailed Description

Analyzes the signal of a nanobeacon in a PMT.

Definition at line 33 of file NBPulse.hh.

Constructor & Destructor Documentation

NBPulse::NBPulse ( )
inline

Default constructor.

Definition at line 66 of file NBPulse.hh.

66 {}
NBPulse::NBPulse ( TH2D *  h)
inline

Constructor.

This constructor initializes the nanobeacon pulse from a 2D histogram containing the time vs ToT hit distribution.

Parameters
htime vs ToT histogram

Definition at line 76 of file NBPulse.hh.

76  :
77  max_ToT (26.4) ,
78  n_bins_right (5) ,
79  n_bins_left (5) ,
81  n_bins_afterpulse (10),
84  n_sigma_peak (5) ,
86  fit_region_min (20.0),
87  fit_region_max (20.0) ,
88  fit_range_min (6.0) ,
89  fit_range_max (15.0) ,
90  landau_m_min (2.0) ,
91  landau_m_max (2.0) ,
92  landau_s_min (1.0) ,
93  landau_s_max (9.0) ,
94  gaus_m_min (2.0) ,
95  gaus_m_max (2.0) ,
96  gaus_s_min (1.0) ,
97  gaus_s_max (9.0) ,
98  l_g_ratio_min (0.0) ,
99  l_g_ratio_max (1.0) ,
100  KS_threshold(1e-3) ,
101  analyzed (false) ,
102  fitted (false) ,
103  is_good (false) ,
104  saturated_hit (false) ,
105  saturated_tot (false) ,
106  weak (false)
107 
108  {
109 
110  time_vs_tot = h ;
111 
112  initialize() ;
113 
114  }
double landau_m_min
Definition: NBPulse.hh:47
double fit_range_min
Definition: NBPulse.hh:47
int n_bins_afterpulse
Definition: NBPulse.hh:45
double landau_s_min
Definition: NBPulse.hh:47
int max_n_bins_baseline
Definition: NBPulse.hh:45
int n_sigma_peak
Definition: NBPulse.hh:45
double gaus_m_max
Definition: NBPulse.hh:47
bool saturated_hit
Definition: NBPulse.hh:55
double max_ToT
Definition: NBPulse.hh:43
void initialize()
Leaves the pulse ready for analysis by projecting the time vs ToT distribution to the time axis and m...
Definition: NBPulse.hh:184
TH2D * time_vs_tot
Definition: NBPulse.hh:35
bool weak
Definition: NBPulse.hh:59
double landau_s_max
Definition: NBPulse.hh:47
double landau_m_max
Definition: NBPulse.hh:47
bool analyzed
Definition: NBPulse.hh:49
double l_g_ratio_max
Definition: NBPulse.hh:47
double fit_region_min
Definition: NBPulse.hh:47
double gaus_s_max
Definition: NBPulse.hh:47
int n_bins_left
Definition: NBPulse.hh:45
int min_n_bins_baseline
Definition: NBPulse.hh:45
double fit_range_max
Definition: NBPulse.hh:47
int n_sigma_before_peak
Definition: NBPulse.hh:45
double l_g_ratio_min
Definition: NBPulse.hh:47
double fit_region_max
Definition: NBPulse.hh:47
bool fitted
Definition: NBPulse.hh:51
double gaus_s_min
Definition: NBPulse.hh:47
double gaus_m_min
Definition: NBPulse.hh:47
double KS_threshold
Definition: NBPulse.hh:47
bool is_good
Definition: NBPulse.hh:53
int n_bins_before_peak
Definition: NBPulse.hh:45
bool saturated_tot
Definition: NBPulse.hh:57
int n_bins_right
Definition: NBPulse.hh:45
NBPulse::~NBPulse ( )
inline

Destructor.

Definition at line 120 of file NBPulse.hh.

120  {
121 
122  delete(htime_full) ;
123 
124  if (analyzed == true){
125 
126  delete(tot_peak) ;
127 
128  delete(htime_peak) ;
129 
130  }
131 
132  }
TH1D * htime_peak
Definition: NBPulse.hh:37
TH1D * tot_peak
Definition: NBPulse.hh:37
bool analyzed
Definition: NBPulse.hh:49
TH1D * htime_full
Definition: NBPulse.hh:37

Member Function Documentation

double NBPulse::Kolmogorov_Smirnov_uniform ( TH1D *  h)
inline

Performs a Kolmogorov-Smirnov compatibility test between the hit time distribution and a uniform distribution.

Parameters
hThe hit time distribution
Returns
The probability that the hit time distribution is compatible with a uniform distribution.

Definition at line 141 of file NBPulse.hh.

141  {
142 
143  TRandom3* r = new TRandom3() ;
144 
145  double low_x = 0 ;
146 
147  double high_x = h->GetBinCenter(h->GetNbinsX() - 1) ;
148 
149  int nentries = h->Integral(1 , h->GetNbinsX() - 1) ;
150 
151  if (h->GetBinCenter(1)>=0){
152 
153  low_x = h->GetBinCenter(1) ;
154 
155  }
156 
157  TH1D* u = (TH1D*)h->Clone("uniform") ;
158 
159  u->Reset() ;
160 
161  for (int i = 0 ; i<nentries ; i++){
162 
163  u->Fill(r->Uniform(low_x, high_x));
164 
165  }
166 
167  double proba = h->KolmogorovTest(u) ;
168 
169  delete(r) ;
170 
171  delete(u) ;
172 
173  return proba ;
174 
175  } ;
void NBPulse::initialize ( )
inline

Leaves the pulse ready for analysis by projecting the time vs ToT distribution to the time axis and moving it out of the Root file directory.

Definition at line 184 of file NBPulse.hh.

184  {
185 
186  htime_full = time_vs_tot->ProjectionX(time_vs_tot->GetName() , 1 , time_vs_tot->GetNbinsY());
187 
188  htime_full->SetDirectory(0) ;
189 
190  RooMsgService::instance().setGlobalKillBelow(RooFit::WARNING);
191 
192  }
#define WARNING(A)
Definition: JMessage.hh:63
TH2D * time_vs_tot
Definition: NBPulse.hh:35
TH1D * htime_full
Definition: NBPulse.hh:37
void NBPulse::analyzeFast ( )
inline

Performs a fast analysis of the hit time distribution.

Definition at line 198 of file NBPulse.hh.

198  {
199 
200  if (analyzed == true) return ;
201 
203 
204  int peakbin = htime_full->GetMaximumBin() ;
205 
206  double height = htime_full->GetBinContent(peakbin) ;
207 
208  int leftbin = peakbin - n_bins_left ;
209 
210  int rightbin = peakbin + n_bins_right ;
211 
212  tpeak = htime_full->GetBinCenter(peakbin) ;
213 
214  tlow = htime_full->GetBinCenter(leftbin) ;
215 
216  thigh = htime_full->GetBinCenter(rightbin) ;
217 
218  htime_peak = (TH1D*)htime_full->Clone();
219 
220  for (int i=0 ; i<htime_peak->GetNbinsX() ; i++){
221 
222  if (i<leftbin || i>rightbin){
223 
224  htime_peak->SetBinContent(i,0);
225 
226  }
227 
228  }
229 
230  htime_peak->SetName((std::string("htime_")+htime_peak->GetName()).c_str());
231 
232  htime_peak->SetDirectory(0) ;
233 
234  tot_peak = (TH1D*)time_vs_tot->ProjectionY((std::string("ToT_") + time_vs_tot->GetName()).c_str() , leftbin , rightbin) ;
235 
236  tot_peak->SetDirectory(0);
237 
238  meanToT = tot_peak->GetMean() ;
239 
240  if (meanToT > max_ToT){
241 
242  saturated_tot = true ;
243 
244  }
245 
246  sigmaToT = tot_peak->GetStdDev() ;
247 
248  double p_KS = Kolmogorov_Smirnov_uniform(htime_full) ;
249 
250  if (baseline.first==0 && baseline.second==0) weak = true ;
251 
252  if (p_KS > KS_threshold) weak = true ;
253 
254  if (height < (baseline.first + n_sigma_peak * baseline.second)) weak = true ;
255 
256  for (int i = 1 ; i<=n_bins_before_peak ; i++){
257 
258  if (htime_full->GetBinContent(peakbin - n_bins_before_peak) < (baseline.first + n_sigma_before_peak * baseline.second)) weak = true ;
259 
260  }
261 
262  if (weak == false && htime_full->GetBinContent(peakbin + n_bins_afterpulse) == 0) saturated_hit = true ;
263 
264  if (weak == false && saturated_hit == false) is_good = true ;
265 
266  nhits = htime_peak->Integral() ;
267 
268  analyzed = true ;
269 
270  }
double meanToT
Definition: NBPulse.hh:43
int n_bins_afterpulse
Definition: NBPulse.hh:45
TH1D * htime_peak
Definition: NBPulse.hh:37
int n_sigma_peak
Definition: NBPulse.hh:45
void computeBaseline(TH1 *h)
Computes baseline of the hit time distribution measured by the PMT.
Definition: NBPulse.hh:279
bool saturated_hit
Definition: NBPulse.hh:55
TH1D * tot_peak
Definition: NBPulse.hh:37
double max_ToT
Definition: NBPulse.hh:43
TH2D * time_vs_tot
Definition: NBPulse.hh:35
bool weak
Definition: NBPulse.hh:59
bool analyzed
Definition: NBPulse.hh:49
TH1D * htime_full
Definition: NBPulse.hh:37
double sigmaToT
Definition: NBPulse.hh:43
int n_bins_left
Definition: NBPulse.hh:45
double Kolmogorov_Smirnov_uniform(TH1D *h)
Performs a Kolmogorov-Smirnov compatibility test between the hit time distribution and a uniform dist...
Definition: NBPulse.hh:141
int n_sigma_before_peak
Definition: NBPulse.hh:45
double thigh
Definition: NBPulse.hh:43
double tpeak
Definition: NBPulse.hh:43
double KS_threshold
Definition: NBPulse.hh:47
bool is_good
Definition: NBPulse.hh:53
pair< double, double > baseline
Definition: NBPulse.hh:41
int n_bins_before_peak
Definition: NBPulse.hh:45
bool saturated_tot
Definition: NBPulse.hh:57
double nhits
Definition: NBPulse.hh:43
double tlow
Definition: NBPulse.hh:43
int n_bins_right
Definition: NBPulse.hh:45
void NBPulse::computeBaseline ( TH1 *  h)
inline

Computes baseline of the hit time distribution measured by the PMT.

For this the first bins of the distribution (where the nanobeacon pulse is absent) are used.

Parameters
hHistogram with the hit time distribution

Definition at line 279 of file NBPulse.hh.

279  {
280 
281  int N = 0;
282 
283  double sum = 0 ;
284 
285  vector<double> bins ;
286 
287  for (int i=1 ; i<max_n_bins_baseline + 1 ; i++){
288 
289  if(h->GetBinContent(i) > 0 ){
290 
291  sum += h->GetBinContent(i);
292 
293  bins.push_back(h->GetBinContent(i));
294 
295  N += 1 ;
296 
297  }
298 
299  }
300 
301  if(N < min_n_bins_baseline){
302 
303  bins.clear();
304 
305  sum = 0 ;
306 
307  N = 0 ;
308 
309  for (int i=1 ; i<max_n_bins_baseline + 1 ; i++){
310 
311  if(h->GetBinContent(h->GetNbinsX() - i)>0){
312 
313  sum += h->GetBinContent(h->GetNbinsX() - i);
314 
315  bins.push_back(h->GetBinContent(h->GetNbinsX() - i));
316 
317  N += 1 ;
318 
319  }
320 
321  }
322 
323  }
324 
325  if(N < min_n_bins_baseline){
326 
327  baseline.first = 0 ;
328 
329  baseline.second = 0 ;
330 
331  return ;
332 
333  }
334 
335  double bg_mean = sum/N ;
336 
337  baseline.first = bg_mean ;
338 
339  double sdev = 0;
340 
341  for (int i=0 ; i<N ; i++){
342 
343  sdev += pow((bins[i] - bg_mean) , 2) ;
344 
345  }
346 
347  baseline.second = sqrt(sdev/N);
348 
349  }
int max_n_bins_baseline
Definition: NBPulse.hh:45
int min_n_bins_baseline
Definition: NBPulse.hh:45
pair< double, double > baseline
Definition: NBPulse.hh:41
void NBPulse::fit ( )
inline

Fits the hit nanobeacon signal to a model composed by a Landau and a Gauss functions.

Definition at line 355 of file NBPulse.hh.

355  {
356 
357  if (fitted==true) return ;
358 
359  if (analyzed==false) analyzeFast() ;
360 
361  w.SetName((std::string("W_") + htime_full->GetName()).c_str());
362 
363  double tpeak = htime_full->GetBinCenter(htime_full->GetMaximumBin()) ;
364 
365  RooRealVar x("time" , "hit time" , tpeak - fit_region_min , tpeak + fit_region_max) ;
366 
367  RooRealVar ml("mean landau" , "mean landau" , tpeak - landau_m_min , tpeak + landau_m_max);
368 
369  RooRealVar sl("sigma landau" , "sigma landau" , landau_s_min , landau_s_max);
370 
371  RooLandau L1("L1" , "Main_Landau" , x , ml , sl) ;
372 
373  RooRealVar mg("mean Gauss" , "mean Gauss" , tpeak - gaus_m_min , tpeak + gaus_m_max);
374 
375  RooRealVar sg("sigma Gauss" , "sigma Gauss" , gaus_s_min , gaus_s_max);
376 
377  RooGaussian G1("G1" , "Main_Gaussian" , x , mg , sg) ;
378 
379  RooRealVar c1("c1" , "c1" , l_g_ratio_min , l_g_ratio_max);
380 
381  RooAddPdf M("Model" , "Model" , RooArgList(G1,L1) , RooArgList(c1));
382 
383  w.import(M) ;
384 
385  RooDataHist Data ("data" , "data" , *w.var("time") , Import(*htime_full));
386 
387  w.pdf("Model")->fitTo(Data , Range(tpeak - fit_range_min , tpeak + fit_range_max) , RooFit::PrintLevel(-1000) , RooFit::PrintEvalErrors(-1000)) ;
388 
389  RooRealVar Tmin ("tmin" , "tmin" , tpeak - fit_range_min) ;
390 
391  RooRealVar Tmax ("tmax" , "tmax" , tpeak + fit_range_max) ;
392 
393  RooRealVar Tpeak ("tpeak" , "tpeak" , tpeak) ;
394 
395  RooRealVar Mpv("mpv" , "mpv" , Tmin.getVal()) ;
396 
397  RooRealVar Diff("peak-mpv" , "peak-mpv" , Tpeak.getVal() - Tmin.getVal()) ;
398 
399  double max_pdfval = 0 ;
400 
401  RooRealVar *var = (RooRealVar*)w.function("Model")->getObservables(RooArgSet(*w.var("time")))->first();
402 
403  for(double x = Tmin.getVal() ; x<Tmax.getVal() ; x+=0.1){
404 
405  var->setVal(x) ;
406 
407  double pdfval = w.function("Model")->getVal(RooArgSet(*var)) ;
408 
409  if (pdfval>max_pdfval){
410 
411  max_pdfval = pdfval ;
412 
413  Mpv.setVal(x) ;
414 
415  }
416 
417  }
418 
419  w.import(Mpv) ;
420 
421  w.import(Tmin) ;
422 
423  w.import(Tmax) ;
424 
425  w.import(Tpeak) ;
426 
427  w.import(Diff) ;
428 
429  fitted = true ;
430 
431  }
double landau_m_min
Definition: NBPulse.hh:47
double fit_range_min
Definition: NBPulse.hh:47
double landau_s_min
Definition: NBPulse.hh:47
double gaus_m_max
Definition: NBPulse.hh:47
double landau_s_max
Definition: NBPulse.hh:47
double landau_m_max
Definition: NBPulse.hh:47
void analyzeFast()
Performs a fast analysis of the hit time distribution.
Definition: NBPulse.hh:198
Double_t G1(const Double_t x)
Integral of method g1.
Definition: JQuantiles.cc:37
bool analyzed
Definition: NBPulse.hh:49
TH1D * htime_full
Definition: NBPulse.hh:37
double l_g_ratio_max
Definition: NBPulse.hh:47
double fit_region_min
Definition: NBPulse.hh:47
double gaus_s_max
Definition: NBPulse.hh:47
double fit_range_max
Definition: NBPulse.hh:47
double l_g_ratio_min
Definition: NBPulse.hh:47
double fit_region_max
Definition: NBPulse.hh:47
bool fitted
Definition: NBPulse.hh:51
double gaus_s_min
Definition: NBPulse.hh:47
TCanvas * c1
Global variables to handle mouse events.
double gaus_m_min
Definition: NBPulse.hh:47
double tpeak
Definition: NBPulse.hh:43
RooWorkspace w
Definition: NBPulse.hh:39
bool NBPulse::IsGood ( )
inline

Check if the peak is good.

Returns
true if the peak is classified as good.

Definition at line 439 of file NBPulse.hh.

439  {
440 
441  return is_good ;
442 
443  }
bool is_good
Definition: NBPulse.hh:53
bool NBPulse::IsFitted ( )
inline

Check if the peak is fitted.

Returns
true if the peak has been fitted.

Definition at line 451 of file NBPulse.hh.

451  {
452 
453  return fitted ;
454 
455  }
bool fitted
Definition: NBPulse.hh:51
bool NBPulse::IsSaturatedHit ( )
inline

Check if the hit time distribution is saturated.

Returns
true if the hit time distributon is classified as saturated.

Definition at line 463 of file NBPulse.hh.

463  {
464 
465  return saturated_hit ;
466 
467  }
bool saturated_hit
Definition: NBPulse.hh:55
bool NBPulse::IsWeak ( )
inline

Check if the hit time distribution is weak.

Returns
true if the hit time distributon is classified as weak.

Definition at line 475 of file NBPulse.hh.

475  {
476 
477  return weak ;
478 
479  }
bool weak
Definition: NBPulse.hh:59
TH1D* NBPulse::getHtime_full ( )
inline

Get the hit time distribution.

Returns
histogram with hit time distribution

Definition at line 488 of file NBPulse.hh.

488  {
489 
490  return htime_full ;
491 
492  }
TH1D * htime_full
Definition: NBPulse.hh:37
TH1D* NBPulse::gettime_peak ( )
inline

Get the hit time distribution.

Returns
histogram with hit time distribution

Definition at line 501 of file NBPulse.hh.

501  {
502 
503  return htime_peak ;
504 
505  }
TH1D * htime_peak
Definition: NBPulse.hh:37
TH1D* NBPulse::gettot_peak ( )
inline

Get the hit time distribution.

Returns
histogram with hit time distribution

Definition at line 514 of file NBPulse.hh.

514  {
515 
516  return tot_peak ;
517 
518  }
TH1D * tot_peak
Definition: NBPulse.hh:37
TH2D* NBPulse::gettime_vs_tot ( )
inline

Get the hit time vs ToT distribution.

Returns
a 2d histogram with hit time vs ToT distribution

Definition at line 528 of file NBPulse.hh.

528  {
529 
530  return time_vs_tot ;
531 
532  }
TH2D * time_vs_tot
Definition: NBPulse.hh:35
RooWorkspace NBPulse::getWorkspace ( )
inline

Get all the information from the fit.

Returns
Roofit Workspace with all the information and parameters from the fit

Definition at line 541 of file NBPulse.hh.

541  {
542 
543  return w ;
544 
545  }
RooWorkspace w
Definition: NBPulse.hh:39
double NBPulse::getPeakTime ( )
inline

Get the arrival time of the nanobeacon pulse.

Returns
pulse arrival time

Definition at line 554 of file NBPulse.hh.

554  {
555 
556  return tpeak ;
557 
558  }
double tpeak
Definition: NBPulse.hh:43
double NBPulse::getMeanToT ( )
inline

Get the mean of the ToT distribution of the hits produced by the nanobeacon.

Returns
mean ToT value

Definition at line 567 of file NBPulse.hh.

567  {
568 
569  return meanToT ;
570 
571  }
double meanToT
Definition: NBPulse.hh:43
double NBPulse::getSigmaToT ( )
inline

Get the standard deviation of the ToT distribution of the hits produced by the nanobeacon.

Returns
standard deviation from ToT distribution

Definition at line 580 of file NBPulse.hh.

580  {
581 
582  return sigmaToT ;
583 
584  }
double sigmaToT
Definition: NBPulse.hh:43
string NBPulse::getName ( )
inline

Get name of the hit time vs ToT histogram.

Returns
The name of the histogram.

Definition at line 593 of file NBPulse.hh.

593  {
594 
595  return time_vs_tot->GetName();
596 
597  }
TH2D * time_vs_tot
Definition: NBPulse.hh:35

Member Data Documentation

TH2D* NBPulse::time_vs_tot
private

Definition at line 35 of file NBPulse.hh.

TH1D* NBPulse::htime_full
private

Definition at line 37 of file NBPulse.hh.

TH1D * NBPulse::htime_peak
private

Definition at line 37 of file NBPulse.hh.

TH1D * NBPulse::tot_peak
private

Definition at line 37 of file NBPulse.hh.

RooWorkspace NBPulse::w
private

Definition at line 39 of file NBPulse.hh.

pair<double , double> NBPulse::baseline
private

Definition at line 41 of file NBPulse.hh.

double NBPulse::tlow
private

Definition at line 43 of file NBPulse.hh.

double NBPulse::thigh
private

Definition at line 43 of file NBPulse.hh.

double NBPulse::tpeak
private

Definition at line 43 of file NBPulse.hh.

double NBPulse::nhits
private

Definition at line 43 of file NBPulse.hh.

double NBPulse::meanToT
private

Definition at line 43 of file NBPulse.hh.

double NBPulse::sigmaToT
private

Definition at line 43 of file NBPulse.hh.

double NBPulse::max_ToT
private

Definition at line 43 of file NBPulse.hh.

int NBPulse::n_bins_right
private

Definition at line 45 of file NBPulse.hh.

int NBPulse::n_bins_left
private

Definition at line 45 of file NBPulse.hh.

int NBPulse::n_bins_before_peak
private

Definition at line 45 of file NBPulse.hh.

int NBPulse::n_bins_afterpulse
private

Definition at line 45 of file NBPulse.hh.

int NBPulse::max_n_bins_baseline
private

Definition at line 45 of file NBPulse.hh.

int NBPulse::min_n_bins_baseline
private

Definition at line 45 of file NBPulse.hh.

int NBPulse::n_sigma_peak
private

Definition at line 45 of file NBPulse.hh.

int NBPulse::n_sigma_before_peak
private

Definition at line 45 of file NBPulse.hh.

double NBPulse::fit_region_min
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::fit_region_max
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::fit_range_min
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::fit_range_max
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::landau_m_min
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::landau_m_max
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::landau_s_min
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::landau_s_max
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::gaus_m_min
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::gaus_m_max
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::gaus_s_min
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::gaus_s_max
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::l_g_ratio_min
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::l_g_ratio_max
private

Definition at line 47 of file NBPulse.hh.

double NBPulse::KS_threshold
private

Definition at line 47 of file NBPulse.hh.

bool NBPulse::analyzed
private

Definition at line 49 of file NBPulse.hh.

bool NBPulse::fitted
private

Definition at line 51 of file NBPulse.hh.

bool NBPulse::is_good
private

Definition at line 53 of file NBPulse.hh.

bool NBPulse::saturated_hit
private

Definition at line 55 of file NBPulse.hh.

bool NBPulse::saturated_tot
private

Definition at line 57 of file NBPulse.hh.

bool NBPulse::weak
private

Definition at line 59 of file NBPulse.hh.


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