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 // Steering class for particle (gamma, hadron) identification and correlation
18 // analysis. It is called by the task class AliAnalysisTaskCaloTrackCorrelation
19 // and it connects the input (ESD/AOD/MonteCarlo) got with AliCaloTrackReader
20 // (produces TClonesArrays of AODs (TParticles in MC case if requested)), with
21 // the analysis classes that derive from AliAnaCaloTrackCorrBaseClass
23 // -- Author: Gustavo Conesa (INFN-LNF, LPSC-Grenoble)
27 // --- ROOT system ---
28 #include "TClonesArray.h"
31 //#include <TObjectTable.h>
33 //---- AliRoot system ----
34 #include "AliAnalysisManager.h"
35 #include "AliInputEventHandler.h"
36 #include "AliESDEvent.h"
37 #include "AliAODEvent.h"
38 #include "AliAnaCaloTrackCorrBaseClass.h"
39 #include "AliAnaCaloTrackCorrMaker.h"
41 ClassImp(AliAnaCaloTrackCorrMaker)
44 //__________________________________________________
45 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker() :
47 fReader(0), fCaloUtils(0),
48 fOutputContainer(new TList ), fAnalysisContainer(new TList ),
49 fMakeHisto(kFALSE), fMakeAOD(kFALSE),
50 fAnaDebug(0), fCuts(new TList),
52 fhNEvents(0), fhNExoticEvents(0),
53 fhNEventsNoTriggerFound(0),
54 fhNPileUpEvents(0), fhNPileUpEventsTriggerBC0(0),
56 fhPileUpClusterMult(0), fhPileUpClusterMultAndSPDPileUp(0),
58 fhCentrality(0), fhEventPlaneAngle(0),
59 fhNMergedFiles(0), fhScaleFactor(0),
60 fhEMCalBCEvent(0), fhEMCalBCEventCut(0),
61 fhTrackBCEvent(0), fhTrackBCEventCut(0),
62 fhPrimaryVertexBC(0), fhTimeStampFraction(0),
63 fhNPileUpVertSPD(0), fhNPileUpVertTracks(0),
65 fhClusterTriggerBC(0), fhClusterTriggerBCExotic(0),
66 fhClusterTriggerBCBadCell(0), fhClusterTriggerBCBadCellExotic(0),
67 fhClusterTriggerBCBadCluster(0), fhClusterTriggerBCBadClusterExotic(0),
68 fhClusterTriggerBCUnMatch(0), fhClusterTriggerBCExoticUnMatch(0),
69 fhClusterTriggerBCBadCellUnMatch(0), fhClusterTriggerBCBadCellExoticUnMatch(0),
70 fhClusterTriggerBCBadClusterUnMatch(0), fhClusterTriggerBCBadClusterExoticUnMatch(0),
71 fhClusterTriggerBCEventBC(0), fhClusterTriggerBCEventBCUnMatch(0),
72 fhClusterTriggerBCExoticEventBC(0), fhClusterTriggerBCExoticEventBCUnMatch(0)
75 if(fAnaDebug > 1 ) printf("*** Analysis Maker Constructor *** \n");
77 //Initialize parameters, pointers and histograms
81 //________________________________________________________________________________________
82 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker(const AliAnaCaloTrackCorrMaker & maker) :
84 fReader(), //(new AliCaloTrackReader(*maker.fReader)),
85 fCaloUtils(),//(new AliCalorimeterUtils(*maker.fCaloUtils)),
86 fOutputContainer(new TList()), fAnalysisContainer(new TList()),
87 fMakeHisto(maker.fMakeHisto), fMakeAOD(maker.fMakeAOD),
88 fAnaDebug(maker.fAnaDebug), fCuts(new TList()),
89 fScaleFactor(maker.fScaleFactor),
90 fhNEvents(maker.fhNEvents),
91 fhNExoticEvents(maker.fhNExoticEvents),
92 fhNEventsNoTriggerFound(maker.fhNEventsNoTriggerFound),
93 fhNPileUpEvents(maker.fhNPileUpEvents),
94 fhNPileUpEventsTriggerBC0(maker.fhNPileUpEventsTriggerBC0),
95 fhZVertex(maker.fhZVertex),
96 fhPileUpClusterMult(maker.fhPileUpClusterMult),
97 fhPileUpClusterMultAndSPDPileUp(maker.fhPileUpClusterMultAndSPDPileUp),
98 fhTrackMult(maker.fhTrackMult),
99 fhCentrality(maker.fhCentrality),
100 fhEventPlaneAngle(maker.fhEventPlaneAngle),
101 fhNMergedFiles(maker.fhNMergedFiles),
102 fhScaleFactor(maker.fhScaleFactor),
103 fhEMCalBCEvent(maker.fhEMCalBCEvent),
104 fhEMCalBCEventCut(maker.fhEMCalBCEventCut),
105 fhTrackBCEvent(maker.fhTrackBCEvent),
106 fhTrackBCEventCut(maker.fhTrackBCEventCut),
107 fhPrimaryVertexBC(maker.fhPrimaryVertexBC),
108 fhTimeStampFraction(maker.fhTimeStampFraction),
109 fhNPileUpVertSPD(maker.fhNPileUpVertSPD),
110 fhNPileUpVertTracks(maker.fhNPileUpVertTracks),
111 fhClusterTriggerBC(maker.fhClusterTriggerBC),
112 fhClusterTriggerBCExotic(maker.fhClusterTriggerBCExotic),
113 fhClusterTriggerBCBadCell(maker.fhClusterTriggerBCBadCell),
114 fhClusterTriggerBCBadCellExotic(maker.fhClusterTriggerBCBadCellExotic),
115 fhClusterTriggerBCBadCluster(maker.fhClusterTriggerBCBadCluster),
116 fhClusterTriggerBCBadClusterExotic(maker.fhClusterTriggerBCBadClusterExotic),
117 fhClusterTriggerBCUnMatch(maker.fhClusterTriggerBCUnMatch),
118 fhClusterTriggerBCExoticUnMatch(maker.fhClusterTriggerBCExoticUnMatch),
119 fhClusterTriggerBCBadCellUnMatch(maker.fhClusterTriggerBCBadCellUnMatch),
120 fhClusterTriggerBCBadCellExoticUnMatch(maker.fhClusterTriggerBCBadCellExoticUnMatch),
121 fhClusterTriggerBCBadClusterUnMatch(maker.fhClusterTriggerBCBadClusterUnMatch),
122 fhClusterTriggerBCBadClusterExoticUnMatch(maker.fhClusterTriggerBCBadClusterExoticUnMatch),
123 fhClusterTriggerBCEventBC(maker.fhClusterTriggerBCEventBC),
124 fhClusterTriggerBCEventBCUnMatch(maker.fhClusterTriggerBCEventBCUnMatch),
125 fhClusterTriggerBCExoticEventBC(maker.fhClusterTriggerBCExoticEventBC),
126 fhClusterTriggerBCExoticEventBCUnMatch(maker.fhClusterTriggerBCExoticEventBCUnMatch)
132 //___________________________________________________
133 AliAnaCaloTrackCorrMaker::~AliAnaCaloTrackCorrMaker()
135 // Remove all owned pointers.
137 // Do not delete it here, already done somewhere else, need to understand where.
138 // if (fOutputContainer) {
139 // fOutputContainer->Clear();
140 // delete fOutputContainer ;
143 if (fAnalysisContainer)
145 fAnalysisContainer->Delete();
146 delete fAnalysisContainer ;
149 if (fReader) delete fReader ;
150 if (fCaloUtils) delete fCaloUtils ;
160 //__________________________________________________________________
161 void AliAnaCaloTrackCorrMaker::AddAnalysis(TObject* ana, Int_t n)
163 // Add analysis depending on AliAnaCaloTrackCorrBaseClass to list
165 if ( fAnalysisContainer)
167 fAnalysisContainer->AddAt(ana,n);
171 printf("AliAnaCaloTrackCorrMaker::AddAnalysis() - AnalysisContainer not initialized\n");
176 //_________________________________________________________
177 TList * AliAnaCaloTrackCorrMaker::FillAndGetAODBranchList()
180 // Get any new output AOD branches from analysis and put them in a list
181 // The list is filled in the maker, and new branch passed to the analysis frame
182 // AliAnalysisTaskCaloTrackCorrelation
184 TList *aodBranchList = fReader->GetAODBranchList() ;
186 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
188 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
189 if(ana->NewOutputAOD()) aodBranchList->Add(ana->GetCreateOutputAODBranch());
192 return aodBranchList ;
196 //_________________________________________________________
197 void AliAnaCaloTrackCorrMaker::FillControlHistograms()
199 // Event control histograms
201 AliVEvent* event = fReader->GetInputEvent();
202 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
203 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
205 fhNEvents ->Fill(0); // Number of events analyzed
207 if( fReader->IsPileUpFromSPD())
208 fhNPileUpEvents->Fill(0.5);
209 //if( event->IsPileupFromSPDInMultBins())
210 // fhNPileUpEvents->Fill(1.5);
211 if( fReader->IsPileUpFromEMCal())
212 fhNPileUpEvents->Fill(2.5);
213 if( fReader->IsPileUpFromSPDOrEMCal() )
214 fhNPileUpEvents->Fill(3.5);
215 if( fReader->IsPileUpFromSPDAndEMCal() )
216 fhNPileUpEvents->Fill(4.5);
217 if( fReader->IsPileUpFromSPDAndNotEMCal() )
218 fhNPileUpEvents->Fill(5.5);
219 if( fReader->IsPileUpFromEMCalAndNotSPD() )
220 fhNPileUpEvents->Fill(6.5);
221 if( fReader->IsPileUpFromNotSPDAndNotEMCal() )
222 fhNPileUpEvents->Fill(7.5);
224 Int_t triggerBC = fReader->GetTriggerClusterBC() ;
225 if( triggerBC == 0 &&
226 !fReader->IsExoticEvent() &&
227 !fReader->IsBadCellTriggerEvent())
229 if( fReader->IsPileUpFromSPD())
230 fhNPileUpEventsTriggerBC0->Fill(0.5);
231 //if( event->IsPileupFromSPDInMultBins())
232 // fhNPileUpEventsTriggerBC0->Fill(1.5);
233 if( fReader->IsPileUpFromEMCal())
234 fhNPileUpEventsTriggerBC0->Fill(2.5);
235 if( fReader->IsPileUpFromSPDOrEMCal() )
236 fhNPileUpEventsTriggerBC0->Fill(3.5);
237 if( fReader->IsPileUpFromSPDAndEMCal() )
238 fhNPileUpEventsTriggerBC0->Fill(4.5);
239 if( fReader->IsPileUpFromSPDAndNotEMCal() )
240 fhNPileUpEventsTriggerBC0->Fill(5.5);
241 if( fReader->IsPileUpFromEMCalAndNotSPD() )
242 fhNPileUpEventsTriggerBC0->Fill(6.5);
243 if( fReader->IsPileUpFromNotSPDAndNotEMCal() )
244 fhNPileUpEventsTriggerBC0->Fill(7.5);
247 if(fReader->IsPileUpFromSPD())
248 fhPileUpClusterMultAndSPDPileUp ->Fill(fReader->GetNPileUpClusters());
250 fhPileUpClusterMult ->Fill(fReader->GetNPileUpClusters ());
251 fhTrackMult ->Fill(fReader->GetTrackMultiplicity());
252 fhCentrality ->Fill(fReader->GetEventCentrality ());
253 fhEventPlaneAngle ->Fill(fReader->GetEventPlaneAngle ());
255 for(Int_t i = 0; i < 19; i++)
257 if(fReader->GetTrackEventBC(i)) fhTrackBCEvent ->Fill(i);
258 if(fReader->GetTrackEventBCcut(i))fhTrackBCEventCut->Fill(i);
259 if(fReader->GetEMCalEventBC(i)) fhEMCalBCEvent ->Fill(i);
260 if(fReader->GetEMCalEventBCcut(i))fhEMCalBCEventCut->Fill(i);
264 event->GetPrimaryVertex()->GetXYZ(v) ;
265 fhZVertex->Fill(v[2]);
267 Int_t bc = fReader->GetVertexBC();
268 if(bc!=AliVTrack::kTOFBCNA)fhPrimaryVertexBC->Fill(bc);
271 // N pile up vertices
272 Int_t nVerticesSPD = -1;
273 Int_t nVerticesTracks = -1;
277 nVerticesSPD = esdevent->GetNumberOfPileupVerticesSPD();
278 nVerticesTracks = esdevent->GetNumberOfPileupVerticesTracks();
283 nVerticesSPD = aodevent->GetNumberOfPileupVerticesSPD();
284 nVerticesTracks = aodevent->GetNumberOfPileupVerticesTracks();
287 fhNPileUpVertSPD ->Fill(nVerticesSPD);
288 fhNPileUpVertTracks->Fill(nVerticesTracks);
291 if(fReader->IsSelectEventTimeStampOn() && esdevent)
293 Int_t timeStamp = esdevent->GetTimeStamp();
294 Float_t timeStampFrac = 1.*(timeStamp-fReader->GetRunTimeStampMin()) /
295 (fReader->GetRunTimeStampMax()-fReader->GetRunTimeStampMin());
297 //printf("stamp %d, min %d, max %d, frac %f\n", timeStamp, fReader->GetRunTimeStampMin(), fReader->GetRunTimeStampMax(), timeStampFrac);
299 fhTimeStampFraction->Fill(timeStampFrac);
303 //_______________________________________________________
304 TList * AliAnaCaloTrackCorrMaker::GetListOfAnalysisCuts()
307 // Get the list of the cuts used for the analysis
308 // The list is filled in the maker, called by the task in LocalInit() and posted there
310 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
312 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
313 TObjString * objstring = ana->GetAnalysisCuts();
315 if(objstring)fCuts->Add(objstring);
322 //___________________________________________________
323 TList *AliAnaCaloTrackCorrMaker::GetOutputContainer()
325 // Fill the output list of histograms during the CreateOutputObjects stage.
327 //Initialize calorimeters geometry pointers
328 //GetCaloUtils()->InitPHOSGeometry();
329 //GetCaloUtils()->InitEMCALGeometry();
331 //General event histograms
333 fhNEvents = new TH1F("hNEvents", "Number of analyzed events" , 1 , 0 , 1 ) ;
334 fhNEvents->SetYTitle("# events");
335 fOutputContainer->Add(fhNEvents);
337 fhNExoticEvents = new TH1F("hNExoticEvents", "Number of analyzed events triggered by exotic cluster" , 1 , 0 , 1 ) ;
338 fhNExoticEvents->SetYTitle("# exotic events");
339 fOutputContainer->Add(fhNExoticEvents);
341 fhNEventsNoTriggerFound = new TH1F("hNEventsNoTriggerFound", "Number of analyzed events triggered but no trigger found" , 1 , 0 , 1 ) ;
342 fhNEventsNoTriggerFound->SetYTitle("# exotic events");
343 fOutputContainer->Add(fhNEventsNoTriggerFound);
347 Float_t minbin =-5.5;
348 Float_t maxbin = 5.5;
349 Int_t labelshift = 6;
351 fhClusterTriggerBCEventBC = new TH2F("hClusterTriggerBCEventBC", "Found trigger BC and Event BC",
352 nbin , minbin ,maxbin,4,0, 4) ;
353 fhClusterTriggerBCEventBC->SetXTitle("cluster trigger BC");
354 for(Int_t i = 1; i < 5; i++)
355 fhClusterTriggerBCEventBC->GetYaxis()->SetBinLabel(i ,Form("BC/4=%d",i));
356 fhClusterTriggerBCEventBC->SetXTitle("cluster trigger BC");
357 for(Int_t i = 1; i < 12; i++)
358 fhClusterTriggerBCEventBC->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
359 fhClusterTriggerBCEventBC->SetYTitle("Event BC%4");
360 fOutputContainer->Add(fhClusterTriggerBCEventBC);
362 fhClusterTriggerBCExoticEventBC = new TH2F("hClusterTriggerBCExoticEventBC", "Found exotic trigger BC and Event BC",
363 nbin , minbin ,maxbin,4,1, 4) ;
364 for(Int_t i = 1; i < 5; i++)
365 fhClusterTriggerBCExoticEventBC->GetYaxis()->SetBinLabel(i ,Form("BC/4=%d",i));
366 fhClusterTriggerBCExoticEventBC->SetXTitle("cluster trigger BC");
367 for(Int_t i = 1; i < 12; i++)
368 fhClusterTriggerBCExoticEventBC->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
369 fhClusterTriggerBCExoticEventBC->SetYTitle("Event BC%4");
370 fOutputContainer->Add(fhClusterTriggerBCExoticEventBC);
372 fhClusterTriggerBCEventBCUnMatch = new TH2F("hClusterTriggerBCEventBCUnMatch", "Found unmatched trigger BC and Event BC",
373 nbin , minbin ,maxbin,4,1, 4) ;
374 for(Int_t i = 1; i < 5; i++)
375 fhClusterTriggerBCEventBCUnMatch->GetYaxis()->SetBinLabel(i ,Form("BC/4=%d",i));
376 fhClusterTriggerBCEventBCUnMatch->SetXTitle("cluster trigger BC");
377 for(Int_t i = 1; i < 12; i++)
378 fhClusterTriggerBCEventBCUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
379 fhClusterTriggerBCEventBCUnMatch->SetYTitle("Event BC%4");
380 fOutputContainer->Add(fhClusterTriggerBCEventBCUnMatch);
382 fhClusterTriggerBCExoticEventBCUnMatch = new TH2F("hClusterTriggerExoticBCEventBCUnMatch", "Found unmatched trigger BC and Event BC",
383 nbin , minbin ,maxbin,4,1, 4) ;
384 for(Int_t i = 1; i < 5; i++)
385 fhClusterTriggerBCExoticEventBCUnMatch->GetYaxis()->SetBinLabel(i ,Form("BC/4=%d",i));
386 fhClusterTriggerBCExoticEventBCUnMatch->SetXTitle("cluster trigger BC");
387 for(Int_t i = 1; i < 12; i++)
388 fhClusterTriggerBCExoticEventBCUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
389 fhClusterTriggerBCExoticEventBCUnMatch->SetYTitle("Event BC%4");
390 fOutputContainer->Add(fhClusterTriggerBCExoticEventBCUnMatch);
392 fhClusterTriggerBC = new TH1F("hClusterTriggerBC",
393 "Number of analyzed events triggered by a cluster in a given BC",
394 nbin , minbin ,maxbin) ;
395 fhClusterTriggerBC->SetYTitle("# events");
396 for(Int_t i = 1; i < 12; i++)
397 fhClusterTriggerBC->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
398 fOutputContainer->Add(fhClusterTriggerBC);
400 fhClusterTriggerBCExotic = new TH1F("hClusterTriggerBCExotic",
401 "Number of analyzed events triggered by a exotic cluster in a given BC",
402 nbin , minbin ,maxbin) ;
403 fhClusterTriggerBCExotic->SetYTitle("# events");
404 for(Int_t i = 1; i < 12; i++)
405 fhClusterTriggerBCExotic->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
406 fOutputContainer->Add(fhClusterTriggerBCExotic);
409 fhClusterTriggerBCBadCell = new TH1F("hClusterTriggerBCBadCell",
410 "Number of analyzed events triggered by a bad cell in a given BC",
411 nbin , minbin ,maxbin) ;
413 fhClusterTriggerBCBadCell->SetYTitle("# events");
414 for(Int_t i = 1; i < 12; i++)
415 fhClusterTriggerBCBadCell->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
416 fOutputContainer->Add(fhClusterTriggerBCBadCell);
418 fhClusterTriggerBCBadCellExotic = new TH1F("hClusterTriggerBCBadCellExotic",
419 "Number of analyzed events triggered by a bad cell & exotic cluster in a given BC",
420 nbin , minbin ,maxbin) ;
421 fhClusterTriggerBCBadCellExotic->SetYTitle("# events");
422 for(Int_t i = 1; i < 12; i++)
423 fhClusterTriggerBCBadCellExotic->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
424 fOutputContainer->Add(fhClusterTriggerBCBadCellExotic);
426 fhClusterTriggerBCBadCluster = new TH1F("hClusterTriggerBCBadCluster",
427 "Number of analyzed events triggered by a bad cluster in a given BC",
428 nbin , minbin ,maxbin) ;
430 fhClusterTriggerBCBadCluster->SetYTitle("# events");
431 for(Int_t i = 1; i < 12; i++)
432 fhClusterTriggerBCBadCluster->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
433 fOutputContainer->Add(fhClusterTriggerBCBadCluster);
436 fhClusterTriggerBCBadClusterExotic = new TH1F("hClusterTriggerBCBadClusterExotic",
437 "Number of analyzed events triggered by a bad cluster & exotic cluster in a given BC",
438 nbin , minbin ,maxbin) ;
440 fhClusterTriggerBCBadClusterExotic->SetYTitle("# events");
441 for(Int_t i = 1; i < 12; i++)
442 fhClusterTriggerBCBadClusterExotic->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
443 fOutputContainer->Add(fhClusterTriggerBCBadClusterExotic);
445 fhClusterTriggerBCUnMatch = new TH1F("hClusterTriggerBCUnMatch",
446 "Number of analyzed events triggered by a cluster (no trigger patch match) in a given BC",
447 nbin , minbin ,maxbin) ;
448 fhClusterTriggerBCUnMatch->SetYTitle("# events");
449 for(Int_t i = 1; i < 12; i++)
450 fhClusterTriggerBCUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
451 fOutputContainer->Add(fhClusterTriggerBCUnMatch);
453 fhClusterTriggerBCExoticUnMatch = new TH1F("hClusterTriggerBCExoticUnMatch",
454 "Number of analyzed events triggered by a exotic cluster (no trigger patch match) in a given BC",
455 nbin , minbin ,maxbin) ;
456 fhClusterTriggerBCExoticUnMatch->SetYTitle("# events");
457 for(Int_t i = 1; i < 12; i++)
458 fhClusterTriggerBCExoticUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
459 fOutputContainer->Add(fhClusterTriggerBCExoticUnMatch);
462 fhClusterTriggerBCBadCellUnMatch = new TH1F("hClusterTriggerBCBadCellUnMatch",
463 "Number of analyzed events triggered by a bad cluster (no trigger patch match) in a given BC",
464 nbin , minbin ,maxbin) ;
465 fhClusterTriggerBCBadCellUnMatch->SetYTitle("# events");
466 for(Int_t i = 1; i < 12; i++)
467 fhClusterTriggerBCBadCellUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
468 fOutputContainer->Add(fhClusterTriggerBCBadCellUnMatch);
471 fhClusterTriggerBCBadCellExoticUnMatch = new TH1F("hClusterTriggerBCBadCellExoticUnMatch",
472 "Number of analyzed events triggered by a bad&exotic cluster (no trigger patch match) in a given BC",
473 nbin , minbin ,maxbin) ;
474 fhClusterTriggerBCBadCellExoticUnMatch->SetYTitle("# events");
475 for(Int_t i = 1; i < 12; i++)
476 fhClusterTriggerBCBadCellExoticUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
477 fOutputContainer->Add(fhClusterTriggerBCBadCellExoticUnMatch);
480 fhClusterTriggerBCBadClusterUnMatch = new TH1F("hClusterTriggerBCBadClusterUnMatch",
481 "Number of analyzed events triggered by a bad cluster (no trigger patch match) in a given BC",
482 nbin , minbin ,maxbin) ;
483 fhClusterTriggerBCBadClusterUnMatch->SetYTitle("# events");
484 for(Int_t i = 1; i < 12; i++)
485 fhClusterTriggerBCBadClusterUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
486 fOutputContainer->Add(fhClusterTriggerBCBadClusterUnMatch);
489 fhClusterTriggerBCBadClusterExoticUnMatch = new TH1F("hClusterTriggerBCBadClusterExoticUnMatch",
490 "Number of analyzed events triggered by a bad&exotic cluster (no trigger patch match) in a given BC",
491 nbin , minbin ,maxbin) ;
492 fhClusterTriggerBCBadClusterExoticUnMatch->SetYTitle("# events");
493 for(Int_t i = 1; i < 12; i++)
494 fhClusterTriggerBCBadClusterExoticUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
495 fOutputContainer->Add(fhClusterTriggerBCBadClusterExoticUnMatch);
498 fhNPileUpEvents = new TH1F("hNPileUpEvents", "Number of events considered as pile-up", 8 , 0 , 8 ) ;
499 fhNPileUpEvents->SetYTitle("# events");
500 fhNPileUpEvents->GetXaxis()->SetBinLabel(1 ,"SPD");
501 fhNPileUpEvents->GetXaxis()->SetBinLabel(2 ,"Multi SPD");
502 fhNPileUpEvents->GetXaxis()->SetBinLabel(3 ,"EMCal");
503 fhNPileUpEvents->GetXaxis()->SetBinLabel(4 ,"EMCal || SPD");
504 fhNPileUpEvents->GetXaxis()->SetBinLabel(5 ,"EMCal && SPD");
505 fhNPileUpEvents->GetXaxis()->SetBinLabel(6 ,"!EMCal && SPD");
506 fhNPileUpEvents->GetXaxis()->SetBinLabel(7 ,"EMCal && !SPD");
507 fhNPileUpEvents->GetXaxis()->SetBinLabel(8 ,"!EMCal && !SPD");
508 fOutputContainer->Add(fhNPileUpEvents);
510 fhNPileUpEventsTriggerBC0 = new TH1F("hNPileUpEventsTriggerBC0","Number of events considered as pile-up, trigger cluster in BC=0", 8 , 0 , 8 ) ;
511 fhNPileUpEventsTriggerBC0->SetYTitle("# events");
512 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(1 ,"SPD");
513 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(2 ,"Multi SPD");
514 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(3 ,"EMCal");
515 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(4 ,"EMCal || SPD");
516 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(5 ,"EMCal && SPD");
517 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(6 ,"!EMCal && SPD");
518 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(7 ,"EMCal && !SPD");
519 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(8 ,"!EMCal && !SPD");
520 fOutputContainer->Add(fhNPileUpEventsTriggerBC0);
523 fhTrackBCEvent = new TH1F("hTrackBCEvent", "Number of events with at least 1 track in a bunch crossing ", 19 , 0 , 19 ) ;
524 fhTrackBCEvent->SetYTitle("# events");
525 fhTrackBCEvent->SetXTitle("Bunch crossing");
526 for(Int_t i = 1; i < 20; i++)
527 fhTrackBCEvent->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
528 fOutputContainer->Add(fhTrackBCEvent);
530 fhTrackBCEventCut = new TH1F("hTrackBCEventCut", "Number of events with at least 1 track in a bunch crossing ", 19 , 0 , 19 ) ;
531 fhTrackBCEventCut->SetYTitle("# events");
532 fhTrackBCEventCut->SetXTitle("Bunch crossing");
533 for(Int_t i = 1; i < 20; i++)
534 fhTrackBCEventCut->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
535 fOutputContainer->Add(fhTrackBCEventCut);
537 fhPrimaryVertexBC = new TH1F("hPrimaryVertexBC", "Number of primary vertex per bunch crossing ", 41 , -20 , 20 ) ;
538 fhPrimaryVertexBC->SetYTitle("# events");
539 fhPrimaryVertexBC->SetXTitle("Bunch crossing");
540 fOutputContainer->Add(fhPrimaryVertexBC);
542 fhEMCalBCEvent = new TH1F("hEMCalBCEvent", "Number of events with at least 1 cluster in a bunch crossing ", 19 , 0 , 19 ) ;
543 fhEMCalBCEvent->SetYTitle("# events");
544 fhEMCalBCEvent->SetXTitle("Bunch crossing");
545 for(Int_t i = 1; i < 20; i++)
546 fhEMCalBCEvent->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
547 fOutputContainer->Add(fhEMCalBCEvent);
549 fhEMCalBCEventCut = new TH1F("hEMCalBCEventCut", "Number of events with at least 1 cluster in a bunch crossing", 19 , 0 , 19 ) ;
550 fhEMCalBCEventCut->SetYTitle("# events");
551 fhEMCalBCEventCut->SetXTitle("Bunch crossing");
552 for(Int_t i = 1; i < 20; i++)
553 fhEMCalBCEventCut->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
554 fOutputContainer->Add(fhEMCalBCEventCut);
556 fhZVertex = new TH1F("hZVertex", " Z vertex distribution" , 200 , -50 , 50 ) ;
557 fhZVertex->SetXTitle("v_{z} (cm)");
558 fOutputContainer->Add(fhZVertex);
560 fhTrackMult = new TH1F("hTrackMult", "Number of tracks per events" , 2000 , 0 , 2000 ) ;
561 fhTrackMult->SetXTitle("# tracks");
562 fOutputContainer->Add(fhTrackMult);
564 fhPileUpClusterMult = new TH1F("hPileUpClusterMult", "Number of clusters per event with large time (|t| > 20 ns)" , 100 , 0 , 100 ) ;
565 fhPileUpClusterMult->SetXTitle("# clusters");
566 fOutputContainer->Add(fhPileUpClusterMult);
568 fhPileUpClusterMultAndSPDPileUp = new TH1F("hPileUpClusterMultAndSPDPileUp", "Number of clusters per event with large time (|t| > 20 ns, events tagged as pile-up by SPD)" , 100 , 0 , 100 ) ;
569 fhPileUpClusterMultAndSPDPileUp->SetXTitle("# clusters");
570 fOutputContainer->Add(fhPileUpClusterMultAndSPDPileUp);
572 fhNPileUpVertSPD = new TH1F ("hNPileUpVertSPD","N pile-up SPD vertex", 50,0,50);
573 fhNPileUpVertSPD->SetYTitle("# vertex ");
574 fOutputContainer->Add(fhNPileUpVertSPD);
576 fhNPileUpVertTracks = new TH1F ("hNPileUpVertTracks","N pile-up Tracks vertex", 50,0,50);
577 fhNPileUpVertTracks->SetYTitle("# vertex ");
578 fOutputContainer->Add(fhNPileUpVertTracks);
580 fhCentrality = new TH1F("hCentrality","Number of events in centrality bin",100,0.,100) ;
581 fhCentrality->SetXTitle("Centrality bin");
582 fOutputContainer->Add(fhCentrality) ;
584 fhEventPlaneAngle=new TH1F("hEventPlaneAngle","Number of events in event plane",100,0.,TMath::Pi()) ;
585 fhEventPlaneAngle->SetXTitle("EP angle (rad)");
586 fOutputContainer->Add(fhEventPlaneAngle) ;
588 if(fReader->IsSelectEventTimeStampOn())
590 fhTimeStampFraction = new TH1F("hTimeStampFraction","Fraction of events within a given time stamp range",150, -1, 2) ;
591 fhTimeStampFraction->SetXTitle("fraction");
592 fOutputContainer->Add(fhTimeStampFraction) ;
597 fhNMergedFiles = new TH1F("hNMergedFiles", "Number of merged output files" , 1 , 0 , 1 ) ;
598 fhNMergedFiles->SetYTitle("# files");
599 fhNMergedFiles->Fill(1); // Fill here with one entry, while merging it will count the rest
600 fOutputContainer->Add(fhNMergedFiles);
602 fhScaleFactor = new TH1F("hScaleFactor", "Number of merged output files" , 1 , 0 , 1 ) ;
603 fhScaleFactor->SetYTitle("scale factor");
604 fhScaleFactor->SetBinContent(1,fScaleFactor); // Fill here
605 fOutputContainer->Add(fhScaleFactor);
608 if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
610 printf("AliAnaCaloTrackCorrMaker::GetOutputContainer() - Analysis job list not initialized!!!\n");
611 return fOutputContainer;
614 const Int_t buffersize = 255;
615 char newname[buffersize];
616 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
619 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
621 if(fMakeHisto) // Analysis with histograms as output on
624 //Fill container with appropriate histograms
625 TList * templist = ana ->GetCreateOutputObjects();
626 templist->SetOwner(kFALSE); //Owner is fOutputContainer.
628 for(Int_t i = 0; i < templist->GetEntries(); i++)
631 //Add only to the histogram name the name of the task
632 if( strcmp((templist->At(i))->ClassName(),"TObjString") )
634 snprintf(newname,buffersize, "%s%s", (ana->GetAddedHistogramsStringToName()).Data(), (templist->At(i))->GetName());
635 //printf("name %s, new name %s\n",(templist->At(i))->GetName(),newname);
636 ((TH1*) templist->At(i))->SetName(newname);
639 //Add histogram to general container
640 fOutputContainer->Add(templist->At(i)) ;
646 }// Analysis with histograms as output on
648 }//Loop on analysis defined
650 return fOutputContainer;
654 //___________________________________
655 void AliAnaCaloTrackCorrMaker::Init()
657 //Init container histograms and other common variables
658 // Fill the output list of histograms during the CreateOutputObjects stage.
662 GetReader()->SetCaloUtils(GetCaloUtils()); // pass the calo utils pointer to the reader
665 if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
667 printf("AliAnaCaloTrackCorrMaker::GetOutputInit() - Analysis job list not initialized!!!\n");
671 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
674 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
676 ana->SetReader(fReader); // Set Reader for each analysis
677 ana->SetCaloUtils(fCaloUtils); // Set CaloUtils for each analysis
681 }//Loop on analysis defined
685 //_____________________________________________
686 void AliAnaCaloTrackCorrMaker::InitParameters()
692 fAnaDebug = 0; // No debugging info displayed by default
696 //______________________________________________________________
697 void AliAnaCaloTrackCorrMaker::Print(const Option_t * opt) const
699 //Print some relevant parameters set for the analysis
704 printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
705 printf("Debug level = %d\n", fAnaDebug ) ;
706 printf("Produce Histo = %d\n", fMakeHisto ) ;
707 printf("Produce AOD = %d\n", fMakeAOD ) ;
708 printf("Number of analysis tasks = %d\n", fAnalysisContainer->GetEntries()) ;
710 if(!strcmp("all",opt))
712 printf("Print analysis Tasks settings :\n") ;
713 for(Int_t iana = 0; iana<fAnalysisContainer->GetEntries(); iana++)
715 ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana))->Print("");
718 printf("Print analysis Reader settings :\n") ;
720 printf("Print analysis Calorimeter Utils settings :\n") ;
721 fCaloUtils->Print("");
727 //_______________________________________________________________________
728 void AliAnaCaloTrackCorrMaker::ProcessEvent(const Int_t iEntry,
729 const char * currentFileName)
731 //Process analysis for this event
733 if(fMakeHisto && !fOutputContainer)
735 printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Histograms not initialized\n");
741 printf("*** AliAnaCaloTrackCorrMaker::ProcessEvent() Event %d *** \n",iEntry);
744 printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Current File Name : %s\n", currentFileName);
745 //printf("fAODBranchList %p, entries %d\n",fAODBranchList,fAODBranchList->GetEntries());
749 //Each event needs an empty branch
750 TList * aodList = fReader->GetAODBranchList();
751 Int_t nAODBranches = aodList->GetEntries();
752 for(Int_t iaod = 0; iaod < nAODBranches; iaod++)
754 TClonesArray *tca = dynamic_cast<TClonesArray*> (aodList->At(iaod));
755 if(tca) tca->Clear("C");
758 //Set geometry matrices before filling arrays, in case recalibration/position calculation etc is needed
759 fCaloUtils->AccessGeometry(fReader->GetInputEvent());
761 //Set the AODB calibration, bad channels etc. parameters at least once
762 fCaloUtils->AccessOADB(fReader->GetInputEvent());
765 //Tell the reader to fill the data in the 3 detector lists
766 Bool_t ok = fReader->FillInputEvent(iEntry, currentFileName);
768 Int_t triggerBC = fReader->GetTriggerClusterBC() ;
769 Bool_t exotic = fReader->IsExoticEvent();
770 Bool_t badCluster = fReader->IsBadCellTriggerEvent();
771 Bool_t badCell = fReader->IsBadMaxCellTriggerEvent();
772 Bool_t triggerMatch= fReader->IsTriggerMatched();
773 Bool_t triggerBCOK = kTRUE;
774 Int_t triggerId = fReader->GetTriggerClusterId() ;
778 //printf("Trigger id %d\n",triggerId);
779 if(triggerId == -2)fhNEventsNoTriggerFound->Fill(0);
780 triggerBCOK = kFALSE;
783 if(exotic) fhNExoticEvents->Fill(0) ;
784 //if(fReader->IsExoticEvent()) printf("Maker: EXOTIC Cluster trigger\n");
788 Int_t eventBC = fReader->GetInputEvent()->GetBunchCrossNumber();
789 if(eventBC%4 < 0 || eventBC%4 > 3 )
790 printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - STRANGE: Trigger BC %d - Event BC %d, modulo4 %d \n",triggerBC,eventBC,eventBC%4);
794 if (!exotic && !badCluster) fhClusterTriggerBC->Fill(triggerBC);
795 else if( exotic && badCluster)
797 fhClusterTriggerBCBadClusterExotic->Fill(triggerBC);
798 if(badCell) fhClusterTriggerBCBadCellExotic->Fill(triggerBC);
800 else if( exotic && !badCluster) fhClusterTriggerBCExotic->Fill(triggerBC);
801 else if( badCluster && !exotic )
803 fhClusterTriggerBCBadCluster ->Fill(triggerBC);
804 if(badCell) fhClusterTriggerBCBadCell->Fill(triggerBC);
807 if(!exotic) fhClusterTriggerBCEventBC ->Fill(triggerBC,eventBC%4);
808 else fhClusterTriggerBCExoticEventBC->Fill(triggerBC,eventBC%4);
812 if (!exotic && !badCluster) fhClusterTriggerBCUnMatch->Fill(triggerBC);
813 else if( exotic && badCluster)
815 fhClusterTriggerBCBadClusterExoticUnMatch->Fill(triggerBC);
816 if(badCell) fhClusterTriggerBCBadCellExoticUnMatch ->Fill(triggerBC);
818 else if( exotic && !badCluster) fhClusterTriggerBCExoticUnMatch->Fill(triggerBC);
819 else if( badCluster && !exotic )
821 fhClusterTriggerBCBadClusterUnMatch->Fill(triggerBC);
822 if(badCell)fhClusterTriggerBCBadCellUnMatch->Fill(triggerBC);
825 if(!exotic) fhClusterTriggerBCEventBCUnMatch ->Fill(triggerBC,eventBC%4);
826 else fhClusterTriggerBCExoticEventBCUnMatch->Fill(triggerBC,eventBC%4);
832 if(fAnaDebug >= 1 )printf("*** Skip event *** %d \n",iEntry);
833 fReader->ResetLists();
837 //Magic line to write events to file
838 if(fReader->WriteDeltaAODToFile())AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kTRUE);
840 //printf(">>>>>>>>>> BEFORE >>>>>>>>>>>\n");
841 //gObjectTable->Print();
843 //Access pointers, and trigger mask check needed in mixing case
844 AliAnalysisManager *manager = AliAnalysisManager::GetAnalysisManager();
845 AliInputEventHandler *inputHandler = dynamic_cast<AliInputEventHandler*>(manager->GetInputEventHandler());
847 UInt_t isMBTrigger = kFALSE;
848 UInt_t isTrigger = kFALSE;
851 isMBTrigger = inputHandler->IsEventSelected() & fReader->GetMixEventTriggerMask();
852 isTrigger = inputHandler->IsEventSelected() & fReader->GetEventTriggerMask();
855 //Loop on analysis algorithms
857 if(fAnaDebug > 0 ) printf("*** Begin analysis *** \n");
859 Int_t nana = fAnalysisContainer->GetEntries() ;
860 for(Int_t iana = 0; iana < nana; iana++)
862 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
864 ana->ConnectInputOutputAODBranches(); //Sets branches for each analysis
866 //Fill pool for mixed event for the analysis that need it
867 if(!fReader->IsEventTriggerAtSEOn() && isMBTrigger)
869 ana->FillEventMixPool();
870 continue; // pool filled do not try to fill AODs or histograms
873 //Make analysis, create aods in aod branch and in some cases fill histograms
874 if(fMakeAOD ) ana->MakeAnalysisFillAOD() ;
876 //Make further analysis with aod branch and fill histograms
877 if(fMakeHisto) ana->MakeAnalysisFillHistograms() ;
881 fReader->ResetLists();
883 // In case of mixing analysis, non triggered events are used,
884 // do not fill control histograms for a non requested triggered event
885 if(!fReader->IsEventTriggerAtSEOn() && !isTrigger)
887 if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis, MB for mixing *** \n");
891 FillControlHistograms();
893 //printf(">>>>>>>>>> AFTER >>>>>>>>>>>\n");
894 //gObjectTable->Print();
896 if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis *** \n");
900 //__________________________________________________________
901 void AliAnaCaloTrackCorrMaker::Terminate(TList * outputList)
903 //Execute Terminate of analysis
904 //Do some final plots.
908 Error("Terminate", "No output list");
912 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
915 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
916 if(ana->MakePlotsOn())ana->Terminate(outputList);
918 }//Loop on analysis defined