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 // Efficiency between different steps of the procedure.
18 // The ouptut of the task is a AliCFContainer from which the efficiencies
20 //-----------------------------------------------------------------------
21 // Author : Marta Verweij - UU
22 //-----------------------------------------------------------------------
25 #ifndef ALIPWG4HIGHPTSPECTRA_CXX
26 #define ALIPWG4HIGHPTSPECTRA_CXX
28 #include "AliPWG4HighPtSpectra.h"
42 #include "AliAnalysisManager.h"
43 #include "AliESDInputHandler.h"
44 #include "AliESDtrack.h"
45 #include "AliESDtrackCuts.h"
46 #include "AliExternalTrackParam.h"
47 #include "AliCentrality.h"
52 #include "TParticle.h"
53 #include "AliMCEvent.h"
54 #include "AliMCEventHandler.h"
55 #include "AliCFContainer.h"
56 #include "AliGenPythiaEventHeader.h"
57 #include "AliGenHijingEventHeader.h"
58 #include "AliGenCocktailEventHeader.h"
59 //#include "$ALICE_ROOT/PWG4/JetTasks/AliAnalysisHelperJetTasks.h"
62 using namespace std; //required for resolving the 'cout' symbol
64 ClassImp(AliPWG4HighPtSpectra)
66 //__________________________________________________________________________
67 AliPWG4HighPtSpectra::AliPWG4HighPtSpectra() : AliAnalysisTask("AliPWG4HighPtSpectra", ""),
75 fTriggerMask(AliVEvent::kMB),
80 fTrackCutsReject(0x0),
81 fbSelectHIJING(kFALSE),
82 fSigmaConstrainedMax(100.),
93 fPtRelUncertainty1PtPrim(0x0),
94 fPtRelUncertainty1PtSec(0x0)
100 //___________________________________________________________________________
101 AliPWG4HighPtSpectra::AliPWG4HighPtSpectra(const Char_t* name) :
102 AliAnalysisTask(name,""),
110 fTriggerMask(AliVEvent::kMB),
115 fTrackCutsReject(0x0),
116 fbSelectHIJING(kFALSE),
117 fSigmaConstrainedMax(100.),
128 fPtRelUncertainty1PtPrim(0x0),
129 fPtRelUncertainty1PtSec(0x0)
132 // Constructor. Initialization of Inputs and Outputs
134 AliDebug(2,Form("AliPWG4HighPtSpectra Calling Constructor"));
135 // Input slot #0 works with a TChain ESD
136 DefineInput(0, TChain::Class());
137 // Output slot #0 writes into a TList
138 DefineOutput(0,TList::Class());
139 // Output slot #1, #2 writes into a AliCFContainer
140 DefineOutput(1,AliCFContainer::Class());
141 DefineOutput(2,AliCFContainer::Class());
142 // Output slot #3 writes into a AliESDtrackCuts
143 DefineOutput(3, AliESDtrackCuts::Class());
144 DefineOutput(4, AliESDtrackCuts::Class());
147 //________________________________________________________________________
148 void AliPWG4HighPtSpectra::LocalInit()
151 // Only called once at beginning
153 PostData(3,fTrackCuts);
156 //________________________________________________________________________
157 void AliPWG4HighPtSpectra::ConnectInputData(Option_t *)
161 AliDebug(2,Form(">> AliPWG4HighPtSpectra::ConnectInputData \n"));
163 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
165 AliDebug(2,Form( "ERROR: Could not read chain from input slot 0 \n"));
169 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
172 AliDebug(2,Form("ERROR: Could not get ESDInputHandler"));
175 fESD = esdH->GetEvent();
177 AliMCEventHandler *eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
179 AliDebug(2,Form( "ERROR: Could not retrieve MC event handler \n"));
182 fMC = eventHandler->MCEvent();
186 //________________________________________________________________________
187 Bool_t AliPWG4HighPtSpectra::SelectEvent() {
189 // Decide if event should be selected for analysis
192 // Checks following requirements:
194 // - trigger info from AliPhysicsSelection
195 // - number of reconstructed tracks > 1
196 // - primary vertex reconstructed
197 // - z-vertex < 10 cm
199 Bool_t selectEvent = kTRUE;
201 //fESD object available?
203 AliDebug(2,Form("ERROR: fInputEvent not available\n"));
204 fNEventReject->Fill("noESD",1);
205 selectEvent = kFALSE;
210 UInt_t isSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
211 if(fTriggerMask != AliVEvent::kAny && !(isSelected&fTriggerMask)) { //Select collison candidates
212 AliDebug(2,Form(" Trigger Selection: event REJECTED ... "));
213 fNEventReject->Fill("Trigger",1);
214 selectEvent = kFALSE;
218 //Check if number of reconstructed tracks is larger than 1
219 if(!fESD->GetNumberOfTracks() || fESD->GetNumberOfTracks()<2) {
220 fNEventReject->Fill("NTracks<2",1);
221 selectEvent = kFALSE;
225 //Check if vertex is reconstructed
226 fVtx = fESD->GetPrimaryVertexSPD();
229 fNEventReject->Fill("noVTX",1);
230 selectEvent = kFALSE;
234 if(!fVtx->GetStatus()) {
235 fNEventReject->Fill("VtxStatus",1);
236 selectEvent = kFALSE;
241 // TString vtxName(fVtx->GetName());
242 if(fVtx->GetNContributors()<2) {
243 fNEventReject->Fill("NCont<2",1);
244 selectEvent = kFALSE;
248 //Check if z-vertex < 10 cm
250 fVtx->GetXYZ(primVtx);
251 if(TMath::Sqrt(primVtx[0]*primVtx[0] + primVtx[1]*primVtx[1])>1. || TMath::Abs(primVtx[2]>10.)){
252 fNEventReject->Fill("ZVTX>10",1);
253 selectEvent = kFALSE;
257 //Centrality selection should only be done in case of PbPb
260 if(fCentClass!=CalculateCentrality(fESD) && fCentClass!=10) {
261 fNEventReject->Fill("cent",1);
262 selectEvent = kFALSE;
266 if(dynamic_cast<AliESDEvent*>(fESD)->GetCentrality()) {
267 cent = dynamic_cast<AliESDEvent*>(fESD)->GetCentrality()->GetCentralityPercentile("V0M");
270 fNEventReject->Fill("cent>90",1);
271 selectEvent = kFALSE;
274 fh1Centrality->Fill(cent);
282 //________________________________________________________________________
283 Int_t AliPWG4HighPtSpectra::CalculateCentrality(AliESDEvent *esd){
289 if(esd->GetCentrality()){
290 cent = esd->GetCentrality()->GetCentralityPercentile("V0M");
303 //_________________________________________________
304 void AliPWG4HighPtSpectra::Exec(Option_t *)
307 // Main loop function
309 AliDebug(2,Form(">> AliPWG4HighPtSpectra::Exec \n"));
312 AliMCEventHandler *eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
314 AliDebug(2,Form( "ERROR: Could not retrieve MC event handler \n"));
317 fMC = eventHandler->MCEvent();
320 // All events without selection
321 fNEventAll->Fill(0.);
325 PostData(0,fHistList);
326 PostData(1,fCFManagerPos->GetParticleContainer());
327 PostData(2,fCFManagerNeg->GetParticleContainer());
334 AliDebug(2,Form("MC particles: %d", fMC->GetNumberOfTracks()));
335 fStack = fMC->Stack(); //Particles Stack
336 AliDebug(2,Form("MC particles stack: %d", fStack->GetNtrack()));
339 Int_t nTracks = fESD->GetNumberOfTracks();
340 AliDebug(2,Form("nTracks %d", nTracks));
343 fNEventReject->Fill("noTrackCuts",1);
345 PostData(0,fHistList);
346 PostData(1,fCFManagerPos->GetParticleContainer());
347 PostData(2,fCFManagerNeg->GetParticleContainer());
351 // Selected events for analysis
352 fNEventSel->Fill(0.);
354 const Int_t nvar = 4;
356 Double_t containerInputRec[nvar] = {0.,0.,0.,0.};
357 Double_t containerInputMC[nvar] = {0.,0.,0.,0.};
358 Double_t containerInputRecMC[nvar] = {0.,0.,0.,0.}; //reconstructed yield as function of MC variable
360 //Now go to rec level
361 for (Int_t iTrack = 0; iTrack<nTracks; iTrack++)
363 //Get track for analysis
364 AliESDtrack *track = 0x0;
365 AliESDtrack *esdtrack = fESD->GetTrack(iTrack);
371 if (!(fTrackCuts->AcceptTrack(esdtrack)))
376 track = AliESDtrackCuts::GetTPCOnlyTrack(fESD,esdtrack->GetID());
377 else if(fTrackType==2 || fTrackType==4) {
378 track = AliESDtrackCuts::GetTPCOnlyTrack(const_cast<AliESDEvent*>(fESD),esdtrack->GetID());
382 AliExternalTrackParam exParam;
383 Bool_t relate = track->RelateToVertexTPC(fVtx,fESD->GetMagneticField(),kVeryBig,&exParam);
385 if(track) delete track;
388 track->Set(exParam.GetX(),exParam.GetAlpha(),exParam.GetParameter(),exParam.GetCovariance());
390 else if(fTrackType==5 || fTrackType==6) {
391 if(fTrackCuts->AcceptTrack(esdtrack)) {
395 if( !(fTrackCutsReject->AcceptTrack(esdtrack)) && fTrackCuts->AcceptTrack(esdtrack) ) {
398 //use TPConly constrained track
399 track = AliESDtrackCuts::GetTPCOnlyTrack(fESD,esdtrack->GetID());
403 AliExternalTrackParam exParam;
404 Bool_t relate = track->RelateToVertexTPC(fVtx,fESD->GetMagneticField(),kVeryBig,&exParam);
406 if(track) delete track;
409 track->Set(exParam.GetX(),exParam.GetAlpha(),exParam.GetParameter(),exParam.GetCovariance());
411 else if(fTrackType==6) {
412 //use global constrained track
413 track = new AliESDtrack(*esdtrack);
414 track->Set(esdtrack->GetConstrainedParam()->GetX(),esdtrack->GetConstrainedParam()->GetAlpha(),esdtrack->GetConstrainedParam()->GetParameter(),esdtrack->GetConstrainedParam()->GetCovariance());
420 else if(fTrackType==7) {
421 //use global constrained track
422 track = new AliESDtrack(*esdtrack);
431 if(fTrackType==2 || fTrackType==4 || fTrackType==5) {
432 //Cut on chi2 of constrained fit
433 if(track->GetConstrainedChi2TPC() > fSigmaConstrainedMax*fSigmaConstrainedMax && fSigmaConstrainedMax>0.) {
434 if(track) delete track;
438 if (!(fTrackCuts->AcceptTrack(track)) && fTrackType!=4 && fTrackType!=5 && fTrackType!=6) {
439 if(fTrackType==1 || fTrackType==2 || fTrackType==7) {
440 if(track) delete track;
446 if(fTrackCutsReject ) {
447 if(fTrackCutsReject->AcceptTrack(track) ) {
448 if(track) delete track;
453 if(esdtrack->GetConstrainedParam())
454 track->Set(esdtrack->GetConstrainedParam()->GetX(),esdtrack->GetConstrainedParam()->GetAlpha(),esdtrack->GetConstrainedParam()->GetParameter(),esdtrack->GetConstrainedParam()->GetCovariance());
458 if(fTrackType==1 || fTrackType==2 || fTrackType==4 || fTrackType==5 || fTrackType==6 || fTrackType==7) {
459 if(track) delete track;
465 containerInputRec[0] = track->Pt();
466 containerInputRec[1] = track->Phi();
467 containerInputRec[2] = track->Eta();
468 containerInputRec[3] = track->GetTPCNclsIter1();
470 if(track->GetSign()>0.) fCFManagerPos->GetParticleContainer()->Fill(containerInputRec,kStepReconstructed);
471 if(track->GetSign()<0.) fCFManagerNeg->GetParticleContainer()->Fill(containerInputRec,kStepReconstructed);
474 //Only fill the MC containers if MC information is available
476 Int_t label = TMath::Abs(track->GetLabel());
477 if(label>fStack->GetNtrack()) {
478 if(fTrackType==1 || fTrackType==2 || fTrackType==7)
482 //Only select particles generated by HIJING if requested
484 if(!IsHIJINGParticle(label)) {
485 if(fTrackType==1 || fTrackType==2 || fTrackType==7)
490 TParticle *particle = fStack->Particle(label) ;
492 if(fTrackType==1 || fTrackType==2 || fTrackType==7)
496 containerInputRecMC[0] = particle->Pt();
497 containerInputRecMC[1] = particle->Phi();
498 containerInputRecMC[2] = particle->Eta();
499 containerInputRecMC[3] = track->GetTPCNclsIter1();
501 //Container with primaries
502 if(fStack->IsPhysicalPrimary(label)) {
503 if(particle->GetPDG()->Charge()>0.) {
504 fCFManagerPos->GetParticleContainer()->Fill(containerInputRecMC,kStepReconstructedMC);
506 if(particle->GetPDG()->Charge()<0.) {
507 fCFManagerNeg->GetParticleContainer()->Fill(containerInputRecMC,kStepReconstructedMC);
510 //Fill pT resolution plots for primaries
511 fPtRelUncertainty1PtPrim->Fill(containerInputRec[0],containerInputRec[0]*TMath::Sqrt(track->GetSigma1Pt2()));
515 //Container with secondaries
516 if (!fStack->IsPhysicalPrimary(label) ) {
517 if(particle->GetPDG()->Charge()>0.) {
518 fCFManagerPos->GetParticleContainer()->Fill(containerInputMC,kStepSecondaries);
520 if(particle->GetPDG()->Charge()<0.) {
521 fCFManagerNeg->GetParticleContainer()->Fill(containerInputMC,kStepSecondaries);
523 //Fill pT resolution plots for primaries
524 fPtRelUncertainty1PtSec->Fill(containerInputRec[0],containerInputRec[0]*TMath::Sqrt(track->GetSigma1Pt2()));
528 if(fTrackType==1 || fTrackType==2 || fTrackType==4 || fTrackType==5 || fTrackType==6 || fTrackType==7) {
529 if(track) delete track;
536 //Fill MC containers if particles are findable
538 for(int iPart = 1; iPart<(fMC->GetNumberOfPrimaries()); iPart++) {
539 AliMCParticle *mcPart = (AliMCParticle*)fMC->GetTrack(iPart);
540 if(!mcPart) continue;
542 //Only select particles generated by HIJING if requested
544 if(!IsHIJINGParticle(iPart))
548 Int_t pdg = mcPart->PdgCode();
550 // select charged pions, protons, kaons , electrons, muons
551 if(TMath::Abs(pdg) == 211 || TMath::Abs(pdg) == 2212 || TMath::Abs(pdg) ==
552 321 || TMath::Abs(pdg) == 11 || TMath::Abs(pdg) == 13){
555 containerInputMC[0] = mcPart->Pt();
556 containerInputMC[1] = mcPart->Phi();
557 containerInputMC[2] = mcPart->Eta();
558 // AliESDtrack *esdtrack = fESD->GetTrack(mcPart->GetLabel());
559 containerInputMC[3] = 159.;
561 if(fStack->IsPhysicalPrimary(iPart)) {
562 if(mcPart->Charge()>0. && fCFManagerPos->CheckParticleCuts(kStepMCAcceptance,mcPart)) fCFManagerPos->GetParticleContainer()->Fill(containerInputMC,kStepMCAcceptance);
563 if(mcPart->Charge()<0. && fCFManagerNeg->CheckParticleCuts(kStepMCAcceptance,mcPart)) fCFManagerNeg->GetParticleContainer()->Fill(containerInputMC,kStepMCAcceptance);
569 PostData(0,fHistList);
570 PostData(1,fCFManagerPos->GetParticleContainer());
571 PostData(2,fCFManagerNeg->GetParticleContainer());
574 //________________________________________________________________________
575 Bool_t AliPWG4HighPtSpectra::PythiaInfoFromFile(const char* currFile,Float_t &fXsec,Float_t &fTrials){
577 // get the cross section and the trails either from pyxsec.root or from pysec_hists.root
578 // This is to called in Notify and should provide the path to the AOD/ESD file
579 // Copied from AliAnalysisTaskJetSpectrum2
582 TString file(currFile);
586 if(file.Contains("root_archive.zip#")){
587 Ssiz_t pos1 = file.Index("root_archive",12,TString::kExact);
588 Ssiz_t pos = file.Index("#",1,pos1,TString::kExact);
589 file.Replace(pos+1,20,"");
592 // not an archive take the basename....
593 file.ReplaceAll(gSystem->BaseName(file.Data()),"");
595 // Printf("%s",file.Data());
597 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
599 // next trial fetch the histgram file
600 fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
602 // not a severe condition but indicates that we have no information
606 // find the tlist we want to be independtent of the name so use the Tkey
607 TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0);
612 TList *list = dynamic_cast<TList*>(key->ReadObj());
617 fXsec = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
618 fTrials = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
621 } // no tree pyxsec.root
623 TTree *xtree = (TTree*)fxsec->Get("Xsection");
629 Double_t xsection = 0;
630 xtree->SetBranchAddress("xsection",&xsection);
631 xtree->SetBranchAddress("ntrials",&ntrials);
639 //________________________________________________________________________
640 Bool_t AliPWG4HighPtSpectra::Notify()
643 // Implemented Notify() to read the cross sections
644 // and number of trials from pyxsec.root
645 // Copied from AliAnalysisTaskJetSpectrum2
648 TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
649 Float_t xsection = 0;
653 TFile *curfile = tree->GetCurrentFile();
655 Error("Notify","No current file");
658 if(!fh1Xsec||!fh1Trials){
659 // Printf("%s%d No Histogram fh1Xsec",(char*)__FILE__,__LINE__);
662 PythiaInfoFromFile(curfile->GetName(),xsection,ftrials);
663 fh1Xsec->Fill("<#sigma>",xsection);
664 fh1Trials->Fill("#sum{ntrials}",ftrials);
669 //________________________________________________________________________
670 AliGenPythiaEventHeader* AliPWG4HighPtSpectra::GetPythiaEventHeader(AliMCEvent *mcEvent){
672 if(!mcEvent)return 0;
673 AliGenEventHeader* genHeader = mcEvent->GenEventHeader();
674 AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
675 if(!pythiaGenHeader){
677 AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(genHeader);
679 if (!genCocktailHeader) {
680 // AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Unknown header type (not Pythia or Cocktail)");
681 // AliWarning(Form("%s %d: Unknown header type (not Pythia or Cocktail)",(char*)__FILE__,__LINE__));
684 TList* headerList = genCocktailHeader->GetHeaders();
685 for (Int_t i=0; i<headerList->GetEntries(); i++) {
686 pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(headerList->At(i));
690 if(!pythiaGenHeader){
691 AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Pythia event header not found");
695 return pythiaGenHeader;
699 //________________________________________________________________________
700 AliGenHijingEventHeader* AliPWG4HighPtSpectra::GetHijingEventHeader(AliMCEvent *mcEvent){
702 if(!mcEvent)return 0;
703 AliGenEventHeader* genHeader = mcEvent->GenEventHeader();
704 AliGenHijingEventHeader* hijingGenHeader = dynamic_cast<AliGenHijingEventHeader*>(genHeader);
705 if(!hijingGenHeader){
707 AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(genHeader);
709 if (!genCocktailHeader) {
710 // AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Unknown header type (not Pythia or Cocktail)");
711 // AliWarning(Form("%s %d: Unknown header type (not Pythia or Cocktail)",(char*)__FILE__,__LINE__));
714 TList* headerList = genCocktailHeader->GetHeaders();
715 for (Int_t i=0; i<headerList->GetEntries(); i++) {
716 hijingGenHeader = dynamic_cast<AliGenHijingEventHeader*>(headerList->At(i));
720 if(!hijingGenHeader){
721 AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Hijing event header not found");
725 return hijingGenHeader;
730 //___________________________________________________________________________
731 void AliPWG4HighPtSpectra::Terminate(Option_t*)
733 // The Terminate() function is the last function to be called during
734 // a query. It always runs on the client, it can be used to present
735 // the results graphically or save the results to file.
739 //___________________________________________________________________________
740 void AliPWG4HighPtSpectra::CreateOutputObjects() {
741 //HERE ONE CAN CREATE OUTPUT OBJECTS, IN PARTICULAR IF THE OBJECT PARAMETERS DON'T NEED
742 //TO BE SET BEFORE THE EXECUTION OF THE TASK
744 AliDebug(2,Form("CreateOutputObjects CreateOutputObjects of task %s", GetName()));
746 Bool_t oldStatus = TH1::AddDirectoryStatus();
747 TH1::AddDirectory(kFALSE);
751 fHistList = new TList();
752 fHistList->SetOwner(kTRUE);
753 fNEventAll = new TH1F("fNEventAll","NEventAll",1,-0.5,0.5);
754 fHistList->Add(fNEventAll);
755 fNEventSel = new TH1F("fNEventSel","NEvent Selected for analysis",1,-0.5,0.5);
756 fHistList->Add(fNEventSel);
758 fNEventReject = new TH1F("fNEventReject","Reason events are rejectected for analysis",20,0,20);
760 fNEventReject->Fill("noESD",0);
761 fNEventReject->Fill("Trigger",0);
762 fNEventReject->Fill("NTracks<2",0);
763 fNEventReject->Fill("noVTX",0);
764 fNEventReject->Fill("VtxStatus",0);
765 fNEventReject->Fill("NCont<2",0);
766 fNEventReject->Fill("ZVTX>10",0);
767 fNEventReject->Fill("cent",0);
768 fNEventReject->Fill("cent>90",0);
769 fHistList->Add(fNEventReject);
771 fh1Centrality = new TH1F("fh1Centrality","fh1Centrality; Centrality %",100,0,100);
772 fHistList->Add(fh1Centrality);
774 fh1Xsec = new TProfile("fh1Xsec","xsec from pyxsec.root",1,0,1);
775 fh1Xsec->GetXaxis()->SetBinLabel(1,"<#sigma>");
776 fHistList->Add(fh1Xsec);
778 fh1Trials = new TH1F("fh1Trials","trials root file",1,0,1);
779 fh1Trials->GetXaxis()->SetBinLabel(1,"#sum{ntrials}");
780 fHistList->Add(fh1Trials);
782 fh1PtHard = new TH1F("fh1PtHard","PYTHIA Pt hard;p_{T,hard}",350,-.5,349.5);
783 fHistList->Add(fh1PtHard);
784 fh1PtHardTrials = new TH1F("fh1PtHardTrials","PYTHIA Pt hard weight with trials;p_{T,hard}",350,-.5,349.5);
785 fHistList->Add(fh1PtHardTrials);
787 Int_t fgkNPtBins = 100;
789 Float_t kMaxPt = 100.;
790 Double_t *binsPt = new Double_t[fgkNPtBins+1];
791 for(Int_t i=0; i<=fgkNPtBins; i++) binsPt[i]=(Double_t)kMinPt + (kMaxPt-kMinPt)/fgkNPtBins*(Double_t)i ;
793 Int_t fgkNRel1PtUncertaintyBins=50;
794 Float_t fgkRel1PtUncertaintyMin = 0.;
795 Float_t fgkRel1PtUncertaintyMax = 1.;
797 Double_t *binsRel1PtUncertainty=new Double_t[fgkNRel1PtUncertaintyBins+1];
798 for(Int_t i=0; i<=fgkNRel1PtUncertaintyBins; i++) binsRel1PtUncertainty[i]=(Double_t)fgkRel1PtUncertaintyMin + (fgkRel1PtUncertaintyMax-fgkRel1PtUncertaintyMin)/fgkNRel1PtUncertaintyBins*(Double_t)i ;
800 fPtRelUncertainty1PtPrim = new TH2F("fPtRelUncertainty1PtPrim","fPtRelUncertainty1PtPrim",fgkNPtBins,binsPt,fgkNRel1PtUncertaintyBins,binsRel1PtUncertainty);
801 fHistList->Add(fPtRelUncertainty1PtPrim);
803 fPtRelUncertainty1PtSec = new TH2F("fPtRelUncertainty1PtSec","fPtRelUncertainty1PtSec",fgkNPtBins,binsPt,fgkNRel1PtUncertaintyBins,binsRel1PtUncertainty);
804 fHistList->Add(fPtRelUncertainty1PtSec);
806 TH1::AddDirectory(oldStatus);
808 PostData(0,fHistList);
809 PostData(1,fCFManagerPos->GetParticleContainer());
810 PostData(2,fCFManagerNeg->GetParticleContainer());
812 if(binsPt) delete [] binsPt;
813 if(binsRel1PtUncertainty) delete [] binsRel1PtUncertainty;
817 //________________________________________________________________________
818 Bool_t AliPWG4HighPtSpectra::IsHIJINGParticle(Int_t label) {
820 // Return kTRUE in case particle is from HIJING event
823 AliGenHijingEventHeader* hijingHeader = GetHijingEventHeader(fMC);
825 Int_t nproduced = hijingHeader->NProduced();
827 TParticle * mom = fStack->Particle(label);
829 Int_t iParent = mom->GetFirstMother();
832 mom = fStack->Particle(iMom);
833 iParent = mom->GetFirstMother();