]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG4/JetTasks/AliPWG4HighPtSpectra.cxx
changed methods in AliAODPWG4Particle (G. Conesa)
[u/mrichter/AliRoot.git] / PWG4 / JetTasks / 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"
44#include "AliESDtrack.h"
45#include "AliESDtrackCuts.h"
46#include "AliExternalTrackParam.h"
2b553e6f 47#include "AliCentrality.h"
65e8ecdd 48
67ebd013 49#include "AliLog.h"
50
fdceab34 51#include "AliStack.h"
52#include "TParticle.h"
fdceab34 53#include "AliMCEvent.h"
54#include "AliMCEventHandler.h"
fdceab34 55#include "AliCFContainer.h"
b4691ee7 56#include "AliGenPythiaEventHeader.h"
57#include "AliGenCocktailEventHeader.h"
67ebd013 58//#include "$ALICE_ROOT/PWG4/JetTasks/AliAnalysisHelperJetTasks.h"
59
60//#include <iostream>
fdceab34 61using namespace std; //required for resolving the 'cout' symbol
62
63ClassImp(AliPWG4HighPtSpectra)
64
65//__________________________________________________________________________
66AliPWG4HighPtSpectra::AliPWG4HighPtSpectra() : AliAnalysisTask("AliPWG4HighPtSpectra", ""),
67 fReadAODData(0),
67ebd013 68 fCFManagerPos(0x0),
69 fCFManagerNeg(0x0),
cce400da 70 fESD(0x0),
71 fMC(0x0),
72 fStack(0x0),
34fc6450 73 fVtx(0x0),
2b553e6f 74 fIsPbPb(0),
75 fCentClass(10),
34fc6450 76 fTrackType(0),
cce400da 77 fTrackCuts(0x0),
2b553e6f 78 fSigmaConstrainedMax(5.),
b4691ee7 79 fAvgTrials(1),
fdceab34 80 fHistList(0),
b5cc0c6d 81 fNEventAll(0),
b4691ee7 82 fNEventSel(0),
cb76764e 83 fNEventReject(0),
2b553e6f 84 fh1Centrality(0x0),
b4691ee7 85 fh1Xsec(0),
86 fh1Trials(0),
87 fh1PtHard(0),
88 fh1PtHardTrials(0)
fdceab34 89{
90 //
91 //Default ctor
92 //
93}
94//___________________________________________________________________________
95AliPWG4HighPtSpectra::AliPWG4HighPtSpectra(const Char_t* name) :
96 AliAnalysisTask(name,""),
97 fReadAODData(0),
67ebd013 98 fCFManagerPos(0x0),
99 fCFManagerNeg(0x0),
cce400da 100 fESD(0x0),
101 fMC(0x0),
102 fStack(0x0),
34fc6450 103 fVtx(0x0),
2b553e6f 104 fIsPbPb(0),
105 fCentClass(10),
34fc6450 106 fTrackType(0),
cce400da 107 fTrackCuts(0x0),
2b553e6f 108 fSigmaConstrainedMax(5.),
b4691ee7 109 fAvgTrials(1),
fdceab34 110 fHistList(0),
b5cc0c6d 111 fNEventAll(0),
b4691ee7 112 fNEventSel(0),
cb76764e 113 fNEventReject(0),
2b553e6f 114 fh1Centrality(0x0),
b4691ee7 115 fh1Xsec(0),
116 fh1Trials(0),
117 fh1PtHard(0),
118 fh1PtHardTrials(0)
fdceab34 119{
120 //
121 // Constructor. Initialization of Inputs and Outputs
122 //
f51451be 123 AliDebug(2,Form("AliPWG4HighPtSpectra Calling Constructor"));
fdceab34 124 // Input slot #0 works with a TChain ESD
125 DefineInput(0, TChain::Class());
e5abcde9 126 // Output slot #0 writes into a TList
fdceab34 127 DefineOutput(0,TList::Class());
e5abcde9 128 // Output slot #1, #2 writes into a AliCFContainer
fdceab34 129 DefineOutput(1,AliCFContainer::Class());
67ebd013 130 DefineOutput(2,AliCFContainer::Class());
e5abcde9 131 // Output slot #3 writes into a AliESDtrackCuts
132 DefineOutput(3, AliESDtrackCuts::Class());
1f329128 133 DefineOutput(4, AliESDtrackCuts::Class());
e5abcde9 134}
135
136//________________________________________________________________________
137void AliPWG4HighPtSpectra::LocalInit()
138{
139 //
140 // Only called once at beginning
141 //
142 PostData(3,fTrackCuts);
fdceab34 143}
144
fdceab34 145//________________________________________________________________________
146void AliPWG4HighPtSpectra::ConnectInputData(Option_t *)
147{
148 // Connect ESD here
149 // Called once
df943115 150 AliDebug(2,Form(">> AliPWG4HighPtSpectra::ConnectInputData \n"));
fdceab34 151
152 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
153 if (!tree) {
cce400da 154 AliDebug(2,Form( "ERROR: Could not read chain from input slot 0 \n"));
155 return;
fdceab34 156 }
cce400da 157
158 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
159
160 if (!esdH) {
161 AliDebug(2,Form("ERROR: Could not get ESDInputHandler"));
162 return;
163 } else
164 fESD = esdH->GetEvent();
67ebd013 165
cce400da 166 AliMCEventHandler *eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
167 if (!eventHandler) {
168 AliDebug(2,Form( "ERROR: Could not retrieve MC event handler \n"));
169 }
170 else
171 fMC = eventHandler->MCEvent();
172
fdceab34 173}
cce400da 174
175//________________________________________________________________________
176Bool_t AliPWG4HighPtSpectra::SelectEvent() {
fdceab34 177 //
cce400da 178 // Decide if event should be selected for analysis
fdceab34 179 //
fdceab34 180
cce400da 181 // Checks following requirements:
182 // - fESD available
183 // - trigger info from AliPhysicsSelection
184 // - number of reconstructed tracks > 1
185 // - primary vertex reconstructed
186 // - z-vertex < 10 cm
187
188 Bool_t selectEvent = kTRUE;
b5cc0c6d 189
cce400da 190 //fESD object available?
fdceab34 191 if (!fESD) {
cce400da 192 AliDebug(2,Form("ERROR: fInputEvent not available\n"));
cb76764e 193 fNEventReject->Fill("noESD",1);
cce400da 194 selectEvent = kFALSE;
195 return selectEvent;
fdceab34 196 }
197
cce400da 198 //Trigger
bfd58011 199 UInt_t isSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
200 if(!(isSelected&AliVEvent::kMB)) { //Select collison candidates
b5cc0c6d 201 AliDebug(2,Form(" Trigger Selection: event REJECTED ... "));
cb76764e 202 fNEventReject->Fill("Trigger",1);
cce400da 203 selectEvent = kFALSE;
204 return selectEvent;
67ebd013 205 }
b5cc0c6d 206
cce400da 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;
211 return selectEvent;
fdceab34 212 }
b4691ee7 213
cce400da 214 //Check if vertex is reconstructed
34fc6450 215 fVtx = fESD->GetPrimaryVertexSPD();
216
217 if(!fVtx) {
218 fNEventReject->Fill("noVTX",1);
219 selectEvent = kFALSE;
220 return selectEvent;
221 }
222
223 if(!fVtx->GetStatus()) {
224 fNEventReject->Fill("VtxStatus",1);
225 selectEvent = kFALSE;
226 return selectEvent;
227 }
228
fdceab34 229 // Need vertex cut
34fc6450 230 // TString vtxName(fVtx->GetName());
231 if(fVtx->GetNContributors()<2) {
232 fNEventReject->Fill("NCont<2",1);
233 selectEvent = kFALSE;
234 return selectEvent;
cd9a6fa2 235 }
cce400da 236
237 //Check if z-vertex < 10 cm
b5cc0c6d 238 double primVtx[3];
34fc6450 239 fVtx->GetXYZ(primVtx);
67ebd013 240 if(TMath::Sqrt(primVtx[0]*primVtx[0] + primVtx[1]*primVtx[1])>1. || TMath::Abs(primVtx[2]>10.)){
cb76764e 241 fNEventReject->Fill("ZVTX>10",1);
cce400da 242 selectEvent = kFALSE;
243 return selectEvent;
67ebd013 244 }
cce400da 245
2b553e6f 246 //Centrality selection should only be done in case of PbPb
247 if(IsPbPb()) {
248 Float_t cent = 0.;
249 if(fCentClass!=CalculateCentrality(fESD) && fCentClass!=10) {
250 fNEventReject->Fill("cent",1);
251 selectEvent = kFALSE;
252 return selectEvent;
253 }
254 else {
255 if(dynamic_cast<AliESDEvent*>(fESD)->GetCentrality()) {
256 cent = dynamic_cast<AliESDEvent*>(fESD)->GetCentrality()->GetCentralityPercentile("V0M");
257 }
258 if(cent>90.) {
259 fNEventReject->Fill("cent>90",1);
260 selectEvent = kFALSE;
261 return selectEvent;
262 }
263 fh1Centrality->Fill(cent);
264 }
265 }
266
cce400da 267 return selectEvent;
268
269}
270
2b553e6f 271//________________________________________________________________________
272Int_t AliPWG4HighPtSpectra::CalculateCentrality(AliESDEvent *esd){
273
274
275 Float_t cent = 999;
276
277 if(esd){
278 if(esd->GetCentrality()){
279 cent = esd->GetCentrality()->GetCentralityPercentile("V0M");
280 }
281 }
282
b1041e3b 283 if(cent<0) return 5;
2b553e6f 284 if(cent>80)return 4;
285 if(cent>50)return 3;
286 if(cent>30)return 2;
287 if(cent>10)return 1;
288 return 0;
289
290}
291
cce400da 292//_________________________________________________
293void AliPWG4HighPtSpectra::Exec(Option_t *)
294{
295 //
296 // Main loop function
297 //
298 AliDebug(2,Form(">> AliPWG4HighPtSpectra::Exec \n"));
299
300 // All events without selection
301 fNEventAll->Fill(0.);
302
303 if(!SelectEvent()) {
cb76764e 304 fNEventReject->Fill("NTracks<2",1);
67ebd013 305 // Post output data
306 PostData(0,fHistList);
307 PostData(1,fCFManagerPos->GetParticleContainer());
308 PostData(2,fCFManagerNeg->GetParticleContainer());
cd9a6fa2 309 return;
310 }
cce400da 311
312 //MCEvent available?
313 //if yes: get stack
314 if(fMC) {
315 AliDebug(2,Form("MC particles: %d", fMC->GetNumberOfTracks()));
316 fStack = fMC->Stack(); //Particles Stack
317 AliDebug(2,Form("MC particles stack: %d", fStack->GetNtrack()));
318 }
319
320 // ---- Get MC Header information (for MC productions in pThard bins) ----
321 Double_t ptHard = 0.;
322 Double_t nTrials = 1; // trials for MC trigger weight for real data
323
324 if(fMC){
325 AliGenPythiaEventHeader* pythiaGenHeader = GetPythiaEventHeader(fMC);
326 if(pythiaGenHeader){
327 nTrials = pythiaGenHeader->Trials();
328 ptHard = pythiaGenHeader->GetPtHard();
329
330 fh1PtHard->Fill(ptHard);
331 fh1PtHardTrials->Fill(ptHard,nTrials);
332
333 fh1Trials->Fill("#sum{ntrials}",fAvgTrials);
334 }
335 }
336
fdceab34 337 Int_t nTracks = fESD->GetNumberOfTracks();
338 AliDebug(2,Form("nTracks %d", nTracks));
fdceab34 339
67ebd013 340 if(!fTrackCuts) {
cb76764e 341 fNEventReject->Fill("noTrackCuts",1);
67ebd013 342 // Post output data
343 PostData(0,fHistList);
344 PostData(1,fCFManagerPos->GetParticleContainer());
345 PostData(2,fCFManagerNeg->GetParticleContainer());
346 return;
347 }
348
b5cc0c6d 349 // Selected events for analysis
350 fNEventSel->Fill(0.);
67ebd013 351
352
ec555a3c 353 Double_t containerInputRec[3] = {0.,0.,0.};
ec555a3c 354 Double_t containerInputMC[3] = {0.,0.,0.};
2b553e6f 355 Double_t containerInputRecMC[3] = {0.,0.,0.}; //reconstructed yield as function of MC variable
ec555a3c 356
fdceab34 357 //Now go to rec level
358 for (Int_t iTrack = 0; iTrack<nTracks; iTrack++)
359 {
2b553e6f 360 //Get track for analysis
361 AliESDtrack *track = 0x0;
362 AliESDtrack *esdtrack = fESD->GetTrack(iTrack);
363 if(!esdtrack) continue;
364
365 if(fTrackType==1)
366 track = AliESDtrackCuts::GetTPCOnlyTrack(fESD,esdtrack->GetID());
367 else if(fTrackType==2) {
368 track = AliESDtrackCuts::GetTPCOnlyTrack(fESD,esdtrack->GetID());
369 if(!track) continue;
370
371 AliExternalTrackParam exParam;
372 Bool_t relate = track->RelateToVertexTPC(fVtx,fESD->GetMagneticField(),kVeryBig,&exParam);
373 if( !relate ) {
374 delete track;
34fc6450 375 continue;
376 }
2b553e6f 377 track->Set(exParam.GetX(),exParam.GetAlpha(),exParam.GetParameter(),exParam.GetCovariance());
bcd60a94 378 }
d3a3f33d 379 else if(fTrackType==7) {
380 //use global constrained track
381 track = esdtrack;
031aac8e 382 // track->Set(esdtrack->GetConstrainedParam()->GetX(),esdtrack->GetConstrainedParam()->GetAlpha(),esdtrack->GetConstrainedParam()->GetParameter(),esdtrack->GetConstrainedParam()->GetCovariance());
d3a3f33d 383
384 }
bcd60a94 385 else
2b553e6f 386 track = esdtrack;
bcd60a94 387
bcd60a94 388 if(!track) {
bcd60a94 389 continue;
390 }
2b553e6f 391
392 if(fTrackType==2) {
393 //Cut on chi2 of constrained fit
394 if(track->GetConstrainedChi2TPC() > fSigmaConstrainedMax*fSigmaConstrainedMax) {
395 delete track;
396 continue;
397 }
398 }
fdceab34 399
34fc6450 400
031aac8e 401
1f329128 402 if (fTrackCuts->AcceptTrack(track)) {
031aac8e 403
404 if(fTrackType==7) {
405 if(esdtrack->GetConstrainedParam())
406 track->Set(esdtrack->GetConstrainedParam()->GetX(),esdtrack->GetConstrainedParam()->GetAlpha(),esdtrack->GetConstrainedParam()->GetParameter(),esdtrack->GetConstrainedParam()->GetCovariance());
407 else
408 continue;
409 }
410
411 //fill the container
412 containerInputRec[0] = track->Pt();
413 containerInputRec[1] = track->Phi();
414 containerInputRec[2] = track->Eta();
415
1f329128 416 if(track->GetSign()>0.) fCFManagerPos->GetParticleContainer()->Fill(containerInputRec,kStepReconstructed);
417 if(track->GetSign()<0.) fCFManagerNeg->GetParticleContainer()->Fill(containerInputRec,kStepReconstructed);
418
67ebd013 419
1f329128 420 //Only fill the MC containers if MC information is available
cce400da 421 if(fMC) {
0cc50ca6 422 Int_t label = TMath::Abs(track->GetLabel());
cce400da 423 TParticle *particle = fStack->Particle(label) ;
bcd60a94 424 if(!particle) {
2b553e6f 425 if(fTrackType==1 || fTrackType==2)
426 delete track;
bcd60a94 427 continue;
428 }
ec555a3c 429 containerInputRecMC[0] = particle->Pt();
430 containerInputRecMC[1] = particle->Phi();
431 containerInputRecMC[2] = particle->Eta();
67ebd013 432
1f329128 433 //Container with primaries
cce400da 434 if(fStack->IsPhysicalPrimary(label)) {
1f329128 435 if(particle->GetPDG()->Charge()>0.) {
ec555a3c 436 fCFManagerPos->GetParticleContainer()->Fill(containerInputRecMC,kStepReconstructedMC);
1f329128 437 }
438 if(particle->GetPDG()->Charge()<0.) {
ec555a3c 439 fCFManagerNeg->GetParticleContainer()->Fill(containerInputRecMC,kStepReconstructedMC);
1f329128 440 }
67ebd013 441 }
442
1f329128 443 //Container with secondaries
cce400da 444 if (!fStack->IsPhysicalPrimary(label) ) {
67ebd013 445 if(particle->GetPDG()->Charge()>0.) {
aa3ba8d2 446 fCFManagerPos->GetParticleContainer()->Fill(containerInputMC,kStepSecondaries);
67ebd013 447 }
448 if(particle->GetPDG()->Charge()<0.) {
aa3ba8d2 449 fCFManagerNeg->GetParticleContainer()->Fill(containerInputMC,kStepSecondaries);
67ebd013 450 }
0cc50ca6 451 }
452 }
67ebd013 453
2b553e6f 454 }//trackCuts global tracks
60829498 455
2b553e6f 456 if(fTrackType==1 || fTrackType==2)
457 delete track;
1f329128 458 }//track loop
67ebd013 459
fdceab34 460
65e8ecdd 461 //Fill MC containters if particles are findable
cce400da 462 if(fMC) {
463 for(int iPart = 1; iPart<(fMC->GetNumberOfPrimaries()); iPart++) {
464 AliMCParticle *mcPart = (AliMCParticle*)fMC->GetTrack(iPart);
465 if(!mcPart) continue;
466 //fill the container
467 containerInputMC[0] = mcPart->Pt();
468 containerInputMC[1] = mcPart->Phi();
469 containerInputMC[2] = mcPart->Eta();
470
471 if(fStack->IsPhysicalPrimary(iPart)) {
472 if(mcPart->Charge()>0. && fCFManagerPos->CheckParticleCuts(kStepMCAcceptance,mcPart)) fCFManagerPos->GetParticleContainer()->Fill(containerInputMC,kStepMCAcceptance);
473 if(mcPart->Charge()<0. && fCFManagerNeg->CheckParticleCuts(kStepMCAcceptance,mcPart)) fCFManagerNeg->GetParticleContainer()->Fill(containerInputMC,kStepMCAcceptance);
0cc50ca6 474 }
cce400da 475 }
0cc50ca6 476 }
67ebd013 477
67ebd013 478 PostData(0,fHistList);
479 PostData(1,fCFManagerPos->GetParticleContainer());
480 PostData(2,fCFManagerNeg->GetParticleContainer());
481
fdceab34 482}
b4691ee7 483//________________________________________________________________________
484Bool_t AliPWG4HighPtSpectra::PythiaInfoFromFile(const char* currFile,Float_t &fXsec,Float_t &fTrials){
485 //
486 // get the cross section and the trails either from pyxsec.root or from pysec_hists.root
487 // This is to called in Notify and should provide the path to the AOD/ESD file
488 // Copied from AliAnalysisTaskJetSpectrum2
489 //
490
491 TString file(currFile);
492 fXsec = 0;
493 fTrials = 1;
494
495 if(file.Contains("root_archive.zip#")){
496 Ssiz_t pos1 = file.Index("root_archive",12,TString::kExact);
497 Ssiz_t pos = file.Index("#",1,pos1,TString::kExact);
498 file.Replace(pos+1,20,"");
499 }
500 else {
501 // not an archive take the basename....
502 file.ReplaceAll(gSystem->BaseName(file.Data()),"");
503 }
cb76764e 504 // Printf("%s",file.Data());
b4691ee7 505
506 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
507 if(!fxsec){
508 // next trial fetch the histgram file
509 fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
510 if(!fxsec){
511 // not a severe condition but inciate that we have no information
512 return kFALSE;
513 }
514 else{
515 // find the tlist we want to be independtent of the name so use the Tkey
516 TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0);
517 if(!key){
518 fxsec->Close();
519 return kFALSE;
520 }
521 TList *list = dynamic_cast<TList*>(key->ReadObj());
522 if(!list){
523 fxsec->Close();
524 return kFALSE;
525 }
526 fXsec = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
527 fTrials = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
528 fxsec->Close();
529 }
530 } // no tree pyxsec.root
531 else {
532 TTree *xtree = (TTree*)fxsec->Get("Xsection");
533 if(!xtree){
534 fxsec->Close();
535 return kFALSE;
536 }
537 UInt_t ntrials = 0;
538 Double_t xsection = 0;
539 xtree->SetBranchAddress("xsection",&xsection);
540 xtree->SetBranchAddress("ntrials",&ntrials);
541 xtree->GetEntry(0);
542 fTrials = ntrials;
543 fXsec = xsection;
544 fxsec->Close();
545 }
546 return kTRUE;
547}
548//________________________________________________________________________
549Bool_t AliPWG4HighPtSpectra::Notify()
550{
551 //
552 // Implemented Notify() to read the cross sections
553 // and number of trials from pyxsec.root
554 // Copied from AliAnalysisTaskJetSpectrum2
555 //
556
557 TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
558 Float_t xsection = 0;
559 Float_t ftrials = 1;
560
561 fAvgTrials = 1;
562 if(tree){
563 TFile *curfile = tree->GetCurrentFile();
564 if (!curfile) {
565 Error("Notify","No current file");
566 return kFALSE;
567 }
568 if(!fh1Xsec||!fh1Trials){
cb76764e 569 // Printf("%s%d No Histogram fh1Xsec",(char*)__FILE__,__LINE__);
b4691ee7 570 return kFALSE;
571 }
572 PythiaInfoFromFile(curfile->GetName(),xsection,ftrials);
573 fh1Xsec->Fill("<#sigma>",xsection);
574 // construct a poor man average trials
575 Float_t nEntries = (Float_t)tree->GetTree()->GetEntries();
576 if(ftrials>=nEntries && nEntries>0.)fAvgTrials = ftrials/nEntries;
577 }
578 return kTRUE;
579}
580
581//________________________________________________________________________
582AliGenPythiaEventHeader* AliPWG4HighPtSpectra::GetPythiaEventHeader(AliMCEvent *mcEvent){
583
584 if(!mcEvent)return 0;
585 AliGenEventHeader* genHeader = mcEvent->GenEventHeader();
586 AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
587 if(!pythiaGenHeader){
588 // cocktail ??
589 AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(genHeader);
590
591 if (!genCocktailHeader) {
cb76764e 592 // AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Unknown header type (not Pythia or Cocktail)");
b4691ee7 593 // AliWarning(Form("%s %d: Unknown header type (not Pythia or Cocktail)",(char*)__FILE__,__LINE__));
594 return 0;
595 }
596 TList* headerList = genCocktailHeader->GetHeaders();
597 for (Int_t i=0; i<headerList->GetEntries(); i++) {
598 pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(headerList->At(i));
599 if (pythiaGenHeader)
600 break;
601 }
602 if(!pythiaGenHeader){
603 AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Pythia event header not found");
604 return 0;
605 }
606 }
607 return pythiaGenHeader;
608
609}
fdceab34 610
611
612//___________________________________________________________________________
613void AliPWG4HighPtSpectra::Terminate(Option_t*)
614{
615 // The Terminate() function is the last function to be called during
616 // a query. It always runs on the client, it can be used to present
617 // the results graphically or save the results to file.
9b58297c 618
fdceab34 619}
620
fdceab34 621//___________________________________________________________________________
622void AliPWG4HighPtSpectra::CreateOutputObjects() {
623 //HERE ONE CAN CREATE OUTPUT OBJECTS, IN PARTICULAR IF THE OBJECT PARAMETERS DON'T NEED
624 //TO BE SET BEFORE THE EXECUTION OF THE TASK
625 //
f51451be 626 AliDebug(2,Form("CreateOutputObjects CreateOutputObjects of task %s", GetName()));
df943115 627
67ebd013 628 Bool_t oldStatus = TH1::AddDirectoryStatus();
629 TH1::AddDirectory(kFALSE);
630
fdceab34 631 //slot #1
b5cc0c6d 632 OpenFile(0);
633 fHistList = new TList();
b4691ee7 634 fHistList->SetOwner(kTRUE);
b5cc0c6d 635 fNEventAll = new TH1F("fNEventAll","NEventAll",1,-0.5,0.5);
636 fHistList->Add(fNEventAll);
637 fNEventSel = new TH1F("fNEventSel","NEvent Selected for analysis",1,-0.5,0.5);
638 fHistList->Add(fNEventSel);
fdceab34 639
cb76764e 640 fNEventReject = new TH1F("fNEventReject","Reason events are rejectected for analysis",20,0,20);
2b553e6f 641 //Set labels
642 fNEventReject->Fill("noESD",0);
643 fNEventReject->Fill("Trigger",0);
644 fNEventReject->Fill("NTracks<2",0);
645 fNEventReject->Fill("noVTX",0);
646 fNEventReject->Fill("VtxStatus",0);
647 fNEventReject->Fill("NCont<2",0);
648 fNEventReject->Fill("ZVTX>10",0);
649 fNEventReject->Fill("cent",0);
650 fNEventReject->Fill("cent>90",0);
cb76764e 651 fHistList->Add(fNEventReject);
652
2b553e6f 653 fh1Centrality = new TH1F("fh1Centrality","fh1Centrality; Centrality %",100,0,100);
654 fHistList->Add(fh1Centrality);
655
b4691ee7 656 fh1Xsec = new TProfile("fh1Xsec","xsec from pyxsec.root",1,0,1);
657 fh1Xsec->GetXaxis()->SetBinLabel(1,"<#sigma>");
658 fHistList->Add(fh1Xsec);
659
660 fh1Trials = new TH1F("fh1Trials","trials root file",1,0,1);
661 fh1Trials->GetXaxis()->SetBinLabel(1,"#sum{ntrials}");
662 fHistList->Add(fh1Trials);
663
664 fh1PtHard = new TH1F("fh1PtHard","PYTHIA Pt hard;p_{T,hard}",350,-.5,349.5);
665 fHistList->Add(fh1PtHard);
666 fh1PtHardTrials = new TH1F("fh1PtHardTrials","PYTHIA Pt hard weight with trials;p_{T,hard}",350,-.5,349.5);
667 fHistList->Add(fh1PtHardTrials);
668
67ebd013 669 TH1::AddDirectory(oldStatus);
670
11245558 671 PostData(0,fHistList);
672 PostData(1,fCFManagerPos->GetParticleContainer());
673 PostData(2,fCFManagerNeg->GetParticleContainer());
674
fdceab34 675}
676
677#endif