AliPhysics  8417398 (8417398)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliNeutralMesonSelection.h
Go to the documentation of this file.
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 //_________________________________________________________________________
18 //_________________________________________________________________________
19 
20 // --- ROOT system ---
21 #include<TObject.h>
22 #include<TArrayD.h>
23 #include<TString.h>
24 
25 class TLorentzVector ;
26 class TList ;
27 class TH2F ;
28 
29 class AliNeutralMesonSelection : public TObject {
30 
31  public:
32 
33  AliNeutralMesonSelection() ; // default ctor
34 
36  virtual ~AliNeutralMesonSelection() { ; }
37 
38  // General
39 
40  TList * GetCreateOutputObjects();
41 
42  void InitParameters();
43 
44  void Print(const Option_t * opt) const;
45 
48 
49  Bool_t SelectPair(TLorentzVector particlei, TLorentzVector particlej, Int_t calo) ;
50 
51  void SetParticle(TString particleName) ; // Do some default settings for "Pi0" or "Eta"
52  TString GetParticle() const { return fParticle ; }
53 
54  Int_t GetDebug() const { return fDebug ; }
55  void SetDebug(Int_t d) { fDebug = d ; }
56 
57  // Asymmetry selection
58 
59  Float_t GetAsymmetryCut() const { return fAsymmetryCut ; }
60  void SetAsymmetryCut(Float_t asy) { fAsymmetryCut = asy ; }
61 
64 
65  //Opening angle selection
66 
67  Double_t GetAngleMaxParam(Int_t i) const { return fAngleMaxParam.At(i) ; }
68  void SetAngleMaxParam(Int_t i, Double_t par) { fAngleMaxParam.AddAt(par,i) ; }
69 
70  void SetShiftMinAngleCut(Float_t a, Float_t b) { fShiftMinAngle[0] = a ;
71  fShiftMinAngle[1] = b ; }
72 
73  void SwitchOnAngleSelection() { fUseAngleCut = kTRUE ; }
74  void SwitchOffAngleSelection() { fUseAngleCut = kFALSE ; }
75 
76  Bool_t IsAngleInWindow(Float_t angle, Float_t e) const ;
77 
78  //Invariant mass selection
79 
80  Double_t GetInvMassMaxCut() const { return fInvMassMaxCut ; }
81  Double_t GetInvMassMinCut() const { return fInvMassMinCut ; }
82 
83  void SetInvMassCutRange(Double_t invmassmin, Double_t invmassmax)
84  { fInvMassMaxCut =invmassmax; fInvMassMinCut =invmassmin ; }
85 
86  void SetSideBandCutRanges( Double_t lmin, Double_t lmax,
87  Double_t rmin, Double_t rmax )
88  { fLeftBandMinCut = lmin ; fLeftBandMaxCut = lmax ;
89  fRightBandMinCut = rmin ; fRightBandMaxCut = rmax ; }
90 
91  void SetInvMassCutMaxParameters(Float_t a, Float_t b, Float_t c)
92  { fInvMassMaxCutParam[0] = a ;
93  fInvMassMaxCutParam[1] = b ;
94  fInvMassMaxCutParam[2] = c ; }
95 
96  Double_t GetMass() const { return fM ; }
97  void SetMass(Double_t m) { fM = m ; }
98 
99 
100  // Decay photon bit setting
101  static const Int_t fgkMaxNDecayBits = 8;
102 
103  enum decayTypes { kPi0 = 0, kEta = 1,
106 
107  UInt_t GetDecayBit() const { return fDecayBit ; }
108 
109  void SetDecayBit(Int_t &tag, UInt_t set) const {
110  // Set bit of type set (decayTypes) in tag
111  tag |= (1<<set) ;
112  }
113 
114  void SetDecayBit(Int_t &tag) const {
115  // Set bit of type set (decayTypes) in tag
116  tag |= (1<<fDecayBit) ;
117  }
118 
119  Bool_t CheckDecayBit(Int_t tag, UInt_t test) const {
120  // Check if in tag the bit test (decayTypes) is set.
121  if (tag & (1<<test) ) return kTRUE ;
122  else return kFALSE ;
123  }
124 
125  Bool_t CheckDecayBit(Int_t tag) const {
126  // Check if in tag the bit test (decayTypes) is set.
127  if (tag & (1<<fDecayBit) ) return kTRUE ;
128  else return kFALSE ;
129  }
130 
131  // Histograms setters and getters
132 
133  virtual void SetHistoERangeAndNBins(Float_t min, Float_t max, Int_t n) {
134  fHistoNEBins = n ;
135  fHistoEMax = max ;
136  fHistoEMin = min ;
137  }
138 
139  Int_t GetHistoNEBins() const { return fHistoNEBins ; }
140  Float_t GetHistoEMin() const { return fHistoEMin ; }
141  Float_t GetHistoEMax() const { return fHistoEMax ; }
142 
143  virtual void SetHistoAngleRangeAndNBins(Float_t min, Float_t max, Int_t n) {
144  fHistoNAngleBins = n ;
145  fHistoAngleMax = max ;
146  fHistoAngleMin = min ;
147  }
148 
149  Int_t GetHistoNAngleBins() const { return fHistoNAngleBins ; }
150  Float_t GetHistoAngleMin() const { return fHistoAngleMin ; }
151  Float_t GetHistoAngleMax() const { return fHistoAngleMax ; }
152 
153  virtual void SetHistoIMRangeAndNBins(Float_t min, Float_t max, Int_t n) {
154  fHistoNIMBins = n ;
155  fHistoIMMax = max ;
156  fHistoIMMin = min ;
157  }
158 
159  Int_t GetHistoNIMBins() const { return fHistoNIMBins ; }
160  Float_t GetHistoIMMin() const { return fHistoIMMin ; }
161  Float_t GetHistoIMMax() const { return fHistoIMMax ; }
162 
163  private:
164 
165  Float_t fAsymmetryCut ;
166 
168 
169  Double_t fM ;
170 
171  Double_t fInvMassMaxCut ;
172 
173  Double_t fInvMassMinCut ;
174 
175  Double_t fInvMassMaxCutParam[3];
176 
177  Double_t fLeftBandMinCut ;
178 
179  Double_t fLeftBandMaxCut ;
180 
181  Double_t fRightBandMinCut ;
182 
183  Double_t fRightBandMaxCut ;
184 
185  TArrayD fAngleMaxParam ;
186 
187  Bool_t fUseAngleCut ;
188 
189  Float_t fShiftMinAngle[2] ;
190 
192 
193  TString fParticle ;
194 
195  UInt_t fDecayBit;
196 
197  Int_t fDebug ;
198 
199  // Histograms
201 
203 
205 
207 
208 
210 
212 
214 
216 
218 
220 
222 
223 
224  // Histograms binning and range
225  Int_t fHistoNEBins ;
226 
227  Float_t fHistoEMax ;
228 
229  Float_t fHistoEMin ;
230 
232 
233  Float_t fHistoAngleMax ;
234 
235  Float_t fHistoAngleMin ;
236 
237 
238  Int_t fHistoNIMBins ;
239 
240  Float_t fHistoIMMax ;
241 
242  Float_t fHistoIMMin ;
243 
246 
249 
251  ClassDef(AliNeutralMesonSelection,8) ;
253 
254 } ;
255 
256 #endif //ALINEUTRALMESONSELECTION_H
257 
258 
259 
Double_t GetAngleMaxParam(Int_t i) const
Float_t fHistoAngleMin
Minimum value of angle histogram range.
Float_t fHistoAngleMax
Maximum value of angle histogram range.
Double_t fLeftBandMaxCut
Side Band selection, max left band cut.
Float_t fAsymmetryCut
Asymmetry cut.
virtual void SetHistoAngleRangeAndNBins(Float_t min, Float_t max, Int_t n)
void InitParameters()
Initialize the parameters of the analysis.
Bool_t AreNeutralMesonSelectionHistosKept() const
TH2F * fhInvMassPairAllCut
! Invariant mass of decay photons, all cuts.
AliNeutralMesonSelection()
Default constructor. Initialize parameters.
Double_t fRightBandMinCut
Side Band selection, min right band cut.
void SetAngleMaxParam(Int_t i, Double_t par)
void SetDecayBit(Int_t &tag) const
void Print(const Option_t *opt) const
Print some relevant parameters set for the analysis.
void SetInvMassCutRange(Double_t invmassmin, Double_t invmassmax)
void KeepNeutralMesonSelectionHistos(Bool_t keep)
Bool_t CheckDecayBit(Int_t tag, UInt_t test) const
Double_t fLeftBandMinCut
Side Band selection, min left band cut.
TH2F * fhAsymmetryOpeningAngleCut
! Asymmetry of decay photons, cut on opening angle.
virtual ~AliNeutralMesonSelection()
Virtual destructor.
TH2F * fhInvMassPairNoCut
! Invariant mass of decay photons, no cuts.
Bool_t fKeepNeutralMesonHistos
Keep neutral meson selection histograms.
Bool_t IsAngleInWindow(Float_t angle, Float_t e) const
Bool_t fUseAngleCut
Select pairs depending on their opening angle.
Double_t fRightBandMaxCut
Side Band selection, max right band cut.
void SetSideBandCutRanges(Double_t lmin, Double_t lmax, Double_t rmin, Double_t rmax)
Int_t fHistoNIMBins
Number of bins in Invariant Mass axis.
Float_t fHistoEMax
Maximum value of pi0 E histogram range.
Double_t fInvMassMaxCutParam[3]
Variable invariant mass max cut, for pi0 in EMCAL.
Double_t fInvMassMaxCut
Invariant Mass cut maximum.
virtual void SetHistoERangeAndNBins(Float_t min, Float_t max, Int_t n)
TArrayD fAngleMaxParam
Maximum opening angle selection parameters.
void SetParticle(TString particleName)
Set some default parameters for selection of pi0 or eta.
Float_t fShiftMinAngle[2]
Correction shift for min angle from true kinematic limit, resolution effects.
Bool_t SelectPair(TLorentzVector particlei, TLorentzVector particlej, Int_t calo)
TH2F * fhAnglePairAllCut
! Aperture angle of decay photons, all cuts.
Double_t fM
Mass of the neutral meson.
Int_t fHistoNEBins
Number of bins in pi0 E axis.
TH2F * fhAnglePairNoCut
! Aperture angle of decay photons, no cuts.
void SetDecayBit(Int_t &tag, UInt_t set) const
TH2F * fhAsymmetryNoCut
! Asymmetry of decay photons, no cuts.
Float_t fHistoEMin
Minimum value of pi0 E histogram range.
Float_t fHistoIMMin
Minimum value of Invariant Mass histogram range.
virtual void SetHistoIMRangeAndNBins(Float_t min, Float_t max, Int_t n)
TH2F * fhAnglePairOpeningAngleCut
! Aperture angle of decay photons, cut on opening angle.
UInt_t fDecayBit
Decay type flag, set while selecting, depending on fParticle and side range. See enum decayTypes for ...
Bool_t CheckDecayBit(Int_t tag) const
Float_t fHistoIMMax
Maximum value of Invariant Mass histogram range.
Class that contains methods to select candidate cluster pairs to neutral meson.
AliNeutralMesonSelection & operator=(const AliNeutralMesonSelection &nm)
Assignment operator not implemented.
TH2F * fhInvMassPairAsymmetryCut
! Invariant mass of decay photons, asymmetry cut.
Bool_t fUseAsymmetryCut
Use the asymmetry cut.
Double_t fInvMassMinCut
Invariant Masscut minimun.
Int_t fHistoNAngleBins
Number of bins in angle axis.
void SetInvMassCutMaxParameters(Float_t a, Float_t b, Float_t c)
TH2F * fhAsymmetryAllCut
! Asymmetry of decay photons, all cuts.
TH2F * fhAnglePairAsymmetryCut
! Aperture angle of decay photons, asymmetry cut.
TString fParticle
Meutral meson name (Pi0, Eta, +SideBand).
TH2F * fhInvMassPairOpeningAngleCut
! Invariant mass of decay photons, cut on opening angle.
void SetShiftMinAngleCut(Float_t a, Float_t b)