1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 //_________________________________________________________________________
17 // Do photon/pi0 analysis for isolation and correlation
18 // at the generator level. Only for kine stack (ESDs)
21 // -- Author: Gustavo Conesa (LPSC-CNRS-Grenoble)
22 //////////////////////////////////////////////////////////////////////////////
24 // --- ROOT system ---
26 #include "TParticle.h"
27 #include "TDatabasePDG.h"
29 //---- ANALYSIS system ----
30 #include "AliAnaGeneratorKine.h"
32 #include "AliGenPythiaEventHeader.h"
34 ClassImp(AliAnaGeneratorKine)
37 //__________________________________________
38 AliAnaGeneratorKine::AliAnaGeneratorKine() :
39 AliAnaCaloTrackCorrBaseClass(),
42 fMinChargedPt(0), fMinNeutralPt(0),
44 fParton2(0), fParton3(0),
45 fParton6(0), fParton7(0),
48 fhPtHard(0), fhPtParton(0), fhPtJet(0),
49 fhPtPartonPtHard(0), fhPtJetPtHard(0), fhPtJetPtParton(0),
50 fhPtPhoton(0), fhPtPi0(0)
54 //Initialize parameters
57 for(Int_t i = 0; i < 4; i++)
59 fhPtPhotonLeading[i] = fhPtPi0Leading[i] = 0;
60 fhPtPhotonLeadingSumPt[i] = fhPtPi0LeadingSumPt[i] = 0;
61 fhPtPhotonLeadingIsolated[i] = fhPtPi0LeadingIsolated[i] = 0;
62 for(Int_t j = 0; j < 2; j++)
64 fhZHardPhoton[j][i] = fhZHardPi0[j][i] = 0;
65 fhZHardPhotonIsolated[j][i] = fhZHardPi0Isolated[j][i] = 0;
66 fhZPartonPhoton[j][i] = fhZPartonPi0[j][i] = 0;
67 fhZPartonPhotonIsolated[j][i] = fhZPartonPi0Isolated[j][i] = 0;
68 fhZJetPhoton[j][i] = fhZJetPi0[j][i] = 0;
69 fhZJetPhotonIsolated[j][i] = fhZJetPi0Isolated[j][i] = 0;
70 fhXEPhoton[j][i] = fhXEPi0[j][i] = 0;
71 fhXEPhotonIsolated[j][i] = fhXEPi0Isolated[j][i] = 0;
72 fhXEUEPhoton[j][i] = fhXEUEPi0[j][i] = 0;
73 fhXEUEPhotonIsolated[j][i] = fhXEUEPi0Isolated[j][i] = 0;
75 fhPtPartonTypeNearPhoton[j][i] = fhPtPartonTypeNearPi0[j][i] = 0;
76 fhPtPartonTypeNearPhotonIsolated[j][i] = fhPtPartonTypeNearPi0Isolated[j][i] = 0;
78 fhPtPartonTypeAwayPhoton[j][i] = fhPtPartonTypeAwayPi0[j][i] = 0;
79 fhPtPartonTypeAwayPhotonIsolated[j][i] = fhPtPartonTypeAwayPi0Isolated[j][i] = 0;
85 //___________________________________________________________________________
86 Bool_t AliAnaGeneratorKine::CorrelateWithPartonOrJet(TLorentzVector trigger,
93 //Correlate trigger with partons or jets, get z
95 if(GetDebug() > 1) printf("AliAnaGeneratorKine::CorrelateWithPartonOrJet() - Start \n");
97 //Get the index of the mother
98 iparton = (fStack->Particle(indexTrig))->GetFirstMother();
99 TParticle * mother = fStack->Particle(iparton);
102 iparton = mother->GetFirstMother();
103 if(iparton < 0) { printf("AliAnaGeneratorKine::CorrelateWithPartonOrJet() - Negative index, skip event\n"); return kFALSE; }
104 mother = fStack->Particle(iparton);
107 //printf("Mother is parton %d with pdg %d\n",imom,mother->GetPdgCode());
111 //printf("This particle is not from hard process - pdg %d, parton index %d\n",pdgTrig, iparton);
115 Float_t ptTrig = trigger.Pt();
116 Float_t partonPt = fParton6->Pt();
117 Float_t jetPt = fJet6.Pt();
120 partonPt = fParton6->Pt();
124 //Get id of parton in near and away side
131 //printf("parton 6 pdg = %d, parton 7 pdg = %d\n",fParton6->GetPdgCode(),fParton7->GetPdgCode());
135 nearPDG = fParton6->GetPdgCode();
136 awayPDG = fParton7->GetPdgCode();
140 nearPDG = fParton7->GetPdgCode();
141 awayPDG = fParton6->GetPdgCode();
144 if (nearPDG == 22) near = 0;
145 else if(nearPDG == 21) near = 1;
148 if (awayPDG == 22) away = 0;
149 else if(awayPDG == 21) away = 1;
152 for( Int_t i = 0; i < 4; i++ )
157 fhPtPartonTypeNearPi0[leading[i]][i]->Fill(ptTrig,near);
158 fhPtPartonTypeAwayPi0[leading[i]][i]->Fill(ptTrig,away);
161 fhPtPartonTypeNearPi0Isolated[leading[i]][i]->Fill(ptTrig,near);
162 fhPtPartonTypeAwayPi0Isolated[leading[i]][i]->Fill(ptTrig,away);
168 fhPtPartonTypeNearPhoton[leading[i]][i]->Fill(ptTrig,near);
169 fhPtPartonTypeAwayPhoton[leading[i]][i]->Fill(ptTrig,away);
172 fhPtPartonTypeNearPhotonIsolated[leading[i]][i]->Fill(ptTrig,near);
173 fhPtPartonTypeAwayPhotonIsolated[leading[i]][i]->Fill(ptTrig,away);
182 Float_t zHard = ptTrig / fPtHard;
183 Float_t zPart = ptTrig / partonPt;
184 Float_t zJet = ptTrig / jetPt;
186 //if(zHard > 1 ) printf("*** Particle energy larger than pT hard z=%f\n",zHard);
188 //printf("Z: hard %2.2f, parton %2.2f, jet %2.2f\n",zHard,zPart,zJet);
190 for( Int_t i = 0; i < 4; i++ )
194 fhZHardPi0[leading[i]][i] ->Fill(ptTrig,zHard);
195 fhZPartonPi0[leading[i]][i]->Fill(ptTrig,zPart);
196 fhZJetPi0[leading[i]][i] ->Fill(ptTrig,zJet );
200 fhZHardPi0Isolated[leading[i]][i] ->Fill(ptTrig,zHard);
201 fhZPartonPi0Isolated[leading[i]][i]->Fill(ptTrig,zPart);
202 fhZJetPi0Isolated[leading[i]][i] ->Fill(ptTrig,zJet);
209 fhZHardPhoton[leading[i]][i] ->Fill(ptTrig,zHard);
210 fhZPartonPhoton[leading[i]][i]->Fill(ptTrig,zPart);
211 fhZJetPhoton[leading[i]][i] ->Fill(ptTrig,zJet );
215 fhZHardPhotonIsolated[leading[i]][i] ->Fill(ptTrig,zHard);
216 fhZPartonPhotonIsolated[leading[i]][i]->Fill(ptTrig,zPart);
217 fhZJetPhotonIsolated[leading[i]][i] ->Fill(ptTrig,zJet);
223 if(GetDebug() > 1) printf("AliAnaGeneratorKine::CorrelateWithPartonOrJet() - End TRUE \n");
229 //____________________________________________________
230 TList * AliAnaGeneratorKine::GetCreateOutputObjects()
232 // Create histograms to be saved in output file
234 TList * outputContainer = new TList() ;
235 outputContainer->SetName("GenKineHistos") ;
237 Int_t nptbins = GetHistogramRanges()->GetHistoPtBins();
238 Float_t ptmax = GetHistogramRanges()->GetHistoPtMax();
239 Float_t ptmin = GetHistogramRanges()->GetHistoPtMin();
241 Int_t nptsumbins = GetHistogramRanges()->GetHistoNPtSumBins();
242 Float_t ptsummax = GetHistogramRanges()->GetHistoPtSumMax();
243 Float_t ptsummin = GetHistogramRanges()->GetHistoPtSumMin();
246 fhPtHard = new TH1F("hPtHard"," pt hard for selected triggers",nptbins,ptmin,ptmax);
247 fhPtHard->SetXTitle("#it{p}_{T}^{hard} (GeV/#it{c})");
248 outputContainer->Add(fhPtHard);
250 fhPtParton = new TH1F("hPtParton"," pt parton for selected triggers",nptbins,ptmin,ptmax);
251 fhPtParton->SetXTitle("#it{p}_{T}^{parton} (GeV/#it{c})");
252 outputContainer->Add(fhPtParton);
254 fhPtJet = new TH1F("hPtJet"," pt jet for selected triggers",nptbins,ptmin,ptmax);
255 fhPtJet->SetXTitle("#it{p}_{T}^{jet} (GeV/#it{c})");
256 outputContainer->Add(fhPtJet);
258 fhPtPartonPtHard = new TH2F("hPtPartonPtHard","parton pt / pt hard for selected triggers",nptbins,ptmin,ptmax,200,0,2);
259 fhPtPartonPtHard->SetXTitle("#it{p}_{T}^{hard} (GeV/#it{c})");
260 fhPtPartonPtHard->SetYTitle("#it{p}_{T}^{parton}/#it{p}_{T}^{hard}");
261 outputContainer->Add(fhPtPartonPtHard);
263 fhPtJetPtHard = new TH2F("hPtJetPtHard","jet pt / pt hard for selected triggers",nptbins,ptmin,ptmax,200,0,2);
264 fhPtJetPtHard->SetXTitle("#it{p}_{T}^{hard} (GeV/#it{c})");
265 fhPtJetPtHard->SetYTitle("#it{p}_{T}^{jet}/#it{p}_{T}^{hard}");
266 outputContainer->Add(fhPtJetPtHard);
268 fhPtJetPtParton = new TH2F("hPtJetPtParton","parton pt / pt hard for selected triggers",nptbins,ptmin,ptmax,200,0,2);
269 fhPtJetPtParton->SetXTitle("#it{p}_{T}^{hard} (GeV/#it{c})");
270 fhPtJetPtParton->SetYTitle("#it{p}_{T}^{jet}/#it{p}_{T}^{parton}");
271 outputContainer->Add(fhPtJetPtParton);
273 fhPtPhoton = new TH1F("hPtPhoton","Input #gamma",nptbins,ptmin,ptmax);
274 fhPtPhoton->SetXTitle("#it{p}_{T} (GeV/#it{c})");
275 outputContainer->Add(fhPtPhoton);
277 fhPtPi0 = new TH1F("hPtPi0","Input #pi^{0}",nptbins,ptmin,ptmax);
278 fhPtPi0->SetXTitle("#it{p}_{T} (GeV/#it{c})");
279 outputContainer->Add(fhPtPi0);
281 TString name [] = {"","_EMC","_Photon","_EMC_Photon"};
282 TString title [] = {"",", neutral in EMCal",", neutral only #gamma-like",", neutral in EMCal and only #gamma-like"};
283 TString leading[] = {"NotLeading","Leading"};
285 for(Int_t i = 0; i < 4; i++)
290 fhPtPhotonLeading[i] = new TH1F(Form("hPtPhotonLeading%s",name[i].Data()),
291 Form("#gamma: Leading of all particles%s",title[i].Data()),
292 nptbins,ptmin,ptmax);
293 fhPtPhotonLeading[i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
294 outputContainer->Add(fhPtPhotonLeading[i]);
296 fhPtPi0Leading[i] = new TH1F(Form("hPtPi0Leading%s",name[i].Data()),
297 Form("#pi^{0}: Leading of all particles%s",title[i].Data()),
298 nptbins,ptmin,ptmax);
299 fhPtPi0Leading[i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
300 outputContainer->Add(fhPtPi0Leading[i]);
302 fhPtPhotonLeadingIsolated[i] = new TH1F(Form("hPtPhotonLeadingIsolated%s",name[i].Data()),
303 Form("#gamma: Leading of all particles%s, isolated",title[i].Data()),
304 nptbins,ptmin,ptmax);
305 fhPtPhotonLeadingIsolated[i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
306 outputContainer->Add(fhPtPhotonLeadingIsolated[i]);
308 fhPtPi0LeadingIsolated[i] = new TH1F(Form("hPtPi0LeadingIsolated%s",name[i].Data()),
309 Form("#pi^{0}: Leading of all particles%s, isolated",title[i].Data()),
310 nptbins,ptmin,ptmax);
311 fhPtPi0LeadingIsolated[i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
312 outputContainer->Add(fhPtPi0LeadingIsolated[i]);
314 fhPtPhotonLeadingSumPt[i] = new TH2F(Form("hPtPhotonLeadingSumPt%s",name[i].Data()),
315 Form("#gamma: Leading of all particles%s",title[i].Data()),
316 nptbins,ptmin,ptmax,nptsumbins,ptsummin,ptsummax);
317 fhPtPhotonLeadingSumPt[i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
318 fhPtPhotonLeadingSumPt[i]->SetYTitle("#Sigma #it{p}_{T} (GeV/#it{c})");
319 outputContainer->Add(fhPtPhotonLeadingSumPt[i]);
321 fhPtPi0LeadingSumPt[i] = new TH2F(Form("hPtPi0LeadingSumPt%s",name[i].Data()),
322 Form("#pi^{0}: Leading of all particles%s",title[i].Data()),
323 nptbins,ptmin,ptmax,nptsumbins,ptsummin,ptsummax);
324 fhPtPi0LeadingSumPt[i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
325 fhPtPi0LeadingSumPt[i]->SetYTitle("#Sigma #it{p}_{T} (GeV/#it{c})");
326 outputContainer->Add(fhPtPi0LeadingSumPt[i]);
329 // Leading or not loop
330 for(Int_t j = 0; j < 2; j++)
334 fhPtPartonTypeNearPhoton[j][i] = new TH2F(Form("hPtPartonTypeNearPhoton%s%s",leading[j].Data(),name[i].Data()),
335 Form("#gamma: %s of all particles%s",leading[j].Data(),title[i].Data()),
336 nptbins,ptmin,ptmax,3,0,3);
337 fhPtPartonTypeNearPhoton[j][i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
338 fhPtPartonTypeNearPhoton[j][i]->SetYTitle("Parton type");
339 fhPtPartonTypeNearPhoton[j][i]->GetYaxis()->SetBinLabel(1,"#gamma");
340 fhPtPartonTypeNearPhoton[j][i]->GetYaxis()->SetBinLabel(2,"g");
341 fhPtPartonTypeNearPhoton[j][i]->GetYaxis()->SetBinLabel(3,"q");
342 outputContainer->Add(fhPtPartonTypeNearPhoton[j][i]);
344 fhPtPartonTypeNearPi0[j][i] = new TH2F(Form("hPtPartonTypeNearPi0%s%s",leading[j].Data(),name[i].Data()),
345 Form("#pi^{0}: %s of all particles%s",leading[j].Data(),title[i].Data()),
346 nptbins,ptmin,ptmax,3,0,3);
347 fhPtPartonTypeNearPi0[j][i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
348 fhPtPartonTypeNearPi0[j][i]->SetYTitle("Parton type");
349 fhPtPartonTypeNearPi0[j][i]->GetYaxis()->SetBinLabel(1,"#gamma");
350 fhPtPartonTypeNearPi0[j][i]->GetYaxis()->SetBinLabel(2,"g");
351 fhPtPartonTypeNearPi0[j][i]->GetYaxis()->SetBinLabel(3,"q");
352 outputContainer->Add(fhPtPartonTypeNearPi0[j][i]);
354 fhPtPartonTypeNearPhotonIsolated[j][i] = new TH2F(Form("hPtPartonTypeNearPhoton%sIsolated%s",leading[j].Data(),name[i].Data()),
355 Form("#gamma: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
356 nptbins,ptmin,ptmax,3,0,3);
357 fhPtPartonTypeNearPhotonIsolated[j][i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
358 fhPtPartonTypeNearPhotonIsolated[j][i]->SetYTitle("Parton type");
359 fhPtPartonTypeNearPhotonIsolated[j][i]->GetYaxis()->SetBinLabel(1,"#gamma");
360 fhPtPartonTypeNearPhotonIsolated[j][i]->GetYaxis()->SetBinLabel(2,"g");
361 fhPtPartonTypeNearPhotonIsolated[j][i]->GetYaxis()->SetBinLabel(3,"q");
362 outputContainer->Add(fhPtPartonTypeNearPhotonIsolated[j][i]);
364 fhPtPartonTypeNearPi0Isolated[j][i] = new TH2F(Form("hPtPartonTypeNearPi0%sIsolated%s",leading[j].Data(),name[i].Data()),
365 Form("#pi^{0}: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
366 nptbins,ptmin,ptmax,3,0,3);
367 fhPtPartonTypeNearPi0Isolated[j][i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
368 fhPtPartonTypeNearPi0Isolated[j][i]->SetYTitle("Parton type");
369 fhPtPartonTypeNearPi0Isolated[j][i]->GetYaxis()->SetBinLabel(1,"#gamma");
370 fhPtPartonTypeNearPi0Isolated[j][i]->GetYaxis()->SetBinLabel(2,"g");
371 fhPtPartonTypeNearPi0Isolated[j][i]->GetYaxis()->SetBinLabel(3,"q");
372 outputContainer->Add(fhPtPartonTypeNearPi0Isolated[j][i]);
377 fhPtPartonTypeAwayPhoton[j][i] = new TH2F(Form("hPtPartonTypeAwayPhoton%s%s",leading[j].Data(),name[i].Data()),
378 Form("#gamma: %s of all particles%s",leading[j].Data(),title[i].Data()),
379 nptbins,ptmin,ptmax,3,0,3);
380 fhPtPartonTypeAwayPhoton[j][i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
381 fhPtPartonTypeAwayPhoton[j][i]->SetYTitle("Parton type");
382 fhPtPartonTypeAwayPhoton[j][i]->GetYaxis()->SetBinLabel(1,"#gamma");
383 fhPtPartonTypeAwayPhoton[j][i]->GetYaxis()->SetBinLabel(2,"g");
384 fhPtPartonTypeAwayPhoton[j][i]->GetYaxis()->SetBinLabel(3,"q");
385 outputContainer->Add(fhPtPartonTypeAwayPhoton[j][i]);
387 fhPtPartonTypeAwayPi0[j][i] = new TH2F(Form("hPtPartonTypeAwayPi0%s%s",leading[j].Data(),name[i].Data()),
388 Form("#pi^{0}: %s of all particles%s",leading[j].Data(),title[i].Data()),
389 nptbins,ptmin,ptmax,3,0,3);
390 fhPtPartonTypeAwayPi0[j][i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
391 fhPtPartonTypeAwayPi0[j][i]->SetYTitle("Parton type");
392 fhPtPartonTypeAwayPi0[j][i]->GetYaxis()->SetBinLabel(1,"#gamma");
393 fhPtPartonTypeAwayPi0[j][i]->GetYaxis()->SetBinLabel(2,"g");
394 fhPtPartonTypeAwayPi0[j][i]->GetYaxis()->SetBinLabel(3,"q");
395 outputContainer->Add(fhPtPartonTypeAwayPi0[j][i]);
397 fhPtPartonTypeAwayPhotonIsolated[j][i] = new TH2F(Form("hPtPartonTypeAwayPhoton%sIsolated%s",leading[j].Data(),name[i].Data()),
398 Form("#gamma: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
399 nptbins,ptmin,ptmax,3,0,3);
400 fhPtPartonTypeAwayPhotonIsolated[j][i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
401 fhPtPartonTypeAwayPhotonIsolated[j][i]->SetYTitle("Parton type");
402 fhPtPartonTypeAwayPhotonIsolated[j][i]->GetYaxis()->SetBinLabel(1,"#gamma");
403 fhPtPartonTypeAwayPhotonIsolated[j][i]->GetYaxis()->SetBinLabel(2,"g");
404 fhPtPartonTypeAwayPhotonIsolated[j][i]->GetYaxis()->SetBinLabel(3,"q");
405 outputContainer->Add(fhPtPartonTypeAwayPhotonIsolated[j][i]);
407 fhPtPartonTypeAwayPi0Isolated[j][i] = new TH2F(Form("hPtPartonTypeAwayPi0%sIsolated%s",leading[j].Data(),name[i].Data()),
408 Form("#pi^{0}: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
409 nptbins,ptmin,ptmax,3,0,3);
410 fhPtPartonTypeAwayPi0Isolated[j][i]->SetXTitle("#it{p}_{T} (GeV/#it{c})");
411 fhPtPartonTypeAwayPi0Isolated[j][i]->SetYTitle("Parton type");
412 fhPtPartonTypeAwayPi0Isolated[j][i]->GetYaxis()->SetBinLabel(1,"#gamma");
413 fhPtPartonTypeAwayPi0Isolated[j][i]->GetYaxis()->SetBinLabel(2,"g");
414 fhPtPartonTypeAwayPi0Isolated[j][i]->GetYaxis()->SetBinLabel(3,"q");
415 outputContainer->Add(fhPtPartonTypeAwayPi0Isolated[j][i]);
419 fhZHardPhoton[j][i] = new TH2F(Form("hZHardPhoton%s%s",leading[j].Data(),name[i].Data()),
420 Form("#it{z}_{Hard} of #gamma: %s of all particles%s",leading[j].Data(),title[i].Data()),
421 nptbins,ptmin,ptmax,200,0,2);
422 fhZHardPhoton[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
423 fhZHardPhoton[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
424 outputContainer->Add(fhZHardPhoton[j][i]);
426 fhZHardPi0[j][i] = new TH2F(Form("hZHardPi0%s%s",leading[j].Data(),name[i].Data()),
427 Form("#it{z}_{Hard} of #pi^{0}: %s of all particles%s",leading[j].Data(),title[i].Data()),
428 nptbins,ptmin,ptmax,200,0,2);
429 fhZHardPi0[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
430 fhZHardPi0[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
431 outputContainer->Add(fhZHardPi0[j][i]);
433 fhZHardPhotonIsolated[j][i] = new TH2F(Form("hZHardPhoton%sIsolated%s",leading[j].Data(),name[i].Data()),
434 Form("#it{z}_{Hard} of #gamma: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
435 nptbins,ptmin,ptmax,200,0,2);
436 fhZHardPhotonIsolated[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
437 fhZHardPhotonIsolated[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
438 outputContainer->Add(fhZHardPhotonIsolated[j][i]);
440 fhZHardPi0Isolated[j][i] = new TH2F(Form("hZHardPi0%sIsolated%s",leading[j].Data(),name[i].Data()),
441 Form("#it{z}_{Hard} of #pi^{0}: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
442 nptbins,ptmin,ptmax,200,0,2);
443 fhZHardPi0Isolated[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
444 fhZHardPi0Isolated[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
445 outputContainer->Add(fhZHardPi0Isolated[j][i]);
449 fhZPartonPhoton[j][i] = new TH2F(Form("hZPartonPhoton%s%s",leading[j].Data(),name[i].Data()),
450 Form("#it{z}_{Parton} of #gamma: %s of all particles%s",leading[j].Data(),title[i].Data()),
451 nptbins,ptmin,ptmax,200,0,2);
452 fhZPartonPhoton[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
453 fhZPartonPhoton[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
454 outputContainer->Add(fhZPartonPhoton[j][i]);
456 fhZPartonPi0[j][i] = new TH2F(Form("hZPartonPi0%s%s",leading[j].Data(),name[i].Data()),
457 Form("#it{z}_{Parton} of #pi^{0}: %s of all particles%s",leading[j].Data(),title[i].Data()),
458 nptbins,ptmin,ptmax,200,0,2);
459 fhZPartonPi0[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
460 fhZPartonPi0[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
461 outputContainer->Add(fhZPartonPi0[j][i]);
463 fhZPartonPhotonIsolated[j][i] = new TH2F(Form("hZPartonPhoton%sIsolated%s",leading[j].Data(),name[i].Data()),
464 Form("#it{z}_{Parton} of #gamma: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
465 nptbins,ptmin,ptmax,200,0,2);
466 fhZPartonPhotonIsolated[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
467 fhZPartonPhotonIsolated[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
468 outputContainer->Add(fhZPartonPhotonIsolated[j][i]);
470 fhZPartonPi0Isolated[j][i] = new TH2F(Form("hZPartonPi0%sIsolated%s",leading[j].Data(),name[i].Data()),
471 Form("#it{z}_{Parton} of #pi^{0}: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
472 nptbins,ptmin,ptmax,200,0,2);
473 fhZPartonPi0Isolated[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
474 fhZPartonPi0Isolated[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
475 outputContainer->Add(fhZPartonPi0Isolated[j][i]);
480 fhZJetPhoton[j][i] = new TH2F(Form("hZJetPhoton%s%s",leading[j].Data(),name[i].Data()),
481 Form("#it{z}_{Jet} of #gamma: %s of all particles%s",leading[j].Data(),title[i].Data()),
482 nptbins,ptmin,ptmax,200,0,2);
483 fhZJetPhoton[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
484 fhZJetPhoton[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
485 outputContainer->Add(fhZJetPhoton[j][i]);
487 fhZJetPi0[j][i] = new TH2F(Form("hZJetPi0%s%s",leading[j].Data(),name[i].Data()),
488 Form("#it{z}_{Jet} of #pi^{0}: %s of all particles%s",leading[j].Data(),title[i].Data()),
489 nptbins,ptmin,ptmax,200,0,2);
490 fhZJetPi0[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
491 fhZJetPi0[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
492 outputContainer->Add(fhZJetPi0[j][i]);
494 fhZJetPhotonIsolated[j][i] = new TH2F(Form("hZJetPhoton%sIsolated%s",leading[j].Data(),name[i].Data()),
495 Form("#it{z}_{Jet} of #gamma: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
496 nptbins,ptmin,ptmax,200,0,2);
497 fhZJetPhotonIsolated[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
498 fhZJetPhotonIsolated[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
499 outputContainer->Add(fhZJetPhotonIsolated[j][i]);
501 fhZJetPi0Isolated[j][i] = new TH2F(Form("hZJetPi0%sIsolated%s",leading[j].Data(),name[i].Data()),
502 Form("#it{z}_{Jet} of #pi^{0}: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
503 nptbins,ptmin,ptmax,200,0,2);
504 fhZJetPi0Isolated[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
505 fhZJetPi0Isolated[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
506 outputContainer->Add(fhZJetPi0Isolated[j][i]);
511 fhXEPhoton[j][i] = new TH2F(Form("hXEPhoton%s%s",leading[j].Data(),name[i].Data()),
512 Form("#it{z}_{Jet} of #gamma: %s of all particles%s",leading[j].Data(),title[i].Data()),
513 nptbins,ptmin,ptmax,200,0,2);
514 fhXEPhoton[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
515 fhXEPhoton[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
516 outputContainer->Add(fhXEPhoton[j][i]);
518 fhXEPi0[j][i] = new TH2F(Form("hXEPi0%s%s",leading[j].Data(),name[i].Data()),
519 Form("#it{z}_{Jet} of #pi^{0}: %s of all particles%s",leading[j].Data(),title[i].Data()),
520 nptbins,ptmin,ptmax,200,0,2);
521 fhXEPi0[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
522 fhXEPi0[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
523 outputContainer->Add(fhXEPi0[j][i]);
525 fhXEPhotonIsolated[j][i] = new TH2F(Form("hXEPhoton%sIsolated%s",leading[j].Data(),name[i].Data()),
526 Form("#it{z}_{Jet} of #gamma: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
527 nptbins,ptmin,ptmax,200,0,2);
528 fhXEPhotonIsolated[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
529 fhXEPhotonIsolated[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
530 outputContainer->Add(fhXEPhotonIsolated[j][i]);
532 fhXEPi0Isolated[j][i] = new TH2F(Form("hXEPi0%sIsolated%s",leading[j].Data(),name[i].Data()),
533 Form("#it{z}_{Jet} of #pi^{0}: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
534 nptbins,ptmin,ptmax,200,0,2);
535 fhXEPi0Isolated[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
536 fhXEPi0Isolated[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
537 outputContainer->Add(fhXEPi0Isolated[j][i]);
542 fhXEUEPhoton[j][i] = new TH2F(Form("hXEUEPhoton%s%s",leading[j].Data(),name[i].Data()),
543 Form("#it{z}_{Jet} of #gamma: %s of all particles%s",leading[j].Data(),title[i].Data()),
544 nptbins,ptmin,ptmax,200,0,2);
545 fhXEUEPhoton[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
546 fhXEUEPhoton[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
547 outputContainer->Add(fhXEUEPhoton[j][i]);
549 fhXEUEPi0[j][i] = new TH2F(Form("hXEUEPi0%s%s",leading[j].Data(),name[i].Data()),
550 Form("#it{z}_{Jet} of #pi^{0}: %s of all particles%s",leading[j].Data(),title[i].Data()),
551 nptbins,ptmin,ptmax,200,0,2);
552 fhXEUEPi0[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
553 fhXEUEPi0[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
554 outputContainer->Add(fhXEUEPi0[j][i]);
556 fhXEUEPhotonIsolated[j][i] = new TH2F(Form("hXEUEPhoton%sIsolated%s",leading[j].Data(),name[i].Data()),
557 Form("#it{z}_{Jet} of #gamma: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
558 nptbins,ptmin,ptmax,200,0,2);
559 fhXEUEPhotonIsolated[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
560 fhXEUEPhotonIsolated[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
561 outputContainer->Add(fhXEUEPhotonIsolated[j][i]);
563 fhXEUEPi0Isolated[j][i] = new TH2F(Form("hXEUEPi0%sIsolated%s",leading[j].Data(),name[i].Data()),
564 Form("#it{z}_{Jet} of #pi^{0}: %s of all particles%s, isolated",leading[j].Data(),title[i].Data()),
565 nptbins,ptmin,ptmax,200,0,2);
566 fhXEUEPi0Isolated[j][i]->SetYTitle("#it{p}_{T}^{particle}/#it{p}_{T}^{hard}");
567 fhXEUEPi0Isolated[j][i]->SetXTitle("#it{p}_{T}^{particle} (GeV/#it{c})");
568 outputContainer->Add(fhXEUEPi0Isolated[j][i]);
572 return outputContainer;
576 //____________________________________________
577 void AliAnaGeneratorKine::GetPartonsAndJets()
579 // Fill data members with partons,jets and generated pt hard
581 if(GetDebug() > 1) printf("AliAnaGeneratorKine::GetPartonsAndJets() - Start \n");
583 fStack = GetMCStack() ;
587 printf("AliAnaGeneratorKine::MakeAnalysisFillHistograms() - No Stack available, STOP\n");
591 fParton2 = fStack->Particle(2) ;
592 fParton3 = fStack->Particle(3) ;
593 fParton6 = fStack->Particle(6) ;
594 fParton7 = fStack->Particle(7) ;
596 Float_t p6phi = fParton6->Phi();
597 if(p6phi < 0) p6phi +=TMath::TwoPi();
598 Float_t p7phi = fParton7->Phi();
599 if(p7phi < 0) p7phi +=TMath::TwoPi();
601 //printf("parton6: pt %2.2f, eta %2.2f, phi %2.2f with pdg %d\n",fParton6->Pt(),fParton6->Eta(),p6phi, fParton6->GetPdgCode());
602 //printf("parton7: pt %2.2f, eta %2.2f, phi %2.2f with pdg %d\n",fParton7->Pt(),fParton7->Eta(),p7phi, fParton7->GetPdgCode());
604 // Get the jets, only for pythia
605 if(!strcmp(GetReader()->GetGenEventHeader()->ClassName(), "AliGenPythiaEventHeader"))
607 AliGenPythiaEventHeader* pygeh= (AliGenPythiaEventHeader*) GetReader()->GetGenEventHeader();
609 fPtHard = pygeh->GetPtHard();
611 //printf("pt Hard %2.2f\n",fPtHard);
613 const Int_t nTriggerJets = pygeh->NTriggerJets();
615 Float_t tmpjet[]={0,0,0,0};
617 // select the closest jet to parton
621 for(Int_t ijet = 0; ijet< nTriggerJets; ijet++)
623 pygeh->TriggerJet(ijet, tmpjet);
625 TLorentzVector jet(tmpjet[0],tmpjet[1],tmpjet[2],tmpjet[3]);
626 Float_t jphi = jet.Phi();
627 if(jphi < 0) jphi +=TMath::TwoPi();
629 Double_t radius6 = GetIsolationCut()->Radius(fParton6->Eta(), p6phi, jet.Eta() , jphi) ;
630 Double_t radius7 = GetIsolationCut()->Radius(fParton7->Eta(), p7phi, jet.Eta() , jphi) ;
632 //printf("jet %d: pt %2.2f, eta %2.2f, phi %2.2f, r6 %2.2f, r7 %2.2f\n",ijet,jet.Pt(),jet.Eta(),jphi,radius6, radius7);
648 //printf("jet6: pt %2.2f, eta %2.2f, phi %2.2f\n",fJet6.Pt(),fJet6.Eta(),fJet6.Phi());
649 //printf("jet7: pt %2.2f, eta %2.2f, phi %2.2f\n",fJet7.Pt(),fJet7.Eta(),fJet6.Phi());
653 fhPtHard ->Fill(fPtHard);
654 fhPtJet ->Fill(fJet6.Pt());
655 fhPtJet ->Fill(fJet7.Pt());
656 fhPtParton ->Fill(fParton6->Pt());
657 fhPtParton ->Fill(fParton7->Pt());
659 fhPtPartonPtHard->Fill(fPtHard, fParton6->Pt()/fPtHard);
660 fhPtPartonPtHard->Fill(fPtHard, fParton7->Pt()/fPtHard);
661 fhPtJetPtHard ->Fill(fPtHard, fJet6.Pt()/fPtHard);
662 fhPtJetPtHard ->Fill(fPtHard, fJet7.Pt()/fPtHard);
663 fhPtJetPtParton ->Fill(fPtHard, fJet6.Pt()/fParton6->Pt());
664 fhPtJetPtParton ->Fill(fPtHard, fJet7.Pt()/fParton7->Pt());
666 if(GetDebug() > 1) printf("AliAnaGeneratorKine::GetPartonsAndJets() - End \n");
670 //_____________________________________________________
671 void AliAnaGeneratorKine::GetXE(TLorentzVector trigger,
679 // Calculate the real XE and the UE XE
681 if(GetDebug() > 1) printf("AliAnaGeneratorKine::GetXE() - Start \n");
683 Float_t ptTrig = trigger.Pt();
684 Float_t phiTrig = trigger.Phi();
685 if(phiTrig < 0 ) phiTrig += TMath::TwoPi();
687 TLorentzVector chPartLV;
689 //Loop on primaries, start from position 8, no partons
690 for(Int_t ipr = 8; ipr < fStack->GetNprimary(); ipr ++ )
692 TParticle * particle = fStack->Particle(ipr) ;
694 if(ipr==indexTrig) continue;
696 Int_t pdg = particle->GetPdgCode();
697 Int_t status = particle->GetStatusCode();
699 // Compare trigger with final state particles
700 if( status != 1 ) continue ;
702 Double_t charge = TDatabasePDG::Instance()->GetParticle(pdg)->Charge();
704 if(charge == 0 ) continue; // construct xe only with charged
706 Float_t pt = particle->Pt();
707 Float_t phi = particle->Phi();
708 if(phi < 0 ) phi += TMath::TwoPi();
710 if( pt < fMinChargedPt) continue ;
712 particle->Momentum(chPartLV);
713 Bool_t inTPC = GetFiducialCut()->IsInFiducialCut(chPartLV,"CTS") ;
717 Float_t xe = -pt/ptTrig*TMath::Cos(phi-phiTrig);
719 // ---------------------------------------------------
720 // Get the index of the mother, get from what parton
721 Int_t ipartonAway = particle->GetFirstMother();
722 if(ipartonAway < 0) return;
724 TParticle * mother = fStack->Particle(ipartonAway);
725 while (ipartonAway > 7)
727 ipartonAway = mother->GetFirstMother();
728 if(ipartonAway < 0) break;
729 mother = fStack->Particle(ipartonAway);
732 //-----------------------------------------
733 // Get XE of particles belonging to the jet
734 // on the opposite side of the trigger
735 if((ipartonAway==6 || ipartonAway==7) && iparton!=ipartonAway)
737 for( Int_t i = 0; i < 4; i++ )
741 fhXEPi0[leading[i]][i] ->Fill(ptTrig,xe);
745 fhXEPi0Isolated[leading[i]][i] ->Fill(ptTrig,xe);
751 fhXEPhoton[leading[i]][i] ->Fill(ptTrig,xe);
755 fhXEPhotonIsolated[leading[i]][i] ->Fill(ptTrig,xe);
762 //----------------------------------------------------------
763 // Get the XE from particles not attached to any of the jets
764 if(ipartonAway!=6 && ipartonAway!=7)
766 for( Int_t i = 0; i < 4; i++ )
770 fhXEUEPi0[leading[i]][i] ->Fill(ptTrig,xe);
774 fhXEUEPi0Isolated[leading[i]][i] ->Fill(ptTrig,xe);
780 fhXEUEPhoton[leading[i]][i] ->Fill(ptTrig,xe);
784 fhXEUEPhotonIsolated[leading[i]][i] ->Fill(ptTrig,xe);
793 if(GetDebug() > 1) printf("AliAnaGeneratorKine::GetPartonsAndJets() - End \n");
797 //________________________________________
798 void AliAnaGeneratorKine::InitParameters()
801 //Initialize the parameters of the analysis.
802 AddToHistogramsName("AnaGenKine_");
804 fTriggerDetector = "EMCAL";
811 //_____________________________________________________________________
812 void AliAnaGeneratorKine::IsLeadingAndIsolated(TLorentzVector trigger,
818 // Check if the trigger is the leading particle and if it is isolated
819 // In case of neutral particles check all neutral or neutral in EMCAL acceptance
821 if(GetDebug() > 1) printf("AliAnaGeneratorKine::GetIsLeadingAndIsolated() - Start \n");
823 Float_t ptMaxCharged = 0; // all charged
824 Float_t ptMaxNeutral = 0; // all neutral
825 Float_t ptMaxNeutEMCAL = 0; // for neutral, select them in EMCAL acceptance
826 Float_t ptMaxNeutPhot = 0; // for neutral, take only photons
827 Float_t ptMaxNeutEMCALPhot = 0; // for neutral, take only photons in EMCAL acceptance
839 Float_t ptTrig = trigger.Pt();
840 Float_t etaTrig = trigger.Eta();
841 Float_t phiTrig = trigger.Phi();
842 if(phiTrig < 0 ) phiTrig += TMath::TwoPi();
844 // Minimum track or cluster energy
847 Float_t ptThresIC = GetIsolationCut()->GetPtThreshold();
848 Float_t sumThresIC = GetIsolationCut()->GetPtThreshold();
849 Float_t rThresIC = GetIsolationCut()->GetConeSize();
850 Float_t isoMethod = GetIsolationCut()->GetICMethod();
854 Int_t nICNeutral = 0;
855 Int_t nICNeutEMCAL = 0;
856 Int_t nICNeutPhot = 0;
857 Int_t nICNeutEMCALPhot = 0;
861 Float_t sumNePtPhot = 0;
862 Float_t sumNePtEMC = 0;
863 Float_t sumNePtEMCPhot = 0;
866 //Loop on primaries, start from position 8, no partons
867 for(Int_t ipr = 8; ipr < fStack->GetNprimary(); ipr ++ )
869 if(ipr == indexTrig) continue;
870 TParticle * particle = fStack->Particle(ipr) ;
872 Int_t imother= particle->GetFirstMother();
873 //printf("Leading ipr %d - mother %d\n",ipr, imother);
875 if(imother==indexTrig) continue ;
877 Int_t pdg = particle->GetPdgCode();
878 Int_t status = particle->GetStatusCode();
880 // Compare trigger with final state particles
881 if( status != 1) continue ;
883 // Select all particles in at least the TPC acceptance
884 Bool_t inTPC = GetFiducialCut()->IsInFiducialCut(trigger,"CTS") ;
887 Float_t pt = particle->Pt();
888 Float_t eta = particle->Eta();
889 Float_t phi = particle->Phi();
890 if(phi < 0 ) phi += TMath::TwoPi();
892 Double_t charge = TDatabasePDG::Instance()->GetParticle(pdg)->Charge();
895 Double_t radius = GetIsolationCut()->Radius(etaTrig, phiTrig, eta , phi) ;
899 if(pt < fMinNeutralPt) continue ;
901 if( ptMaxNeutral < pt ) ptMaxNeutral = pt;
903 if( radius < rThresIC )
905 if( pt > ptThresIC ) nICNeutral++ ;
909 Bool_t phPDG = kFALSE;
910 if(pdg==22 || pdg==111) phPDG = kTRUE;
912 //if(pt > ptTrig) printf(" --- pdg %d, phPDG %d pT %2.2f, pTtrig %2.2f, eta %2.2f, phi %2.2f\n",pdg,phPDG,pt,ptTrig,particle->Eta(), particle->Phi()*TMath::RadToDeg());
915 if( ptMaxNeutPhot < pt) ptMaxNeutPhot = pt;
917 if( radius < rThresIC )
919 if(pt > ptThresIC) nICNeutPhot++ ;
924 //Calorimeter acceptance
925 Bool_t inCalo = GetFiducialCut()->IsInFiducialCut(trigger,GetCalorimeter()) ;
926 if(!inCalo) continue;
928 if( ptMaxNeutEMCAL < pt ) ptMaxNeutEMCAL = pt;
929 if( radius < rThresIC )
931 if( pt > ptThresIC ) nICNeutEMCAL++ ;
937 if( ptMaxNeutEMCALPhot < pt ) ptMaxNeutEMCALPhot = pt;
938 if( radius < rThresIC )
940 if (pt > ptThresIC) nICNeutEMCALPhot++ ;
941 sumNePtEMCPhot += pt;
947 if( pt < fMinChargedPt) continue ;
949 if( ptMaxCharged < pt ) ptMaxCharged = pt;
951 if( radius < rThresIC )
953 //printf("UE track? pTtrig %2.2f, pt %2.2f, etaTrig %2.2f, eta %2.2f, phiTrig %2.2f, phi %2.2f, radius %2.2f\n",
954 // ptTrig, pt,etaTrig, eta, phiTrig*TMath::RadToDeg(), phi*TMath::RadToDeg(),radius);
955 if( pt > ptThresIC ) nICTrack++ ;
963 if(ptTrig > ptMaxCharged)
965 //printf("pt charged %2.2f, pt neutral %2.2f, pt neutral emcal %2.2f, pt photon %2.2f, pt photon emcal %2.2f\n",
966 // ptMaxCharged, ptMaxNeutral, ptMaxNeutEMCAL,ptMaxNeutPhot, ptMaxNeutEMCALPhot);
967 if(ptTrig > ptMaxNeutral ) leading[0] = kTRUE ;
968 if(ptTrig > ptMaxNeutEMCAL ) leading[1] = kTRUE ;
969 if(ptTrig > ptMaxNeutPhot ) leading[2] = kTRUE ;
970 if(ptTrig > ptMaxNeutEMCALPhot) leading[3] = kTRUE ;
973 //printf("N in cone over threshold: tracks %d, neutral %d, neutral emcal %d, photon %d, photon emcal %d\n",
974 // nICTrack, nICNeutral ,nICNeutEMCAL,nICNeutPhot, nICNeutEMCALPhot);
977 // Isolation decision
978 if( isoMethod == AliIsolationCut::kPtThresIC )
982 if(nICNeutral == 0 ) isolated[0] = kTRUE ;
983 if(nICNeutEMCAL == 0 ) isolated[1] = kTRUE ;
984 if(nICNeutPhot == 0 ) isolated[2] = kTRUE ;
985 if(nICNeutEMCALPhot == 0 ) isolated[3] = kTRUE ;
988 else if( isoMethod == AliIsolationCut::kSumPtIC )
990 if(sumChPt + sumNePt < sumThresIC ) isolated[0] = kTRUE ;
991 if(sumChPt + sumNePtEMC < sumThresIC ) isolated[1] = kTRUE ;
992 if(sumChPt + sumNePtPhot < sumThresIC ) isolated[2] = kTRUE ;
993 if(sumChPt + sumNePtEMCPhot < sumThresIC ) isolated[3] = kTRUE ;
997 //----------------------------------------------------
998 // Fill histograms if conditions apply for all 4 cases
999 for( Int_t i = 0; i < 4; i++ )
1005 fhPtPi0Leading[i]->Fill(ptTrig);
1007 if (i == 0) fhPtPi0LeadingSumPt[i]->Fill(ptTrig, sumChPt + sumNePt);
1008 else if(i == 1) fhPtPi0LeadingSumPt[i]->Fill(ptTrig, sumChPt + sumNePtEMC);
1009 else if(i == 2) fhPtPi0LeadingSumPt[i]->Fill(ptTrig, sumChPt + sumNePtPhot);
1010 else if(i == 3) fhPtPi0LeadingSumPt[i]->Fill(ptTrig, sumChPt + sumNePtEMCPhot);
1012 if(isolated[i]) fhPtPi0LeadingIsolated[i]->Fill(ptTrig);
1015 else if(pdgTrig==22)
1019 fhPtPhotonLeading[i]->Fill(ptTrig);
1021 if (i == 0) fhPtPhotonLeadingSumPt[i]->Fill(ptTrig, sumChPt + sumNePt);
1022 else if(i == 1) fhPtPhotonLeadingSumPt[i]->Fill(ptTrig, sumChPt + sumNePtEMC);
1023 else if(i == 2) fhPtPhotonLeadingSumPt[i]->Fill(ptTrig, sumChPt + sumNePtPhot);
1024 else if(i == 3) fhPtPhotonLeadingSumPt[i]->Fill(ptTrig, sumChPt + sumNePtEMCPhot);
1026 if(isolated[i]) fhPtPhotonLeadingIsolated[i]->Fill(ptTrig);
1029 } // conditions loop
1031 if(GetDebug() > 1) printf("AliAnaGeneratorKine::IsLeadingAndIsolated() - End \n");
1035 //_____________________________________________________
1036 void AliAnaGeneratorKine::MakeAnalysisFillHistograms()
1038 //Particle-Parton Correlation Analysis, fill histograms
1040 if(GetDebug() > 1) printf("AliAnaGeneratorKine::MakeAnalysisFillHistograms() - Start \n");
1042 TLorentzVector trigger;
1044 GetPartonsAndJets();
1046 for(Int_t ipr = 0; ipr < fStack->GetNprimary(); ipr ++ )
1048 TParticle * particle = fStack->Particle(ipr) ;
1050 Int_t pdgTrig = particle->GetPdgCode();
1051 Int_t statusTrig = particle->GetStatusCode();
1052 Int_t imother = particle->GetFirstMother();
1053 Float_t ptTrig = particle->Pt();
1055 // Select final state photons (prompt, fragmentation) or pi0s
1057 //Check the origin of the photon, accept if prompt or initial/final state radiation
1058 if(pdgTrig == 22 && statusTrig == 1)
1060 if(imother > 8) continue;
1062 // If not photon, trigger on pi0
1063 else if(pdgTrig != 111) continue;
1066 // Acceptance and kinematical cuts
1067 if( ptTrig < GetMinPt() ) continue ;
1069 // Recover the kinematics:
1070 particle->Momentum(trigger);
1072 Bool_t in = GetFiducialCutForTrigger()->IsInFiducialCut(trigger,fTriggerDetector) ;
1073 if(! in ) continue ;
1075 if( GetDebug() > 2) printf("Select trigger particle %d: pdg %d status %d, mother index %d, pT %2.2f, eta %2.2f, phi %2.2f \n",
1076 ipr, pdgTrig, statusTrig, imother, ptTrig, particle->Eta(), particle->Phi()*TMath::RadToDeg());
1080 // printf("\t pi0 daughters %d, %d\n", particle->GetDaughter(0), particle->GetDaughter(1));
1083 if (pdgTrig==22 ) fhPtPhoton->Fill(ptTrig);
1084 else if(pdgTrig==111) fhPtPi0 ->Fill(ptTrig);
1086 // Check if it is leading
1088 Bool_t isolated[4] ;
1090 IsLeadingAndIsolated(trigger, ipr, pdgTrig, leading, isolated);
1093 Int_t ok = CorrelateWithPartonOrJet(trigger, ipr, pdgTrig, leading, isolated, iparton);
1096 GetXE(trigger,ipr,pdgTrig,leading,isolated,iparton) ;
1100 if(GetDebug() > 1) printf("AliAnaGeneratorKine::MakeAnalysisFillHistograms() - End fill histograms \n");