]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/PartCorrBase/AliNeutralMesonSelection.h
Centrality update (Alberica)
[u/mrichter/AliRoot.git] / PWG4 / PartCorrBase / AliNeutralMesonSelection.h
1 #ifndef ALINEUTRALMESONSELECTION_H
2 #define ALINEUTRALMESONSELECTION_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice     */
5 /* $Id: AliNeutralMesonSelection.h 27413 2008-07-18 13:28:12Z gconesab $ */
6
7 //_________________________________________________________________________
8 // Class that contains methods to select candidate pairs to neutral meson 
9 // 2 main selections, invariant mass around pi0 (also any other mass),
10 // apperture angle to distinguish from combinatorial.
11 // There is a 3rd cut based on the gamma correlation on phi or pt.
12 //-- Author: Gustavo Conesa (INFN-LNF)
13
14 // --- ROOT system ---
15 #include<TObject.h>
16 #include<TArrayD.h>
17
18 class TLorentzVector ;
19 class TList ;
20 class TH2F ;
21
22 //--- ANALYSIS system ---
23
24 class AliNeutralMesonSelection : public TObject {
25   
26  public: 
27   AliNeutralMesonSelection() ; // default ctor
28   virtual ~AliNeutralMesonSelection() { ; } //virtual dtor  
29  private:
30   AliNeutralMesonSelection(const AliNeutralMesonSelection & g) ; // cpy ctor
31   AliNeutralMesonSelection & operator = (const AliNeutralMesonSelection & g) ;//cpy assignment
32   
33  public:
34   // General
35
36   TList *  GetCreateOutputObjects();
37   
38   void     InitParameters();    
39   
40   void     Print(const Option_t * opt) const;
41
42   Bool_t   AreNeutralMesonSelectionHistosKept()   const { return fKeepNeutralMesonHistos ; }
43   void     KeepNeutralMesonSelectionHistos(Bool_t keep) { fKeepNeutralMesonHistos = keep ; }
44   
45   Bool_t   SelectPair(TLorentzVector particlei,  TLorentzVector particlej, TString calo)  ;
46   
47   void     SetParticle(TString particleName) ;  // Do some default settings for "Pi0" or "Eta"
48   
49   // Asymmetry selection
50     
51   Float_t  GetAsymmetryCut()                     const { return fAsymmetryCut            ; }
52   void     SetAsymmetryCut(Float_t asy)                { fAsymmetryCut = asy             ; }    
53   
54   void     SwitchOnAsymmetryCut()                      { fUseAsymmetryCut = kTRUE        ; }
55   void     SwitchOffAsymmetryCut()                     { fUseAsymmetryCut = kFALSE       ; }
56   
57   //Opening angle selection 
58   
59   Double_t GetAngleMaxParam(Int_t i)              const { return fAngleMaxParam.At(i)    ; }
60   void     SetAngleMaxParam(Int_t i, Double_t par)      { fAngleMaxParam.AddAt(par,i)    ; }
61   
62   void     SetShiftMinAngleCut(Float_t a, Float_t b)    { fShiftMinAngle[0] = a          ; 
63                                                           fShiftMinAngle[1] = b          ; }
64   
65   void     SwitchOnAngleSelection()                     { fUseAngleCut = kTRUE           ; }
66   void     SwitchOffAngleSelection()                    { fUseAngleCut = kFALSE          ; }
67
68   Bool_t   IsAngleInWindow(const Float_t angle, const Float_t e) const ;
69   
70   //Invariant mass selection
71   
72   Double_t GetInvMassMaxCut()                     const { return fInvMassMaxCut          ; }
73   Double_t GetInvMassMinCut()                     const { return fInvMassMinCut          ; }
74   
75   void     SetInvMassCutRange(Double_t invmassmin, Double_t invmassmax)
76             { fInvMassMaxCut =invmassmax;  fInvMassMinCut =invmassmin                    ; }    
77   
78   void     SetInvMassCutMaxParameters(Float_t a, Float_t b, Float_t c)
79             { fInvMassMaxCutParam[0] = a ; 
80               fInvMassMaxCutParam[1] = b ; 
81               fInvMassMaxCutParam[2] = c ; }    
82
83   Double_t GetMass()                              const { return fM                      ; }
84   void     SetMass(Double_t m)                          { fM = m                         ; }
85       
86   //Histogrammes setters and getters
87   
88   virtual void SetHistoERangeAndNBins(Float_t min, Float_t max, Int_t n) {
89     fHistoNEBins = n ;
90     fHistoEMax = max ;
91     fHistoEMin = min ;
92   }
93   
94   Int_t   GetHistoNEBins()     const { return fHistoNEBins    ; }
95   Float_t GetHistoEMin()       const { return fHistoEMin      ; }
96   Float_t GetHistoEMax()       const { return fHistoEMax      ; }
97     
98   virtual void SetHistoAngleRangeAndNBins(Float_t min, Float_t max, Int_t n) {
99     fHistoNAngleBins = n ;
100     fHistoAngleMax = max ;
101     fHistoAngleMin = min ;
102   }
103   
104   Int_t   GetHistoNAngleBins() const { return fHistoNAngleBins ; }
105   Float_t GetHistoAngleMin()   const { return fHistoAngleMin   ; }
106   Float_t GetHistoAngleMax()   const { return fHistoAngleMax   ; }
107   
108   virtual void SetHistoIMRangeAndNBins(Float_t min, Float_t max, Int_t n) {
109     fHistoNIMBins = n ;
110     fHistoIMMax = max ;
111     fHistoIMMin = min ;
112   }
113   
114   Int_t   GetHistoNIMBins()    const { return fHistoNIMBins    ; }
115   Float_t GetHistoIMMin()      const { return fHistoIMMin      ; }
116   Float_t GetHistoIMMax()      const { return fHistoIMMax      ; }
117   
118   
119  private:
120   
121   Float_t  fAsymmetryCut  ;               // Asymmetry cut
122   Bool_t   fUseAsymmetryCut;              // Use the asymmetry cut
123   
124   Double_t fM ;                           // Mass of the neutral meson
125   Double_t fInvMassMaxCut ;               // Invariant Mass cut maximum
126   Double_t fInvMassMinCut ;               // Invariant Masscut minimun
127   Double_t fInvMassMaxCutParam[3];        // Variable invariant mass max cut, for pi0 in EMCAL
128   
129   TArrayD  fAngleMaxParam ;               // Max opening angle selection parameters
130   Bool_t   fUseAngleCut   ;               // Select pairs depending on their opening angle
131   Float_t  fShiftMinAngle[2] ;            // Correction shift for min angle from true kinematic limit, resolution effects
132   
133   Bool_t   fKeepNeutralMesonHistos ;      // Keep neutral meson selection histograms
134
135   //Histograms
136   TH2F *   fhAnglePairNoCut ;             //! Aperture angle of decay photons, no cuts
137   TH2F *   fhAnglePairOpeningAngleCut ;   //! Aperture angle of decay photons, cut on opening angle
138   TH2F *   fhAnglePairAsymmetryCut ;      //! Aperture angle of decay photons, asymmetry cut
139   TH2F *   fhAnglePairAllCut ;            //! Aperture angle of decay photons, all cuts
140   
141   TH2F *   fhInvMassPairNoCut ;           //! Invariant mass of decay photons, no cuts
142   TH2F *   fhInvMassPairOpeningAngleCut ; //! Invariant mass of decay photons, cut on opening angle
143   TH2F *   fhInvMassPairAsymmetryCut ;    //! Invariant mass of decay photons, asymmetry cut  
144   TH2F *   fhInvMassPairAllCut ;          //! Invariant mass of decay photons, all cuts
145
146   TH2F *   fhAsymmetryNoCut ;             //! Asymmetry of decay photons, no cuts
147   TH2F *   fhAsymmetryOpeningAngleCut ;   //! Asymmetry of decay photons, cut on opening angle
148   TH2F *   fhAsymmetryAllCut ;            //! Asymmetry of decay photons, all cuts
149   
150   //Histograms binning and range    
151   Int_t    fHistoNEBins ;                 // Number of bins in pi0 E axis
152   Float_t  fHistoEMax ;                   // Maximum value of pi0 E histogram range
153   Float_t  fHistoEMin ;                   // Minimum value of pi0 E histogram range
154   
155   Int_t    fHistoNAngleBins ;             // Number of bins in angle axis
156   Float_t  fHistoAngleMax ;               // Maximum value of angle histogram range
157   Float_t  fHistoAngleMin ;               // Minimum value of angle histogram range
158   
159   Int_t    fHistoNIMBins ;                // Number of bins in Invariant Mass axis
160   Float_t  fHistoIMMax ;                  // Maximum value of Invariant Mass histogram range
161   Float_t  fHistoIMMin ;                  // Minimum value of Invariant Mass histogram range
162   
163   ClassDef(AliNeutralMesonSelection,6)
164     
165 } ;
166
167 #endif //ALINEUTRALMESONSELECTION_H
168
169
170