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 $ */
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)
14 // --- ROOT system ---
18 class TLorentzVector ;
22 //--- ANALYSIS system ---
24 class AliNeutralMesonSelection : public TObject {
27 AliNeutralMesonSelection() ; // default ctor
28 virtual ~AliNeutralMesonSelection() { ; } //virtual dtor
32 TList * GetCreateOutputObjects();
34 void InitParameters();
36 void Print(const Option_t * opt) const;
38 Bool_t AreNeutralMesonSelectionHistosKept() const { return fKeepNeutralMesonHistos ; }
39 void KeepNeutralMesonSelectionHistos(Bool_t keep) { fKeepNeutralMesonHistos = keep ; }
41 Bool_t SelectPair(TLorentzVector particlei, TLorentzVector particlej, TString calo) ;
43 void SetParticle(TString particleName) ; // Do some default settings for "Pi0" or "Eta"
45 // Asymmetry selection
47 Float_t GetAsymmetryCut() const { return fAsymmetryCut ; }
48 void SetAsymmetryCut(Float_t asy) { fAsymmetryCut = asy ; }
50 void SwitchOnAsymmetryCut() { fUseAsymmetryCut = kTRUE ; }
51 void SwitchOffAsymmetryCut() { fUseAsymmetryCut = kFALSE ; }
53 //Opening angle selection
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) ; }
58 void SetShiftMinAngleCut(Float_t a, Float_t b) { fShiftMinAngle[0] = a ;
59 fShiftMinAngle[1] = b ; }
61 void SwitchOnAngleSelection() { fUseAngleCut = kTRUE ; }
62 void SwitchOffAngleSelection() { fUseAngleCut = kFALSE ; }
64 Bool_t IsAngleInWindow(const Float_t angle, const Float_t e) const ;
66 //Invariant mass selection
68 Double_t GetInvMassMaxCut() const { return fInvMassMaxCut ; }
69 Double_t GetInvMassMinCut() const { return fInvMassMinCut ; }
71 void SetInvMassCutRange(Double_t invmassmin, Double_t invmassmax)
72 { fInvMassMaxCut =invmassmax; fInvMassMinCut =invmassmin ; }
74 void SetInvMassCutMaxParameters(Float_t a, Float_t b, Float_t c)
75 { fInvMassMaxCutParam[0] = a ;
76 fInvMassMaxCutParam[1] = b ;
77 fInvMassMaxCutParam[2] = c ; }
79 Double_t GetMass() const { return fM ; }
80 void SetMass(Double_t m) { fM = m ; }
82 //Histogrammes setters and getters
84 virtual void SetHistoERangeAndNBins(Float_t min, Float_t max, Int_t n) {
90 Int_t GetHistoNEBins() const { return fHistoNEBins ; }
91 Float_t GetHistoEMin() const { return fHistoEMin ; }
92 Float_t GetHistoEMax() const { return fHistoEMax ; }
94 virtual void SetHistoAngleRangeAndNBins(Float_t min, Float_t max, Int_t n) {
95 fHistoNAngleBins = n ;
96 fHistoAngleMax = max ;
97 fHistoAngleMin = min ;
100 Int_t GetHistoNAngleBins() const { return fHistoNAngleBins ; }
101 Float_t GetHistoAngleMin() const { return fHistoAngleMin ; }
102 Float_t GetHistoAngleMax() const { return fHistoAngleMax ; }
104 virtual void SetHistoIMRangeAndNBins(Float_t min, Float_t max, Int_t n) {
110 Int_t GetHistoNIMBins() const { return fHistoNIMBins ; }
111 Float_t GetHistoIMMin() const { return fHistoIMMin ; }
112 Float_t GetHistoIMMax() const { return fHistoIMMax ; }
117 Float_t fAsymmetryCut ; // Asymmetry cut
118 Bool_t fUseAsymmetryCut; // Use the asymmetry cut
120 Double_t fM ; // Mass of the neutral meson
121 Double_t fInvMassMaxCut ; // Invariant Mass cut maximum
122 Double_t fInvMassMinCut ; // Invariant Masscut minimun
123 Double_t fInvMassMaxCutParam[3]; // Variable invariant mass max cut, for pi0 in EMCAL
125 TArrayD fAngleMaxParam ; // Max opening angle selection parameters
126 Bool_t fUseAngleCut ; // Select pairs depending on their opening angle
127 Float_t fShiftMinAngle[2] ; // Correction shift for min angle from true kinematic limit, resolution effects
129 Bool_t fKeepNeutralMesonHistos ; // Keep neutral meson selection histograms
132 TH2F * fhAnglePairNoCut ; //! Aperture angle of decay photons, no cuts
133 TH2F * fhAnglePairOpeningAngleCut ; //! Aperture angle of decay photons, cut on opening angle
134 TH2F * fhAnglePairAsymmetryCut ; //! Aperture angle of decay photons, asymmetry cut
135 TH2F * fhAnglePairAllCut ; //! Aperture angle of decay photons, all cuts
137 TH2F * fhInvMassPairNoCut ; //! Invariant mass of decay photons, no cuts
138 TH2F * fhInvMassPairOpeningAngleCut ; //! Invariant mass of decay photons, cut on opening angle
139 TH2F * fhInvMassPairAsymmetryCut ; //! Invariant mass of decay photons, asymmetry cut
140 TH2F * fhInvMassPairAllCut ; //! Invariant mass of decay photons, all cuts
142 TH2F * fhAsymmetryNoCut ; //! Asymmetry of decay photons, no cuts
143 TH2F * fhAsymmetryOpeningAngleCut ; //! Asymmetry of decay photons, cut on opening angle
144 TH2F * fhAsymmetryAllCut ; //! Asymmetry of decay photons, all cuts
146 //Histograms binning and range
147 Int_t fHistoNEBins ; // Number of bins in pi0 E axis
148 Float_t fHistoEMax ; // Maximum value of pi0 E histogram range
149 Float_t fHistoEMin ; // Minimum value of pi0 E histogram range
151 Int_t fHistoNAngleBins ; // Number of bins in angle axis
152 Float_t fHistoAngleMax ; // Maximum value of angle histogram range
153 Float_t fHistoAngleMin ; // Minimum value of angle histogram range
155 Int_t fHistoNIMBins ; // Number of bins in Invariant Mass axis
156 Float_t fHistoIMMax ; // Maximum value of Invariant Mass histogram range
157 Float_t fHistoIMMin ; // Minimum value of Invariant Mass histogram range
159 AliNeutralMesonSelection(const AliNeutralMesonSelection & g) ; // cpy ctor
160 AliNeutralMesonSelection & operator = (const AliNeutralMesonSelection & g) ; // cpy assignment
162 ClassDef(AliNeutralMesonSelection,6)
166 #endif //ALINEUTRALMESONSELECTION_H