]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/CaloTrackCorrBase/AliNeutralMesonSelection.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWG / CaloTrackCorrBase / 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
6 //_________________________________________________________________________
7 // Class that contains methods to select candidate pairs to neutral meson 
8 // 2 main selections, invariant mass around pi0 (also any other mass),
9 // apperture angle to distinguish from combinatorial.
10 // There is a 3rd cut based on the gamma correlation on phi or pt.
11 //-- Author: Gustavo Conesa (INFN-LNF)
12
13 // --- ROOT system ---
14 #include<TObject.h>
15 #include<TArrayD.h>
16 #include<TString.h>
17
18 class TLorentzVector ;
19 class TList   ;
20 class TH2F    ;
21 class TString ;
22
23 class AliNeutralMesonSelection : public TObject {
24   
25  public: 
26   AliNeutralMesonSelection() ; // default ctor
27   virtual ~AliNeutralMesonSelection() { ; } //virtual dtor  
28
29   // General
30
31   TList *  GetCreateOutputObjects();
32   
33   void     InitParameters();    
34   
35   void     Print(const Option_t * opt) const;
36
37   Bool_t   AreNeutralMesonSelectionHistosKept()   const { return fKeepNeutralMesonHistos ; }
38   void     KeepNeutralMesonSelectionHistos(Bool_t keep) { fKeepNeutralMesonHistos = keep ; }
39   
40   Bool_t   SelectPair(TLorentzVector particlei,  TLorentzVector particlej, TString calo)  ;
41   
42   void     SetParticle(TString particleName) ;  // Do some default settings for "Pi0" or "Eta"
43   TString  GetParticle()                          const { return fParticle               ; }
44   
45   // Asymmetry selection
46     
47   Float_t  GetAsymmetryCut()                      const { return fAsymmetryCut           ; }
48   void     SetAsymmetryCut(Float_t asy)                 { fAsymmetryCut = asy            ; }
49   
50   void     SwitchOnAsymmetryCut()                       { fUseAsymmetryCut = kTRUE       ; }
51   void     SwitchOffAsymmetryCut()                      { fUseAsymmetryCut = kFALSE      ; }
52   
53   //Opening angle selection 
54   
55   Double_t GetAngleMaxParam(Int_t i)              const { return fAngleMaxParam.At(i)    ; }
56   void     SetAngleMaxParam(Int_t i, Double_t par)      { fAngleMaxParam.AddAt(par,i)    ; }
57   
58   void     SetShiftMinAngleCut(Float_t a, Float_t b)    { fShiftMinAngle[0] = a          ;
59                                                           fShiftMinAngle[1] = b          ; }
60   
61   void     SwitchOnAngleSelection()                     { fUseAngleCut = kTRUE           ; }
62   void     SwitchOffAngleSelection()                    { fUseAngleCut = kFALSE          ; }
63
64   Bool_t   IsAngleInWindow(Float_t angle, Float_t e) const ;
65   
66   //Invariant mass selection
67   
68   Double_t GetInvMassMaxCut()                     const { return fInvMassMaxCut          ; }
69   Double_t GetInvMassMinCut()                     const { return fInvMassMinCut          ; }
70   
71   void     SetInvMassCutRange(Double_t invmassmin, Double_t invmassmax)
72             { fInvMassMaxCut =invmassmax;  fInvMassMinCut =invmassmin                    ; }    
73   
74   void     SetSideBandCutRanges( Double_t lmin, Double_t lmax, 
75                                  Double_t rmin, Double_t rmax )
76             { fLeftBandMinCut  = lmin ; fLeftBandMaxCut  = lmax ; 
77               fRightBandMinCut = rmin ; fRightBandMaxCut = rmax ; }     
78   
79   void     SetInvMassCutMaxParameters(Float_t a, Float_t b, Float_t c)
80             { fInvMassMaxCutParam[0] = a ; 
81               fInvMassMaxCutParam[1] = b ; 
82               fInvMassMaxCutParam[2] = c ; }    
83
84   Double_t GetMass()                              const { return fM                      ; }
85   void     SetMass(Double_t m)                          { fM = m                         ; }
86   
87   
88   // Decay photon bit setting
89   
90   enum decayTypes { kPi0, kEta, kPi0Side, kEtaSide} ;
91
92   UInt_t    GetDecayBit()                         const { return fDecayBit               ; }
93   
94   void    SetDecayBit(Int_t &tag, UInt_t set) const {
95     // Set bit of type set (decayTypes) in tag
96     tag |= (1<<set) ;
97   }
98   
99   void    SetDecayBit(Int_t &tag) const {
100     // Set bit of type set (decayTypes) in tag
101     tag |= (1<<fDecayBit) ;
102   }
103   
104   Bool_t  CheckDecayBit(Int_t tag, UInt_t test) const {
105     // Check if in tag the bit test (decayTypes) is set.
106     if (tag & (1<<test) ) return  kTRUE ;
107     else return kFALSE ;
108   }
109
110   Bool_t  CheckDecayBit(Int_t tag) const {
111     // Check if in tag the bit test (decayTypes) is set.
112     if (tag & (1<<fDecayBit) ) return  kTRUE ;
113     else return kFALSE ;
114   }
115   
116   // Histograms setters and getters
117   
118   virtual void SetHistoERangeAndNBins(Float_t min, Float_t max, Int_t n) {
119     fHistoNEBins = n ;
120     fHistoEMax = max ;
121     fHistoEMin = min ;
122   }
123   
124   Int_t   GetHistoNEBins()     const { return fHistoNEBins    ; }
125   Float_t GetHistoEMin()       const { return fHistoEMin      ; }
126   Float_t GetHistoEMax()       const { return fHistoEMax      ; }
127     
128   virtual void SetHistoAngleRangeAndNBins(Float_t min, Float_t max, Int_t n) {
129     fHistoNAngleBins = n ;
130     fHistoAngleMax = max ;
131     fHistoAngleMin = min ;
132   }
133   
134   Int_t   GetHistoNAngleBins() const { return fHistoNAngleBins ; }
135   Float_t GetHistoAngleMin()   const { return fHistoAngleMin   ; }
136   Float_t GetHistoAngleMax()   const { return fHistoAngleMax   ; }
137   
138   virtual void SetHistoIMRangeAndNBins(Float_t min, Float_t max, Int_t n) {
139     fHistoNIMBins = n ;
140     fHistoIMMax = max ;
141     fHistoIMMin = min ;
142   }
143   
144   Int_t   GetHistoNIMBins()    const { return fHistoNIMBins    ; }
145   Float_t GetHistoIMMin()      const { return fHistoIMMin      ; }
146   Float_t GetHistoIMMax()      const { return fHistoIMMax      ; }
147   
148  private:
149   
150   Float_t  fAsymmetryCut  ;               // Asymmetry cut
151   Bool_t   fUseAsymmetryCut;              // Use the asymmetry cut
152
153   Double_t fM ;                           // Mass of the neutral meson
154   Double_t fInvMassMaxCut ;               // Invariant Mass cut maximum
155   Double_t fInvMassMinCut ;               // Invariant Masscut minimun
156   Double_t fInvMassMaxCutParam[3];        // Variable invariant mass max cut, for pi0 in EMCAL
157   
158   Double_t fLeftBandMinCut  ;             // Side Band selection, min left  band cut
159   Double_t fLeftBandMaxCut  ;             // Side Band selection, max left  band cut
160   Double_t fRightBandMinCut ;             // Side Band selection, min right band cut
161   Double_t fRightBandMaxCut ;             // Side Band selection, max right band cut
162   
163   TArrayD  fAngleMaxParam ;               // Max opening angle selection parameters
164   Bool_t   fUseAngleCut   ;               // Select pairs depending on their opening angle
165   Float_t  fShiftMinAngle[2] ;            // Correction shift for min angle from true kinematic limit, resolution effects
166   
167   Bool_t   fKeepNeutralMesonHistos ;      // Keep neutral meson selection histograms
168   
169   TString  fParticle ;                    // neutral meson name (Pi0, Eta, +SideBand)
170   UInt_t   fDecayBit;                     // Decay type flag, set when fParticle is set
171
172   //Histograms
173   TH2F *   fhAnglePairNoCut ;             //! Aperture angle of decay photons, no cuts
174   TH2F *   fhAnglePairOpeningAngleCut ;   //! Aperture angle of decay photons, cut on opening angle
175   TH2F *   fhAnglePairAsymmetryCut ;      //! Aperture angle of decay photons, asymmetry cut
176   TH2F *   fhAnglePairAllCut ;            //! Aperture angle of decay photons, all cuts
177   
178   TH2F *   fhInvMassPairNoCut ;           //! Invariant mass of decay photons, no cuts
179   TH2F *   fhInvMassPairOpeningAngleCut ; //! Invariant mass of decay photons, cut on opening angle
180   TH2F *   fhInvMassPairAsymmetryCut ;    //! Invariant mass of decay photons, asymmetry cut  
181   TH2F *   fhInvMassPairAllCut ;          //! Invariant mass of decay photons, all cuts
182
183   TH2F *   fhAsymmetryNoCut ;             //! Asymmetry of decay photons, no cuts
184   TH2F *   fhAsymmetryOpeningAngleCut ;   //! Asymmetry of decay photons, cut on opening angle
185   TH2F *   fhAsymmetryAllCut ;            //! Asymmetry of decay photons, all cuts
186   
187   //Histograms binning and range    
188   Int_t    fHistoNEBins ;                 // Number of bins in pi0 E axis
189   Float_t  fHistoEMax ;                   // Maximum value of pi0 E histogram range
190   Float_t  fHistoEMin ;                   // Minimum value of pi0 E histogram range
191   
192   Int_t    fHistoNAngleBins ;             // Number of bins in angle axis
193   Float_t  fHistoAngleMax ;               // Maximum value of angle histogram range
194   Float_t  fHistoAngleMin ;               // Minimum value of angle histogram range
195   
196   Int_t    fHistoNIMBins ;                // Number of bins in Invariant Mass axis
197   Float_t  fHistoIMMax ;                  // Maximum value of Invariant Mass histogram range
198   Float_t  fHistoIMMin ;                  // Minimum value of Invariant Mass histogram range  
199   
200   AliNeutralMesonSelection(              const AliNeutralMesonSelection & g) ; // cpy ctor
201   AliNeutralMesonSelection & operator = (const AliNeutralMesonSelection & g) ; // cpy assignment
202   
203   ClassDef(AliNeutralMesonSelection,8)
204     
205 } ;
206
207 #endif //ALINEUTRALMESONSELECTION_H
208
209
210