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),
55 fhXVertex(0), fhYVertex(0), fhZVertex(0),
56 fhXVertexExotic(0), fhYVertexExotic(0), fhZVertexExotic(0),
57 fhPileUpClusterMult(0), fhPileUpClusterMultAndSPDPileUp(0),
59 fhCentrality(0), fhEventPlaneAngle(0),
60 fhNMergedFiles(0), fhScaleFactor(0),
61 fhEMCalBCEvent(0), fhEMCalBCEventCut(0),
62 fhTrackBCEvent(0), fhTrackBCEventCut(0),
63 fhPrimaryVertexBC(0), fhTimeStampFraction(0),
64 fhNPileUpVertSPD(0), fhNPileUpVertTracks(0),
66 fhClusterTriggerBC(0), fhClusterTriggerBCExotic(0),
67 fhClusterTriggerBCBadCell(0), fhClusterTriggerBCBadCellExotic(0),
68 fhClusterTriggerBCBadCluster(0), fhClusterTriggerBCBadClusterExotic(0),
69 fhClusterTriggerBCUnMatch(0), fhClusterTriggerBCExoticUnMatch(0),
70 fhClusterTriggerBCBadCellUnMatch(0), fhClusterTriggerBCBadCellExoticUnMatch(0),
71 fhClusterTriggerBCBadClusterUnMatch(0), fhClusterTriggerBCBadClusterExoticUnMatch(0),
72 fhClusterTriggerBCEventBC(0), fhClusterTriggerBCEventBCUnMatch(0),
73 fhClusterTriggerBCExoticEventBC(0), fhClusterTriggerBCExoticEventBCUnMatch(0)
76 if(fAnaDebug > 1 ) printf("*** Analysis Maker Constructor *** \n");
78 //Initialize parameters, pointers and histograms
82 //________________________________________________________________________________________
83 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker(const AliAnaCaloTrackCorrMaker & maker) :
85 fReader(), //(new AliCaloTrackReader(*maker.fReader)),
86 fCaloUtils(),//(new AliCalorimeterUtils(*maker.fCaloUtils)),
87 fOutputContainer(new TList()), fAnalysisContainer(new TList()),
88 fMakeHisto(maker.fMakeHisto), fMakeAOD(maker.fMakeAOD),
89 fAnaDebug(maker.fAnaDebug), fCuts(new TList()),
90 fScaleFactor(maker.fScaleFactor),
91 fhNEvents(maker.fhNEvents),
92 fhNExoticEvents(maker.fhNExoticEvents),
93 fhNEventsNoTriggerFound(maker.fhNEventsNoTriggerFound),
94 fhNPileUpEvents(maker.fhNPileUpEvents),
95 fhNPileUpEventsTriggerBC0(maker.fhNPileUpEventsTriggerBC0),
96 fhXVertex(maker.fhXVertex),
97 fhYVertex(maker.fhYVertex),
98 fhZVertex(maker.fhZVertex),
99 fhXVertexExotic(maker.fhXVertexExotic),
100 fhYVertexExotic(maker.fhYVertexExotic),
101 fhZVertexExotic(maker.fhZVertexExotic),
102 fhPileUpClusterMult(maker.fhPileUpClusterMult),
103 fhPileUpClusterMultAndSPDPileUp(maker.fhPileUpClusterMultAndSPDPileUp),
104 fhTrackMult(maker.fhTrackMult),
105 fhCentrality(maker.fhCentrality),
106 fhEventPlaneAngle(maker.fhEventPlaneAngle),
107 fhNMergedFiles(maker.fhNMergedFiles),
108 fhScaleFactor(maker.fhScaleFactor),
109 fhEMCalBCEvent(maker.fhEMCalBCEvent),
110 fhEMCalBCEventCut(maker.fhEMCalBCEventCut),
111 fhTrackBCEvent(maker.fhTrackBCEvent),
112 fhTrackBCEventCut(maker.fhTrackBCEventCut),
113 fhPrimaryVertexBC(maker.fhPrimaryVertexBC),
114 fhTimeStampFraction(maker.fhTimeStampFraction),
115 fhNPileUpVertSPD(maker.fhNPileUpVertSPD),
116 fhNPileUpVertTracks(maker.fhNPileUpVertTracks),
117 fhClusterTriggerBC(maker.fhClusterTriggerBC),
118 fhClusterTriggerBCExotic(maker.fhClusterTriggerBCExotic),
119 fhClusterTriggerBCBadCell(maker.fhClusterTriggerBCBadCell),
120 fhClusterTriggerBCBadCellExotic(maker.fhClusterTriggerBCBadCellExotic),
121 fhClusterTriggerBCBadCluster(maker.fhClusterTriggerBCBadCluster),
122 fhClusterTriggerBCBadClusterExotic(maker.fhClusterTriggerBCBadClusterExotic),
123 fhClusterTriggerBCUnMatch(maker.fhClusterTriggerBCUnMatch),
124 fhClusterTriggerBCExoticUnMatch(maker.fhClusterTriggerBCExoticUnMatch),
125 fhClusterTriggerBCBadCellUnMatch(maker.fhClusterTriggerBCBadCellUnMatch),
126 fhClusterTriggerBCBadCellExoticUnMatch(maker.fhClusterTriggerBCBadCellExoticUnMatch),
127 fhClusterTriggerBCBadClusterUnMatch(maker.fhClusterTriggerBCBadClusterUnMatch),
128 fhClusterTriggerBCBadClusterExoticUnMatch(maker.fhClusterTriggerBCBadClusterExoticUnMatch),
129 fhClusterTriggerBCEventBC(maker.fhClusterTriggerBCEventBC),
130 fhClusterTriggerBCEventBCUnMatch(maker.fhClusterTriggerBCEventBCUnMatch),
131 fhClusterTriggerBCExoticEventBC(maker.fhClusterTriggerBCExoticEventBC),
132 fhClusterTriggerBCExoticEventBCUnMatch(maker.fhClusterTriggerBCExoticEventBCUnMatch)
138 //___________________________________________________
139 AliAnaCaloTrackCorrMaker::~AliAnaCaloTrackCorrMaker()
141 // Remove all owned pointers.
143 // Do not delete it here, already done somewhere else, need to understand where.
144 // if (fOutputContainer) {
145 // fOutputContainer->Clear();
146 // delete fOutputContainer ;
149 if (fAnalysisContainer)
151 fAnalysisContainer->Delete();
152 delete fAnalysisContainer ;
155 if (fReader) delete fReader ;
156 if (fCaloUtils) delete fCaloUtils ;
166 //__________________________________________________________________
167 void AliAnaCaloTrackCorrMaker::AddAnalysis(TObject* ana, Int_t n)
169 // Add analysis depending on AliAnaCaloTrackCorrBaseClass to list
171 if ( fAnalysisContainer)
173 fAnalysisContainer->AddAt(ana,n);
177 printf("AliAnaCaloTrackCorrMaker::AddAnalysis() - AnalysisContainer not initialized\n");
182 //_________________________________________________________
183 TList * AliAnaCaloTrackCorrMaker::FillAndGetAODBranchList()
186 // Get any new output AOD branches from analysis and put them in a list
187 // The list is filled in the maker, and new branch passed to the analysis frame
188 // AliAnalysisTaskCaloTrackCorrelation
190 TList *aodBranchList = fReader->GetAODBranchList() ;
192 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
194 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
195 if(ana->NewOutputAOD()) aodBranchList->Add(ana->GetCreateOutputAODBranch());
198 return aodBranchList ;
202 //_________________________________________________________
203 void AliAnaCaloTrackCorrMaker::FillControlHistograms()
205 // Event control histograms
207 AliVEvent* event = fReader->GetInputEvent();
208 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
209 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
211 fhNEvents ->Fill(0); // Number of events analyzed
213 if( fReader->IsPileUpFromSPD())
214 fhNPileUpEvents->Fill(0.5);
215 //if( event->IsPileupFromSPDInMultBins())
216 // fhNPileUpEvents->Fill(1.5);
217 if( fReader->IsPileUpFromEMCal())
218 fhNPileUpEvents->Fill(2.5);
219 if( fReader->IsPileUpFromSPDOrEMCal() )
220 fhNPileUpEvents->Fill(3.5);
221 if( fReader->IsPileUpFromSPDAndEMCal() )
222 fhNPileUpEvents->Fill(4.5);
223 if( fReader->IsPileUpFromSPDAndNotEMCal() )
224 fhNPileUpEvents->Fill(5.5);
225 if( fReader->IsPileUpFromEMCalAndNotSPD() )
226 fhNPileUpEvents->Fill(6.5);
227 if( fReader->IsPileUpFromNotSPDAndNotEMCal() )
228 fhNPileUpEvents->Fill(7.5);
230 Int_t triggerBC = fReader->GetTriggerClusterBC() ;
231 if( triggerBC == 0 &&
232 !fReader->IsExoticEvent() &&
233 !fReader->IsBadCellTriggerEvent())
235 if( fReader->IsPileUpFromSPD())
236 fhNPileUpEventsTriggerBC0->Fill(0.5);
237 //if( event->IsPileupFromSPDInMultBins())
238 // fhNPileUpEventsTriggerBC0->Fill(1.5);
239 if( fReader->IsPileUpFromEMCal())
240 fhNPileUpEventsTriggerBC0->Fill(2.5);
241 if( fReader->IsPileUpFromSPDOrEMCal() )
242 fhNPileUpEventsTriggerBC0->Fill(3.5);
243 if( fReader->IsPileUpFromSPDAndEMCal() )
244 fhNPileUpEventsTriggerBC0->Fill(4.5);
245 if( fReader->IsPileUpFromSPDAndNotEMCal() )
246 fhNPileUpEventsTriggerBC0->Fill(5.5);
247 if( fReader->IsPileUpFromEMCalAndNotSPD() )
248 fhNPileUpEventsTriggerBC0->Fill(6.5);
249 if( fReader->IsPileUpFromNotSPDAndNotEMCal() )
250 fhNPileUpEventsTriggerBC0->Fill(7.5);
253 if(fReader->IsPileUpFromSPD())
254 fhPileUpClusterMultAndSPDPileUp ->Fill(fReader->GetNPileUpClusters());
256 fhPileUpClusterMult ->Fill(fReader->GetNPileUpClusters ());
257 fhTrackMult ->Fill(fReader->GetTrackMultiplicity());
258 fhCentrality ->Fill(fReader->GetEventCentrality ());
259 fhEventPlaneAngle ->Fill(fReader->GetEventPlaneAngle ());
261 for(Int_t i = 0; i < 19; i++)
263 if(fReader->GetTrackEventBC(i)) fhTrackBCEvent ->Fill(i);
264 if(fReader->GetTrackEventBCcut(i))fhTrackBCEventCut->Fill(i);
265 if(fReader->GetEMCalEventBC(i)) fhEMCalBCEvent ->Fill(i);
266 if(fReader->GetEMCalEventBCcut(i))fhEMCalBCEventCut->Fill(i);
270 event->GetPrimaryVertex()->GetXYZ(v) ;
271 fhXVertex->Fill(v[0]);
272 fhYVertex->Fill(v[1]);
273 fhZVertex->Fill(v[2]);
275 Int_t bc = fReader->GetVertexBC();
276 if(bc!=AliVTrack::kTOFBCNA)fhPrimaryVertexBC->Fill(bc);
279 // N pile up vertices
280 Int_t nVerticesSPD = -1;
281 Int_t nVerticesTracks = -1;
285 nVerticesSPD = esdevent->GetNumberOfPileupVerticesSPD();
286 nVerticesTracks = esdevent->GetNumberOfPileupVerticesTracks();
291 nVerticesSPD = aodevent->GetNumberOfPileupVerticesSPD();
292 nVerticesTracks = aodevent->GetNumberOfPileupVerticesTracks();
295 fhNPileUpVertSPD ->Fill(nVerticesSPD);
296 fhNPileUpVertTracks->Fill(nVerticesTracks);
299 if(fReader->IsSelectEventTimeStampOn() && esdevent)
301 Int_t timeStamp = esdevent->GetTimeStamp();
302 Float_t timeStampFrac = 1.*(timeStamp-fReader->GetRunTimeStampMin()) /
303 (fReader->GetRunTimeStampMax()-fReader->GetRunTimeStampMin());
305 //printf("stamp %d, min %d, max %d, frac %f\n", timeStamp, fReader->GetRunTimeStampMin(), fReader->GetRunTimeStampMax(), timeStampFrac);
307 fhTimeStampFraction->Fill(timeStampFrac);
311 //_______________________________________________________
312 TList * AliAnaCaloTrackCorrMaker::GetListOfAnalysisCuts()
315 // Get the list of the cuts used for the analysis
316 // The list is filled in the maker, called by the task in LocalInit() and posted there
318 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
320 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
321 TObjString * objstring = ana->GetAnalysisCuts();
323 if(objstring)fCuts->Add(objstring);
330 //___________________________________________________
331 TList *AliAnaCaloTrackCorrMaker::GetOutputContainer()
333 // Fill the output list of histograms during the CreateOutputObjects stage.
335 //Initialize calorimeters geometry pointers
336 //GetCaloUtils()->InitPHOSGeometry();
337 //GetCaloUtils()->InitEMCALGeometry();
339 //General event histograms
341 fhNEvents = new TH1F("hNEvents", "Number of analyzed events" , 1 , 0 , 1 ) ;
342 fhNEvents->SetYTitle("# events");
343 fOutputContainer->Add(fhNEvents);
345 fhNExoticEvents = new TH1F("hNExoticEvents", "Number of analyzed events triggered by exotic cluster" , 1 , 0 , 1 ) ;
346 fhNExoticEvents->SetYTitle("# exotic events");
347 fOutputContainer->Add(fhNExoticEvents);
349 fhNEventsNoTriggerFound = new TH1F("hNEventsNoTriggerFound", "Number of analyzed events triggered but no trigger found" , 1 , 0 , 1 ) ;
350 fhNEventsNoTriggerFound->SetYTitle("# exotic events");
351 fOutputContainer->Add(fhNEventsNoTriggerFound);
355 Float_t minbin =-5.5;
356 Float_t maxbin = 5.5;
357 Int_t labelshift = 6;
359 fhClusterTriggerBCEventBC = new TH2F("hClusterTriggerBCEventBC", "Found trigger BC and Event BC",
360 nbin , minbin ,maxbin,4,0, 4) ;
361 fhClusterTriggerBCEventBC->SetXTitle("cluster trigger BC");
362 for(Int_t i = 1; i < 5; i++)
363 fhClusterTriggerBCEventBC->GetYaxis()->SetBinLabel(i ,Form("BC/4=%d",i));
364 fhClusterTriggerBCEventBC->SetXTitle("cluster trigger BC");
365 for(Int_t i = 1; i < 12; i++)
366 fhClusterTriggerBCEventBC->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
367 fhClusterTriggerBCEventBC->SetYTitle("Event BC%4");
368 fOutputContainer->Add(fhClusterTriggerBCEventBC);
370 fhClusterTriggerBCExoticEventBC = new TH2F("hClusterTriggerBCExoticEventBC", "Found exotic trigger BC and Event BC",
371 nbin , minbin ,maxbin,4,1, 4) ;
372 for(Int_t i = 1; i < 5; i++)
373 fhClusterTriggerBCExoticEventBC->GetYaxis()->SetBinLabel(i ,Form("BC/4=%d",i));
374 fhClusterTriggerBCExoticEventBC->SetXTitle("cluster trigger BC");
375 for(Int_t i = 1; i < 12; i++)
376 fhClusterTriggerBCExoticEventBC->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
377 fhClusterTriggerBCExoticEventBC->SetYTitle("Event BC%4");
378 fOutputContainer->Add(fhClusterTriggerBCExoticEventBC);
380 fhClusterTriggerBCEventBCUnMatch = new TH2F("hClusterTriggerBCEventBCUnMatch", "Found unmatched trigger BC and Event BC",
381 nbin , minbin ,maxbin,4,1, 4) ;
382 for(Int_t i = 1; i < 5; i++)
383 fhClusterTriggerBCEventBCUnMatch->GetYaxis()->SetBinLabel(i ,Form("BC/4=%d",i));
384 fhClusterTriggerBCEventBCUnMatch->SetXTitle("cluster trigger BC");
385 for(Int_t i = 1; i < 12; i++)
386 fhClusterTriggerBCEventBCUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
387 fhClusterTriggerBCEventBCUnMatch->SetYTitle("Event BC%4");
388 fOutputContainer->Add(fhClusterTriggerBCEventBCUnMatch);
390 fhClusterTriggerBCExoticEventBCUnMatch = new TH2F("hClusterTriggerExoticBCEventBCUnMatch", "Found unmatched trigger BC and Event BC",
391 nbin , minbin ,maxbin,4,1, 4) ;
392 for(Int_t i = 1; i < 5; i++)
393 fhClusterTriggerBCExoticEventBCUnMatch->GetYaxis()->SetBinLabel(i ,Form("BC/4=%d",i));
394 fhClusterTriggerBCExoticEventBCUnMatch->SetXTitle("cluster trigger BC");
395 for(Int_t i = 1; i < 12; i++)
396 fhClusterTriggerBCExoticEventBCUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
397 fhClusterTriggerBCExoticEventBCUnMatch->SetYTitle("Event BC%4");
398 fOutputContainer->Add(fhClusterTriggerBCExoticEventBCUnMatch);
400 fhClusterTriggerBC = new TH1F("hClusterTriggerBC",
401 "Number of analyzed events triggered by a cluster in a given BC",
402 nbin , minbin ,maxbin) ;
403 fhClusterTriggerBC->SetYTitle("# events");
404 for(Int_t i = 1; i < 12; i++)
405 fhClusterTriggerBC->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
406 fOutputContainer->Add(fhClusterTriggerBC);
408 fhClusterTriggerBCExotic = new TH1F("hClusterTriggerBCExotic",
409 "Number of analyzed events triggered by a exotic cluster in a given BC",
410 nbin , minbin ,maxbin) ;
411 fhClusterTriggerBCExotic->SetYTitle("# events");
412 for(Int_t i = 1; i < 12; i++)
413 fhClusterTriggerBCExotic->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
414 fOutputContainer->Add(fhClusterTriggerBCExotic);
417 fhClusterTriggerBCBadCell = new TH1F("hClusterTriggerBCBadCell",
418 "Number of analyzed events triggered by a bad cell in a given BC",
419 nbin , minbin ,maxbin) ;
421 fhClusterTriggerBCBadCell->SetYTitle("# events");
422 for(Int_t i = 1; i < 12; i++)
423 fhClusterTriggerBCBadCell->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
424 fOutputContainer->Add(fhClusterTriggerBCBadCell);
426 fhClusterTriggerBCBadCellExotic = new TH1F("hClusterTriggerBCBadCellExotic",
427 "Number of analyzed events triggered by a bad cell & exotic cluster in a given BC",
428 nbin , minbin ,maxbin) ;
429 fhClusterTriggerBCBadCellExotic->SetYTitle("# events");
430 for(Int_t i = 1; i < 12; i++)
431 fhClusterTriggerBCBadCellExotic->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
432 fOutputContainer->Add(fhClusterTriggerBCBadCellExotic);
434 fhClusterTriggerBCBadCluster = new TH1F("hClusterTriggerBCBadCluster",
435 "Number of analyzed events triggered by a bad cluster in a given BC",
436 nbin , minbin ,maxbin) ;
438 fhClusterTriggerBCBadCluster->SetYTitle("# events");
439 for(Int_t i = 1; i < 12; i++)
440 fhClusterTriggerBCBadCluster->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
441 fOutputContainer->Add(fhClusterTriggerBCBadCluster);
444 fhClusterTriggerBCBadClusterExotic = new TH1F("hClusterTriggerBCBadClusterExotic",
445 "Number of analyzed events triggered by a bad cluster & exotic cluster in a given BC",
446 nbin , minbin ,maxbin) ;
448 fhClusterTriggerBCBadClusterExotic->SetYTitle("# events");
449 for(Int_t i = 1; i < 12; i++)
450 fhClusterTriggerBCBadClusterExotic->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
451 fOutputContainer->Add(fhClusterTriggerBCBadClusterExotic);
453 fhClusterTriggerBCUnMatch = new TH1F("hClusterTriggerBCUnMatch",
454 "Number of analyzed events triggered by a cluster (no trigger patch match) in a given BC",
455 nbin , minbin ,maxbin) ;
456 fhClusterTriggerBCUnMatch->SetYTitle("# events");
457 for(Int_t i = 1; i < 12; i++)
458 fhClusterTriggerBCUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
459 fOutputContainer->Add(fhClusterTriggerBCUnMatch);
461 fhClusterTriggerBCExoticUnMatch = new TH1F("hClusterTriggerBCExoticUnMatch",
462 "Number of analyzed events triggered by a exotic cluster (no trigger patch match) in a given BC",
463 nbin , minbin ,maxbin) ;
464 fhClusterTriggerBCExoticUnMatch->SetYTitle("# events");
465 for(Int_t i = 1; i < 12; i++)
466 fhClusterTriggerBCExoticUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
467 fOutputContainer->Add(fhClusterTriggerBCExoticUnMatch);
470 fhClusterTriggerBCBadCellUnMatch = new TH1F("hClusterTriggerBCBadCellUnMatch",
471 "Number of analyzed events triggered by a bad cluster (no trigger patch match) in a given BC",
472 nbin , minbin ,maxbin) ;
473 fhClusterTriggerBCBadCellUnMatch->SetYTitle("# events");
474 for(Int_t i = 1; i < 12; i++)
475 fhClusterTriggerBCBadCellUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
476 fOutputContainer->Add(fhClusterTriggerBCBadCellUnMatch);
479 fhClusterTriggerBCBadCellExoticUnMatch = new TH1F("hClusterTriggerBCBadCellExoticUnMatch",
480 "Number of analyzed events triggered by a bad&exotic cluster (no trigger patch match) in a given BC",
481 nbin , minbin ,maxbin) ;
482 fhClusterTriggerBCBadCellExoticUnMatch->SetYTitle("# events");
483 for(Int_t i = 1; i < 12; i++)
484 fhClusterTriggerBCBadCellExoticUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
485 fOutputContainer->Add(fhClusterTriggerBCBadCellExoticUnMatch);
488 fhClusterTriggerBCBadClusterUnMatch = new TH1F("hClusterTriggerBCBadClusterUnMatch",
489 "Number of analyzed events triggered by a bad cluster (no trigger patch match) in a given BC",
490 nbin , minbin ,maxbin) ;
491 fhClusterTriggerBCBadClusterUnMatch->SetYTitle("# events");
492 for(Int_t i = 1; i < 12; i++)
493 fhClusterTriggerBCBadClusterUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
494 fOutputContainer->Add(fhClusterTriggerBCBadClusterUnMatch);
497 fhClusterTriggerBCBadClusterExoticUnMatch = new TH1F("hClusterTriggerBCBadClusterExoticUnMatch",
498 "Number of analyzed events triggered by a bad&exotic cluster (no trigger patch match) in a given BC",
499 nbin , minbin ,maxbin) ;
500 fhClusterTriggerBCBadClusterExoticUnMatch->SetYTitle("# events");
501 for(Int_t i = 1; i < 12; i++)
502 fhClusterTriggerBCBadClusterExoticUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
503 fOutputContainer->Add(fhClusterTriggerBCBadClusterExoticUnMatch);
506 fhNPileUpEvents = new TH1F("hNPileUpEvents", "Number of events considered as pile-up", 8 , 0 , 8 ) ;
507 fhNPileUpEvents->SetYTitle("# events");
508 fhNPileUpEvents->GetXaxis()->SetBinLabel(1 ,"SPD");
509 fhNPileUpEvents->GetXaxis()->SetBinLabel(2 ,"Multi SPD");
510 fhNPileUpEvents->GetXaxis()->SetBinLabel(3 ,"EMCal");
511 fhNPileUpEvents->GetXaxis()->SetBinLabel(4 ,"EMCal || SPD");
512 fhNPileUpEvents->GetXaxis()->SetBinLabel(5 ,"EMCal && SPD");
513 fhNPileUpEvents->GetXaxis()->SetBinLabel(6 ,"!EMCal && SPD");
514 fhNPileUpEvents->GetXaxis()->SetBinLabel(7 ,"EMCal && !SPD");
515 fhNPileUpEvents->GetXaxis()->SetBinLabel(8 ,"!EMCal && !SPD");
516 fOutputContainer->Add(fhNPileUpEvents);
518 fhNPileUpEventsTriggerBC0 = new TH1F("hNPileUpEventsTriggerBC0","Number of events considered as pile-up, trigger cluster in BC=0", 8 , 0 , 8 ) ;
519 fhNPileUpEventsTriggerBC0->SetYTitle("# events");
520 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(1 ,"SPD");
521 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(2 ,"Multi SPD");
522 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(3 ,"EMCal");
523 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(4 ,"EMCal || SPD");
524 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(5 ,"EMCal && SPD");
525 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(6 ,"!EMCal && SPD");
526 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(7 ,"EMCal && !SPD");
527 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(8 ,"!EMCal && !SPD");
528 fOutputContainer->Add(fhNPileUpEventsTriggerBC0);
531 fhTrackBCEvent = new TH1F("hTrackBCEvent", "Number of events with at least 1 track in a bunch crossing ", 19 , 0 , 19 ) ;
532 fhTrackBCEvent->SetYTitle("# events");
533 fhTrackBCEvent->SetXTitle("Bunch crossing");
534 for(Int_t i = 1; i < 20; i++)
535 fhTrackBCEvent->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
536 fOutputContainer->Add(fhTrackBCEvent);
538 fhTrackBCEventCut = new TH1F("hTrackBCEventCut", "Number of events with at least 1 track in a bunch crossing ", 19 , 0 , 19 ) ;
539 fhTrackBCEventCut->SetYTitle("# events");
540 fhTrackBCEventCut->SetXTitle("Bunch crossing");
541 for(Int_t i = 1; i < 20; i++)
542 fhTrackBCEventCut->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
543 fOutputContainer->Add(fhTrackBCEventCut);
545 fhPrimaryVertexBC = new TH1F("hPrimaryVertexBC", "Number of primary vertex per bunch crossing ", 41 , -20 , 20 ) ;
546 fhPrimaryVertexBC->SetYTitle("# events");
547 fhPrimaryVertexBC->SetXTitle("Bunch crossing");
548 fOutputContainer->Add(fhPrimaryVertexBC);
550 fhEMCalBCEvent = new TH1F("hEMCalBCEvent", "Number of events with at least 1 cluster in a bunch crossing ", 19 , 0 , 19 ) ;
551 fhEMCalBCEvent->SetYTitle("# events");
552 fhEMCalBCEvent->SetXTitle("Bunch crossing");
553 for(Int_t i = 1; i < 20; i++)
554 fhEMCalBCEvent->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
555 fOutputContainer->Add(fhEMCalBCEvent);
557 fhEMCalBCEventCut = new TH1F("hEMCalBCEventCut", "Number of events with at least 1 cluster in a bunch crossing", 19 , 0 , 19 ) ;
558 fhEMCalBCEventCut->SetYTitle("# events");
559 fhEMCalBCEventCut->SetXTitle("Bunch crossing");
560 for(Int_t i = 1; i < 20; i++)
561 fhEMCalBCEventCut->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
562 fOutputContainer->Add(fhEMCalBCEventCut);
564 fhXVertex = new TH1F("hXVertex", " X vertex distribution" , 200 , -4 , 4 ) ;
565 fhXVertex->SetXTitle("v_{x} (cm)");
566 fOutputContainer->Add(fhXVertex);
568 fhYVertex = new TH1F("hYVertex", " Y vertex distribution" , 200 , -4 , 4 ) ;
569 fhYVertex->SetXTitle("v_{y} (cm)");
570 fOutputContainer->Add(fhYVertex);
572 fhZVertex = new TH1F("hZVertex", " Z vertex distribution" , 200 , -50 , 50 ) ;
573 fhZVertex->SetXTitle("v_{z} (cm)");
574 fOutputContainer->Add(fhZVertex);
576 fhXVertexExotic = new TH1F("hXVertexExotic", " X vertex distribution in exotic events" , 200 , -4 , 4 ) ;
577 fhXVertexExotic->SetXTitle("v_{x} (cm)");
578 fOutputContainer->Add(fhXVertexExotic);
580 fhYVertexExotic = new TH1F("hYVertexExotic", " Y vertex distribution in exotic events" , 200 , -4 , 4 ) ;
581 fhYVertexExotic->SetXTitle("v_{y} (cm)");
582 fOutputContainer->Add(fhYVertexExotic);
584 fhZVertexExotic = new TH1F("hZVertexExotic", " Z vertex distribution in exotic events" , 200 , -50 , 50 ) ;
585 fhZVertexExotic->SetXTitle("v_{z} (cm)");
586 fOutputContainer->Add(fhZVertexExotic);
588 fhTrackMult = new TH1F("hTrackMult", "Number of tracks per events" , 2000 , 0 , 2000 ) ;
589 fhTrackMult->SetXTitle("# tracks");
590 fOutputContainer->Add(fhTrackMult);
592 fhPileUpClusterMult = new TH1F("hPileUpClusterMult", "Number of clusters per event with large time (|t| > 20 ns)" , 100 , 0 , 100 ) ;
593 fhPileUpClusterMult->SetXTitle("# clusters");
594 fOutputContainer->Add(fhPileUpClusterMult);
596 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 ) ;
597 fhPileUpClusterMultAndSPDPileUp->SetXTitle("# clusters");
598 fOutputContainer->Add(fhPileUpClusterMultAndSPDPileUp);
600 fhNPileUpVertSPD = new TH1F ("hNPileUpVertSPD","N pile-up SPD vertex", 50,0,50);
601 fhNPileUpVertSPD->SetYTitle("# vertex ");
602 fOutputContainer->Add(fhNPileUpVertSPD);
604 fhNPileUpVertTracks = new TH1F ("hNPileUpVertTracks","N pile-up Tracks vertex", 50,0,50);
605 fhNPileUpVertTracks->SetYTitle("# vertex ");
606 fOutputContainer->Add(fhNPileUpVertTracks);
608 fhCentrality = new TH1F("hCentrality","Number of events in centrality bin",100,0.,100) ;
609 fhCentrality->SetXTitle("Centrality bin");
610 fOutputContainer->Add(fhCentrality) ;
612 fhEventPlaneAngle=new TH1F("hEventPlaneAngle","Number of events in event plane",100,0.,TMath::Pi()) ;
613 fhEventPlaneAngle->SetXTitle("EP angle (rad)");
614 fOutputContainer->Add(fhEventPlaneAngle) ;
616 if(fReader->IsSelectEventTimeStampOn())
618 fhTimeStampFraction = new TH1F("hTimeStampFraction","Fraction of events within a given time stamp range",150, -1, 2) ;
619 fhTimeStampFraction->SetXTitle("fraction");
620 fOutputContainer->Add(fhTimeStampFraction) ;
625 fhNMergedFiles = new TH1F("hNMergedFiles", "Number of merged output files" , 1 , 0 , 1 ) ;
626 fhNMergedFiles->SetYTitle("# files");
627 fhNMergedFiles->Fill(1); // Fill here with one entry, while merging it will count the rest
628 fOutputContainer->Add(fhNMergedFiles);
630 fhScaleFactor = new TH1F("hScaleFactor", "Number of merged output files" , 1 , 0 , 1 ) ;
631 fhScaleFactor->SetYTitle("scale factor");
632 fhScaleFactor->SetBinContent(1,fScaleFactor); // Fill here
633 fOutputContainer->Add(fhScaleFactor);
636 if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
638 printf("AliAnaCaloTrackCorrMaker::GetOutputContainer() - Analysis job list not initialized!!!\n");
639 return fOutputContainer;
642 const Int_t buffersize = 255;
643 char newname[buffersize];
644 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
647 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
649 if(fMakeHisto) // Analysis with histograms as output on
652 //Fill container with appropriate histograms
653 TList * templist = ana ->GetCreateOutputObjects();
654 templist->SetOwner(kFALSE); //Owner is fOutputContainer.
656 for(Int_t i = 0; i < templist->GetEntries(); i++)
659 //Add only to the histogram name the name of the task
660 if( strcmp((templist->At(i))->ClassName(),"TObjString") )
662 snprintf(newname,buffersize, "%s%s", (ana->GetAddedHistogramsStringToName()).Data(), (templist->At(i))->GetName());
663 //printf("name %s, new name %s\n",(templist->At(i))->GetName(),newname);
664 ((TH1*) templist->At(i))->SetName(newname);
667 //Add histogram to general container
668 fOutputContainer->Add(templist->At(i)) ;
674 }// Analysis with histograms as output on
676 }//Loop on analysis defined
678 return fOutputContainer;
682 //___________________________________
683 void AliAnaCaloTrackCorrMaker::Init()
685 //Init container histograms and other common variables
686 // Fill the output list of histograms during the CreateOutputObjects stage.
690 GetReader()->SetCaloUtils(GetCaloUtils()); // pass the calo utils pointer to the reader
693 if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
695 printf("AliAnaCaloTrackCorrMaker::GetOutputInit() - Analysis job list not initialized!!!\n");
699 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
702 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
704 ana->SetReader(fReader); // Set Reader for each analysis
705 ana->SetCaloUtils(fCaloUtils); // Set CaloUtils for each analysis
709 }//Loop on analysis defined
713 //_____________________________________________
714 void AliAnaCaloTrackCorrMaker::InitParameters()
720 fAnaDebug = 0; // No debugging info displayed by default
724 //______________________________________________________________
725 void AliAnaCaloTrackCorrMaker::Print(const Option_t * opt) const
727 //Print some relevant parameters set for the analysis
732 printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
733 printf("Debug level = %d\n", fAnaDebug ) ;
734 printf("Produce Histo = %d\n", fMakeHisto ) ;
735 printf("Produce AOD = %d\n", fMakeAOD ) ;
736 printf("Number of analysis tasks = %d\n", fAnalysisContainer->GetEntries()) ;
738 if(!strcmp("all",opt))
740 printf("Print analysis Tasks settings :\n") ;
741 for(Int_t iana = 0; iana<fAnalysisContainer->GetEntries(); iana++)
743 ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana))->Print("");
746 printf("Print analysis Reader settings :\n") ;
748 printf("Print analysis Calorimeter Utils settings :\n") ;
749 fCaloUtils->Print("");
755 //_______________________________________________________________________
756 void AliAnaCaloTrackCorrMaker::ProcessEvent(const Int_t iEntry,
757 const char * currentFileName)
759 //Process analysis for this event
761 if(fMakeHisto && !fOutputContainer)
763 printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Histograms not initialized\n");
769 printf("*** AliAnaCaloTrackCorrMaker::ProcessEvent() Event %d *** \n",iEntry);
772 printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Current File Name : %s\n", currentFileName);
773 //printf("fAODBranchList %p, entries %d\n",fAODBranchList,fAODBranchList->GetEntries());
777 //Each event needs an empty branch
778 TList * aodList = fReader->GetAODBranchList();
779 Int_t nAODBranches = aodList->GetEntries();
780 for(Int_t iaod = 0; iaod < nAODBranches; iaod++)
782 TClonesArray *tca = dynamic_cast<TClonesArray*> (aodList->At(iaod));
783 if(tca) tca->Clear("C");
786 //Set geometry matrices before filling arrays, in case recalibration/position calculation etc is needed
787 fCaloUtils->AccessGeometry(fReader->GetInputEvent());
789 //Set the AODB calibration, bad channels etc. parameters at least once
790 fCaloUtils->AccessOADB(fReader->GetInputEvent());
793 //Tell the reader to fill the data in the 3 detector lists
794 Bool_t ok = fReader->FillInputEvent(iEntry, currentFileName);
796 Int_t triggerBC = fReader->GetTriggerClusterBC() ;
797 Bool_t exotic = fReader->IsExoticEvent();
798 Bool_t badCluster = fReader->IsBadCellTriggerEvent();
799 Bool_t badCell = fReader->IsBadMaxCellTriggerEvent();
800 Bool_t triggerMatch= fReader->IsTriggerMatched();
801 Bool_t triggerBCOK = kTRUE;
802 Int_t triggerId = fReader->GetTriggerClusterId() ;
806 //printf("Trigger id %d\n",triggerId);
807 if(triggerId == -2)fhNEventsNoTriggerFound->Fill(0);
808 triggerBCOK = kFALSE;
813 fhNExoticEvents->Fill(0) ;
815 fReader->GetInputEvent()->GetPrimaryVertex()->GetXYZ(v) ;
816 fhXVertexExotic->Fill(v[0]);
817 fhYVertexExotic->Fill(v[1]);
818 fhZVertexExotic->Fill(v[2]);
820 //if(fReader->IsExoticEvent()) printf("Maker: EXOTIC Cluster trigger\n");
824 Int_t eventBC = fReader->GetInputEvent()->GetBunchCrossNumber();
825 if(eventBC%4 < 0 || eventBC%4 > 3 )
826 printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - STRANGE: Trigger BC %d - Event BC %d, modulo4 %d \n",triggerBC,eventBC,eventBC%4);
830 if (!exotic && !badCluster) fhClusterTriggerBC->Fill(triggerBC);
831 else if( exotic && badCluster)
833 fhClusterTriggerBCBadClusterExotic->Fill(triggerBC);
834 if(badCell) fhClusterTriggerBCBadCellExotic->Fill(triggerBC);
836 else if( exotic && !badCluster) fhClusterTriggerBCExotic->Fill(triggerBC);
837 else if( badCluster && !exotic )
839 fhClusterTriggerBCBadCluster ->Fill(triggerBC);
840 if(badCell) fhClusterTriggerBCBadCell->Fill(triggerBC);
843 if(!exotic) fhClusterTriggerBCEventBC ->Fill(triggerBC,eventBC%4);
844 else fhClusterTriggerBCExoticEventBC->Fill(triggerBC,eventBC%4);
848 if (!exotic && !badCluster) fhClusterTriggerBCUnMatch->Fill(triggerBC);
849 else if( exotic && badCluster)
851 fhClusterTriggerBCBadClusterExoticUnMatch->Fill(triggerBC);
852 if(badCell) fhClusterTriggerBCBadCellExoticUnMatch ->Fill(triggerBC);
854 else if( exotic && !badCluster) fhClusterTriggerBCExoticUnMatch->Fill(triggerBC);
855 else if( badCluster && !exotic )
857 fhClusterTriggerBCBadClusterUnMatch->Fill(triggerBC);
858 if(badCell)fhClusterTriggerBCBadCellUnMatch->Fill(triggerBC);
861 if(!exotic) fhClusterTriggerBCEventBCUnMatch ->Fill(triggerBC,eventBC%4);
862 else fhClusterTriggerBCExoticEventBCUnMatch->Fill(triggerBC,eventBC%4);
868 if(fAnaDebug >= 1 )printf("*** Skip event *** %d \n",iEntry);
869 fReader->ResetLists();
873 //Magic line to write events to file
874 if(fReader->WriteDeltaAODToFile())AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kTRUE);
876 //printf(">>>>>>>>>> BEFORE >>>>>>>>>>>\n");
877 //gObjectTable->Print();
879 //Access pointers, and trigger mask check needed in mixing case
880 AliAnalysisManager *manager = AliAnalysisManager::GetAnalysisManager();
881 AliInputEventHandler *inputHandler = dynamic_cast<AliInputEventHandler*>(manager->GetInputEventHandler());
883 UInt_t isMBTrigger = kFALSE;
884 UInt_t isTrigger = kFALSE;
887 isMBTrigger = inputHandler->IsEventSelected() & fReader->GetMixEventTriggerMask();
888 isTrigger = inputHandler->IsEventSelected() & fReader->GetEventTriggerMask();
891 //Loop on analysis algorithms
893 if(fAnaDebug > 0 ) printf("*** Begin analysis *** \n");
895 Int_t nana = fAnalysisContainer->GetEntries() ;
896 for(Int_t iana = 0; iana < nana; iana++)
898 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
900 ana->ConnectInputOutputAODBranches(); //Sets branches for each analysis
902 //Fill pool for mixed event for the analysis that need it
903 if(!fReader->IsEventTriggerAtSEOn() && isMBTrigger)
905 ana->FillEventMixPool();
906 continue; // pool filled do not try to fill AODs or histograms
909 //Make analysis, create aods in aod branch and in some cases fill histograms
910 if(fMakeAOD ) ana->MakeAnalysisFillAOD() ;
912 //Make further analysis with aod branch and fill histograms
913 if(fMakeHisto) ana->MakeAnalysisFillHistograms() ;
917 fReader->ResetLists();
919 // In case of mixing analysis, non triggered events are used,
920 // do not fill control histograms for a non requested triggered event
921 if(!fReader->IsEventTriggerAtSEOn() && !isTrigger)
923 if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis, MB for mixing *** \n");
927 FillControlHistograms();
929 //printf(">>>>>>>>>> AFTER >>>>>>>>>>>\n");
930 //gObjectTable->Print();
932 if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis *** \n");
936 //__________________________________________________________
937 void AliAnaCaloTrackCorrMaker::Terminate(TList * outputList)
939 //Execute Terminate of analysis
940 //Do some final plots.
944 Error("Terminate", "No output list");
948 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
951 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
952 if(ana->MakePlotsOn())ana->Terminate(outputList);
954 }//Loop on analysis defined