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),
64 fhClusterTriggerBC(0), fhClusterTriggerBCExotic(0),
65 fhClusterTriggerBCBad(0), fhClusterTriggerBCBadExotic(0),
66 fhClusterTriggerBCUnMatch(0), fhClusterTriggerBCExoticUnMatch(0),
67 fhClusterTriggerBCBadUnMatch(0), fhClusterTriggerBCBadExoticUnMatch(0)
70 if(fAnaDebug > 1 ) printf("*** Analysis Maker Constructor *** \n");
72 //Initialize parameters, pointers and histograms
76 //________________________________________________________________________________________
77 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker(const AliAnaCaloTrackCorrMaker & maker) :
79 fReader(), //(new AliCaloTrackReader(*maker.fReader)),
80 fCaloUtils(),//(new AliCalorimeterUtils(*maker.fCaloUtils)),
81 fOutputContainer(new TList()), fAnalysisContainer(new TList()),
82 fMakeHisto(maker.fMakeHisto), fMakeAOD(maker.fMakeAOD),
83 fAnaDebug(maker.fAnaDebug), fCuts(new TList()),
84 fScaleFactor(maker.fScaleFactor),
85 fhNEvents(maker.fhNEvents),
86 fhNExoticEvents(maker.fhNExoticEvents),
87 fhNEventsNoTriggerFound(maker.fhNEventsNoTriggerFound),
88 fhNPileUpEvents(maker.fhNPileUpEvents),
89 fhNPileUpEventsTriggerBC0(maker.fhNPileUpEventsTriggerBC0),
90 fhZVertex(maker.fhZVertex),
91 fhPileUpClusterMult(maker.fhPileUpClusterMult),
92 fhPileUpClusterMultAndSPDPileUp(maker.fhPileUpClusterMultAndSPDPileUp),
93 fhTrackMult(maker.fhTrackMult),
94 fhCentrality(maker.fhCentrality),
95 fhEventPlaneAngle(maker.fhEventPlaneAngle),
96 fhNMergedFiles(maker.fhNMergedFiles),
97 fhScaleFactor(maker.fhScaleFactor),
98 fhEMCalBCEvent(maker.fhEMCalBCEvent),
99 fhEMCalBCEventCut(maker.fhEMCalBCEventCut),
100 fhTrackBCEvent(maker.fhTrackBCEvent),
101 fhTrackBCEventCut(maker.fhTrackBCEventCut),
102 fhPrimaryVertexBC(maker.fhPrimaryVertexBC),
103 fhTimeStampFraction(maker.fhTimeStampFraction),
104 fhNPileUpVertSPD(maker.fhNPileUpVertSPD),
105 fhNPileUpVertTracks(maker.fhNPileUpVertTracks),
106 fhClusterTriggerBC(maker.fhClusterTriggerBC),
107 fhClusterTriggerBCExotic(maker.fhClusterTriggerBCExotic),
108 fhClusterTriggerBCBad(maker.fhClusterTriggerBCBad),
109 fhClusterTriggerBCBadExotic(maker.fhClusterTriggerBCBadExotic),
110 fhClusterTriggerBCUnMatch(maker.fhClusterTriggerBCUnMatch),
111 fhClusterTriggerBCExoticUnMatch(maker.fhClusterTriggerBCExoticUnMatch),
112 fhClusterTriggerBCBadUnMatch(maker.fhClusterTriggerBCBadUnMatch),
113 fhClusterTriggerBCBadExoticUnMatch(maker.fhClusterTriggerBCBadExoticUnMatch)
118 //___________________________________________________
119 AliAnaCaloTrackCorrMaker::~AliAnaCaloTrackCorrMaker()
121 // Remove all owned pointers.
123 // Do not delete it here, already done somewhere else, need to understand where.
124 // if (fOutputContainer) {
125 // fOutputContainer->Clear();
126 // delete fOutputContainer ;
129 if (fAnalysisContainer)
131 fAnalysisContainer->Delete();
132 delete fAnalysisContainer ;
135 if (fReader) delete fReader ;
136 if (fCaloUtils) delete fCaloUtils ;
146 //__________________________________________________________________
147 void AliAnaCaloTrackCorrMaker::AddAnalysis(TObject* ana, Int_t n)
149 // Add analysis depending on AliAnaCaloTrackCorrBaseClass to list
151 if ( fAnalysisContainer)
153 fAnalysisContainer->AddAt(ana,n);
157 printf("AliAnaCaloTrackCorrMaker::AddAnalysis() - AnalysisContainer not initialized\n");
162 //_________________________________________________________
163 TList * AliAnaCaloTrackCorrMaker::FillAndGetAODBranchList()
166 // Get any new output AOD branches from analysis and put them in a list
167 // The list is filled in the maker, and new branch passed to the analysis frame
168 // AliAnalysisTaskCaloTrackCorrelation
170 TList *aodBranchList = fReader->GetAODBranchList() ;
172 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
174 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
175 if(ana->NewOutputAOD()) aodBranchList->Add(ana->GetCreateOutputAODBranch());
178 return aodBranchList ;
182 //_________________________________________________________
183 void AliAnaCaloTrackCorrMaker::FillControlHistograms()
185 // Event control histograms
187 AliVEvent* event = fReader->GetInputEvent();
188 AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (event);
189 AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (event);
191 fhNEvents ->Fill(0); // Number of events analyzed
193 if( fReader->IsPileUpFromSPD())
194 fhNPileUpEvents->Fill(0.5);
195 //if( event->IsPileupFromSPDInMultBins())
196 // fhNPileUpEvents->Fill(1.5);
197 if( fReader->IsPileUpFromEMCal())
198 fhNPileUpEvents->Fill(2.5);
199 if( fReader->IsPileUpFromSPDOrEMCal() )
200 fhNPileUpEvents->Fill(3.5);
201 if( fReader->IsPileUpFromSPDAndEMCal() )
202 fhNPileUpEvents->Fill(4.5);
203 if( fReader->IsPileUpFromSPDAndNotEMCal() )
204 fhNPileUpEvents->Fill(5.5);
205 if( fReader->IsPileUpFromEMCalAndNotSPD() )
206 fhNPileUpEvents->Fill(6.5);
207 if( fReader->IsPileUpFromNotSPDAndNotEMCal() )
208 fhNPileUpEvents->Fill(7.5);
210 Int_t triggerBC = fReader->GetTriggerClusterBC() ;
211 if( triggerBC == 0 &&
212 !fReader->IsExoticEvent() &&
213 !fReader->IsBadCellTriggerEvent())
215 if( fReader->IsPileUpFromSPD())
216 fhNPileUpEventsTriggerBC0->Fill(0.5);
217 //if( event->IsPileupFromSPDInMultBins())
218 // fhNPileUpEventsTriggerBC0->Fill(1.5);
219 if( fReader->IsPileUpFromEMCal())
220 fhNPileUpEventsTriggerBC0->Fill(2.5);
221 if( fReader->IsPileUpFromSPDOrEMCal() )
222 fhNPileUpEventsTriggerBC0->Fill(3.5);
223 if( fReader->IsPileUpFromSPDAndEMCal() )
224 fhNPileUpEventsTriggerBC0->Fill(4.5);
225 if( fReader->IsPileUpFromSPDAndNotEMCal() )
226 fhNPileUpEventsTriggerBC0->Fill(5.5);
227 if( fReader->IsPileUpFromEMCalAndNotSPD() )
228 fhNPileUpEventsTriggerBC0->Fill(6.5);
229 if( fReader->IsPileUpFromNotSPDAndNotEMCal() )
230 fhNPileUpEventsTriggerBC0->Fill(7.5);
233 if(fReader->IsPileUpFromSPD())
234 fhPileUpClusterMultAndSPDPileUp ->Fill(fReader->GetNPileUpClusters());
236 fhPileUpClusterMult ->Fill(fReader->GetNPileUpClusters ());
237 fhTrackMult ->Fill(fReader->GetTrackMultiplicity());
238 fhCentrality ->Fill(fReader->GetEventCentrality ());
239 fhEventPlaneAngle ->Fill(fReader->GetEventPlaneAngle ());
241 for(Int_t i = 0; i < 19; i++)
243 if(fReader->GetTrackEventBC(i)) fhTrackBCEvent ->Fill(i);
244 if(fReader->GetTrackEventBCcut(i))fhTrackBCEventCut->Fill(i);
245 if(fReader->GetEMCalEventBC(i)) fhEMCalBCEvent ->Fill(i);
246 if(fReader->GetEMCalEventBCcut(i))fhEMCalBCEventCut->Fill(i);
250 event->GetPrimaryVertex()->GetXYZ(v) ;
251 fhZVertex->Fill(v[2]);
253 Int_t bc = fReader->GetVertexBC();
254 if(bc!=AliVTrack::kTOFBCNA)fhPrimaryVertexBC->Fill(bc);
257 // N pile up vertices
258 Int_t nVerticesSPD = -1;
259 Int_t nVerticesTracks = -1;
263 nVerticesSPD = esdevent->GetNumberOfPileupVerticesSPD();
264 nVerticesTracks = esdevent->GetNumberOfPileupVerticesTracks();
269 nVerticesSPD = aodevent->GetNumberOfPileupVerticesSPD();
270 nVerticesTracks = aodevent->GetNumberOfPileupVerticesTracks();
273 fhNPileUpVertSPD ->Fill(nVerticesSPD);
274 fhNPileUpVertTracks->Fill(nVerticesTracks);
277 if(fReader->IsSelectEventTimeStampOn() && esdevent)
279 Int_t timeStamp = esdevent->GetTimeStamp();
280 Float_t timeStampFrac = 1.*(timeStamp-fReader->GetRunTimeStampMin()) /
281 (fReader->GetRunTimeStampMax()-fReader->GetRunTimeStampMin());
283 //printf("stamp %d, min %d, max %d, frac %f\n", timeStamp, fReader->GetRunTimeStampMin(), fReader->GetRunTimeStampMax(), timeStampFrac);
285 fhTimeStampFraction->Fill(timeStampFrac);
289 //_______________________________________________________
290 TList * AliAnaCaloTrackCorrMaker::GetListOfAnalysisCuts()
293 // Get the list of the cuts used for the analysis
294 // The list is filled in the maker, called by the task in LocalInit() and posted there
296 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
298 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
299 TObjString * objstring = ana->GetAnalysisCuts();
301 if(objstring)fCuts->Add(objstring);
308 //___________________________________________________
309 TList *AliAnaCaloTrackCorrMaker::GetOutputContainer()
311 // Fill the output list of histograms during the CreateOutputObjects stage.
313 //Initialize calorimeters geometry pointers
314 //GetCaloUtils()->InitPHOSGeometry();
315 //GetCaloUtils()->InitEMCALGeometry();
317 //General event histograms
319 fhNEvents = new TH1F("hNEvents", "Number of analyzed events" , 1 , 0 , 1 ) ;
320 fhNEvents->SetYTitle("# events");
321 fOutputContainer->Add(fhNEvents);
323 fhNExoticEvents = new TH1F("hNExoticEvents", "Number of analyzed events triggered by exotic cluster" , 1 , 0 , 1 ) ;
324 fhNExoticEvents->SetYTitle("# exotic events");
325 fOutputContainer->Add(fhNExoticEvents);
327 fhNEventsNoTriggerFound = new TH1F("hNEventsNoTriggerFound", "Number of analyzed events triggered but no trigger found" , 1 , 0 , 1 ) ;
328 fhNEventsNoTriggerFound->SetYTitle("# exotic events");
329 fOutputContainer->Add(fhNEventsNoTriggerFound);
332 Float_t minbin =-5.5;
333 Float_t maxbin = 5.5;
334 Int_t labelshift = 6;
335 fhClusterTriggerBC = new TH1F("hClusterTriggerBC",
336 "Number of analyzed events triggered by a cluster in a given BC",
337 nbin , minbin ,maxbin) ;
338 fhClusterTriggerBC->SetYTitle("# events");
339 for(Int_t i = 1; i < 12; i++)
340 fhClusterTriggerBC->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
341 fOutputContainer->Add(fhClusterTriggerBC);
343 fhClusterTriggerBCExotic = new TH1F("hClusterTriggerBCExotic",
344 "Number of analyzed events triggered by a exotic cluster in a given BC",
345 nbin , minbin ,maxbin) ;
346 fhClusterTriggerBCExotic->SetYTitle("# events");
347 for(Int_t i = 1; i < 12; i++)
348 fhClusterTriggerBCExotic->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
349 fOutputContainer->Add(fhClusterTriggerBCExotic);
351 fhClusterTriggerBCBad = new TH1F("hClusterTriggerBCBad",
352 "Number of analyzed events triggered by a bad cluster in a given BC",
353 nbin , minbin ,maxbin) ;
354 fhClusterTriggerBCBad->SetYTitle("# events");
355 for(Int_t i = 1; i < 12; i++)
356 fhClusterTriggerBCBad->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
357 fOutputContainer->Add(fhClusterTriggerBCBad);
360 fhClusterTriggerBCBadExotic = new TH1F("hClusterTriggerBCBadExotic",
361 "Number of analyzed events triggered by a bad&exotic cluster in a given BC",
362 nbin , minbin ,maxbin) ;
363 fhClusterTriggerBCBadExotic->SetYTitle("# events");
364 for(Int_t i = 1; i < 12; i++)
365 fhClusterTriggerBCBadExotic->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
366 fOutputContainer->Add(fhClusterTriggerBCBadExotic);
368 fhClusterTriggerBCUnMatch = new TH1F("hClusterTriggerBCUnMatch",
369 "Number of analyzed events triggered by a cluster (no trigger patch match) in a given BC",
370 nbin , minbin ,maxbin) ;
371 fhClusterTriggerBCUnMatch->SetYTitle("# events");
372 for(Int_t i = 1; i < 12; i++)
373 fhClusterTriggerBCUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
374 fOutputContainer->Add(fhClusterTriggerBCUnMatch);
376 fhClusterTriggerBCExoticUnMatch = new TH1F("hClusterTriggerBCExoticUnMatch",
377 "Number of analyzed events triggered by a exotic cluster (no trigger patch match) in a given BC",
378 nbin , minbin ,maxbin) ;
379 fhClusterTriggerBCExoticUnMatch->SetYTitle("# events");
380 for(Int_t i = 1; i < 12; i++)
381 fhClusterTriggerBCExoticUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
382 fOutputContainer->Add(fhClusterTriggerBCExoticUnMatch);
384 fhClusterTriggerBCBadUnMatch = new TH1F("hClusterTriggerBCBadUnMatch",
385 "Number of analyzed events triggered by a bad cluster (no trigger patch match) in a given BC",
386 nbin , minbin ,maxbin) ;
387 fhClusterTriggerBCBadUnMatch->SetYTitle("# events");
388 for(Int_t i = 1; i < 12; i++)
389 fhClusterTriggerBCBadUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
390 fOutputContainer->Add(fhClusterTriggerBCBadUnMatch);
393 fhClusterTriggerBCBadExoticUnMatch = new TH1F("hClusterTriggerBCBadExoticUnMatch",
394 "Number of analyzed events triggered by a bad&exotic cluster (no trigger patch match) in a given BC",
395 nbin , minbin ,maxbin) ;
396 fhClusterTriggerBCBadExoticUnMatch->SetYTitle("# events");
397 for(Int_t i = 1; i < 12; i++)
398 fhClusterTriggerBCBadExoticUnMatch->GetXaxis()->SetBinLabel(i ,Form("BC%d",i-labelshift));
399 fOutputContainer->Add(fhClusterTriggerBCBadExoticUnMatch);
402 fhNPileUpEvents = new TH1F("hNPileUpEvents", "Number of events considered as pile-up", 8 , 0 , 8 ) ;
403 fhNPileUpEvents->SetYTitle("# events");
404 fhNPileUpEvents->GetXaxis()->SetBinLabel(1 ,"SPD");
405 fhNPileUpEvents->GetXaxis()->SetBinLabel(2 ,"Multi SPD");
406 fhNPileUpEvents->GetXaxis()->SetBinLabel(3 ,"EMCal");
407 fhNPileUpEvents->GetXaxis()->SetBinLabel(4 ,"EMCal || SPD");
408 fhNPileUpEvents->GetXaxis()->SetBinLabel(5 ,"EMCal && SPD");
409 fhNPileUpEvents->GetXaxis()->SetBinLabel(6 ,"!EMCal && SPD");
410 fhNPileUpEvents->GetXaxis()->SetBinLabel(7 ,"EMCal && !SPD");
411 fhNPileUpEvents->GetXaxis()->SetBinLabel(8 ,"!EMCal && !SPD");
412 fOutputContainer->Add(fhNPileUpEvents);
414 fhNPileUpEventsTriggerBC0 = new TH1F("hNPileUpEventsTriggerBC0","Number of events considered as pile-up, trigger cluster in BC=0", 8 , 0 , 8 ) ;
415 fhNPileUpEventsTriggerBC0->SetYTitle("# events");
416 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(1 ,"SPD");
417 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(2 ,"Multi SPD");
418 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(3 ,"EMCal");
419 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(4 ,"EMCal || SPD");
420 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(5 ,"EMCal && SPD");
421 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(6 ,"!EMCal && SPD");
422 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(7 ,"EMCal && !SPD");
423 fhNPileUpEventsTriggerBC0->GetXaxis()->SetBinLabel(8 ,"!EMCal && !SPD");
424 fOutputContainer->Add(fhNPileUpEventsTriggerBC0);
427 fhTrackBCEvent = new TH1F("hTrackBCEvent", "Number of events with at least 1 track in a bunch crossing ", 19 , 0 , 19 ) ;
428 fhTrackBCEvent->SetYTitle("# events");
429 fhTrackBCEvent->SetXTitle("Bunch crossing");
430 for(Int_t i = 1; i < 20; i++)
431 fhTrackBCEvent->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
432 fOutputContainer->Add(fhTrackBCEvent);
434 fhTrackBCEventCut = new TH1F("hTrackBCEventCut", "Number of events with at least 1 track in a bunch crossing ", 19 , 0 , 19 ) ;
435 fhTrackBCEventCut->SetYTitle("# events");
436 fhTrackBCEventCut->SetXTitle("Bunch crossing");
437 for(Int_t i = 1; i < 20; i++)
438 fhTrackBCEventCut->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
439 fOutputContainer->Add(fhTrackBCEventCut);
441 fhPrimaryVertexBC = new TH1F("hPrimaryVertexBC", "Number of primary vertex per bunch crossing ", 41 , -20 , 20 ) ;
442 fhPrimaryVertexBC->SetYTitle("# events");
443 fhPrimaryVertexBC->SetXTitle("Bunch crossing");
444 fOutputContainer->Add(fhPrimaryVertexBC);
446 fhEMCalBCEvent = new TH1F("hEMCalBCEvent", "Number of events with at least 1 cluster in a bunch crossing ", 19 , 0 , 19 ) ;
447 fhEMCalBCEvent->SetYTitle("# events");
448 fhEMCalBCEvent->SetXTitle("Bunch crossing");
449 for(Int_t i = 1; i < 20; i++)
450 fhEMCalBCEvent->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
451 fOutputContainer->Add(fhEMCalBCEvent);
453 fhEMCalBCEventCut = new TH1F("hEMCalBCEventCut", "Number of events with at least 1 cluster in a bunch crossing", 19 , 0 , 19 ) ;
454 fhEMCalBCEventCut->SetYTitle("# events");
455 fhEMCalBCEventCut->SetXTitle("Bunch crossing");
456 for(Int_t i = 1; i < 20; i++)
457 fhEMCalBCEventCut->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
458 fOutputContainer->Add(fhEMCalBCEventCut);
460 fhZVertex = new TH1F("hZVertex", " Z vertex distribution" , 200 , -50 , 50 ) ;
461 fhZVertex->SetXTitle("v_{z} (cm)");
462 fOutputContainer->Add(fhZVertex);
464 fhTrackMult = new TH1F("hTrackMult", "Number of tracks per events" , 2000 , 0 , 2000 ) ;
465 fhTrackMult->SetXTitle("# tracks");
466 fOutputContainer->Add(fhTrackMult);
468 fhPileUpClusterMult = new TH1F("hPileUpClusterMult", "Number of clusters per event with large time (|t| > 20 ns)" , 100 , 0 , 100 ) ;
469 fhPileUpClusterMult->SetXTitle("# clusters");
470 fOutputContainer->Add(fhPileUpClusterMult);
472 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 ) ;
473 fhPileUpClusterMultAndSPDPileUp->SetXTitle("# clusters");
474 fOutputContainer->Add(fhPileUpClusterMultAndSPDPileUp);
476 fhNPileUpVertSPD = new TH1F ("hNPileUpVertSPD","N pile-up SPD vertex", 50,0,50);
477 fhNPileUpVertSPD->SetYTitle("# vertex ");
478 fOutputContainer->Add(fhNPileUpVertSPD);
480 fhNPileUpVertTracks = new TH1F ("hNPileUpVertTracks","N pile-up Tracks vertex", 50,0,50);
481 fhNPileUpVertTracks->SetYTitle("# vertex ");
482 fOutputContainer->Add(fhNPileUpVertTracks);
484 fhCentrality = new TH1F("hCentrality","Number of events in centrality bin",100,0.,100) ;
485 fhCentrality->SetXTitle("Centrality bin");
486 fOutputContainer->Add(fhCentrality) ;
488 fhEventPlaneAngle=new TH1F("hEventPlaneAngle","Number of events in event plane",100,0.,TMath::Pi()) ;
489 fhEventPlaneAngle->SetXTitle("EP angle (rad)");
490 fOutputContainer->Add(fhEventPlaneAngle) ;
492 if(fReader->IsSelectEventTimeStampOn())
494 fhTimeStampFraction = new TH1F("hTimeStampFraction","Fraction of events within a given time stamp range",150, -1, 2) ;
495 fhTimeStampFraction->SetXTitle("fraction");
496 fOutputContainer->Add(fhTimeStampFraction) ;
501 fhNMergedFiles = new TH1F("hNMergedFiles", "Number of merged output files" , 1 , 0 , 1 ) ;
502 fhNMergedFiles->SetYTitle("# files");
503 fhNMergedFiles->Fill(1); // Fill here with one entry, while merging it will count the rest
504 fOutputContainer->Add(fhNMergedFiles);
506 fhScaleFactor = new TH1F("hScaleFactor", "Number of merged output files" , 1 , 0 , 1 ) ;
507 fhScaleFactor->SetYTitle("scale factor");
508 fhScaleFactor->SetBinContent(1,fScaleFactor); // Fill here
509 fOutputContainer->Add(fhScaleFactor);
512 if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
514 printf("AliAnaCaloTrackCorrMaker::GetOutputContainer() - Analysis job list not initialized!!!\n");
515 return fOutputContainer;
518 const Int_t buffersize = 255;
519 char newname[buffersize];
520 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
523 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
525 if(fMakeHisto) // Analysis with histograms as output on
528 //Fill container with appropriate histograms
529 TList * templist = ana ->GetCreateOutputObjects();
530 templist->SetOwner(kFALSE); //Owner is fOutputContainer.
532 for(Int_t i = 0; i < templist->GetEntries(); i++)
535 //Add only to the histogram name the name of the task
536 if( strcmp((templist->At(i))->ClassName(),"TObjString") )
538 snprintf(newname,buffersize, "%s%s", (ana->GetAddedHistogramsStringToName()).Data(), (templist->At(i))->GetName());
539 //printf("name %s, new name %s\n",(templist->At(i))->GetName(),newname);
540 ((TH1*) templist->At(i))->SetName(newname);
543 //Add histogram to general container
544 fOutputContainer->Add(templist->At(i)) ;
550 }// Analysis with histograms as output on
552 }//Loop on analysis defined
554 return fOutputContainer;
558 //___________________________________
559 void AliAnaCaloTrackCorrMaker::Init()
561 //Init container histograms and other common variables
562 // Fill the output list of histograms during the CreateOutputObjects stage.
566 GetReader()->SetCaloUtils(GetCaloUtils()); // pass the calo utils pointer to the reader
569 if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
571 printf("AliAnaCaloTrackCorrMaker::GetOutputInit() - Analysis job list not initialized!!!\n");
575 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
578 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
580 ana->SetReader(fReader); // Set Reader for each analysis
581 ana->SetCaloUtils(fCaloUtils); // Set CaloUtils for each analysis
585 }//Loop on analysis defined
589 //_____________________________________________
590 void AliAnaCaloTrackCorrMaker::InitParameters()
596 fAnaDebug = 0; // No debugging info displayed by default
600 //______________________________________________________________
601 void AliAnaCaloTrackCorrMaker::Print(const Option_t * opt) const
603 //Print some relevant parameters set for the analysis
608 printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
609 printf("Debug level = %d\n", fAnaDebug ) ;
610 printf("Produce Histo = %d\n", fMakeHisto ) ;
611 printf("Produce AOD = %d\n", fMakeAOD ) ;
612 printf("Number of analysis tasks = %d\n", fAnalysisContainer->GetEntries()) ;
614 if(!strcmp("all",opt))
616 printf("Print analysis Tasks settings :\n") ;
617 for(Int_t iana = 0; iana<fAnalysisContainer->GetEntries(); iana++)
619 ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana))->Print("");
622 printf("Print analysis Reader settings :\n") ;
624 printf("Print analysis Calorimeter Utils settings :\n") ;
625 fCaloUtils->Print("");
631 //_______________________________________________________________________
632 void AliAnaCaloTrackCorrMaker::ProcessEvent(const Int_t iEntry,
633 const char * currentFileName)
635 //Process analysis for this event
637 if(fMakeHisto && !fOutputContainer)
639 printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Histograms not initialized\n");
645 printf("*** AliAnaCaloTrackCorrMaker::ProcessEvent() Event %d *** \n",iEntry);
648 printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Current File Name : %s\n", currentFileName);
649 //printf("fAODBranchList %p, entries %d\n",fAODBranchList,fAODBranchList->GetEntries());
653 //Each event needs an empty branch
654 TList * aodList = fReader->GetAODBranchList();
655 Int_t nAODBranches = aodList->GetEntries();
656 for(Int_t iaod = 0; iaod < nAODBranches; iaod++)
658 TClonesArray *tca = dynamic_cast<TClonesArray*> (aodList->At(iaod));
659 if(tca) tca->Clear("C");
662 //Set geometry matrices before filling arrays, in case recalibration/position calculation etc is needed
663 fCaloUtils->AccessGeometry(fReader->GetInputEvent());
665 //Set the AODB calibration, bad channels etc. parameters at least once
666 fCaloUtils->AccessOADB(fReader->GetInputEvent());
669 //Tell the reader to fill the data in the 3 detector lists
670 Bool_t ok = fReader->FillInputEvent(iEntry, currentFileName);
672 Int_t triggerBC = fReader->GetTriggerClusterBC() ;
673 Bool_t exotic = fReader->IsExoticEvent();
674 Bool_t badCell = fReader->IsBadCellTriggerEvent();
675 Bool_t triggerMatch= fReader->IsTriggerMatched();
676 Bool_t triggerBCOK = kTRUE;
677 Int_t triggerId = fReader->GetTriggerClusterId() ;
680 //printf("Trigger id %d\n",triggerId);
681 if(triggerId == -2)fhNEventsNoTriggerFound->Fill(0);
682 triggerBCOK = kFALSE;
685 if(exotic) fhNExoticEvents->Fill(0) ;
686 //if(fReader->IsExoticEvent()) printf("Maker: EXOTIC Cluster trigger\n");
692 if (!exotic && !badCell) fhClusterTriggerBC ->Fill(triggerBC);
693 else if( exotic && badCell) fhClusterTriggerBCBadExotic->Fill(triggerBC);
694 else if( exotic && !badCell) fhClusterTriggerBCExotic ->Fill(triggerBC);
695 else if( badCell && !exotic ) fhClusterTriggerBCBad ->Fill(triggerBC);
699 if (!exotic && !badCell) fhClusterTriggerBCUnMatch ->Fill(triggerBC);
700 else if( exotic && badCell) fhClusterTriggerBCBadExoticUnMatch->Fill(triggerBC);
701 else if( exotic && !badCell) fhClusterTriggerBCExoticUnMatch ->Fill(triggerBC);
702 else if( badCell && !exotic ) fhClusterTriggerBCBadUnMatch ->Fill(triggerBC);
706 if(!ok && triggerBC > -9999) printf("Maker: Cluster trigger BC = %d\n",triggerBC);
710 if(fAnaDebug >= 1 )printf("*** Skip event *** %d \n",iEntry);
711 fReader->ResetLists();
715 //Magic line to write events to file
716 if(fReader->WriteDeltaAODToFile())AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kTRUE);
718 //printf(">>>>>>>>>> BEFORE >>>>>>>>>>>\n");
719 //gObjectTable->Print();
721 //Access pointers, and trigger mask check needed in mixing case
722 AliAnalysisManager *manager = AliAnalysisManager::GetAnalysisManager();
723 AliInputEventHandler *inputHandler = dynamic_cast<AliInputEventHandler*>(manager->GetInputEventHandler());
725 UInt_t isMBTrigger = kFALSE;
726 UInt_t isTrigger = kFALSE;
729 isMBTrigger = inputHandler->IsEventSelected() & fReader->GetMixEventTriggerMask();
730 isTrigger = inputHandler->IsEventSelected() & fReader->GetEventTriggerMask();
733 //Loop on analysis algorithms
735 if(fAnaDebug > 0 ) printf("*** Begin analysis *** \n");
737 Int_t nana = fAnalysisContainer->GetEntries() ;
738 for(Int_t iana = 0; iana < nana; iana++)
740 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
742 ana->ConnectInputOutputAODBranches(); //Sets branches for each analysis
744 //Fill pool for mixed event for the analysis that need it
745 if(!fReader->IsEventTriggerAtSEOn() && isMBTrigger)
747 ana->FillEventMixPool();
748 continue; // pool filled do not try to fill AODs or histograms
751 //Make analysis, create aods in aod branch and in some cases fill histograms
752 if(fMakeAOD ) ana->MakeAnalysisFillAOD() ;
754 //Make further analysis with aod branch and fill histograms
755 if(fMakeHisto) ana->MakeAnalysisFillHistograms() ;
759 fReader->ResetLists();
761 // In case of mixing analysis, non triggered events are used,
762 // do not fill control histograms for a non requested triggered event
763 if(!fReader->IsEventTriggerAtSEOn() && !isTrigger)
765 if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis, MB for mixing *** \n");
769 FillControlHistograms();
771 //printf(">>>>>>>>>> AFTER >>>>>>>>>>>\n");
772 //gObjectTable->Print();
774 if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis *** \n");
778 //__________________________________________________________
779 void AliAnaCaloTrackCorrMaker::Terminate(TList * outputList)
781 //Execute Terminate of analysis
782 //Do some final plots.
786 Error("Terminate", "No output list");
790 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
793 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
794 if(ana->MakePlotsOn())ana->Terminate(outputList);
796 }//Loop on analysis defined