]>
Commit | Line | Data |
---|---|---|
1c5acb87 | 1 | #ifndef ALIANAPARTICLEJETLEADINGCONECORRELATION_H |
2 | #define ALIANAPARTICLEJETLEADINGCONECORRELATION_H | |
3 | /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
4 | * See cxx source for full Copyright notice */ | |
5 | /* $Id: $ */ | |
6 | ||
7 | //_________________________________________________________________________ | |
8 | // Class that contains the algorithm for the reconstruction of jet, cone around leading particle | |
9 | // The seed is a backward particle (direct photon) | |
78219bac | 10 | // 1) Take the a trigger particle found stored in AliAODPWG4ParticleCorrelation, |
1c5acb87 | 11 | // 2) Search for the highest pt leading particle opposite to the trigger within a phi, pt window |
12 | // 3) Take all particles around leading in a cone R with pt larger than threshold and construct the jet | |
13 | // | |
14 | // Class created from old AliPHOSGammaJet | |
15 | // (see AliRoot versions previous Release 4-09) | |
16 | // | |
17 | //-- Author: Gustavo Conesa (INFN-LNF) | |
18 | ||
19 | // --- ROOT system --- | |
20 | class TH2F; | |
21 | ||
22 | //---- Analysis system ---- | |
23 | class AliAODTrack; | |
0ae57829 | 24 | class AliVCluster; |
1c5acb87 | 25 | class AliCaloTrackReader; |
26 | class AliNeutralMesonSelection; | |
1c5acb87 | 27 | |
28 | #include "AliAnaPartCorrBaseClass.h" | |
29 | ||
30 | class AliAnaParticleJetLeadingConeCorrelation : public AliAnaPartCorrBaseClass { | |
31 | ||
32 | public: | |
1c5acb87 | 33 | AliAnaParticleJetLeadingConeCorrelation() ; // default ctor |
78219bac | 34 | virtual ~AliAnaParticleJetLeadingConeCorrelation() ; //virtual dtor |
35 | private: | |
1c5acb87 | 36 | AliAnaParticleJetLeadingConeCorrelation(const AliAnaParticleJetLeadingConeCorrelation & g) ; // cpy ctor |
37 | AliAnaParticleJetLeadingConeCorrelation & operator = (const AliAnaParticleJetLeadingConeCorrelation & g) ;//cpy assignment | |
78219bac | 38 | |
39 | public: | |
40 | ||
1c5acb87 | 41 | TList * GetCreateOutputObjects(); |
42 | ||
43 | void InitParameters(); | |
44 | ||
45 | void Print(const Option_t * opt) const; | |
46 | ||
47 | Bool_t AreJetsRecalculated() const {return fReMakeJet ; } | |
48 | void SwitchOnJetsRecalculation(){fReMakeJet = kTRUE; } | |
49 | void SwitchOffJetsRecalculation(){fReMakeJet = kFALSE; } | |
50 | ||
51 | Bool_t AreJetsOnlyInCTS() const {return fJetsOnlyInCTS ; } | |
52 | void SwitchOnJetsOnlyInCTS(){fJetsOnlyInCTS = kTRUE; } | |
53 | void SwitchOffJetsOnlyInCTS(){fJetsOnlyInCTS = kFALSE; } | |
54 | ||
55 | Bool_t AreSeveralConeAndPtCuts() const {return fSeveralConeAndPtCuts ; } | |
56 | void SwitchOnSeveralConeAndPtCuts(){fSeveralConeAndPtCuts = kTRUE ;} | |
57 | void SwitchOffSeveralConeAndPtCuts(){fSeveralConeAndPtCuts = kFALSE ;} | |
58 | ||
59 | Bool_t IsPbPb() const {return fPbPb ; } | |
60 | void SetppCollisions(){fPbPb = kFALSE; } | |
61 | void SetPbPbCollisions(){fPbPb = kTRUE; } | |
62 | ||
63 | Double_t GetDeltaPhiMaxCut() const {return fDeltaPhiMaxCut ; } | |
64 | Double_t GetDeltaPhiMinCut() const {return fDeltaPhiMinCut ; } | |
65 | Double_t GetLeadingRatioMaxCut() const {return fLeadingRatioMaxCut ; } | |
66 | Double_t GetLeadingRatioMinCut() const {return fLeadingRatioMinCut ; } | |
67 | ||
68 | Double_t GetPtTriggerSelectionCut() const {return fPtTriggerSelectionCut ; } | |
69 | Double_t GetJetRatioMaxCut() const {return fJetRatioMaxCut ; } | |
70 | Double_t GetJetRatioMinCut() const {return fJetRatioMinCut ; } | |
71 | ||
72 | void SetPtTriggerSelectionCut(Double_t cut){fPtTriggerSelectionCut = cut; } | |
73 | void SetJetSelectionMode(UInt_t select){ fSelect= select ; } | |
74 | ||
75 | Int_t GetJetNCones() const {return fJetNCone ; } | |
76 | Int_t GetJetNPtThres() const {return fJetNPt ; } | |
77 | Float_t GetJetCone() const {return fJetCone ; } | |
78 | Float_t GetJetPtThreshold() const {return fJetPtThreshold ; } | |
79 | Float_t GetJetPtThresPbPb() const {return fJetPtThresPbPb ; } | |
80 | Float_t GetJetCones(Int_t i) const {return fJetCones[i] ; } | |
81 | Float_t GetJetPtThreshold(Int_t i) const {return fJetPtThres[i] ; } | |
82 | TString GetJetConeName(Int_t i) const {return fJetNameCones[i] ; } | |
83 | TString GetJetPtThresName(Int_t i) const {return fJetNamePtThres[i] ; } | |
84 | ||
85 | ||
86 | void SetDeltaPhiCutRange(Double_t phimin, Double_t phimax) | |
87 | {fDeltaPhiMaxCut =phimax; fDeltaPhiMinCut =phimin;} | |
88 | void SetLeadingRatioCutRange(Double_t ratiomin, Double_t ratiomax) | |
89 | {fLeadingRatioMaxCut =ratiomax; fLeadingRatioMinCut = ratiomin ; } | |
90 | ||
91 | void SetJetNCones(Int_t n){fJetNCone = n ; } | |
92 | void SetJetNPtThresholds(Int_t n){fJetNPt = n ; } | |
93 | void SetJetCones(Int_t i, Float_t cone, TString sc) {fJetCones[i] = cone ; fJetNameCones[i] = sc; }; | |
94 | void SetCone(Float_t cone) {fJetCone = cone; } | |
95 | void SetJetPtThreshold(Float_t pt){fJetPtThreshold = pt; }; | |
96 | void SetJetPtThresPbPb(Float_t pt){fJetPtThresPbPb = pt; }; | |
97 | void SetJetPtThresholds(Int_t i,Float_t pt, TString spt){fJetPtThres[i] = pt ; fJetNamePtThres[i] = spt; }; | |
98 | ||
99 | void SetJetRatioCutRange(Double_t ratiomin, Double_t ratiomax) | |
100 | {fJetRatioMaxCut =ratiomax; fJetRatioMinCut = ratiomin ; } | |
101 | void SetJetCTSRatioCutRange(Double_t ratiomin, Double_t ratiomax) | |
102 | {fJetCTSRatioMaxCut =ratiomax; fJetCTSRatioMinCut = ratiomin ; } | |
103 | ||
477d6cee | 104 | Bool_t OnlyIsolated() const {return fSelectIsolated ; } |
105 | void SelectIsolated(Bool_t select) {fSelectIsolated = select ; } | |
106 | ||
1c5acb87 | 107 | private: |
108 | ||
109 | Double_t CalculateJetRatioLimit(const Double_t ptTrig, const Double_t *param, const Double_t *x) const ; | |
110 | ||
111 | void FillJetHistos(AliAODPWG4ParticleCorrelation * particle, const TLorentzVector leading, const TLorentzVector jet, const TString type, const TString lastname); | |
112 | ||
113 | TList * GetOutputContainer() const {return fOutCont; } | |
114 | ||
115 | Bool_t IsJetSelected(const Double_t ptTrig, const Double_t ptjet) const ; | |
116 | Bool_t IsParticleInJetCone(const Double_t eta, Double_t phi, const Double_t etal, Double_t phil) const ; | |
117 | ||
233e0df8 | 118 | void GetLeadingCharge(AliAODPWG4ParticleCorrelation* const particle, TLorentzVector & pLeading) const ; |
decca433 | 119 | void GetLeadingPi0 (AliAODPWG4ParticleCorrelation* const particle, TLorentzVector & pLeading) ; |
120 | Bool_t GetLeadingParticle(AliAODPWG4ParticleCorrelation *particle, TLorentzVector & pLeading) ; | |
1c5acb87 | 121 | |
122 | void MakeAnalysisFillAOD(); | |
123 | void MakeAnalysisFillHistograms(); | |
124 | void MakeAODJet(AliAODPWG4ParticleCorrelation * particle, const TLorentzVector pLeading) const ; | |
125 | void MakeJetFromAOD(AliAODPWG4ParticleCorrelation * particle, const TLorentzVector pLeading, | |
126 | TLorentzVector & jet, TLorentzVector & bkg) const ; | |
127 | ||
0ae57829 | 128 | Bool_t SelectCluster(AliVCluster * calo, Double_t *vertex, TLorentzVector & mom, Int_t & pdg) ; |
1c5acb87 | 129 | |
130 | private: | |
131 | ||
132 | Bool_t fJetsOnlyInCTS ; // Jets measured only in TPC+ITS. | |
133 | Bool_t fPbPb; // PbPb event | |
134 | Bool_t fSeveralConeAndPtCuts; // To play with the jet cone size and pt th. | |
135 | Bool_t fReMakeJet ; //Re make the jet reconstruction from AODParticleCorrelation input | |
136 | ||
137 | //Leading particle selection parameters | |
138 | Double_t fDeltaPhiMaxCut ; // Minimum Delta Phi Gamma-Leading | |
139 | Double_t fDeltaPhiMinCut ; // Maximum Delta Phi Gamma-Leading | |
140 | Double_t fLeadingRatioMaxCut ; // Leading /gamma Ratio cut maximum | |
141 | Double_t fLeadingRatioMinCut ; // Leading/gamma Ratio cut minimum | |
142 | ||
143 | //Jet selection parameters | |
144 | //Fixed cuts (old) | |
145 | Double_t fJetCTSRatioMaxCut ; // Jet(CTS) /gamma Ratio cut maximum | |
146 | Double_t fJetCTSRatioMinCut ; // Jet(CTS) /gamma Ratio cut maximum | |
147 | Double_t fJetRatioMaxCut ; // Jet(EMCAL+CTS)/gamma Ratio cut maximum | |
148 | Double_t fJetRatioMinCut ; // Jet(EMCAL+CTS)/gamma Ratio cut minimum | |
149 | ||
150 | //Cuts depending on jet pt | |
151 | Double_t fJetE1[2]; //Rec. jet energy parameters | |
152 | Double_t fJetE2[2]; //Rec. jet energy parameters | |
153 | Double_t fJetSigma1[2];//Rec. sigma of jet energy parameters | |
154 | Double_t fJetSigma2[2];//Rec. sigma of jet energy parameters | |
155 | Double_t fBkgMean[6]; //Background mean energy | |
156 | Double_t fBkgRMS[6]; //Background RMS | |
157 | Double_t fJetXMin1[6]; //X Factor to set jet min limit for pp | |
158 | Double_t fJetXMin2[6]; //X Factor to set jet min limit for PbPb | |
159 | Double_t fJetXMax1[6]; //X Factor to set jet max limit for pp | |
160 | Double_t fJetXMax2[6]; //X Factor to set jet max limit for PbPb | |
161 | ||
162 | Int_t fJetNCone ; // Number of jet cones sizes, maximum 5 | |
163 | Int_t fJetNPt ; // Number of jet particle pT threshold, maximum 5 | |
164 | Double_t fJetCone ; // Jet cone sizes under study (!fSeveralConeAndPtCuts) | |
165 | Double_t fJetCones[5]; // Jet cone sizes under study (fSeveralConeAndPtCuts) | |
166 | TString fJetNameCones[5]; // String name of cone to append to histos | |
167 | Double_t fJetPtThreshold; // Jet pT threshold under study(!fSeveralConeAndPtCuts) | |
168 | Double_t fJetPtThresPbPb; // Jet pT threshold under study(!fSeveralConeAndPtCuts) | |
169 | Double_t fJetPtThres[5]; // Jet pT threshold under study(fSeveralConeAndPtCuts) | |
170 | TString fJetNamePtThres[5]; // String name of pt th to append to histos | |
171 | Double_t fPtTriggerSelectionCut; // Jet pt to change to low pt jets analysis | |
172 | UInt_t fSelect ; //kTRUE: Selects all jets, no limits. | |
477d6cee | 173 | Bool_t fSelectIsolated ; // Select only trigger particles isolated |
174 | ||
1c5acb87 | 175 | //Histograms |
176 | //Leading particle distributions | |
177 | TList * fOutCont ; //! Container for histograms | |
178 | ||
179 | TH2F * fhChargedLeadingPt ; //! Pt(Pt trigger) distribution of charged hadrons | |
180 | TH2F * fhChargedLeadingPhi ; //! Phi(Pt trigger) distribution of charged hadrons | |
181 | TH2F * fhChargedLeadingEta ; //! Eta(Pt trigger) distribution of charged hadrons | |
182 | TH2F * fhChargedLeadingDeltaPt ; //! Difference of charged hadron and trigger pT as function of trigger p | |
183 | TH2F * fhChargedLeadingDeltaPhi ; //! Difference of charged hadron and trigger phi as function of trigger pT | |
184 | TH2F * fhChargedLeadingDeltaEta ; //! Difference of charged particle and trigger eta as function of trigger pT | |
185 | TH2F * fhChargedLeadingRatioPt ; //! Ratio of Pt leading charge and trigger | |
186 | ||
187 | TH2F * fhNeutralLeadingPt ; //! Pt(Pt trigger) distribution of neutral hadrons | |
188 | TH2F * fhNeutralLeadingPhi ; //! Phi(Pt trigger) distribution of neutral hadrons | |
189 | TH2F * fhNeutralLeadingEta ; //! Eta(Pt trigger) distribution of neutral hadrons | |
190 | TH2F * fhNeutralLeadingDeltaPt ; //! Difference of neutral hadron and trigger pT as function of trigger pT | |
191 | TH2F * fhNeutralLeadingDeltaPhi ; //! Difference of neutral hadron and trigger phi as function of trigger pT | |
192 | TH2F * fhNeutralLeadingDeltaEta ; //! Difference of charged particle and trigger eta as function of trigger pT | |
193 | TH2F * fhNeutralLeadingRatioPt ; //! Ratio of Pt leading neutral and trigger | |
9415d854 | 194 | |
195 | TH2F * fhChargedLeadingXi ; //! Ln (pt leading charge / pt trigger) | |
196 | TH2F * fhNeutralLeadingXi ; //! Ln (pt leading neutral / pt trigger) | |
dde5a268 | 197 | |
198 | TH2F * fhChargedLeadingDeltaPhiRatioPt30 ; //! Difference of charged hadron and trigger phi as function of pT leading / trigger pT, pT Trigger > 30 GeV | |
199 | TH2F * fhNeutralLeadingDeltaPhiRatioPt30 ; //! Difference of neutral hadron and trigger phi as function of pT leading / trigger pT, pT Trigger > 30 GeV | |
200 | TH2F * fhChargedLeadingDeltaPhiRatioPt50 ; //! Difference of charged hadron and trigger phi as function of pT leading / trigger pT, pT Trigger > 50 GeV | |
201 | TH2F * fhNeutralLeadingDeltaPhiRatioPt50 ; //! Difference of neutral hadron and trigger phi as function of pT leading / trigger pT, pT Trigger > 50 GeV | |
9415d854 | 202 | |
1c5acb87 | 203 | // Jet distributions |
204 | // Fixed cone and pt threshold | |
205 | TH2F * fhJetPt ; //! leading pt jet vs pt trigger | |
206 | TH2F * fhJetRatioPt ; //! Ratio of pt jet and pt trigger | |
207 | TH2F * fhJetDeltaPhi ; //! Delta phi jet-trigger | |
208 | TH2F * fhJetDeltaEta ; //! Delta eta jet-trigger | |
209 | TH2F * fhJetLeadingRatioPt ; //! Ratio of pt leading and pt jet | |
210 | TH2F * fhJetLeadingDeltaPhi ; //! Delta phi jet-leading | |
211 | TH2F * fhJetLeadingDeltaEta ; //! Delta eta jet-leading | |
212 | TH2F * fhJetFFz; //! Accepted reconstructed jet fragmentation function, z=ptjet/pttrig | |
213 | TH2F * fhJetFFxi; //! Accepted reconstructed jet fragmentation function, xsi = ln(pttrig/ptjet) | |
214 | TH2F * fhJetFFpt; //! Jet particle pt distribution in cone | |
215 | TH2F * fhJetNTracksInCone ; //! jet multiplicity in cone | |
216 | ||
217 | TH2F * fhBkgPt ; //! leading pt bakground vs pt trigger | |
218 | TH2F * fhBkgRatioPt ; //! Ratio of pt background and pt trigger | |
219 | TH2F * fhBkgDeltaPhi ; //! Delta phi background-trigger | |
220 | TH2F * fhBkgDeltaEta ; //! Delta eta background-trigger | |
221 | TH2F * fhBkgLeadingRatioPt ; //! Ratio of pt leading and pt background | |
222 | TH2F * fhBkgLeadingDeltaPhi ; //! Delta phi background-leading | |
223 | TH2F * fhBkgLeadingDeltaEta ; //! Delta eta background-leading | |
224 | TH2F * fhBkgFFz; //! Accepted reconstructed background fragmentation function, z=ptjet/pttrig | |
225 | TH2F * fhBkgFFxi; //! Accepted reconstructed background fragmentation function, xsi = ln(pttrig/ptjet) | |
226 | TH2F * fhBkgFFpt; //! Background particle pt distribution in cone | |
227 | TH2F * fhBkgNTracksInCone ; //! Background multiplicity in cone | |
228 | ||
229 | // Variable cone and pt threshold | |
230 | ||
231 | TH2F * fhJetPts[5][5]; //! leading pt jet vs pt trigger | |
232 | TH2F * fhJetRatioPts[5][5]; //! Ratio of pt jet and pt trigger | |
233 | TH2F * fhJetDeltaPhis[5][5]; //! Delta phi jet-trigger | |
234 | TH2F * fhJetDeltaEtas[5][5]; //! Delta eta jet-trigger | |
235 | TH2F * fhJetLeadingRatioPts[5][5]; //! Ratio of pt leading and pt jet | |
236 | TH2F * fhJetLeadingDeltaPhis[5][5]; //! Delta phi jet-leading | |
237 | TH2F * fhJetLeadingDeltaEtas[5][5]; //! Delta eta jet-leading | |
238 | TH2F * fhJetFFzs[5][5]; //! Accepted reconstructed jet fragmentation function, z=ptjet/pttrig | |
239 | TH2F * fhJetFFxis[5][5]; //! Accepted reconstructed jet fragmentation function, xsi = ln(pttrig/ptjet) | |
240 | TH2F * fhJetFFpts[5][5]; //! Jet particle pt distribution in cone | |
241 | TH2F * fhJetNTracksInCones[5][5]; //! jet multiplicity in cone | |
242 | ||
243 | TH2F * fhBkgPts[5][5]; //! leading pt bakground vs pt trigger | |
244 | TH2F * fhBkgRatioPts[5][5]; //! Ratio of pt background and pt trigger | |
245 | TH2F * fhBkgDeltaPhis[5][5]; //! Delta phi background-trigger | |
246 | TH2F * fhBkgDeltaEtas[5][5]; //! Delta eta background-trigger | |
247 | TH2F * fhBkgLeadingRatioPts[5][5]; //! Ratio of pt leading and pt background | |
248 | TH2F * fhBkgLeadingDeltaPhis[5][5]; //! Delta phi background-leading | |
249 | TH2F * fhBkgLeadingDeltaEtas[5][5]; //! Delta eta background-leading | |
250 | TH2F * fhBkgFFzs[5][5]; //! Accepted reconstructed background fragmentation function, z=ptjet/pttrig | |
251 | TH2F * fhBkgFFxis[5][5]; //! Accepted reconstructed background fragmentation function, xsi = ln(pttrig/ptjet) | |
252 | TH2F * fhBkgFFpts[5][5]; //! Background particle pt distribution in cone | |
253 | TH2F * fhBkgNTracksInCones[5][5]; //! Background multiplicity in cone | |
254 | ||
255 | ||
256 | ClassDef(AliAnaParticleJetLeadingConeCorrelation,1) | |
257 | } ; | |
258 | ||
259 | ||
260 | #endif //ALIANAPARTICLEJETLEADINGCONECORRELATION_H | |
261 | ||
262 | ||
263 |