FEDRA emulsion software from the OPERA Collaboration
EdbPeakProb Class Reference

prob the peak probability (2-dim only for the moment) More...

#include <EdbPeakProb.h>

Inheritance diagram for EdbPeakProb:
Collaboration diagram for EdbPeakProb:

Public Member Functions

 EdbPeakProb ()
 
Float_t Prob ()
 
Float_t Probability2D (TH2F *hd, float &xpeak, float &ypeak)
 
virtual ~EdbPeakProb ()
 

Public Attributes

TH1F * eHbin
 service spectrum histo More...
 
TH2F * eHD
 input histogram to be analysed More...
 
TH2F * eHpeak
 service peak histo More...
 
Int_t ePeakNmax
 apriori limits in the number of entries to accept the found peak More...
 
Int_t ePeakNmin
 apriori limits in the number of entries to accept the found peak More...
 
Float_t ePeakRMS
 found peak parameters More...
 
Float_t ePeakRMSmax
 apriori limits in RMS to accept the found peak More...
 
Float_t ePeakRMSmin
 apriori limits in RMS to accept the found peak More...
 
Float_t ePeakX
 
Float_t ePeakY
 
Float_t eProb
 the probability the found peak is not the BG fluctuation More...
 
Int_t eVerbosity
 0-no any message, 1-print, 2-plot; default=2 More...
 

Detailed Description

prob the peak probability (2-dim only for the moment)

Constructor & Destructor Documentation

◆ EdbPeakProb()

EdbPeakProb::EdbPeakProb ( )

◆ ~EdbPeakProb()

EdbPeakProb::~EdbPeakProb ( )
virtual
41{
42 SafeDelete(eHbin);
43 SafeDelete(eHpeak);
44}
TH2F * eHpeak
service peak histo
Definition: EdbPeakProb.h:22
TH1F * eHbin
service spectrum histo
Definition: EdbPeakProb.h:21

Member Function Documentation

◆ Prob()

Float_t EdbPeakProb::Prob ( )
48{
49 float prob = Probability2D( eHD, ePeakX, ePeakY);
50 return prob;
51}
Float_t Probability2D(TH2F *hd, float &xpeak, float &ypeak)
Definition: EdbPeakProb.cxx:54
TH2F * eHD
input histogram to be analysed
Definition: EdbPeakProb.h:19
Float_t ePeakY
Definition: EdbPeakProb.h:16
Float_t ePeakX
Definition: EdbPeakProb.h:16

◆ Probability2D()

float EdbPeakProb::Probability2D ( TH2F *  hd,
float &  xpeak,
float &  ypeak 
)

estimate the probability that the highest bin of the histogram do not belong to the background
Method: estimate the deviation of the tail of the coincidence distribution from exponential shape

55{
58
59 float prob=0;
60
61 float maxbin = hd->GetMaximum();
62 eHbin = new TH1F("hbin","hbin", (int)maxbin+1, 0, maxbin+1);
63 int nx = hd->GetNbinsX();
64 int ny = hd->GetNbinsY();
65 for(int ix=1; ix<=nx; ix++)
66 for(int iy=1; iy<=ny; iy++)
67 eHbin->Fill( hd->GetBinContent(ix,iy) );
68
69 float tailbin = eHbin->GetMean()+eHbin->GetRMS();
70 eHbin->Fit("expo","QR0","",tailbin,maxbin);
71 float a = eHbin->GetFunction("expo")->GetParameter(0);
72 float b = eHbin->GetFunction("expo")->GetParameter(1);
73
74 float lowlim = (int)((-a)/b+0.99999);
75 if(lowlim<1) lowlim=1;
76
77 float real = eHbin->GetBinContent((int)maxbin+1);
78 if( lowlim>maxbin ) prob=0;
79 else prob = 1.- eHbin->GetFunction("expo")->Integral(maxbin,2*(maxbin+1)+10) / real;
80
81 // calculate mean and rms for bins>lowlim to check the peak integrity
82
83 eHpeak = (TH2F*)hd->Clone("hdpeak");
84 eHpeak->Reset();
85
86 float bin;
87 for(int ix=1; ix<=nx; ix++)
88 for(int iy=1; iy<=ny; iy++) {
89 bin = hd->GetBinContent(hd->GetBin(ix, iy));
90 if(bin<lowlim) continue;
91 eHpeak->Fill( hd->GetXaxis()->GetBinCenter( ix ),
92 hd->GetYaxis()->GetBinCenter( iy ),
93 bin );
94 }
95 Log(3, "AlignmentProbability", "peak: ( %f %f ) rms: ( %f %f )",
96 eHpeak->GetMean(1), eHpeak->GetMean(2), eHpeak->GetRMS(1),eHpeak->GetRMS(2));
97
98 if(eVerbosity>1) {
99 TCanvas *c = new TCanvas("aa","aa");
100 c->Divide(2,2);
101 c->cd(1); hd->Draw("colz");
102 c->GetPad(2)->SetLogy();
103 c->cd(2); eHbin->Draw(); eHbin->GetFunction("expo")->Draw("same");
104 //c->cd(3); hd->SetMinimum(lowlim); hd->Draw();
105 c->cd(3);
106 eHpeak->Draw("colz");
107 } else {
108 SafeDelete(eHbin);
109 SafeDelete(eHpeak);
110 }
111
112 Log(2,"Probability2D","%f BGlim=%6.1f, maxbin=%6.1f",prob,lowlim,maxbin);
113 return prob;
114}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
void a()
Definition: check_aligned.C:59
Int_t eVerbosity
0-no any message, 1-print, 2-plot; default=2
Definition: EdbPeakProb.h:24
new TCanvas()

Member Data Documentation

◆ eHbin

TH1F* EdbPeakProb::eHbin

service spectrum histo

◆ eHD

TH2F* EdbPeakProb::eHD

input histogram to be analysed

◆ eHpeak

TH2F* EdbPeakProb::eHpeak

service peak histo

◆ ePeakNmax

Int_t EdbPeakProb::ePeakNmax

apriori limits in the number of entries to accept the found peak

◆ ePeakNmin

Int_t EdbPeakProb::ePeakNmin

apriori limits in the number of entries to accept the found peak

◆ ePeakRMS

Float_t EdbPeakProb::ePeakRMS

found peak parameters

◆ ePeakRMSmax

Float_t EdbPeakProb::ePeakRMSmax

apriori limits in RMS to accept the found peak

◆ ePeakRMSmin

Float_t EdbPeakProb::ePeakRMSmin

apriori limits in RMS to accept the found peak

◆ ePeakX

Float_t EdbPeakProb::ePeakX

◆ ePeakY

Float_t EdbPeakProb::ePeakY

◆ eProb

Float_t EdbPeakProb::eProb

the probability the found peak is not the BG fluctuation

◆ eVerbosity

Int_t EdbPeakProb::eVerbosity

0-no any message, 1-print, 2-plot; default=2


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