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 "AliAnaCaloTrackCorrBaseClass.h"
37 #include "AliAnaCaloTrackCorrMaker.h"
39 ClassImp(AliAnaCaloTrackCorrMaker)
42 //__________________________________________________
43 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker() :
45 fReader(0), fCaloUtils(0),
46 fOutputContainer(new TList ), fAnalysisContainer(new TList ),
47 fMakeHisto(kFALSE), fMakeAOD(kFALSE),
48 fAnaDebug(0), fCuts(new TList),
50 fhNEvents(0), fhZVertex(0),
51 fhTrackMult(0), fhCentrality(0),
53 fhNMergedFiles(0), fhScaleFactor(0)
56 if(fAnaDebug > 1 ) printf("*** Analysis Maker Constructor *** \n");
58 //Initialize parameters, pointers and histograms
62 //________________________________________________________________________________________
63 AliAnaCaloTrackCorrMaker::AliAnaCaloTrackCorrMaker(const AliAnaCaloTrackCorrMaker & maker) :
65 fReader(), //(new AliCaloTrackReader(*maker.fReader)),
66 fCaloUtils(),//(new AliCalorimeterUtils(*maker.fCaloUtils)),
67 fOutputContainer(new TList()), fAnalysisContainer(new TList()),
68 fMakeHisto(maker.fMakeHisto), fMakeAOD(maker.fMakeAOD),
69 fAnaDebug(maker.fAnaDebug), fCuts(new TList()),
70 fScaleFactor(maker.fScaleFactor),
71 fhNEvents(maker.fhNEvents), fhZVertex(maker.fhZVertex),
72 fhTrackMult(maker.fhTrackMult),fhCentrality(maker.fhCentrality),
73 fhEventPlaneAngle(maker.fhEventPlaneAngle),
74 fhNMergedFiles(maker.fhNMergedFiles),
75 fhScaleFactor(maker.fhScaleFactor)
80 //___________________________________________________
81 AliAnaCaloTrackCorrMaker::~AliAnaCaloTrackCorrMaker()
83 // Remove all owned pointers.
85 // Do not delete it here, already done somewhere else, need to understand where.
86 // if (fOutputContainer) {
87 // fOutputContainer->Clear();
88 // delete fOutputContainer ;
91 if (fAnalysisContainer)
93 fAnalysisContainer->Delete();
94 delete fAnalysisContainer ;
97 if (fReader) delete fReader ;
98 if (fCaloUtils) delete fCaloUtils ;
108 //__________________________________________________________________
109 void AliAnaCaloTrackCorrMaker::AddAnalysis(TObject* ana, Int_t n)
111 // Add analysis depending on AliAnaCaloTrackCorrBaseClass to list
113 if ( fAnalysisContainer)
115 fAnalysisContainer->AddAt(ana,n);
119 printf("AliAnaCaloTrackCorrMaker::AddAnalysis() - AnalysisContainer not initialized\n");
124 //_________________________________________________________
125 TList * AliAnaCaloTrackCorrMaker::FillAndGetAODBranchList()
128 // Get any new output AOD branches from analysis and put them in a list
129 // The list is filled in the maker, and new branch passed to the analysis frame
130 // AliAnalysisTaskCaloTrackCorrelation
132 TList *aodBranchList = fReader->GetAODBranchList() ;
134 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
136 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
137 if(ana->NewOutputAOD()) aodBranchList->Add(ana->GetCreateOutputAODBranch());
140 return aodBranchList ;
144 //_______________________________________________________
145 TList * AliAnaCaloTrackCorrMaker::GetListOfAnalysisCuts()
148 // Get the list of the cuts used for the analysis
149 // The list is filled in the maker, called by the task in LocalInit() and posted there
151 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
153 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
154 TObjString * objstring = ana->GetAnalysisCuts();
156 if(objstring)fCuts->Add(objstring);
163 //___________________________________________________
164 TList *AliAnaCaloTrackCorrMaker::GetOutputContainer()
166 // Fill the output list of histograms during the CreateOutputObjects stage.
168 //Initialize calorimeters geometry pointers
169 //GetCaloUtils()->InitPHOSGeometry();
170 //GetCaloUtils()->InitEMCALGeometry();
172 //General event histograms
174 fhNEvents = new TH1I("hNEvents", "Number of analyzed events" , 1 , 0 , 1 ) ;
175 fhNEvents->SetYTitle("# events");
176 fOutputContainer->Add(fhNEvents);
178 fhZVertex = new TH1F("hZVertex", " Z vertex distribution" , 200 , -50 , 50 ) ;
179 fhZVertex->SetXTitle("v_{z} (cm)");
180 fOutputContainer->Add(fhZVertex);
182 fhTrackMult = new TH1I("hTrackMult", "Number of tracks per events" , 2000 , 0 , 2000 ) ;
183 fhTrackMult->SetXTitle("# tracks");
184 fOutputContainer->Add(fhTrackMult);
186 fhCentrality = new TH1F("hCentrality","Number of events in centrality bin",100,0.,100) ;
187 fhCentrality->SetXTitle("Centrality bin");
188 fOutputContainer->Add(fhCentrality) ;
190 fhEventPlaneAngle=new TH1F("hEventPlaneAngle","Number of events in event plane",100,0.,TMath::Pi()) ;
191 fhEventPlaneAngle->SetXTitle("EP angle (rad)");
192 fOutputContainer->Add(fhEventPlaneAngle) ;
196 fhNMergedFiles = new TH1I("hNMergedFiles", "Number of merged output files" , 1 , 0 , 1 ) ;
197 fhNMergedFiles->SetYTitle("# files");
198 fhNMergedFiles->Fill(1); // Fill here with one entry, while merging it will count the rest
199 fOutputContainer->Add(fhNMergedFiles);
201 fhScaleFactor = new TH1F("hScaleFactor", "Number of merged output files" , 1 , 0 , 1 ) ;
202 fhScaleFactor->SetYTitle("scale factor");
203 fhScaleFactor->SetBinContent(1,fScaleFactor); // Fill here
204 fOutputContainer->Add(fhScaleFactor);
207 if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
209 printf("AliAnaCaloTrackCorrMaker::GetOutputContainer() - Analysis job list not initialized!!!\n");
210 return fOutputContainer;
213 const Int_t buffersize = 255;
214 char newname[buffersize];
215 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
218 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
220 if(fMakeHisto) // Analysis with histograms as output on
223 //Fill container with appropriate histograms
224 TList * templist = ana ->GetCreateOutputObjects();
225 templist->SetOwner(kFALSE); //Owner is fOutputContainer.
227 for(Int_t i = 0; i < templist->GetEntries(); i++)
230 //Add only to the histogram name the name of the task
231 if( strcmp((templist->At(i))->ClassName(),"TObjString") )
233 snprintf(newname,buffersize, "%s%s", (ana->GetAddedHistogramsStringToName()).Data(), (templist->At(i))->GetName());
234 ((TH1*) templist->At(i))->SetName(newname);
237 //Add histogram to general container
238 fOutputContainer->Add(templist->At(i)) ;
244 }// Analysis with histograms as output on
246 }//Loop on analysis defined
248 return fOutputContainer;
252 //___________________________________
253 void AliAnaCaloTrackCorrMaker::Init()
255 //Init container histograms and other common variables
256 // Fill the output list of histograms during the CreateOutputObjects stage.
260 GetReader()->SetCaloUtils(GetCaloUtils()); // pass the calo utils pointer to the reader
263 if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0)
265 printf("AliAnaCaloTrackCorrMaker::GetOutputInit() - Analysis job list not initialized!!!\n");
269 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
272 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
274 ana->SetReader(fReader); // Set Reader for each analysis
275 ana->SetCaloUtils(fCaloUtils); // Set CaloUtils for each analysis
279 }//Loop on analysis defined
283 //_____________________________________________
284 void AliAnaCaloTrackCorrMaker::InitParameters()
290 fAnaDebug = 0; // No debugging info displayed by default
294 //______________________________________________________________
295 void AliAnaCaloTrackCorrMaker::Print(const Option_t * opt) const
297 //Print some relevant parameters set for the analysis
302 printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
303 printf("Debug level = %d\n", fAnaDebug ) ;
304 printf("Produce Histo = %d\n", fMakeHisto ) ;
305 printf("Produce AOD = %d\n", fMakeAOD ) ;
306 printf("Number of analysis tasks = %d\n", fAnalysisContainer->GetEntries()) ;
308 if(!strcmp("all",opt))
310 printf("Print analysis Tasks settings :\n") ;
311 for(Int_t iana = 0; iana<fAnalysisContainer->GetEntries(); iana++)
313 ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana))->Print("");
316 printf("Print analysis Reader settings :\n") ;
318 printf("Print analysis Calorimeter Utils settings :\n") ;
319 fCaloUtils->Print("");
325 //_______________________________________________________________________
326 void AliAnaCaloTrackCorrMaker::ProcessEvent(const Int_t iEntry,
327 const char * currentFileName)
329 //Process analysis for this event
331 if(fMakeHisto && !fOutputContainer)
333 printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Histograms not initialized\n");
339 printf("*** AliAnaCaloTrackCorrMaker::ProcessEvent() Event %d *** \n",iEntry);
342 printf("AliAnaCaloTrackCorrMaker::ProcessEvent() - Current File Name : %s\n", currentFileName);
343 //printf("fAODBranchList %p, entries %d\n",fAODBranchList,fAODBranchList->GetEntries());
347 //Each event needs an empty branch
348 TList * aodList = fReader->GetAODBranchList();
349 Int_t nAODBranches = aodList->GetEntries();
350 for(Int_t iaod = 0; iaod < nAODBranches; iaod++)
352 TClonesArray *tca = dynamic_cast<TClonesArray*> (aodList->At(iaod));
353 if(tca) tca->Clear("C");
356 //Set geometry matrices before filling arrays, in case recalibration/position calculation etc is needed
357 fCaloUtils->AccessGeometry(fReader->GetInputEvent());
359 //Set the AODB calibration, bad channels etc. parameters at least once
360 fCaloUtils->AccessOADB(fReader->GetInputEvent());
363 //Tell the reader to fill the data in the 3 detector lists
364 Bool_t ok = fReader->FillInputEvent(iEntry, currentFileName);
367 if(fAnaDebug >= 1 )printf("*** Skip event *** %d \n",iEntry);
371 //Magic line to write events to file
372 if(fReader->WriteDeltaAODToFile())AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kTRUE);
374 //printf(">>>>>>>>>> BEFORE >>>>>>>>>>>\n");
375 //gObjectTable->Print();
377 //Loop on analysis algorithms
379 if(fAnaDebug > 0 ) printf("*** Begin analysis *** \n");
381 Int_t nana = fAnalysisContainer->GetEntries() ;
382 for(Int_t iana = 0; iana < nana; iana++)
384 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
386 ana->ConnectInputOutputAODBranches(); //Sets branches for each analysis
387 //Make analysis, create aods in aod branch or AODCaloClusters
388 if(fMakeAOD ) ana->MakeAnalysisFillAOD() ;
389 //Make further analysis with aod branch and fill histograms
390 if(fMakeHisto) ana->MakeAnalysisFillHistograms() ;
394 fReader->ResetLists();
396 // In case of mixing analysis, non triggered events are used,
397 // do not fill control histograms for a non requested triggered event
398 if(!fReader->IsEventTriggerAtSEOn())
400 AliAnalysisManager *manager = AliAnalysisManager::GetAnalysisManager();
401 AliInputEventHandler *inputHandler = dynamic_cast<AliInputEventHandler*>(manager->GetInputEventHandler());
403 if(!inputHandler) return ;
405 UInt_t isTrigger = inputHandler->IsEventSelected() & fReader->GetEventTriggerMask();
408 if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis, MB for mixing *** \n");
414 // Event control histograms
416 fhNEvents ->Fill(0); // Number of events analyzed
417 fhTrackMult ->Fill(fReader->GetTrackMultiplicity());
418 fhCentrality ->Fill(fReader->GetEventCentrality ());
419 fhEventPlaneAngle->Fill(fReader->GetEventPlaneAngle ());
423 fReader->GetInputEvent()->GetPrimaryVertex()->GetXYZ(v) ;
424 fhZVertex->Fill(v[2]);
427 //printf(">>>>>>>>>> AFTER >>>>>>>>>>>\n");
428 //gObjectTable->Print();
430 if(fAnaDebug > 0 ) printf("AliAnaCaloTrackMaker::ProcessEvent() - *** End analysis *** \n");
434 //__________________________________________________________
435 void AliAnaCaloTrackCorrMaker::Terminate(TList * outputList)
437 //Execute Terminate of analysis
438 //Do some final plots.
442 Error("Terminate", "No output list");
446 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++)
449 AliAnaCaloTrackCorrBaseClass * ana = ((AliAnaCaloTrackCorrBaseClass *) fAnalysisContainer->At(iana)) ;
450 if(ana->MakePlotsOn())ana->Terminate(outputList);
452 }//Loop on analysis defined