]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/dielectron/AliDielectronSignalBase.h
Signal extraction and error calc. updates (Ionut), updated configs (Jens)
[u/mrichter/AliRoot.git] / PWG3 / 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     kFitted = 0,
36     kLikeSign,
37     kEventMixing
38   };
39
40   AliDielectronSignalBase();
41   AliDielectronSignalBase(const char*name, const char* title);
42   
43   virtual ~AliDielectronSignalBase();
44   
45   void SetIntegralRange(Double_t min, Double_t max) {fIntMin=min;fIntMax=max;}
46   void SetFitRange(Double_t min, Double_t max) {fFitMin=min; fFitMax=max;}
47   void SetRebin(Int_t factor) {fRebin = factor;}
48   void SetMethod(EBackgroundMethod method) {fMethod = method;}
49
50   const TVectorD& GetValues() const {return fValues;}
51   const TVectorD& GetErrors() const {return fErrors;}
52
53   Double_t GetIntegralMin()          const { return fIntMin; }
54   Double_t GetIntegralMax()          const { return fIntMax; }
55   Double_t GetSignal()               const { return fValues(0);}
56   Double_t GetSignalError()          const { return fErrors(0);}
57   Double_t GetBackground()           const { return fValues(1);}
58   Double_t GetBackgroundError()      const { return fErrors(1);}
59   Double_t GetSignificance()         const { return fValues(2);}
60   Double_t GetSignificanceError()    const { return fErrors(2);}
61   Double_t GetSB()                   const { return fValues(3);}
62   Double_t GetSBError()              const { return fErrors(3);}
63   Double_t GetMass()                 const { return fValues(4);}
64   Double_t GetMassError()            const { return fErrors(4);}
65   Double_t GetMassWidth()            const { return fValues(5);}
66   Double_t GetMassWidthError()       const { return fErrors(5);}
67
68   TH1F* GetSignalHistogram()        {return fHistSignal;}
69   TH1F* GetBackgroundHistogram()    {return fHistBackground;}
70
71   virtual void Print(Option_t *option="") const;
72
73   /**
74   This function needs to be implemented by the signal extraction classes.
75   Here all the work should be done.
76   
77   The signal extraction is done on the mass spectra.
78   The TObjArray should contain the Inv. Mass spectra of the 10 possible combinations
79      for single and mixed events defined in AliDielectron.cxx
80   */
81   virtual void Process(TObjArray * const /*arrhist*/) = 0;
82     
83 protected: 
84
85   TH1F *fHistSignal;                  // histogram of pure signal
86   TH1F *fHistBackground;              // histogram of background (fitted=0, like-sign=1, event mixing=2)
87   TH1F *fHistDataPM;                  // histogram of selected +- pair candidates
88   TH1F *fHistDataPP;                  // histogram of selected ++ pair candidates
89   TH1F *fHistDataMM;                  // histogram of selected -- pair candidates
90
91   TVectorD fValues;                   // values
92   TVectorD fErrors;                   // value errors
93
94   Double_t fIntMin;                   // signal extraction range min
95   Double_t fIntMax;                   // signal extraction range max
96   Double_t fFitMin;                   // fit range lowest inv. mass
97   Double_t fFitMax;                   // fit range highest inv. mass
98
99   Int_t fRebin;                       // histogram rebin factor
100   EBackgroundMethod fMethod;          // method for background substraction
101
102   Bool_t fProcessed;                  // flag
103
104   void SetSignificanceAndSOB();       // calculate the significance and S/B
105   TPaveText* DrawStats(Double_t x1=0., Double_t y1=0., Double_t x2=0., Double_t y2=0.);
106
107   AliDielectronSignalBase(const AliDielectronSignalBase &c);
108   AliDielectronSignalBase &operator=(const AliDielectronSignalBase &c);
109
110   ClassDef(AliDielectronSignalBase,3) // Dielectron SignalBase
111 };
112
113 inline void AliDielectronSignalBase::SetSignificanceAndSOB()
114 {
115   //
116   // Calculate S/B and significance
117   //
118   // Signal/Background
119   fValues(3) = (fValues(1)>0 ? fValues(0)/fValues(1) : 0);
120   Float_t epsSig = (fValues(0)>0 ? fErrors(0)/fValues(0) : 0);
121   Float_t epsBknd = (fValues(1)>0 ? fErrors(1)/fValues(1) : 0);
122   fErrors(3) = fValues(3)*TMath::Sqrt(epsSig*epsSig + epsBknd*epsBknd);
123   // Significance
124   fValues(2) = ((fValues(0)+fValues(1))>0 ? fValues(0)/TMath::Sqrt(fValues(0)+fValues(1)) : 0);
125   Float_t s = fValues(0); Float_t b = fValues(1);
126   fErrors(2) = ((s+b)>0 ? TMath::Sqrt((s*(s+2*b)*(s+2*b)+b*s*s)/(4*TMath::Power(s+b,3))) : 0);
127 }
128
129 #endif