1 #ifndef ALIDIELECTRONSIGNALBASE_H
2 #define ALIDIELECTRONSIGNALBASE_H
4 /* Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5 * See cxx source for full Copyright notice */
7 //#############################################################
9 //# Class AliDielectronSignalBase #
10 //# Manage Cuts on the legs of the pair #
13 //# Anton Andronic, GSI / A.Andronic@gsi.de #
14 //# Ionut C. Arsene, GSI / I.C.Arsene@gsi.de #
15 //# Julian Book, Uni Ffm / Julian.Book@cern.ch #
16 //# Frederick Kramer, Uni Ffm, / Frederick.Kramer@cern.ch #
17 //# Magnus Mager, CERN / Magnus.Mager@cern.ch #
18 //# WooJin J. Park, GSI / W.J.Park@gsi.de #
19 //# Jens Wiechula, Uni HD / Jens.Wiechula@cern.ch #
21 //#############################################################
32 class AliDielectronSignalBase : public TNamed {
34 enum EBackgroundMethod {
45 AliDielectronSignalBase();
46 AliDielectronSignalBase(const char*name, const char* title);
48 virtual ~AliDielectronSignalBase();
50 void SetIntegralRange(Double_t min, Double_t max) {fIntMin=min;fIntMax=max;}
51 void SetFitRange(Double_t min, Double_t max) {fFitMin=min; fFitMax=max;}
52 void SetRebin(Int_t factor) {fRebin = factor;}
53 void SetMethod(EBackgroundMethod method) {fMethod = method;}
54 Int_t GetMethod() const { return (Int_t)fMethod; }
56 const TVectorD& GetValues() const {return fValues;}
57 const TVectorD& GetErrors() const {return fErrors;}
59 Double_t GetIntegralMin() const { return fIntMin; }
60 Double_t GetIntegralMax() const { return fIntMax; }
61 Double_t GetSignal() const { return fValues(0);}
62 Double_t GetSignalError() const { return fErrors(0);}
63 Double_t GetBackground() const { return fValues(1);}
64 Double_t GetBackgroundError() const { return fErrors(1);}
65 Double_t GetSignificance() const { return fValues(2);}
66 Double_t GetSignificanceError() const { return fErrors(2);}
67 Double_t GetSB() const { return fValues(3);}
68 Double_t GetSBError() const { return fErrors(3);}
69 Double_t GetMass() const { return fValues(4);}
70 Double_t GetMassError() const { return fErrors(4);}
71 Double_t GetMassWidth() const { return fValues(5);}
72 Double_t GetMassWidthError() const { return fErrors(5);}
74 TH1* GetSignalHistogram() const {return fHistSignal;}
75 TH1* GetBackgroundHistogram() const {return fHistBackground;}
76 TH1* GetUnlikeSignHistogram() const {return fHistDataPM;}
77 TH1* GetRfactorHistogram() const {return fHistRfactor;}
79 void SetScaleRawToBackground(Double_t intMin, Double_t intMax) { fScaleMin=intMin; fScaleMax=intMax; }
80 Double_t GetScaleMin() const { return fScaleMin; }
81 Double_t GetScaleMax() const { return fScaleMax; }
83 Double_t GetScaleFactor() const { return fScaleFactor; }
85 void SetMixingCorrection(Bool_t mixcorr=kTRUE) { fMixingCorr=mixcorr; }
87 static Double_t ScaleHistograms(TH1* histRaw, TH1* histBackground, Double_t intMin, Double_t intMax);
89 virtual void Print(Option_t *option="") const;
92 This function needs to be implemented by the signal extraction classes.
93 Here all the work should be done.
95 The signal extraction is done on the mass spectra.
96 The TObjArray should contain the Inv. Mass spectra of the 10 possible combinations
97 for single and mixed events defined in AliDielectron.cxx
99 virtual void Process(TObjArray * const /*arrhist*/) = 0;
103 TH1 *fHistSignal; // histogram of pure signal
104 TH1 *fHistBackground; // histogram of background (fitted=0, like-sign=1, event mixing=2)
105 TH1 *fHistDataPM; // histogram of selected +- pair candidates
106 TH1 *fHistDataPP; // histogram of selected ++ pair candidates
107 TH1 *fHistDataMM; // histogram of selected -- pair candidates
108 TH1 *fHistDataME; // histogram of selected +- pair candidates from mixed event
109 TH1 *fHistRfactor; // histogram of R factors
111 TVectorD fValues; // values
112 TVectorD fErrors; // value errors
114 Double_t fIntMin; // signal extraction range min
115 Double_t fIntMax; // signal extraction range max
116 Double_t fFitMin; // fit range lowest inv. mass
117 Double_t fFitMax; // fit range highest inv. mass
119 Int_t fRebin; // histogram rebin factor
120 EBackgroundMethod fMethod; // method for background substraction
121 Double_t fScaleMin; // min for scaling of raw and background histogram
122 Double_t fScaleMax; // max for scaling of raw and background histogram
123 Double_t fScaleFactor; // scale factor of raw to background histogram scaling
124 Bool_t fMixingCorr; // switch for bin by bin correction with R factor
126 Bool_t fProcessed; // flag
128 void SetSignificanceAndSOB(); // calculate the significance and S/B
129 TPaveText* DrawStats(Double_t x1=0., Double_t y1=0., Double_t x2=0., Double_t y2=0.);
131 AliDielectronSignalBase(const AliDielectronSignalBase &c);
132 AliDielectronSignalBase &operator=(const AliDielectronSignalBase &c);
134 ClassDef(AliDielectronSignalBase,4) // Dielectron SignalBase
137 inline void AliDielectronSignalBase::SetSignificanceAndSOB()
140 // Calculate S/B and significance
143 fValues(3) = (fValues(1)>0 ? fValues(0)/fValues(1) : 0);
144 Float_t epsSig = (fValues(0)>0 ? fErrors(0)/fValues(0) : 0);
145 Float_t epsBknd = (fValues(1)>0 ? fErrors(1)/fValues(1) : 0);
146 fErrors(3) = fValues(3)*TMath::Sqrt(epsSig*epsSig + epsBknd*epsBknd);
148 fValues(2) = ((fValues(0)+fValues(1))>0 ? fValues(0)/TMath::Sqrt(fValues(0)+fValues(1)) : 0);
149 Float_t s = (fValues(0)>0?fValues(0):0); Float_t b = fValues(1);
150 fErrors(2) = ((s+b)>0 ? TMath::Sqrt((s*(s+2*b)*(s+2*b)+b*s*s)/(4*TMath::Power(s+b,3))) : 0);