]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGJE/AliPWG4HighPtSpectra.cxx
Adding cuts on vertex contributors (for pA and PileUp rejection) to AliFemto (ESD...
[u/mrichter/AliRoot.git] / PWGJE / AliPWG4HighPtSpectra.cxx
CommitLineData
fdceab34 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16//-----------------------------------------------------------------------
2b553e6f 17// Efficiency between different steps of the procedure.
fdceab34 18// The ouptut of the task is a AliCFContainer from which the efficiencies
19// can be calculated
20//-----------------------------------------------------------------------
21// Author : Marta Verweij - UU
22//-----------------------------------------------------------------------
23
24
67ebd013 25#ifndef ALIPWG4HIGHPTSPECTRA_CXX
26#define ALIPWG4HIGHPTSPECTRA_CXX
fdceab34 27
28#include "AliPWG4HighPtSpectra.h"
29
67ebd013 30#include "TVector3.h"
31#include <iostream>
32#include "TH1.h"
33#include "TH2.h"
34#include "TH3.h"
b4691ee7 35#include "TProfile.h"
67ebd013 36#include "TList.h"
37#include "TChain.h"
b4691ee7 38#include "TKey.h"
39#include "TSystem.h"
40#include "TFile.h"
67ebd013 41
42#include "AliAnalysisManager.h"
43#include "AliESDInputHandler.h"
9477b466 44#include "AliAODInputHandler.h"
67ebd013 45#include "AliESDtrack.h"
9477b466 46#include "AliAODTrack.h"
47#include "AliAODMCParticle.h"
67ebd013 48#include "AliESDtrackCuts.h"
49#include "AliExternalTrackParam.h"
2b553e6f 50#include "AliCentrality.h"
65e8ecdd 51
67ebd013 52#include "AliLog.h"
53
fdceab34 54#include "AliStack.h"
55#include "TParticle.h"
fdceab34 56#include "AliMCEvent.h"
57#include "AliMCEventHandler.h"
fdceab34 58#include "AliCFContainer.h"
b4691ee7 59#include "AliGenPythiaEventHeader.h"
03372fd1 60#include "AliGenHijingEventHeader.h"
b4691ee7 61#include "AliGenCocktailEventHeader.h"
3abf7501 62#include "AliAODMCHeader.h"
9477b466 63#include "AliESDVertex.h"
67ebd013 64//#include "$ALICE_ROOT/PWG4/JetTasks/AliAnalysisHelperJetTasks.h"
65
66//#include <iostream>
fdceab34 67using namespace std; //required for resolving the 'cout' symbol
68
69ClassImp(AliPWG4HighPtSpectra)
70
71//__________________________________________________________________________
72AliPWG4HighPtSpectra::AliPWG4HighPtSpectra() : AliAnalysisTask("AliPWG4HighPtSpectra", ""),
73 fReadAODData(0),
9477b466 74 fNoPythiaInfo(0),
67ebd013 75 fCFManagerPos(0x0),
76 fCFManagerNeg(0x0),
cce400da 77 fESD(0x0),
9477b466 78 fAOD(0x0),
cce400da 79 fMC(0x0),
80 fStack(0x0),
3abf7501 81 fArrayMCAOD(0x0),
34fc6450 82 fVtx(0x0),
5d87a047 83 fTriggerMask(AliVEvent::kMB),
2b553e6f 84 fIsPbPb(0),
85 fCentClass(10),
34fc6450 86 fTrackType(0),
cce400da 87 fTrackCuts(0x0),
327d12da 88 fTrackCutsReject(0x0),
9477b466 89 fFilterMask(0),
03372fd1 90 fbSelectHIJING(kFALSE),
dae7dd67 91 fSigmaConstrainedMax(100.),
b4691ee7 92 fAvgTrials(1),
fdceab34 93 fHistList(0),
b5cc0c6d 94 fNEventAll(0),
b4691ee7 95 fNEventSel(0),
cb76764e 96 fNEventReject(0),
2b553e6f 97 fh1Centrality(0x0),
b4691ee7 98 fh1Xsec(0),
99 fh1Trials(0),
100 fh1PtHard(0),
a337a5a9 101 fh1PtHardTrials(0),
102 fPtRelUncertainty1PtPrim(0x0),
103 fPtRelUncertainty1PtSec(0x0)
fdceab34 104{
105 //
106 //Default ctor
107 //
108}
9477b466 109
fdceab34 110//___________________________________________________________________________
111AliPWG4HighPtSpectra::AliPWG4HighPtSpectra(const Char_t* name) :
112 AliAnalysisTask(name,""),
113 fReadAODData(0),
9477b466 114 fNoPythiaInfo(0),
67ebd013 115 fCFManagerPos(0x0),
116 fCFManagerNeg(0x0),
cce400da 117 fESD(0x0),
9477b466 118 fAOD(0x0),
cce400da 119 fMC(0x0),
120 fStack(0x0),
3abf7501 121 fArrayMCAOD(0x0),
34fc6450 122 fVtx(0x0),
5d87a047 123 fTriggerMask(AliVEvent::kMB),
2b553e6f 124 fIsPbPb(0),
125 fCentClass(10),
34fc6450 126 fTrackType(0),
cce400da 127 fTrackCuts(0x0),
327d12da 128 fTrackCutsReject(0x0),
9477b466 129 fFilterMask(0),
03372fd1 130 fbSelectHIJING(kFALSE),
dae7dd67 131 fSigmaConstrainedMax(100.),
b4691ee7 132 fAvgTrials(1),
fdceab34 133 fHistList(0),
b5cc0c6d 134 fNEventAll(0),
b4691ee7 135 fNEventSel(0),
cb76764e 136 fNEventReject(0),
2b553e6f 137 fh1Centrality(0x0),
b4691ee7 138 fh1Xsec(0),
139 fh1Trials(0),
140 fh1PtHard(0),
a337a5a9 141 fh1PtHardTrials(0),
142 fPtRelUncertainty1PtPrim(0x0),
143 fPtRelUncertainty1PtSec(0x0)
fdceab34 144{
145 //
146 // Constructor. Initialization of Inputs and Outputs
147 //
f51451be 148 AliDebug(2,Form("AliPWG4HighPtSpectra Calling Constructor"));
fdceab34 149 // Input slot #0 works with a TChain ESD
150 DefineInput(0, TChain::Class());
e5abcde9 151 // Output slot #0 writes into a TList
fdceab34 152 DefineOutput(0,TList::Class());
e5abcde9 153 // Output slot #1, #2 writes into a AliCFContainer
fdceab34 154 DefineOutput(1,AliCFContainer::Class());
67ebd013 155 DefineOutput(2,AliCFContainer::Class());
e5abcde9 156 // Output slot #3 writes into a AliESDtrackCuts
157 DefineOutput(3, AliESDtrackCuts::Class());
1f329128 158 DefineOutput(4, AliESDtrackCuts::Class());
e5abcde9 159}
160
161//________________________________________________________________________
162void AliPWG4HighPtSpectra::LocalInit()
163{
164 //
165 // Only called once at beginning
166 //
167 PostData(3,fTrackCuts);
fdceab34 168}
169
fdceab34 170//________________________________________________________________________
171void AliPWG4HighPtSpectra::ConnectInputData(Option_t *)
172{
173 // Connect ESD here
174 // Called once
df943115 175 AliDebug(2,Form(">> AliPWG4HighPtSpectra::ConnectInputData \n"));
fdceab34 176
177 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
178 if (!tree) {
cce400da 179 AliDebug(2,Form( "ERROR: Could not read chain from input slot 0 \n"));
180 return;
fdceab34 181 }
cce400da 182
9477b466 183 if(fReadAODData) {
184 AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
185
186 if (!aodH) {
187 AliDebug(2,Form("ERROR: Could not get AODInputHandler"));
188 return;
189 } else
190 fAOD = aodH->GetEvent();
191 } else {
192 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
67ebd013 193
9477b466 194 if (!esdH) {
195 AliDebug(2,Form("ERROR: Could not get ESDInputHandler"));
196 return;
197 } else
198 fESD = esdH->GetEvent();
199 }
200
cce400da 201 AliMCEventHandler *eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
202 if (!eventHandler) {
203 AliDebug(2,Form( "ERROR: Could not retrieve MC event handler \n"));
204 }
205 else
206 fMC = eventHandler->MCEvent();
207
fdceab34 208}
cce400da 209
210//________________________________________________________________________
9477b466 211Bool_t AliPWG4HighPtSpectra::SelectEvent()
212{
fdceab34 213 //
cce400da 214 // Decide if event should be selected for analysis
fdceab34 215 //
fdceab34 216
cce400da 217 // Checks following requirements:
218 // - fESD available
219 // - trigger info from AliPhysicsSelection
220 // - number of reconstructed tracks > 1
221 // - primary vertex reconstructed
222 // - z-vertex < 10 cm
223
224 Bool_t selectEvent = kTRUE;
b5cc0c6d 225
9477b466 226 //fESD or fAOD object available?
227 if (!fESD && !fAOD) {
cce400da 228 AliDebug(2,Form("ERROR: fInputEvent not available\n"));
9477b466 229 fNEventReject->Fill("noESD/AOD",1);
cce400da 230 selectEvent = kFALSE;
231 return selectEvent;
fdceab34 232 }
233
cce400da 234 //Trigger
3abf7501 235 UInt_t isSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
5d87a047 236 if(fTriggerMask != AliVEvent::kAny && !(isSelected&fTriggerMask)) { //Select collison candidates
b5cc0c6d 237 AliDebug(2,Form(" Trigger Selection: event REJECTED ... "));
cb76764e 238 fNEventReject->Fill("Trigger",1);
cce400da 239 selectEvent = kFALSE;
240 return selectEvent;
67ebd013 241 }
b5cc0c6d 242
cce400da 243 //Check if number of reconstructed tracks is larger than 1
3abf7501 244 if(!fReadAODData) {
9477b466 245 if(!fESD->GetNumberOfTracks() || fESD->GetNumberOfTracks()<2) {
246 fNEventReject->Fill("NTracks<2",1);
247 selectEvent = kFALSE;
248 return selectEvent;
249 }
fdceab34 250 }
b4691ee7 251
cce400da 252 //Check if vertex is reconstructed
9477b466 253 if(fESD)
254 fVtx = fESD->GetPrimaryVertexSPD();
255 if(fAOD)
3abf7501 256 fVtx = fAOD->GetPrimaryVertex();
34fc6450 257
258 if(!fVtx) {
259 fNEventReject->Fill("noVTX",1);
260 selectEvent = kFALSE;
261 return selectEvent;
262 }
263
3abf7501 264 if(!fReadAODData){
9477b466 265 const AliESDVertex* esdVtx = fESD->GetPrimaryVertexSPD();
266 if(!esdVtx->GetStatus()) {
267 fNEventReject->Fill("VtxStatus",1);
268 selectEvent = kFALSE;
269 return selectEvent;
270 }
34fc6450 271 }
272
fdceab34 273 // Need vertex cut
34fc6450 274 // TString vtxName(fVtx->GetName());
275 if(fVtx->GetNContributors()<2) {
276 fNEventReject->Fill("NCont<2",1);
277 selectEvent = kFALSE;
278 return selectEvent;
cd9a6fa2 279 }
cce400da 280
281 //Check if z-vertex < 10 cm
b5cc0c6d 282 double primVtx[3];
34fc6450 283 fVtx->GetXYZ(primVtx);
67ebd013 284 if(TMath::Sqrt(primVtx[0]*primVtx[0] + primVtx[1]*primVtx[1])>1. || TMath::Abs(primVtx[2]>10.)){
cb76764e 285 fNEventReject->Fill("ZVTX>10",1);
cce400da 286 selectEvent = kFALSE;
287 return selectEvent;
67ebd013 288 }
cce400da 289
2b553e6f 290 //Centrality selection should only be done in case of PbPb
291 if(IsPbPb()) {
292 Float_t cent = 0.;
9477b466 293 Int_t calccent = 0;
294 if(fReadAODData)
295 calccent = CalculateCentrality(fAOD);
296 else
297 calccent = CalculateCentrality(fESD);
298 if(fCentClass!=calccent && fCentClass!=10) {
2b553e6f 299 fNEventReject->Fill("cent",1);
300 selectEvent = kFALSE;
301 return selectEvent;
302 }
303 else {
9477b466 304 if(fReadAODData) {
305 if(dynamic_cast<AliAODEvent*>(fAOD)->GetCentrality()) {
306 cent = dynamic_cast<AliAODEvent*>(fAOD)->GetCentrality()->GetCentralityPercentile("V0M");
307 }
2b553e6f 308 }
9477b466 309 else {
310 if(dynamic_cast<AliESDEvent*>(fESD)->GetCentrality()) {
311 cent = dynamic_cast<AliESDEvent*>(fESD)->GetCentrality()->GetCentralityPercentile("V0M");
312 }
313 }
314
2b553e6f 315 if(cent>90.) {
316 fNEventReject->Fill("cent>90",1);
317 selectEvent = kFALSE;
9477b466 318 return selectEvent;
2b553e6f 319 }
320 fh1Centrality->Fill(cent);
321 }
322 }
323
cce400da 324 return selectEvent;
325
326}
327
2b553e6f 328//________________________________________________________________________
9477b466 329Int_t AliPWG4HighPtSpectra::CalculateCentrality(AliVEvent *event)
330{
2b553e6f 331 Float_t cent = 999;
332
9477b466 333 if(event){
334 if(event->GetCentrality()){
335 cent = event->GetCentrality()->GetCentralityPercentile("V0M");
2b553e6f 336 }
337 }
338
9477b466 339 if(cent<0) return 5;
2b553e6f 340 if(cent>80)return 4;
341 if(cent>50)return 3;
342 if(cent>30)return 2;
343 if(cent>10)return 1;
344 return 0;
2b553e6f 345}
346
cce400da 347//_________________________________________________
348void AliPWG4HighPtSpectra::Exec(Option_t *)
349{
350 //
351 // Main loop function
352 //
353 AliDebug(2,Form(">> AliPWG4HighPtSpectra::Exec \n"));
354
9477b466 355 AliVEvent* event;
356 if(fReadAODData)
357 event = fAOD;
358 else {
359 event = fESD;
360 if(!fMC) {
361 AliMCEventHandler *eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
362 if (!eventHandler) {
363 AliDebug(2,Form( "ERROR: Could not retrieve MC event handler \n"));
364 }
365 else
366 fMC = eventHandler->MCEvent();
0714a353 367 }
0714a353 368 }
369
cce400da 370 // All events without selection
371 fNEventAll->Fill(0.);
372
373 if(!SelectEvent()) {
67ebd013 374 // Post output data
375 PostData(0,fHistList);
376 PostData(1,fCFManagerPos->GetParticleContainer());
377 PostData(2,fCFManagerNeg->GetParticleContainer());
cd9a6fa2 378 return;
379 }
cce400da 380
9477b466 381 if(fReadAODData){
382 //Open the MC particles in the case of AODanalysis
3abf7501 383 fArrayMCAOD = (TClonesArray*) fAOD->GetList()->FindObject(AliAODMCParticle::StdBranchName());
9477b466 384 }
385 else {
386 //In case of ESD event: MCEvent available? if yes: get MC particles stack
387 if(fMC) {
388 AliDebug(2,Form("MC particles: %d", fMC->GetNumberOfTracks()));
389 fStack = fMC->Stack(); //Particles Stack
390 AliDebug(2,Form("MC particles stack: %d", fStack->GetNtrack()));
391 }
cce400da 392 }
393
9477b466 394 Int_t nTracks = event->GetNumberOfTracks();
fdceab34 395 AliDebug(2,Form("nTracks %d", nTracks));
fdceab34 396
9477b466 397 if((!fTrackCuts && !fReadAODData) || (fFilterMask==0 && fReadAODData)) { //if the TrackCuts are missing in ESD analysis or the FilterBit is missing in AOD analysis
cb76764e 398 fNEventReject->Fill("noTrackCuts",1);
67ebd013 399 // Post output data
400 PostData(0,fHistList);
401 PostData(1,fCFManagerPos->GetParticleContainer());
402 PostData(2,fCFManagerNeg->GetParticleContainer());
403 return;
404 }
405
b5cc0c6d 406 // Selected events for analysis
407 fNEventSel->Fill(0.);
67ebd013 408
dae7dd67 409 const Int_t nvar = 4;
67ebd013 410
dae7dd67 411 Double_t containerInputRec[nvar] = {0.,0.,0.,0.};
412 Double_t containerInputMC[nvar] = {0.,0.,0.,0.};
413 Double_t containerInputRecMC[nvar] = {0.,0.,0.,0.}; //reconstructed yield as function of MC variable
ec555a3c 414
fdceab34 415 //Now go to rec level
9477b466 416 if(fReadAODData) {//AOD analysis
9477b466 417 for (Int_t iTrack = 0; iTrack<nTracks; iTrack++)
418 {
419 //Get track for analysis
420 AliVTrack *track = 0x0;
421
422 AliAODTrack *aodtrack = fAOD->GetTrack(iTrack);
423 if(!aodtrack)
424 continue;
3abf7501 425 if( !aodtrack->TestFilterBit(fFilterMask) )
9477b466 426 continue;
427 else {
3abf7501 428 track = static_cast<AliVTrack*>(aodtrack);
9477b466 429 }
430 //fill the container
431 containerInputRec[0] = track->Pt();
432 containerInputRec[1] = track->Phi();
433 containerInputRec[2] = track->Eta();
3abf7501 434 containerInputRec[3] = track->GetTPCNcls();
9477b466 435
436 if(track->Charge()>0.) {
3abf7501 437 fCFManagerPos->GetParticleContainer()->Fill(containerInputRec,kStepReconstructed);
9477b466 438 }
439 if(track->Charge()<0.) {
440 fCFManagerNeg->GetParticleContainer()->Fill(containerInputRec,kStepReconstructed);
9477b466 441 }
442
3abf7501 443 if(fArrayMCAOD) {
9477b466 444 Int_t label = TMath::Abs(track->GetLabel());
3abf7501 445 if(label>fArrayMCAOD->GetEntries()) {
9477b466 446 continue;
447 }
3abf7501 448 AliAODMCParticle *particle = (AliAODMCParticle*) fArrayMCAOD->At(label);
9477b466 449 if(!particle) {
9477b466 450 continue;
451 }
452 //Only select particles generated by HIJING if requested
453 if(fbSelectHIJING) {
454 if(!IsHIJINGParticle(label)) {
9477b466 455 continue;
456 }
457 }
458 containerInputRecMC[0] = particle->Pt();
459 containerInputRecMC[1] = particle->Phi();
460 containerInputRecMC[2] = particle->Eta();
461 containerInputRecMC[3] = track->GetTPCNcls();
462
463 //Container with primaries
464 if(particle->IsPhysicalPrimary()) {
9477b466 465 if(particle->Charge()>0.) {
466 fCFManagerPos->GetParticleContainer()->Fill(containerInputRecMC,kStepReconstructedMC);
9477b466 467 }
468 if(particle->Charge()<0.) {
469 fCFManagerNeg->GetParticleContainer()->Fill(containerInputRecMC,kStepReconstructedMC);
9477b466 470 }
471 //Fill pT resolution plots for primaries
3abf7501 472 //fPtRelUncertainty1PtPrim->Fill(containerInputRec[0],containerInputRec[0]*TMath::Sqrt(track->GetSigma1Pt2())); //This has not been implemented in AOD analysis, since they are also produced by the AddTaskPWG4HighPtTrackQA.C macro
9477b466 473 }
474
475 //Container with secondaries
476 if (!particle->IsPhysicalPrimary() ) {
477 if(particle->Charge()>0.) {
478 fCFManagerPos->GetParticleContainer()->Fill(containerInputRecMC,kStepSecondaries);
9477b466 479 }
480 if(particle->Charge()<0.) {
481 fCFManagerNeg->GetParticleContainer()->Fill(containerInputRecMC,kStepSecondaries);
9477b466 482 }
483 //Fill pT resolution plots for primaries
3abf7501 484 //fPtRelUncertainty1PtSec->Fill(containerInputRec[0],containerInputRec[0]*TMath::Sqrt(track->GetSigma1Pt2())); //This has not been implemented in AOD analysis, since they are also produced by the AddTaskPWG4HighPtTrackQA.C macro
9477b466 485 }
486 }
487 }//track loop
488
489 //Fill MC containers if particles are findable
3abf7501 490 if(fArrayMCAOD) {
491 int noPart = fArrayMCAOD->GetEntriesFast();
492 for(int iPart = 1; iPart<noPart; iPart++) {
493 AliAODMCParticle *mcPart = (AliAODMCParticle*) fArrayMCAOD->At(iPart);
9477b466 494 if(!mcPart) continue;
495
496 //Only select particles generated by HIJING if requested
497 if(fbSelectHIJING) {
498 if(!IsHIJINGParticle(iPart))
499 continue;
500 }
501
502 Int_t pdg = mcPart->PdgCode();
503
504 // select charged pions, protons, kaons , electrons, muons
505 if(TMath::Abs(pdg) == 211 || TMath::Abs(pdg) == 2212 || TMath::Abs(pdg) ==
506 321 || TMath::Abs(pdg) == 11 || TMath::Abs(pdg) == 13){
507
508 //fill the container
509 containerInputMC[0] = mcPart->Pt();
510 containerInputMC[1] = mcPart->Phi();
511 containerInputMC[2] = mcPart->Eta();
512 // AliESDtrack *esdtrack = fESD->GetTrack(mcPart->GetLabel());
513 containerInputMC[3] = 159.;
514
515 if(mcPart->IsPhysicalPrimary()) {
516 if(mcPart->Charge()>0. && fCFManagerPos->CheckParticleCuts(kStepMCAcceptance,mcPart)) {
517 fCFManagerPos->GetParticleContainer()->Fill(containerInputMC,kStepMCAcceptance);
9477b466 518 }
519 if(mcPart->Charge()<0. && fCFManagerNeg->CheckParticleCuts(kStepMCAcceptance,mcPart)) {
520 fCFManagerNeg->GetParticleContainer()->Fill(containerInputMC,kStepMCAcceptance);
9477b466 521 }
522 }
523 }
524 }
525 }
526
527 }//end of AOD analysis
528 else {//ESD analysis
529 for (Int_t iTrack = 0; iTrack<nTracks; iTrack++)
530 {
2b553e6f 531 //Get track for analysis
532 AliESDtrack *track = 0x0;
533 AliESDtrack *esdtrack = fESD->GetTrack(iTrack);
2461511f 534 if(!esdtrack)
535 continue;
536
2b553e6f 537
2461511f 538 if(fTrackType==4) {
539 if (!(fTrackCuts->AcceptTrack(esdtrack)))
540 continue;
8da84c2c 541 }
2b553e6f 542
2461511f 543 if(fTrackType==1)
544 track = AliESDtrackCuts::GetTPCOnlyTrack(fESD,esdtrack->GetID());
545 else if(fTrackType==2 || fTrackType==4) {
546 track = AliESDtrackCuts::GetTPCOnlyTrack(const_cast<AliESDEvent*>(fESD),esdtrack->GetID());
547 if(!track)
548 continue;
549
2b553e6f 550 AliExternalTrackParam exParam;
9477b466 551 Bool_t relate = track->RelateToVertexTPC(fESD->GetPrimaryVertexSPD(),fESD->GetMagneticField(),kVeryBig,&exParam);
2b553e6f 552 if( !relate ) {
2461511f 553 if(track) delete track;
34fc6450 554 continue;
555 }
2b553e6f 556 track->Set(exParam.GetX(),exParam.GetAlpha(),exParam.GetParameter(),exParam.GetCovariance());
bcd60a94 557 }
2461511f 558 else if(fTrackType==5 || fTrackType==6) {
559 if(fTrackCuts->AcceptTrack(esdtrack)) {
560 continue;
561 }
562 else {
563 if( !(fTrackCutsReject->AcceptTrack(esdtrack)) && fTrackCuts->AcceptTrack(esdtrack) ) {
564
565 if(fTrackType==5) {
566 //use TPConly constrained track
567 track = AliESDtrackCuts::GetTPCOnlyTrack(fESD,esdtrack->GetID());
568 if(!track)
569 continue;
570
571 AliExternalTrackParam exParam;
9477b466 572 Bool_t relate = track->RelateToVertexTPC(fESD->GetPrimaryVertexSPD(),fESD->GetMagneticField(),kVeryBig,&exParam);
2461511f 573 if( !relate ) {
574 if(track) delete track;
575 continue;
576 }
577 track->Set(exParam.GetX(),exParam.GetAlpha(),exParam.GetParameter(),exParam.GetCovariance());
578 }
579 else if(fTrackType==6) {
580 //use global constrained track
581 track = new AliESDtrack(*esdtrack);
582 track->Set(esdtrack->GetConstrainedParam()->GetX(),esdtrack->GetConstrainedParam()->GetAlpha(),esdtrack->GetConstrainedParam()->GetParameter(),esdtrack->GetConstrainedParam()->GetCovariance());
2461511f 583 }
584 }
585 }
586 }
d3a3f33d 587 else if(fTrackType==7) {
588 //use global constrained track
42881dab 589 track = new AliESDtrack(*esdtrack);
d3a3f33d 590 }
bcd60a94 591 else
2b553e6f 592 track = esdtrack;
bcd60a94 593
2461511f 594 if(!track)
595 continue;
596
597
598 if(fTrackType==2 || fTrackType==4 || fTrackType==5) {
2b553e6f 599 //Cut on chi2 of constrained fit
2461511f 600 if(track->GetConstrainedChi2TPC() > fSigmaConstrainedMax*fSigmaConstrainedMax && fSigmaConstrainedMax>0.) {
601 if(track) delete track;
2b553e6f 602 continue;
603 }
604 }
2461511f 605 if (!(fTrackCuts->AcceptTrack(track)) && fTrackType!=4 && fTrackType!=5 && fTrackType!=6) {
606 if(fTrackType==1 || fTrackType==2 || fTrackType==7) {
607 if(track) delete track;
608 }
609 continue;
610 }
fdceab34 611
2461511f 612 if(fTrackType==7) {
613 if(fTrackCutsReject ) {
614 if(fTrackCutsReject->AcceptTrack(track) ) {
615 if(track) delete track;
616 continue;
327d12da 617 }
031aac8e 618 }
2461511f 619 if(esdtrack->GetConstrainedParam())
620 track->Set(esdtrack->GetConstrainedParam()->GetX(),esdtrack->GetConstrainedParam()->GetAlpha(),esdtrack->GetConstrainedParam()->GetParameter(),esdtrack->GetConstrainedParam()->GetCovariance());
621 }
031aac8e 622
2461511f 623 if(!track) {
624 if(fTrackType==1 || fTrackType==2 || fTrackType==4 || fTrackType==5 || fTrackType==6 || fTrackType==7) {
625 if(track) delete track;
626 }
627 continue;
628 }
629
630 //fill the container
631 containerInputRec[0] = track->Pt();
632 containerInputRec[1] = track->Phi();
633 containerInputRec[2] = track->Eta();
634 containerInputRec[3] = track->GetTPCNclsIter1();
031aac8e 635
2461511f 636 if(track->GetSign()>0.) fCFManagerPos->GetParticleContainer()->Fill(containerInputRec,kStepReconstructed);
637 if(track->GetSign()<0.) fCFManagerNeg->GetParticleContainer()->Fill(containerInputRec,kStepReconstructed);
1f329128 638
67ebd013 639
2461511f 640 //Only fill the MC containers if MC information is available
641 if(fMC) {
642 Int_t label = TMath::Abs(track->GetLabel());
643 if(label>fStack->GetNtrack()) {
644 if(fTrackType==1 || fTrackType==2 || fTrackType==7)
645 delete track;
646 continue;
647 }
648 //Only select particles generated by HIJING if requested
649 if(fbSelectHIJING) {
650 if(!IsHIJINGParticle(label)) {
5d87a047 651 if(fTrackType==1 || fTrackType==2 || fTrackType==7)
652 delete track;
653 continue;
654 }
2461511f 655 }
656 TParticle *particle = fStack->Particle(label) ;
657 if(!particle) {
658 if(fTrackType==1 || fTrackType==2 || fTrackType==7)
659 delete track;
660 continue;
661 }
662 containerInputRecMC[0] = particle->Pt();
663 containerInputRecMC[1] = particle->Phi();
664 containerInputRecMC[2] = particle->Eta();
665 containerInputRecMC[3] = track->GetTPCNclsIter1();
666
667 //Container with primaries
668 if(fStack->IsPhysicalPrimary(label)) {
669 if(particle->GetPDG()->Charge()>0.) {
670 fCFManagerPos->GetParticleContainer()->Fill(containerInputRecMC,kStepReconstructedMC);
03372fd1 671 }
2461511f 672 if(particle->GetPDG()->Charge()<0.) {
673 fCFManagerNeg->GetParticleContainer()->Fill(containerInputRecMC,kStepReconstructedMC);
bcd60a94 674 }
2461511f 675 //Fill pT resolution plots for primaries
676 fPtRelUncertainty1PtPrim->Fill(containerInputRec[0],containerInputRec[0]*TMath::Sqrt(track->GetSigma1Pt2()));
2461511f 677 }
67ebd013 678
2461511f 679 //Container with secondaries
680 if (!fStack->IsPhysicalPrimary(label) ) {
681 if(particle->GetPDG()->Charge()>0.) {
682 fCFManagerPos->GetParticleContainer()->Fill(containerInputMC,kStepSecondaries);
683 }
684 if(particle->GetPDG()->Charge()<0.) {
685 fCFManagerNeg->GetParticleContainer()->Fill(containerInputMC,kStepSecondaries);
0cc50ca6 686 }
2461511f 687 //Fill pT resolution plots for primaries
688 fPtRelUncertainty1PtSec->Fill(containerInputRec[0],containerInputRec[0]*TMath::Sqrt(track->GetSigma1Pt2()));
0cc50ca6 689 }
2461511f 690 }
60829498 691
2461511f 692 if(fTrackType==1 || fTrackType==2 || fTrackType==4 || fTrackType==5 || fTrackType==6 || fTrackType==7) {
42881dab 693 if(track) delete track;
694 }
1f329128 695 }//track loop
03372fd1 696
9477b466 697 //Fill MC containers if particles are findable
698 if(fMC) {
699 for(int iPart = 1; iPart<(fMC->GetNumberOfPrimaries()); iPart++) {
700 AliMCParticle *mcPart = (AliMCParticle*)fMC->GetTrack(iPart);
701 if(!mcPart) continue;
702
703 //Only select particles generated by HIJING if requested
704 if(fbSelectHIJING) {
705 if(!IsHIJINGParticle(iPart))
706 continue;
707 }
708
709 Int_t pdg = mcPart->PdgCode();
710
711 // select charged pions, protons, kaons , electrons, muons
712 if(TMath::Abs(pdg) == 211 || TMath::Abs(pdg) == 2212 || TMath::Abs(pdg) ==
713 321 || TMath::Abs(pdg) == 11 || TMath::Abs(pdg) == 13){
714
715 //fill the container
716 containerInputMC[0] = mcPart->Pt();
717 containerInputMC[1] = mcPart->Phi();
718 containerInputMC[2] = mcPart->Eta();
719 // AliESDtrack *esdtrack = fESD->GetTrack(mcPart->GetLabel());
720 containerInputMC[3] = 159.;
721
722 if(fStack->IsPhysicalPrimary(iPart)) {
723 if(mcPart->Charge()>0. && fCFManagerPos->CheckParticleCuts(kStepMCAcceptance,mcPart)) fCFManagerPos->GetParticleContainer()->Fill(containerInputMC,kStepMCAcceptance);
724 if(mcPart->Charge()<0. && fCFManagerNeg->CheckParticleCuts(kStepMCAcceptance,mcPart)) fCFManagerNeg->GetParticleContainer()->Fill(containerInputMC,kStepMCAcceptance);
725 }
726 }
0cc50ca6 727 }
cce400da 728 }
9477b466 729 }//end of ESD analysis
730
67ebd013 731 PostData(0,fHistList);
732 PostData(1,fCFManagerPos->GetParticleContainer());
733 PostData(2,fCFManagerNeg->GetParticleContainer());
734
fdceab34 735}
9477b466 736
b4691ee7 737//________________________________________________________________________
9477b466 738Bool_t AliPWG4HighPtSpectra::PythiaInfoFromFile(const char* currFile,Float_t &fXsec,Float_t &fTrials)
739{
b4691ee7 740 //
741 // get the cross section and the trails either from pyxsec.root or from pysec_hists.root
742 // This is to called in Notify and should provide the path to the AOD/ESD file
743 // Copied from AliAnalysisTaskJetSpectrum2
744 //
745
746 TString file(currFile);
747 fXsec = 0;
748 fTrials = 1;
749
750 if(file.Contains("root_archive.zip#")){
751 Ssiz_t pos1 = file.Index("root_archive",12,TString::kExact);
752 Ssiz_t pos = file.Index("#",1,pos1,TString::kExact);
753 file.Replace(pos+1,20,"");
754 }
755 else {
756 // not an archive take the basename....
757 file.ReplaceAll(gSystem->BaseName(file.Data()),"");
758 }
cb76764e 759 // Printf("%s",file.Data());
b4691ee7 760
761 TFile *fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec.root")); // problem that we cannot really test the existance of a file in a archive so we have to lvie with open error message from root
762 if(!fxsec){
763 // next trial fetch the histgram file
764 fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
765 if(!fxsec){
2461511f 766 // not a severe condition but indicates that we have no information
b4691ee7 767 return kFALSE;
768 }
769 else{
770 // find the tlist we want to be independtent of the name so use the Tkey
771 TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0);
772 if(!key){
773 fxsec->Close();
774 return kFALSE;
775 }
776 TList *list = dynamic_cast<TList*>(key->ReadObj());
777 if(!list){
778 fxsec->Close();
779 return kFALSE;
780 }
781 fXsec = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
782 fTrials = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
783 fxsec->Close();
784 }
785 } // no tree pyxsec.root
786 else {
787 TTree *xtree = (TTree*)fxsec->Get("Xsection");
788 if(!xtree){
789 fxsec->Close();
790 return kFALSE;
791 }
792 UInt_t ntrials = 0;
793 Double_t xsection = 0;
794 xtree->SetBranchAddress("xsection",&xsection);
795 xtree->SetBranchAddress("ntrials",&ntrials);
796 xtree->GetEntry(0);
797 fTrials = ntrials;
798 fXsec = xsection;
799 fxsec->Close();
800 }
801 return kTRUE;
802}
9477b466 803
b4691ee7 804//________________________________________________________________________
805Bool_t AliPWG4HighPtSpectra::Notify()
806{
807 //
808 // Implemented Notify() to read the cross sections
809 // and number of trials from pyxsec.root
810 // Copied from AliAnalysisTaskJetSpectrum2
811 //
812
9477b466 813 if(fNoPythiaInfo)
814 return kTRUE;
815
b4691ee7 816 TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
817 Float_t xsection = 0;
818 Float_t ftrials = 1;
819
b4691ee7 820 if(tree){
821 TFile *curfile = tree->GetCurrentFile();
822 if (!curfile) {
823 Error("Notify","No current file");
824 return kFALSE;
825 }
826 if(!fh1Xsec||!fh1Trials){
cb76764e 827 // Printf("%s%d No Histogram fh1Xsec",(char*)__FILE__,__LINE__);
b4691ee7 828 return kFALSE;
829 }
adf3920d 830 PythiaInfoFromFile(curfile->GetName(),xsection,ftrials);
b4691ee7 831 fh1Xsec->Fill("<#sigma>",xsection);
adf3920d 832 fh1Trials->Fill("#sum{ntrials}",ftrials);
833 }
b4691ee7 834 return kTRUE;
835}
836
837//________________________________________________________________________
9477b466 838AliGenPythiaEventHeader* AliPWG4HighPtSpectra::GetPythiaEventHeader(AliMCEvent *mcEvent)
839{
b4691ee7 840
841 if(!mcEvent)return 0;
842 AliGenEventHeader* genHeader = mcEvent->GenEventHeader();
843 AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
844 if(!pythiaGenHeader){
845 // cocktail ??
846 AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(genHeader);
847
848 if (!genCocktailHeader) {
cb76764e 849 // AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Unknown header type (not Pythia or Cocktail)");
b4691ee7 850 // AliWarning(Form("%s %d: Unknown header type (not Pythia or Cocktail)",(char*)__FILE__,__LINE__));
851 return 0;
852 }
853 TList* headerList = genCocktailHeader->GetHeaders();
854 for (Int_t i=0; i<headerList->GetEntries(); i++) {
855 pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(headerList->At(i));
856 if (pythiaGenHeader)
857 break;
858 }
859 if(!pythiaGenHeader){
860 AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Pythia event header not found");
861 return 0;
862 }
863 }
864 return pythiaGenHeader;
865
866}
fdceab34 867
03372fd1 868//________________________________________________________________________
9477b466 869AliGenHijingEventHeader* AliPWG4HighPtSpectra::GetHijingEventHeader(AliMCEvent *mcEvent)
870{
03372fd1 871
872 if(!mcEvent)return 0;
873 AliGenEventHeader* genHeader = mcEvent->GenEventHeader();
874 AliGenHijingEventHeader* hijingGenHeader = dynamic_cast<AliGenHijingEventHeader*>(genHeader);
875 if(!hijingGenHeader){
876 // cocktail ??
877 AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(genHeader);
878
879 if (!genCocktailHeader) {
3abf7501 880 AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Unknown header type (not Pythia or Cocktail)");
03372fd1 881 // AliWarning(Form("%s %d: Unknown header type (not Pythia or Cocktail)",(char*)__FILE__,__LINE__));
882 return 0;
883 }
884 TList* headerList = genCocktailHeader->GetHeaders();
885 for (Int_t i=0; i<headerList->GetEntries(); i++) {
886 hijingGenHeader = dynamic_cast<AliGenHijingEventHeader*>(headerList->At(i));
887 if (hijingGenHeader)
888 break;
889 }
890 if(!hijingGenHeader){
891 AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Hijing event header not found");
892 return 0;
893 }
894 }
895 return hijingGenHeader;
896
897}
898
3abf7501 899//________________________________________________________________________
900AliGenHijingEventHeader* AliPWG4HighPtSpectra::GetHijingEventHeaderAOD(AliAODEvent *aodEvent)
901{
902 AliGenHijingEventHeader* hijingGenHeader = 0x0;
903 if(aodEvent) {
904 AliAODMCHeader* mcHeader = (AliAODMCHeader*) aodEvent->GetList()->FindObject(AliAODMCHeader::StdBranchName());
905 hijingGenHeader = (AliGenHijingEventHeader*) mcHeader->GetCocktailHeaders()->FindObject("Hijing");
906 }
907 return hijingGenHeader;
908}
909
fdceab34 910//___________________________________________________________________________
911void AliPWG4HighPtSpectra::Terminate(Option_t*)
912{
913 // The Terminate() function is the last function to be called during
914 // a query. It always runs on the client, it can be used to present
915 // the results graphically or save the results to file.
9b58297c 916
fdceab34 917}
918
fdceab34 919//___________________________________________________________________________
9477b466 920void AliPWG4HighPtSpectra::CreateOutputObjects()
921{
fdceab34 922 //HERE ONE CAN CREATE OUTPUT OBJECTS, IN PARTICULAR IF THE OBJECT PARAMETERS DON'T NEED
923 //TO BE SET BEFORE THE EXECUTION OF THE TASK
924 //
f51451be 925 AliDebug(2,Form("CreateOutputObjects CreateOutputObjects of task %s", GetName()));
df943115 926
67ebd013 927 Bool_t oldStatus = TH1::AddDirectoryStatus();
928 TH1::AddDirectory(kFALSE);
929
fdceab34 930 //slot #1
b5cc0c6d 931 OpenFile(0);
932 fHistList = new TList();
b4691ee7 933 fHistList->SetOwner(kTRUE);
b5cc0c6d 934 fNEventAll = new TH1F("fNEventAll","NEventAll",1,-0.5,0.5);
935 fHistList->Add(fNEventAll);
936 fNEventSel = new TH1F("fNEventSel","NEvent Selected for analysis",1,-0.5,0.5);
937 fHistList->Add(fNEventSel);
fdceab34 938
cb76764e 939 fNEventReject = new TH1F("fNEventReject","Reason events are rejectected for analysis",20,0,20);
2b553e6f 940 //Set labels
9477b466 941 fNEventReject->Fill("noESD/AOD",0);
2b553e6f 942 fNEventReject->Fill("Trigger",0);
943 fNEventReject->Fill("NTracks<2",0);
944 fNEventReject->Fill("noVTX",0);
945 fNEventReject->Fill("VtxStatus",0);
946 fNEventReject->Fill("NCont<2",0);
947 fNEventReject->Fill("ZVTX>10",0);
948 fNEventReject->Fill("cent",0);
949 fNEventReject->Fill("cent>90",0);
cb76764e 950 fHistList->Add(fNEventReject);
951
2b553e6f 952 fh1Centrality = new TH1F("fh1Centrality","fh1Centrality; Centrality %",100,0,100);
953 fHistList->Add(fh1Centrality);
954
b4691ee7 955 fh1Xsec = new TProfile("fh1Xsec","xsec from pyxsec.root",1,0,1);
956 fh1Xsec->GetXaxis()->SetBinLabel(1,"<#sigma>");
957 fHistList->Add(fh1Xsec);
958
959 fh1Trials = new TH1F("fh1Trials","trials root file",1,0,1);
960 fh1Trials->GetXaxis()->SetBinLabel(1,"#sum{ntrials}");
961 fHistList->Add(fh1Trials);
962
963 fh1PtHard = new TH1F("fh1PtHard","PYTHIA Pt hard;p_{T,hard}",350,-.5,349.5);
964 fHistList->Add(fh1PtHard);
965 fh1PtHardTrials = new TH1F("fh1PtHardTrials","PYTHIA Pt hard weight with trials;p_{T,hard}",350,-.5,349.5);
966 fHistList->Add(fh1PtHardTrials);
967
a337a5a9 968 Int_t fgkNPtBins = 100;
969 Float_t kMinPt = 0.;
970 Float_t kMaxPt = 100.;
971 Double_t *binsPt = new Double_t[fgkNPtBins+1];
972 for(Int_t i=0; i<=fgkNPtBins; i++) binsPt[i]=(Double_t)kMinPt + (kMaxPt-kMinPt)/fgkNPtBins*(Double_t)i ;
973
974 Int_t fgkNRel1PtUncertaintyBins=50;
975 Float_t fgkRel1PtUncertaintyMin = 0.;
976 Float_t fgkRel1PtUncertaintyMax = 1.;
977
978 Double_t *binsRel1PtUncertainty=new Double_t[fgkNRel1PtUncertaintyBins+1];
979 for(Int_t i=0; i<=fgkNRel1PtUncertaintyBins; i++) binsRel1PtUncertainty[i]=(Double_t)fgkRel1PtUncertaintyMin + (fgkRel1PtUncertaintyMax-fgkRel1PtUncertaintyMin)/fgkNRel1PtUncertaintyBins*(Double_t)i ;
980
981 fPtRelUncertainty1PtPrim = new TH2F("fPtRelUncertainty1PtPrim","fPtRelUncertainty1PtPrim",fgkNPtBins,binsPt,fgkNRel1PtUncertaintyBins,binsRel1PtUncertainty);
982 fHistList->Add(fPtRelUncertainty1PtPrim);
983
984 fPtRelUncertainty1PtSec = new TH2F("fPtRelUncertainty1PtSec","fPtRelUncertainty1PtSec",fgkNPtBins,binsPt,fgkNRel1PtUncertaintyBins,binsRel1PtUncertainty);
985 fHistList->Add(fPtRelUncertainty1PtSec);
986
67ebd013 987 TH1::AddDirectory(oldStatus);
988
f6d670e1 989 OpenFile(1);
990 OpenFile(2);
991
11245558 992 PostData(0,fHistList);
993 PostData(1,fCFManagerPos->GetParticleContainer());
994 PostData(2,fCFManagerNeg->GetParticleContainer());
995
a337a5a9 996 if(binsPt) delete [] binsPt;
997 if(binsRel1PtUncertainty) delete [] binsRel1PtUncertainty;
998
fdceab34 999}
1000
03372fd1 1001//________________________________________________________________________
3abf7501 1002Bool_t AliPWG4HighPtSpectra::IsHIJINGParticle(Int_t label)
9477b466 1003{
03372fd1 1004 //
1005 // Return kTRUE in case particle is from HIJING event
1006 //
1007
3abf7501 1008 AliGenHijingEventHeader* hijingHeader;
1009 if(fReadAODData)
1010 hijingHeader = GetHijingEventHeaderAOD(fAOD);
1011 else
1012 hijingHeader = GetHijingEventHeader(fMC);
03372fd1 1013
1014 Int_t nproduced = hijingHeader->NProduced();
1015
3abf7501 1016 if(fReadAODData) {
1017 AliAODMCParticle* mom = (AliAODMCParticle*) fArrayMCAOD->At(label);
1018 Int_t iMom = label;
1019 Int_t iParent = mom->GetMother();
1020 while(iParent!=-1){
1021 iMom = iParent;
1022 mom = (AliAODMCParticle*) fArrayMCAOD->At(iMom);
1023 iParent = mom->GetMother();
1024 }
03372fd1 1025
3abf7501 1026 if(iMom<nproduced) {
1027 return kTRUE;
1028 } else {
1029 return kFALSE;
1030 }
1031 }
1032 else { //ESD analysis
1033 TParticle * mom = fStack->Particle(label);
1034 Int_t iMom = label;
1035 Int_t iParent = mom->GetFirstMother();
1036 while(iParent!=-1){
1037 iMom = iParent;
1038 mom = fStack->Particle(iMom);
1039 iParent = mom->GetFirstMother();
1040 }
1041
1042 if(iMom<nproduced) {
1043 return kTRUE;
1044 } else {
1045 return kFALSE;
1046 }
03372fd1 1047 }
03372fd1 1048}
1049
fdceab34 1050#endif