]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGDQ/dielectron/AliDielectronSignalBase.h
o small update for signal extraction
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / AliDielectronSignalBase.h
1 #ifndef ALIDIELECTRONSIGNALBASE_H
2 #define ALIDIELECTRONSIGNALBASE_H
3
4 /* Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 //#############################################################
8 //#                                                           # 
9 //#         Class AliDielectronSignalBase                       #
10 //#         Manage Cuts on the legs of the pair               #
11 //#                                                           #
12 //#  Authors:                                                 #
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      #
20 //#                                                           #
21 //#############################################################
22
23
24 #include <TNamed.h>
25 #include <TVectorT.h>
26 #include <TMath.h>
27
28 class TObjArray;
29 class TPaveText;
30 class TH1F;
31
32 class AliDielectronSignalBase : public TNamed {
33 public:
34   enum EBackgroundMethod {
35     kFittedMC = 0,
36     kFitted,
37     kLikeSign,
38     kLikeSignArithm,
39     kLikeSignRcorr,
40     kLikeSignArithmRcorr,
41     kEventMixing,
42     kRotation
43   };
44
45   AliDielectronSignalBase();
46   AliDielectronSignalBase(const char*name, const char* title);
47   
48   virtual ~AliDielectronSignalBase();
49   
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; }
55
56   const TVectorD& GetValues() const {return fValues;}
57   const TVectorD& GetErrors() const {return fErrors;}
58
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);}
73  
74   TH1* GetSignalHistogram()      const {return fHistSignal;}
75   TH1* GetBackgroundHistogram()  const {return fHistBackground;}
76   TH1* GetUnlikeSignHistogram()  const {return fHistDataPM;}
77   TH1* GetRfactorHistogram()     const {return fHistRfactor;}
78   
79   void SetScaleRawToBackground(Double_t intMin, Double_t intMax) { fScaleMin=intMin; fScaleMax=intMax; }
80   Double_t GetScaleFactor() const { return fScaleFactor; }
81
82   void SetMixingCorrection(Bool_t mixcorr=kTRUE) { fMixingCorr=mixcorr; }
83   
84   static Double_t ScaleHistograms(TH1* histRaw, TH1* histBackground, Double_t intMin, Double_t intMax);
85   
86   virtual void Print(Option_t *option="") const;
87
88   /**
89   This function needs to be implemented by the signal extraction classes.
90   Here all the work should be done.
91   
92   The signal extraction is done on the mass spectra.
93   The TObjArray should contain the Inv. Mass spectra of the 10 possible combinations
94      for single and mixed events defined in AliDielectron.cxx
95   */
96   virtual void Process(TObjArray * const /*arrhist*/) = 0;
97     
98 protected: 
99
100   TH1 *fHistSignal;                  // histogram of pure signal
101   TH1 *fHistBackground;              // histogram of background (fitted=0, like-sign=1, event mixing=2)
102   TH1 *fHistDataPM;                  // histogram of selected +- pair candidates
103   TH1 *fHistDataPP;                  // histogram of selected ++ pair candidates
104   TH1 *fHistDataMM;                  // histogram of selected -- pair candidates
105   TH1 *fHistDataME;                  // histogram of selected +- pair candidates from mixed event
106   TH1 *fHistRfactor;                 // histogram of R factors
107   
108   TVectorD fValues;                  // values
109   TVectorD fErrors;                  // value errors
110
111   Double_t fIntMin;                  // signal extraction range min
112   Double_t fIntMax;                  // signal extraction range max
113   Double_t fFitMin;                  // fit range lowest inv. mass
114   Double_t fFitMax;                  // fit range highest inv. mass
115
116   Int_t fRebin;                      // histogram rebin factor
117   EBackgroundMethod fMethod;         // method for background substraction
118   Double_t fScaleMin;                // min for scaling of raw and background histogram
119   Double_t fScaleMax;                // max for scaling of raw and background histogram
120   Double_t fScaleFactor;             // scale factor of raw to background histogram scaling
121   Bool_t fMixingCorr;                // switch for bin by bin correction with R factor
122   
123   Bool_t fProcessed;                 // flag
124   
125   void SetSignificanceAndSOB();      // calculate the significance and S/B
126   TPaveText* DrawStats(Double_t x1=0., Double_t y1=0., Double_t x2=0., Double_t y2=0.);
127
128   AliDielectronSignalBase(const AliDielectronSignalBase &c);
129   AliDielectronSignalBase &operator=(const AliDielectronSignalBase &c);
130
131   ClassDef(AliDielectronSignalBase,4) // Dielectron SignalBase
132 };
133
134 inline void AliDielectronSignalBase::SetSignificanceAndSOB()
135 {
136   //
137   // Calculate S/B and significance
138   //
139   // Signal/Background
140   fValues(3) = (fValues(1)>0 ? fValues(0)/fValues(1) : 0);
141   Float_t epsSig = (fValues(0)>0 ? fErrors(0)/fValues(0) : 0);
142   Float_t epsBknd = (fValues(1)>0 ? fErrors(1)/fValues(1) : 0);
143   fErrors(3) = fValues(3)*TMath::Sqrt(epsSig*epsSig + epsBknd*epsBknd);
144   // Significance
145   fValues(2) = ((fValues(0)+fValues(1))>0 ? fValues(0)/TMath::Sqrt(fValues(0)+fValues(1)) : 0);
146   Float_t s = fValues(0); Float_t b = fValues(1);
147   fErrors(2) = ((s+b)>0 ? TMath::Sqrt((s*(s+2*b)*(s+2*b)+b*s*s)/(4*TMath::Power(s+b,3))) : 0);
148 }
149
150 #endif