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 // Count events with different selections
19 // It produces a histogram with the number of events with 9 bins:
20 // 0: all events (that passed the physics selection if it was on)
21 // 1: same but cross check that event pointer did exist
22 // 2: passes vertex cut
23 // 3: passes track number cut, tracks for eta < 0.8
29 // 9: not pileup from SPD
34 // 14: 10 && 2 && 3 && 5
38 // Author: Gustavo Conesa Balbastre (LPSC)
40 //_________________________________________________________________________
46 #include <TProfile2D.h>
48 #include <TClonesArray.h>
49 #include <TGeoGlobalMagField.h>
50 #include "AliAODHeader.h"
51 #include "AliTriggerAnalysis.h"
52 #include "AliESDEvent.h"
53 #include "AliAODEvent.h"
54 #include "AliESDtrackCuts.h"
55 #include "AliAnalysisManager.h"
56 #include "AliInputEventHandler.h"
58 #include "AliAnalysisTaskCounter.h"
59 ClassImp(AliAnalysisTaskCounter)
61 //______________________________________________________________
62 AliAnalysisTaskCounter::AliAnalysisTaskCounter(const char *name)
63 : AliAnalysisTaskSE(name),
64 fAcceptFastCluster(kTRUE),
66 fTrackMultEtaCut(0.8),
68 fOutputContainer(0x0),
69 fESDtrackCuts(AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()),
70 fTriggerAnalysis (new AliTriggerAnalysis),
73 fhXVertex(0), fhYVertex(0), fhZVertex(0),
74 fhXGoodVertex(0),fhYGoodVertex(0),fhZGoodVertex(0),
75 fhCentrality(0), fhEventPlaneAngle(0),
76 fh1Xsec(0), fh1Trials(0)
79 DefineOutput(1, TList::Class());
82 //______________________________________________
83 AliAnalysisTaskCounter::AliAnalysisTaskCounter()
84 : AliAnalysisTaskSE("DefaultAnalysis_AliAnalysisTaskCounter"),
85 fAcceptFastCluster(kTRUE),
87 fTrackMultEtaCut(0.8),
89 fOutputContainer(0x0),
90 fESDtrackCuts(AliESDtrackCuts::GetStandardITSTPCTrackCuts2010()),
91 fTriggerAnalysis (new AliTriggerAnalysis),
94 fhXVertex(0), fhYVertex(0), fhZVertex(0),
95 fhXGoodVertex(0),fhYGoodVertex(0),fhZGoodVertex(0),
96 fhCentrality(0), fhEventPlaneAngle(0),
97 fh1Xsec(0), fh1Trials(0)
100 DefineOutput(1, TList::Class());
103 //__________________________________________________
104 AliAnalysisTaskCounter::~AliAnalysisTaskCounter()
108 if (AliAnalysisManager::GetAnalysisManager()->IsProofMode()) return;
112 fOutputContainer->Delete() ;
113 delete fOutputContainer ;
116 if(fESDtrackCuts) delete fESDtrackCuts;
117 if(fTriggerAnalysis) delete fTriggerAnalysis;
122 //____________________________________________________
123 void AliAnalysisTaskCounter::UserCreateOutputObjects()
127 fOutputContainer = new TList();
129 fh1Xsec = new TH1F("hXsec","xsec from pyxsec.root",1,0,1);
130 fh1Xsec->GetXaxis()->SetBinLabel(1,"<#sigma>");
131 fOutputContainer->Add(fh1Xsec);
133 fh1Trials = new TH1F("hTrials","trials root file",1,0,1);
134 fh1Trials->GetXaxis()->SetBinLabel(1,"#sum{ntrials}");
135 fOutputContainer->Add(fh1Trials);
137 fhZVertex = new TH1F("hZVertex", " Z vertex distribution" , 200 , -50 , 50 ) ;
138 fhZVertex->SetXTitle("v_{z} (cm)");
139 fOutputContainer->Add(fhZVertex);
141 fhZGoodVertex = new TH1F("hZGoodVertex", " Good Z vertex distribution" , 200 , -50 , 50 ) ;
142 fhZGoodVertex->SetXTitle("v_{z} (cm)");
143 fOutputContainer->Add(fhZGoodVertex);
145 fhXVertex = new TH1F("hXVertex", " X vertex distribution" , 200 , -2 , 2 ) ;
146 fhXVertex->SetXTitle("v_{x} (cm)");
147 fOutputContainer->Add(fhXVertex);
149 fhXGoodVertex = new TH1F("hXGoodVertex", " Good X vertex distribution" , 200 , -2 , 2 ) ;
150 fhXGoodVertex->SetXTitle("v_{x} (cm)");
151 fOutputContainer->Add(fhXGoodVertex);
153 fhYVertex = new TH1F("hYVertex", " Y vertex distribution" , 200 , -2 , 2 ) ;
154 fhYVertex->SetXTitle("v_{y} (cm)");
155 fOutputContainer->Add(fhYVertex);
157 fhYGoodVertex = new TH1F("hYGoodVertex", " Good Y vertex distribution" , 200 , -2 , 2 ) ;
158 fhYGoodVertex->SetXTitle("v_{y} (cm)");
159 fOutputContainer->Add(fhYGoodVertex);
161 fhCentrality = new TH1F("hCentrality","Number of events in centrality bin, |vz|<10 cm, method <V0M> ",100,0.,100.) ;
162 fhCentrality->SetXTitle("Centrality bin");
163 fOutputContainer->Add(fhCentrality) ;
165 fhEventPlaneAngle=new TH1F("hEventPlaneAngle","Number of events in event plane, |vz|<10 cm, method <V0> ",100,0.,TMath::Pi()) ;
166 fhEventPlaneAngle->SetXTitle("EP angle (rad)");
167 fOutputContainer->Add(fhEventPlaneAngle) ;
169 fhNEvents = new TH1I("hNEvents", "Number of analyzed events", 21, 0, 21) ;
170 fhNEvents->SetXTitle("Selection");
171 fhNEvents->SetYTitle("# events");
172 fhNEvents->GetXaxis()->SetBinLabel(1 ,"1 = PS");
173 fhNEvents->GetXaxis()->SetBinLabel(2 ,"2 = 1 & ESD");
174 fhNEvents->GetXaxis()->SetBinLabel(3 ,"3 = 2 & |Z|<10");
175 fhNEvents->GetXaxis()->SetBinLabel(4 ,"4 = 2 & |eta|<0.8");
176 fhNEvents->GetXaxis()->SetBinLabel(5 ,"5 = 3 & 4");
177 fhNEvents->GetXaxis()->SetBinLabel(6 ,"6 = 2 & V0AND");
178 fhNEvents->GetXaxis()->SetBinLabel(7 ,"7 = 3 & 6");
179 fhNEvents->GetXaxis()->SetBinLabel(8 ,"8 = 4 & 6");
180 fhNEvents->GetXaxis()->SetBinLabel(9 ,"9 = 5 & 6");
181 fhNEvents->GetXaxis()->SetBinLabel(10,"10 = 2 & not pileup");
182 fhNEvents->GetXaxis()->SetBinLabel(11,"11 = 2 & good vertex");
183 fhNEvents->GetXaxis()->SetBinLabel(12,"12 = 3 & 11");
184 fhNEvents->GetXaxis()->SetBinLabel(13,"13 = 4 & 11");
185 fhNEvents->GetXaxis()->SetBinLabel(14,"14 = 6 & 11");
186 fhNEvents->GetXaxis()->SetBinLabel(15,"15 = 9 & 11");
187 fhNEvents->GetXaxis()->SetBinLabel(16,"16 = 10 & 11");
188 fhNEvents->GetXaxis()->SetBinLabel(17,"17 = 6 & 10");
189 fhNEvents->GetXaxis()->SetBinLabel(17,"17 = 1 & |Z|<50");
190 fhNEvents->GetXaxis()->SetBinLabel(18,"18 = Reject EMCAL 1");
191 fhNEvents->GetXaxis()->SetBinLabel(19,"19 = 18 & 2");
192 fhNEvents->GetXaxis()->SetBinLabel(20,"20 = Reject EMCAL 2");
193 fhNEvents->GetXaxis()->SetBinLabel(21,"20 = 20 & 2");
195 fOutputContainer->Add(fhNEvents);
197 fOutputContainer->SetOwner(kTRUE);
199 PostData(1,fOutputContainer);
204 //_______________________________________________
205 void AliAnalysisTaskCounter::UserExec(Option_t *)
208 // Called for each event
210 //printf("___ Event __ %d __\n",(Int_t)Entry());
214 fhNEvents->Fill(0.5);
216 AliVEvent * event = InputEvent();
217 AliESDEvent * esdevent = dynamic_cast<AliESDEvent*> (event);
218 AliAODEvent * aodevent = dynamic_cast<AliAODEvent*> (event);
220 // Init mag field for tracks in case of ESDs, needed, not clear why
221 if (!TGeoGlobalMagField::Instance()->GetField() && esdevent) esdevent->InitMagneticField();
223 TString triggerclasses = event->GetFiredTriggerClasses();
225 //printf("Trigger class fired: %s \n",event->GetFiredTriggerClasses().Data());
227 if (triggerclasses.Contains("FAST") && !triggerclasses.Contains("ALL") && !fAcceptFastCluster)
229 //printf("Do not count events from fast cluster, trigger name %s\n",triggerclasses.Data());
233 fhNEvents->Fill(1.5);
236 Bool_t bSelectVZ = kFALSE;
237 Bool_t bV0AND = kFALSE;
238 Bool_t bPileup = kFALSE;
239 Bool_t bGoodV = kFALSE;
240 Bool_t bSelectTrack = kFALSE;
243 //---------------------------------
244 //Get the primary vertex, cut on Z
245 //---------------------------------
247 event->GetPrimaryVertex()->GetXYZ(v) ;
248 fhXVertex->Fill(v[0]);
249 fhYVertex->Fill(v[1]);
250 fhZVertex->Fill(v[2]);
252 if(TMath::Abs(v[2]) < fZVertexCut)
255 fhNEvents->Fill(2.5);
257 //else printf("Vertex out %f \n",v[2]);
259 //--------------------------------------------------
260 //Count tracks, cut on number of tracks in eta < 0.8
261 //--------------------------------------------------
262 Int_t nTracks = event->GetNumberOfTracks() ;
263 for (Int_t itrack = 0; itrack < nTracks; itrack++)
264 {////////////// track loop
265 AliVTrack * track = (AliVTrack*)event->GetTrack(itrack) ; // retrieve track from esd
268 if(esdevent && !fESDtrackCuts->AcceptTrack((AliESDtrack*)track)) continue;
271 if(aodevent && !((AliAODTrack*)track)->IsHybridGlobalConstrainedGlobal()) continue ;
273 //Do not count tracks out of acceptance cut
274 if(TMath::Abs(track->Eta())< fTrackMultEtaCut) trackMult++;
277 //printf("AliAnalysisTaskCounter::UserExec() - Track Mult %d \n",trackMult);
279 //--------------------------------------------------
280 // At least one track
281 //--------------------------------------------------
284 bSelectTrack = kTRUE;
285 fhNEvents->Fill(3.5);
286 if(bSelectVZ) fhNEvents->Fill(4.5);
289 //---------------------------------
291 //---------------------------------
293 if(esdevent) bV0AND = fTriggerAnalysis->IsOfflineTriggerFired(esdevent, AliTriggerAnalysis::kV0AND);
297 fhNEvents->Fill(5.5);
298 if (bSelectVZ) fhNEvents->Fill(6.5);
299 if (bSelectTrack) fhNEvents->Fill(7.5);
300 if (bSelectVZ && bSelectTrack) fhNEvents->Fill(8.5);
303 //---------------------------------
305 //---------------------------------
306 bPileup = event->IsPileupFromSPD(3, 0.8, 3., 2., 5.); //Default values, if not it does not compile
307 //bPileup = event->IsPileupFromSPD();
311 fhNEvents->Fill(9.5);
312 if(bV0AND) fhNEvents->Fill(16.5);
315 //---------------------------------
317 //---------------------------------
318 bGoodV = CheckForPrimaryVertex();
320 //Remove events with vertex (0,0,0), bad vertex reconstruction
321 if(TMath::Abs(v[0]) < 1.e-6 &&
322 TMath::Abs(v[1]) < 1.e-6 &&
323 TMath::Abs(v[2]) < 1.e-6) bGoodV = kFALSE;
327 fhXGoodVertex->Fill(v[0]);
328 fhYGoodVertex->Fill(v[1]);
329 fhZGoodVertex->Fill(v[2]);
331 fhNEvents->Fill(10.5);
332 if(bSelectVZ) fhNEvents->Fill(11.5);
333 if(bSelectTrack) fhNEvents->Fill(12.5);
334 if(bV0AND) fhNEvents->Fill(13.5);
335 if(bSelectVZ && bSelectTrack && bV0AND)
336 fhNEvents->Fill(14.5);
337 if(!bPileup) fhNEvents->Fill(15.5);
339 if(TMath::Abs(v[2]) < 10.)
341 if(InputEvent()->GetCentrality())
342 fhCentrality->Fill(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M"));
344 if(InputEvent()->GetEventplane())
346 Float_t ep = InputEvent()->GetEventplane()->GetEventplane("V0", InputEvent());
348 ep+=TMath::Pi()/2.; // put same range as for <Q> method, [0,pi]
350 fhEventPlaneAngle->Fill(ep);
356 //printf("AliAnalysisTaskCounter::UserExec() : z vertex %d, good vertex %d, v0and %d, pile up %d, track mult %d\n ", bSelectVZ, bGoodV, bV0AND, bPileup, trackMult);
358 // Events that could be rejected in EMCAL
359 // LHC11a, SM4 and some SM3 events cut with this
360 Bool_t bEMCALRejected = kFALSE;
361 for (Int_t i = 0; i < InputEvent()->GetNumberOfCaloClusters(); i++)
363 AliVCluster *clus = InputEvent()->GetCaloCluster(i);
365 if ((clus->E() > 500 && clus->GetNCells() > 200 ) || clus->GetNCells() > 200)
368 //printf("Counter: Reject event with cluster: E %f, ncells %d\n",clus->E(),clus->GetNCells());
370 fhNEvents->Fill(17.5);
371 if(bSelectVZ) fhNEvents->Fill(18.5);
372 bEMCALRejected = kTRUE;
378 //LHC11a, 3 last runs, cut with this
381 // Count number of cells in SM3 with energy larger than 0.1, cut on this number
384 for(Int_t icell = 0; icell < event->GetEMCALCells()->GetNumberOfCells(); icell++)
386 if(event->GetEMCALCells()->GetAmplitude(icell) > 0.1 && event->GetEMCALCells()->GetCellNumber(icell)/(24*48)==3) ncellsSM3++;
387 if(event->GetEMCALCells()->GetAmplitude(icell) > 0.1 && event->GetEMCALCells()->GetCellNumber(icell)/(24*48)==4) ncellsSM4++;
391 if(triggerclasses.Contains("EMC")) ncellcut = 35;
393 if( ncellsSM3 >= ncellcut || ncellsSM4 >= 100 )
395 //printf("Counter: reject event with ncells in SM3: ncells %d\n",ncells);
397 fhNEvents->Fill(19.5);
398 if(bSelectVZ) fhNEvents->Fill(20.5);
403 PostData(1,fOutputContainer);
408 //____________________________________________________
409 Bool_t AliAnalysisTaskCounter::CheckForPrimaryVertex()
411 //Check if the vertex was well reconstructed, copy from V0Reader of conversion group
412 //It only works for ESDs
414 AliESDEvent * event = dynamic_cast<AliESDEvent*> (InputEvent());
417 if(event->GetPrimaryVertexTracks()->GetNContributors() > 0)
422 if(event->GetPrimaryVertexTracks()->GetNContributors() < 1)
425 if(event->GetPrimaryVertexSPD()->GetNContributors() > 0)
427 //cout<<"spd vertex type::"<< fESDEvent->GetPrimaryVertex()->GetName() << endl;
431 if(event->GetPrimaryVertexSPD()->GetNContributors() < 1)
433 // cout<<"bad vertex type::"<< fESDEvent->GetPrimaryVertex()->GetName() << endl;
439 //return fInputEvent->GetPrimaryVertex()->GetNContributors()>0;
443 //_____________________________________________________
444 void AliAnalysisTaskCounter::FinishTaskOutput()
446 // Put in the output some event summary histograms
448 AliAnalysisManager *am = AliAnalysisManager::GetAnalysisManager();
449 AliInputEventHandler *inputH = dynamic_cast<AliInputEventHandler*>(am->GetInputEventHandler());
451 TH2F *histStat = dynamic_cast<TH2F*>(inputH->GetStatistics());
452 TH2F *histBin0 = dynamic_cast<TH2F*>(inputH->GetStatistics("BIN0"));
455 fOutputContainer->Add(histStat);
456 else if(DebugLevel() > 1)
457 printf("AliAnalysisTaskCounter::FinishTaskOutput() - Stat histogram not available check, \n if ESDs, that AliPhysicsSelection was on, \n if AODs, if EventStat_temp.root exists \n");
460 fOutputContainer->Add(histBin0);
465 //_____________________________________
466 Bool_t AliAnalysisTaskCounter::Notify()
469 // Implemented Notify() to read the cross sections
470 // and number of trials from pyxsec.root
473 // Fetch the aod also from the input in,
474 // have todo it in notify
476 Float_t xsection = 0;
480 TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
481 if(!tree) return kFALSE;
483 TFile *curfile = tree->GetCurrentFile();
485 if(!curfile) return kFALSE;
487 if(fCurrFileName == curfile->GetName()) return kFALSE;
489 fCurrFileName = TString(curfile->GetName());
491 if(!fh1Xsec||!fh1Trials)
493 Printf("%s%d No Histogram fh1Xsec",(char*)__FILE__,__LINE__);
497 Bool_t ok = PythiaInfoFromFile(fCurrFileName,xsection,trials);
499 if(!ok) return kFALSE;
501 fh1Xsec->Fill("<#sigma>",xsection);
503 // construct a poor man average trials
504 Float_t nEntries = (Float_t)tree->GetTree()->GetEntries();
506 if(trials >= nEntries && nEntries > 0.) fAvgTrials = trials/nEntries;
508 fh1Trials->Fill("#sum{ntrials}",trials);
510 printf("AliAnalysisTaskCounter::Notify() - xs %f, trial %f, avg trials %f\n",xsection,trials, fAvgTrials);
512 if(fDebug) Printf("Reading File %s",fInputHandler->GetTree()->GetCurrentFile()->GetName());
517 //_____________________________________________________________________________________________________
518 Bool_t AliAnalysisTaskCounter::PythiaInfoFromFile(TString file,Float_t & xsec,Float_t & trials)
521 // get the cross section and the trails either from pyxsec.root or from pysec_hists.root
522 // This is to called in Notify and should provide the path to the AOD/ESD file
527 if(file.Contains("root_archive.zip#"))
529 Ssiz_t pos1 = file.Index("root_archive",12,0,TString::kExact);
530 Ssiz_t pos = file.Index("#",1,pos1,TString::kExact);
531 Ssiz_t pos2 = file.Index(".root",5,TString::kExact);
532 file.Replace(pos+1,pos2-pos1,"");
536 // not an archive take the basename....
537 file.ReplaceAll(gSystem->BaseName(file.Data()),"");
540 //Printf("%s",file.Data());
542 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
545 // next trial fetch the histgram file
546 fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
549 // not a severe condition but inciate that we have no information
554 // find the tlist we want to be independtent of the name so use the Tkey
555 TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0);
562 TList *list = dynamic_cast<TList*>(key->ReadObj());
569 xsec = ((TProfile*)list->FindObject("h1Xsec")) ->GetBinContent(1);
570 trials = ((TH1F*) list->FindObject("h1Trials"))->GetBinContent(1);
573 } // no tree pyxsec.root
576 TTree *xtree = (TTree*)fxsec->Get("Xsection");
584 Double_t xsection = 0;
585 xtree->SetBranchAddress("xsection",&xsection);
586 xtree->SetBranchAddress("ntrials",&ntrials);