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 "AliGenCocktailEventHeader.h"
58 //#include "$ALICE_ROOT/PWG4/JetTasks/AliAnalysisHelperJetTasks.h"
61 using namespace std; //required for resolving the 'cout' symbol
63 ClassImp(AliPWG4HighPtSpectra)
65 //__________________________________________________________________________
66 AliPWG4HighPtSpectra::AliPWG4HighPtSpectra() : AliAnalysisTask("AliPWG4HighPtSpectra", ""),
78 fSigmaConstrainedMax(5.),
94 //___________________________________________________________________________
95 AliPWG4HighPtSpectra::AliPWG4HighPtSpectra(const Char_t* name) :
96 AliAnalysisTask(name,""),
108 fSigmaConstrainedMax(5.),
121 // Constructor. Initialization of Inputs and Outputs
123 AliDebug(2,Form("AliPWG4HighPtSpectra Calling Constructor"));
124 // Input slot #0 works with a TChain ESD
125 DefineInput(0, TChain::Class());
126 // Output slot #0 writes into a TList
127 DefineOutput(0,TList::Class());
128 // Output slot #1, #2 writes into a AliCFContainer
129 DefineOutput(1,AliCFContainer::Class());
130 DefineOutput(2,AliCFContainer::Class());
131 // Output slot #3 writes into a AliESDtrackCuts
132 DefineOutput(3, AliESDtrackCuts::Class());
133 DefineOutput(4, AliESDtrackCuts::Class());
136 //________________________________________________________________________
137 void AliPWG4HighPtSpectra::LocalInit()
140 // Only called once at beginning
142 PostData(3,fTrackCuts);
145 //________________________________________________________________________
146 void AliPWG4HighPtSpectra::ConnectInputData(Option_t *)
150 AliDebug(2,Form(">> AliPWG4HighPtSpectra::ConnectInputData \n"));
152 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
154 AliDebug(2,Form( "ERROR: Could not read chain from input slot 0 \n"));
158 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
161 AliDebug(2,Form("ERROR: Could not get ESDInputHandler"));
164 fESD = esdH->GetEvent();
166 AliMCEventHandler *eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
168 AliDebug(2,Form( "ERROR: Could not retrieve MC event handler \n"));
171 fMC = eventHandler->MCEvent();
175 //________________________________________________________________________
176 Bool_t AliPWG4HighPtSpectra::SelectEvent() {
178 // Decide if event should be selected for analysis
181 // Checks following requirements:
183 // - trigger info from AliPhysicsSelection
184 // - number of reconstructed tracks > 1
185 // - primary vertex reconstructed
186 // - z-vertex < 10 cm
188 Bool_t selectEvent = kTRUE;
190 //fESD object available?
192 AliDebug(2,Form("ERROR: fInputEvent not available\n"));
193 fNEventReject->Fill("noESD",1);
194 selectEvent = kFALSE;
199 UInt_t isSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
200 if(!(isSelected&AliVEvent::kMB)) { //Select collison candidates
201 AliDebug(2,Form(" Trigger Selection: event REJECTED ... "));
202 fNEventReject->Fill("Trigger",1);
203 selectEvent = kFALSE;
207 //Check if number of reconstructed tracks is larger than 1
208 if(!fESD->GetNumberOfTracks() || fESD->GetNumberOfTracks()<2) {
209 fNEventReject->Fill("NTracks<2",1);
210 selectEvent = kFALSE;
214 //Check if vertex is reconstructed
215 fVtx = fESD->GetPrimaryVertexSPD();
218 fNEventReject->Fill("noVTX",1);
219 selectEvent = kFALSE;
223 if(!fVtx->GetStatus()) {
224 fNEventReject->Fill("VtxStatus",1);
225 selectEvent = kFALSE;
230 // TString vtxName(fVtx->GetName());
231 if(fVtx->GetNContributors()<2) {
232 fNEventReject->Fill("NCont<2",1);
233 selectEvent = kFALSE;
237 //Check if z-vertex < 10 cm
239 fVtx->GetXYZ(primVtx);
240 if(TMath::Sqrt(primVtx[0]*primVtx[0] + primVtx[1]*primVtx[1])>1. || TMath::Abs(primVtx[2]>10.)){
241 fNEventReject->Fill("ZVTX>10",1);
242 selectEvent = kFALSE;
246 //Centrality selection should only be done in case of PbPb
249 if(fCentClass!=CalculateCentrality(fESD) && fCentClass!=10) {
250 fNEventReject->Fill("cent",1);
251 selectEvent = kFALSE;
255 if(dynamic_cast<AliESDEvent*>(fESD)->GetCentrality()) {
256 cent = dynamic_cast<AliESDEvent*>(fESD)->GetCentrality()->GetCentralityPercentile("V0M");
259 fNEventReject->Fill("cent>90",1);
260 selectEvent = kFALSE;
263 fh1Centrality->Fill(cent);
271 //________________________________________________________________________
272 Int_t AliPWG4HighPtSpectra::CalculateCentrality(AliESDEvent *esd){
278 if(esd->GetCentrality()){
279 cent = esd->GetCentrality()->GetCentralityPercentile("V0M");
291 //_________________________________________________
292 void AliPWG4HighPtSpectra::Exec(Option_t *)
295 // Main loop function
297 AliDebug(2,Form(">> AliPWG4HighPtSpectra::Exec \n"));
299 // All events without selection
300 fNEventAll->Fill(0.);
303 fNEventReject->Fill("NTracks<2",1);
305 PostData(0,fHistList);
306 PostData(1,fCFManagerPos->GetParticleContainer());
307 PostData(2,fCFManagerNeg->GetParticleContainer());
314 AliDebug(2,Form("MC particles: %d", fMC->GetNumberOfTracks()));
315 fStack = fMC->Stack(); //Particles Stack
316 AliDebug(2,Form("MC particles stack: %d", fStack->GetNtrack()));
319 // ---- Get MC Header information (for MC productions in pThard bins) ----
320 Double_t ptHard = 0.;
321 Double_t nTrials = 1; // trials for MC trigger weight for real data
324 AliGenPythiaEventHeader* pythiaGenHeader = GetPythiaEventHeader(fMC);
326 nTrials = pythiaGenHeader->Trials();
327 ptHard = pythiaGenHeader->GetPtHard();
329 fh1PtHard->Fill(ptHard);
330 fh1PtHardTrials->Fill(ptHard,nTrials);
332 fh1Trials->Fill("#sum{ntrials}",fAvgTrials);
336 Int_t nTracks = fESD->GetNumberOfTracks();
337 AliDebug(2,Form("nTracks %d", nTracks));
340 fNEventReject->Fill("noTrackCuts",1);
342 PostData(0,fHistList);
343 PostData(1,fCFManagerPos->GetParticleContainer());
344 PostData(2,fCFManagerNeg->GetParticleContainer());
348 // Selected events for analysis
349 fNEventSel->Fill(0.);
352 Double_t containerInputRec[3] = {0.,0.,0.};
353 Double_t containerInputMC[3] = {0.,0.,0.};
354 Double_t containerInputRecMC[3] = {0.,0.,0.}; //reconstructed yield as function of MC variable
356 //Now go to rec level
357 for (Int_t iTrack = 0; iTrack<nTracks; iTrack++)
359 //Get track for analysis
360 AliESDtrack *track = 0x0;
361 AliESDtrack *esdtrack = fESD->GetTrack(iTrack);
362 if(!esdtrack) continue;
365 track = AliESDtrackCuts::GetTPCOnlyTrack(fESD,esdtrack->GetID());
366 else if(fTrackType==2) {
367 track = AliESDtrackCuts::GetTPCOnlyTrack(fESD,esdtrack->GetID());
370 AliExternalTrackParam exParam;
371 Bool_t relate = track->RelateToVertexTPC(fVtx,fESD->GetMagneticField(),kVeryBig,&exParam);
376 track->Set(exParam.GetX(),exParam.GetAlpha(),exParam.GetParameter(),exParam.GetCovariance());
382 if(fTrackType==1 || fTrackType==2) delete track;
387 //Cut on chi2 of constrained fit
388 if(track->GetConstrainedChi2TPC() > fSigmaConstrainedMax*fSigmaConstrainedMax) {
396 containerInputRec[0] = track->Pt();
397 containerInputRec[1] = track->Phi();
398 containerInputRec[2] = track->Eta();
400 if (fTrackCuts->AcceptTrack(track)) {
401 if(track->GetSign()>0.) fCFManagerPos->GetParticleContainer()->Fill(containerInputRec,kStepReconstructed);
402 if(track->GetSign()<0.) fCFManagerNeg->GetParticleContainer()->Fill(containerInputRec,kStepReconstructed);
405 //Only fill the MC containers if MC information is available
407 Int_t label = TMath::Abs(track->GetLabel());
408 TParticle *particle = fStack->Particle(label) ;
410 if(fTrackType==1 || fTrackType==2)
414 containerInputRecMC[0] = particle->Pt();
415 containerInputRecMC[1] = particle->Phi();
416 containerInputRecMC[2] = particle->Eta();
418 //Container with primaries
419 if(fStack->IsPhysicalPrimary(label)) {
420 if(particle->GetPDG()->Charge()>0.) {
421 fCFManagerPos->GetParticleContainer()->Fill(containerInputRecMC,kStepReconstructedMC);
423 if(particle->GetPDG()->Charge()<0.) {
424 fCFManagerNeg->GetParticleContainer()->Fill(containerInputRecMC,kStepReconstructedMC);
428 //Container with secondaries
429 if (!fStack->IsPhysicalPrimary(label) ) {
430 if(particle->GetPDG()->Charge()>0.) {
431 fCFManagerPos->GetParticleContainer()->Fill(containerInputMC,kStepSecondaries);
433 if(particle->GetPDG()->Charge()<0.) {
434 fCFManagerNeg->GetParticleContainer()->Fill(containerInputMC,kStepSecondaries);
439 }//trackCuts global tracks
441 if(fTrackType==1 || fTrackType==2)
446 //Fill MC containters if particles are findable
448 for(int iPart = 1; iPart<(fMC->GetNumberOfPrimaries()); iPart++) {
449 AliMCParticle *mcPart = (AliMCParticle*)fMC->GetTrack(iPart);
450 if(!mcPart) continue;
452 containerInputMC[0] = mcPart->Pt();
453 containerInputMC[1] = mcPart->Phi();
454 containerInputMC[2] = mcPart->Eta();
456 if(fStack->IsPhysicalPrimary(iPart)) {
457 if(mcPart->Charge()>0. && fCFManagerPos->CheckParticleCuts(kStepMCAcceptance,mcPart)) fCFManagerPos->GetParticleContainer()->Fill(containerInputMC,kStepMCAcceptance);
458 if(mcPart->Charge()<0. && fCFManagerNeg->CheckParticleCuts(kStepMCAcceptance,mcPart)) fCFManagerNeg->GetParticleContainer()->Fill(containerInputMC,kStepMCAcceptance);
463 PostData(0,fHistList);
464 PostData(1,fCFManagerPos->GetParticleContainer());
465 PostData(2,fCFManagerNeg->GetParticleContainer());
468 //________________________________________________________________________
469 Bool_t AliPWG4HighPtSpectra::PythiaInfoFromFile(const char* currFile,Float_t &fXsec,Float_t &fTrials){
471 // get the cross section and the trails either from pyxsec.root or from pysec_hists.root
472 // This is to called in Notify and should provide the path to the AOD/ESD file
473 // Copied from AliAnalysisTaskJetSpectrum2
476 TString file(currFile);
480 if(file.Contains("root_archive.zip#")){
481 Ssiz_t pos1 = file.Index("root_archive",12,TString::kExact);
482 Ssiz_t pos = file.Index("#",1,pos1,TString::kExact);
483 file.Replace(pos+1,20,"");
486 // not an archive take the basename....
487 file.ReplaceAll(gSystem->BaseName(file.Data()),"");
489 // Printf("%s",file.Data());
491 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
493 // next trial fetch the histgram file
494 fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
496 // not a severe condition but inciate that we have no information
500 // find the tlist we want to be independtent of the name so use the Tkey
501 TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0);
506 TList *list = dynamic_cast<TList*>(key->ReadObj());
511 fXsec = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
512 fTrials = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
515 } // no tree pyxsec.root
517 TTree *xtree = (TTree*)fxsec->Get("Xsection");
523 Double_t xsection = 0;
524 xtree->SetBranchAddress("xsection",&xsection);
525 xtree->SetBranchAddress("ntrials",&ntrials);
533 //________________________________________________________________________
534 Bool_t AliPWG4HighPtSpectra::Notify()
537 // Implemented Notify() to read the cross sections
538 // and number of trials from pyxsec.root
539 // Copied from AliAnalysisTaskJetSpectrum2
542 TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
543 Float_t xsection = 0;
548 TFile *curfile = tree->GetCurrentFile();
550 Error("Notify","No current file");
553 if(!fh1Xsec||!fh1Trials){
554 // Printf("%s%d No Histogram fh1Xsec",(char*)__FILE__,__LINE__);
557 PythiaInfoFromFile(curfile->GetName(),xsection,ftrials);
558 fh1Xsec->Fill("<#sigma>",xsection);
559 // construct a poor man average trials
560 Float_t nEntries = (Float_t)tree->GetTree()->GetEntries();
561 if(ftrials>=nEntries && nEntries>0.)fAvgTrials = ftrials/nEntries;
566 //________________________________________________________________________
567 AliGenPythiaEventHeader* AliPWG4HighPtSpectra::GetPythiaEventHeader(AliMCEvent *mcEvent){
569 if(!mcEvent)return 0;
570 AliGenEventHeader* genHeader = mcEvent->GenEventHeader();
571 AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
572 if(!pythiaGenHeader){
574 AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(genHeader);
576 if (!genCocktailHeader) {
577 // AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Unknown header type (not Pythia or Cocktail)");
578 // AliWarning(Form("%s %d: Unknown header type (not Pythia or Cocktail)",(char*)__FILE__,__LINE__));
581 TList* headerList = genCocktailHeader->GetHeaders();
582 for (Int_t i=0; i<headerList->GetEntries(); i++) {
583 pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(headerList->At(i));
587 if(!pythiaGenHeader){
588 AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Pythia event header not found");
592 return pythiaGenHeader;
597 //___________________________________________________________________________
598 void AliPWG4HighPtSpectra::Terminate(Option_t*)
600 // The Terminate() function is the last function to be called during
601 // a query. It always runs on the client, it can be used to present
602 // the results graphically or save the results to file.
606 //___________________________________________________________________________
607 void AliPWG4HighPtSpectra::CreateOutputObjects() {
608 //HERE ONE CAN CREATE OUTPUT OBJECTS, IN PARTICULAR IF THE OBJECT PARAMETERS DON'T NEED
609 //TO BE SET BEFORE THE EXECUTION OF THE TASK
611 AliDebug(2,Form("CreateOutputObjects CreateOutputObjects of task %s", GetName()));
613 Bool_t oldStatus = TH1::AddDirectoryStatus();
614 TH1::AddDirectory(kFALSE);
618 fHistList = new TList();
619 fHistList->SetOwner(kTRUE);
620 fNEventAll = new TH1F("fNEventAll","NEventAll",1,-0.5,0.5);
621 fHistList->Add(fNEventAll);
622 fNEventSel = new TH1F("fNEventSel","NEvent Selected for analysis",1,-0.5,0.5);
623 fHistList->Add(fNEventSel);
625 fNEventReject = new TH1F("fNEventReject","Reason events are rejectected for analysis",20,0,20);
627 fNEventReject->Fill("noESD",0);
628 fNEventReject->Fill("Trigger",0);
629 fNEventReject->Fill("NTracks<2",0);
630 fNEventReject->Fill("noVTX",0);
631 fNEventReject->Fill("VtxStatus",0);
632 fNEventReject->Fill("NCont<2",0);
633 fNEventReject->Fill("ZVTX>10",0);
634 fNEventReject->Fill("cent",0);
635 fNEventReject->Fill("cent>90",0);
636 fHistList->Add(fNEventReject);
638 fh1Centrality = new TH1F("fh1Centrality","fh1Centrality; Centrality %",100,0,100);
639 fHistList->Add(fh1Centrality);
641 fh1Xsec = new TProfile("fh1Xsec","xsec from pyxsec.root",1,0,1);
642 fh1Xsec->GetXaxis()->SetBinLabel(1,"<#sigma>");
643 fHistList->Add(fh1Xsec);
645 fh1Trials = new TH1F("fh1Trials","trials root file",1,0,1);
646 fh1Trials->GetXaxis()->SetBinLabel(1,"#sum{ntrials}");
647 fHistList->Add(fh1Trials);
649 fh1PtHard = new TH1F("fh1PtHard","PYTHIA Pt hard;p_{T,hard}",350,-.5,349.5);
650 fHistList->Add(fh1PtHard);
651 fh1PtHardTrials = new TH1F("fh1PtHardTrials","PYTHIA Pt hard weight with trials;p_{T,hard}",350,-.5,349.5);
652 fHistList->Add(fh1PtHardTrials);
654 TH1::AddDirectory(oldStatus);
656 PostData(0,fHistList);
657 PostData(1,fCFManagerPos->GetParticleContainer());
658 PostData(2,fCFManagerNeg->GetParticleContainer());