]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG4/PartCorrBase/AliAnaPartCorrMaker.cxx
Add calorimeter utility class for geometry access, cell indexing bad channels rejecti...
[u/mrichter/AliRoot.git] / PWG4 / PartCorrBase / AliAnaPartCorrMaker.cxx
CommitLineData
1c5acb87 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15/* $Id: $ */
16
17//_________________________________________________________________________
18// Steering class for particle (gamma, hadron) identification and correlation analysis
19// It is called by the task class AliAnalysisTaskParticleCorrelation and it connects the input
20// (ESD/AOD/MonteCarlo) got with AliCaloTrackReader (produces TClonesArrays of AODs
21// (TParticles in MC case if requested)), with the
22// analysis classes that derive from AliAnaPartCorrBaseClass
23//
24// -- Author: Gustavo Conesa (INFN-LNF)
25
fbeaf916 26#include <cstdlib>
27
1c5acb87 28// --- ROOT system ---
477d6cee 29#include "TClonesArray.h"
477d6cee 30#include "TList.h"
031ac63e 31#include "TH1I.h"
1c5acb87 32//#include "Riostream.h"
1d4cea01 33//#include <TObjectTable.h>
1c5acb87 34
35//---- AliRoot system ----
36#include "AliAnaPartCorrBaseClass.h"
37#include "AliAnaPartCorrMaker.h"
38#include "AliCaloTrackReader.h"
1c5acb87 39
40
41ClassImp(AliAnaPartCorrMaker)
42
43
44//____________________________________________________________________________
45AliAnaPartCorrMaker::AliAnaPartCorrMaker() :
46TObject(),
47fOutputContainer(new TList ), fAnalysisContainer(new TList ),
48fMakeHisto(0), fMakeAOD(0), fAnaDebug(0),
765d44e7 49fReader(0x0), fCaloUtils(0x0), fAODBranchList(new TList ),
031ac63e 50fhNEvents(0x0)
1c5acb87 51{
477d6cee 52 //Default Ctor
53 if(fAnaDebug > 1 ) printf("*** Analysis Maker Constructor *** \n");
54
55 //Initialize parameters, pointers and histograms
56 if(!fReader)
765d44e7 57 fReader = new AliCaloTrackReader();
58 if(!fCaloUtils)
59 fCaloUtils = new AliCalorimeterUtils();
60
477d6cee 61 InitParameters();
1c5acb87 62}
63
64//____________________________________________________________________________
d151d2ee 65AliAnaPartCorrMaker::AliAnaPartCorrMaker(const AliAnaPartCorrMaker & maker) :
1c5acb87 66TObject(),
d151d2ee 67fOutputContainer(new TList()), fAnalysisContainer(new TList()),
68fMakeHisto(maker.fMakeHisto), fMakeAOD(maker.fMakeAOD), fAnaDebug(maker.fAnaDebug),
765d44e7 69fReader(new AliCaloTrackReader(*maker.fReader)),
70fCaloUtils(new AliCalorimeterUtils(*maker.fCaloUtils)),
71fAODBranchList(new TList()),
031ac63e 72fhNEvents(maker.fhNEvents)
1c5acb87 73{
477d6cee 74 // cpy ctor
1c5acb87 75
76}
77
78//_________________________________________________________________________
d151d2ee 79//AliAnaPartCorrMaker & AliAnaPartCorrMaker::operator = (const AliAnaPartCorrMaker & source)
80//{
81// // assignment operator
82//
83// if(this == &source)return *this;
84// ((TObject *)this)->operator=(source);
85//
86// delete fOutputContainer; fOutputContainer = source.fOutputContainer ;
87// delete fAnalysisContainer; fAnalysisContainer = source.fAnalysisContainer ;
88// fAnaDebug = source.fAnaDebug;
89// fMakeHisto = source.fMakeHisto;
90// fMakeAOD = source.fMakeAOD;
91//
92// delete fReader ; fReader = new AliCaloTrackReader(*source.fReader) ;
93// delete fAODBranchList; fAODBranchList = source.fAODBranchList;
94//
95// return *this;
96//
97//}
1c5acb87 98
99//____________________________________________________________________________
100AliAnaPartCorrMaker::~AliAnaPartCorrMaker()
101{
477d6cee 102 // Remove all pointers.
103
104 // Protection added in case of NULL pointers (MG)
105 if (fOutputContainer) {
106 fOutputContainer->Clear();
107 delete fOutputContainer ;
108 }
109
110 if (fAnalysisContainer) {
111 fAnalysisContainer->Clear();
112 delete fAnalysisContainer ;
113 }
114
115 if (fReader) delete fReader ;
116
765d44e7 117 if (fCaloUtils) delete fCaloUtils ;
118
119
477d6cee 120 if(fAODBranchList){
1c5acb87 121// for(Int_t iaod = 0; iaod < fAODBranchList->GetEntries(); iaod++)
122// fAODBranchList->At(iaod)->Clear();
123
477d6cee 124 fAODBranchList->Clear();
125 delete fAODBranchList ;
126 }
127
1c5acb87 128}
129
130//________________________________________________________________________
131TList * AliAnaPartCorrMaker::GetAODBranchList()
132{
133
134// Get any new output AOD branches from analysis and put them in a list
135// The list is filled in the maker, and new branch passed to the analysis frame
136// AliAnalysisTaskPartCorr
137
477d6cee 138 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++){
139
140 AliAnaPartCorrBaseClass * ana = ((AliAnaPartCorrBaseClass *) fAnalysisContainer->At(iana)) ;
141 if(ana->NewOutputAOD()) fAODBranchList->Add(ana->GetCreateOutputAODBranch());
142 }
143
144 return fAODBranchList ;
145
1c5acb87 146}
147
2e557d1c 148//________________________________________________________________________
149TList *AliAnaPartCorrMaker::GetOutputContainer()
150{
151// Fill the output list of histograms during the CreateOutputObjects stage.
477d6cee 152 if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0){
a8a55c9d 153 printf("AliAnaPartCorrMaker::GetOutputContainer() - Analysis job list not initialized!!!\n");
154 //abort();
477d6cee 155 }
a3aebfff 156
765d44e7 157 //Initialize the geometry pointers
158 //GetCaloUtils()->InitPHOSGeometry();
159 //GetCaloUtils()->InitEMCALGeometry();
160
a3aebfff 161 char newname[128];
477d6cee 162 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++){
163 AliAnaPartCorrBaseClass * ana = ((AliAnaPartCorrBaseClass *) fAnalysisContainer->At(iana)) ;
164 if(fMakeHisto){// Analysis with histograms as output on
165 //Fill container with appropriate histograms
1cd71065 166 TList * templist = ana ->GetCreateOutputObjects();
a33b76ed 167 templist->SetOwner(kFALSE); //Owner is fOutputContainer.
a3aebfff 168 for(Int_t i = 0; i < templist->GetEntries(); i++){
169
170 //Add only to the histogram name the name of the task
171 if( strcmp((templist->At(i))->ClassName(),"TObjString") ) {
172 sprintf(newname,"%s%s", (ana->GetAddedHistogramsStringToName()).Data(), (templist->At(i))->GetName());
173 ((TH1*) templist->At(i))->SetName(newname);
174 }
175 //Add histogram to general container
477d6cee 176 fOutputContainer->Add(templist->At(i)) ;
a3aebfff 177 }
a33b76ed 178 delete templist;
477d6cee 179 }// Analysis with histograms as output on
180 }//Loop on analysis defined
a3aebfff 181
031ac63e 182
183 fhNEvents = new TH1I("hNEvents", "Number of analyzed events" , 1 , 0 , 1 ) ;
184 fOutputContainer->Add(fhNEvents);
185
477d6cee 186 return fOutputContainer;
a3aebfff 187
2e557d1c 188}
189
1c5acb87 190//________________________________________________________________________
191void AliAnaPartCorrMaker::Init()
192{
477d6cee 193 //Init container histograms and other common variables
194 // Fill the output list of histograms during the CreateOutputObjects stage.
195
196 if(!fAnalysisContainer || fAnalysisContainer->GetEntries()==0){
a8a55c9d 197 printf("AliAnaPartCorrMaker::GetOutputInit() - Analysis job list not initialized!!!\n");
198 //abort();
477d6cee 199 }
591cc579 200
765d44e7 201 //Initialize the geometry pointers
202 GetCaloUtils()->InitPHOSGeometry();
203 GetCaloUtils()->InitEMCALGeometry();
204
591cc579 205 //Initialize reader
206 fReader->Init();
765d44e7 207 fReader->SetCaloUtils(fCaloUtils); // pass the calo utils pointer to the reader
591cc579 208
765d44e7 209 //fCaloUtils->Init();
477d6cee 210 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++){
211
212 AliAnaPartCorrBaseClass * ana = ((AliAnaPartCorrBaseClass *) fAnalysisContainer->At(iana)) ;
213 ana->SetReader(fReader); //SetReader for each analysis
765d44e7 214 ana->SetCaloUtils(fCaloUtils); //Set CaloUtils for each analysis
215
477d6cee 216 ana->Init();
217
218 }//Loop on analysis defined
1c5acb87 219}
220
221//____________________________________________________________________________
222void AliAnaPartCorrMaker::InitParameters()
477d6cee 223{
224 //Init data members
225
226 fMakeHisto = kTRUE;
1cd71065 227 fMakeAOD = kTRUE;
228 fAnaDebug = 0; // No debugging info displayed by default
591cc579 229
1c5acb87 230}
231
232//__________________________________________________________________
233void AliAnaPartCorrMaker::Print(const Option_t * opt) const
477d6cee 234{
235 //Print some relevant parameters set for the analysis
236
237 if(! opt)
238 return;
239
240 printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
241 printf("Debug level = %d\n", fAnaDebug) ;
242 printf("Produce Histo = %d\n", fMakeHisto) ;
243 printf("Produce AOD = %d\n", fMakeAOD) ;
1cd71065 244 printf("Number of analysis tasks = %d\n", fAnalysisContainer->GetEntries()) ;
245 if(!strcmp("all",opt)){
246 printf("Print analysis Tasks settings :\n") ;
247 for(Int_t iana = 0; iana<fAnalysisContainer->GetEntries(); iana++){
248 ((AliAnaPartCorrBaseClass *) fAnalysisContainer->At(iana))->Print("");
249 }
250
251 printf("Print analysis Reader settings :\n") ;
252 fReader->Print("");
765d44e7 253 printf("Print analysis Calorimeter Utils settings :\n") ;
254 fCaloUtils->Print("");
1cd71065 255 }
1c5acb87 256}
257
258
259//____________________________________________________________________________
29b2ceec 260void AliAnaPartCorrMaker::ProcessEvent(const Int_t iEntry, const char * currentFileName){
477d6cee 261 //Process analysis for this event
262
263 if(fMakeHisto && !fOutputContainer){
264 printf("AliAnaPartCorrMaker::ProcessEvent() - Histograms not initialized\n");
265 abort();
266 }
a79a2424 267
902aa95c 268 if(fAnaDebug >= 0 ){
a79a2424 269 printf("*** Event %d *** \n",iEntry);
902aa95c 270 if(fAnaDebug > 1 ) {
271 printf("AliAnaPartCorrMaker::ProcessEvent() - Current File Name : %s\n", currentFileName);
272 //printf("fAODBranchList %p, entries %d\n",fAODBranchList,fAODBranchList->GetEntries());
273 }
a79a2424 274 }
902aa95c 275 //Each event needs an empty branch
276 Int_t nAODBranches = fAODBranchList->GetEntries();
4a745797 277 for(Int_t iaod = 0; iaod < nAODBranches; iaod++){
278 //fAODBranchList->At(iaod)->Clear();
279 TClonesArray *tca = dynamic_cast<TClonesArray*> (fAODBranchList->At(iaod));
280 if(tca) tca->Delete();
281 }
902aa95c 282
477d6cee 283 //Tell the reader to fill the data in the 3 detector lists
29b2ceec 284 Bool_t ok = fReader->FillInputEvent(iEntry, currentFileName);
285 if(!ok){
72d2488e 286 if(fAnaDebug >= 1 )printf("*** Skip event *** %d \n",iEntry);
29b2ceec 287 return ;
288 }
902aa95c 289
765d44e7 290 fCaloUtils->SetGeometryTransformationMatrices(fReader->GetInputEvent());
291
1d4cea01 292 //printf(">>>>>>>>>> BEFORE >>>>>>>>>>>\n");
293 //gObjectTable->Print();
477d6cee 294 //Loop on analysis algorithms
295 if(fAnaDebug > 0 ) printf("*** Begin analysis *** \n");
296 Int_t nana = fAnalysisContainer->GetEntries() ;
297 for(Int_t iana = 0; iana < nana; iana++){
298 AliAnaPartCorrBaseClass * ana = ((AliAnaPartCorrBaseClass *) fAnalysisContainer->At(iana)) ;
1e86c71e 299
300 ana->ConnectInputOutputAODBranches(); //Sets branches for each analysis
477d6cee 301 //Make analysis, create aods in aod branch or AODCaloClusters
1e86c71e 302 if(fMakeAOD) ana->MakeAnalysisFillAOD() ;
477d6cee 303 //Make further analysis with aod branch and fill histograms
304 if(fMakeHisto) ana->MakeAnalysisFillHistograms() ;
305
306 }
031ac63e 307
308 fhNEvents->Fill(0); //Event analyzed
309
477d6cee 310 fReader->ResetLists();
311
1d4cea01 312 //printf(">>>>>>>>>> AFTER >>>>>>>>>>>\n");
313 //gObjectTable->Print();
5efec477 314
477d6cee 315 if(fAnaDebug > 0 ) printf("*** End analysis *** \n");
316
1c5acb87 317}
6639984f 318
319//________________________________________________________________________
a5cc4f03 320void AliAnaPartCorrMaker::Terminate(TList * outputList)
6639984f 321{
477d6cee 322 //Execute Terminate of analysis
323 //Do some final plots.
324
a5cc4f03 325 if (!outputList) {
326 Error("Terminate", "No output list");
327 return;
328 }
329
477d6cee 330 for(Int_t iana = 0; iana < fAnalysisContainer->GetEntries(); iana++){
331
332 AliAnaPartCorrBaseClass * ana = ((AliAnaPartCorrBaseClass *) fAnalysisContainer->At(iana)) ;
a5cc4f03 333 ana->Terminate(outputList);
477d6cee 334
335 }//Loop on analysis defined
6639984f 336}