]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/EMCALTasks/AliAnalysisTaskEMCALClusterize.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGGA / EMCALTasks / AliAnalysisTaskEMCALClusterize.cxx
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
16 //_________________________________________________________________________
17 // This analysis provides a new list of clusters to be used in other analysis
18 //
19 // Author: Gustavo Conesa Balbastre,
20 //         Adapted from analysis class from Deepa Thomas
21 //
22 // $Id$
23 //_________________________________________________________________________
24
25 // --- Root ---
26 #include <TString.h>
27 #include <TRefArray.h>
28 #include <TClonesArray.h>
29 #include <TTree.h>
30 #include <TGeoManager.h>
31 #include <TROOT.h>
32 #include <TInterpreter.h>
33 #include <TFile.h>
34
35 // --- AliRoot Analysis Steering
36 #include "AliAnalysisTask.h"
37 #include "AliAnalysisManager.h"
38 #include "AliESDEvent.h"
39 #include "AliGeomManager.h"
40 #include "AliVCaloCells.h"
41 #include "AliAODCaloCluster.h"
42 #include "AliCDBManager.h"
43 #include "AliCDBStorage.h"
44 #include "AliCDBEntry.h"
45 #include "AliLog.h"
46 #include "AliVEventHandler.h"
47 #include "AliAODInputHandler.h"
48 #include "AliOADBContainer.h"
49 #include "AliAODMCParticle.h"
50
51 // --- EMCAL
52 #include "AliEMCALAfterBurnerUF.h"
53 #include "AliEMCALGeometry.h"
54 #include "AliEMCALClusterizerNxN.h"
55 #include "AliEMCALClusterizerv1.h"
56 #include "AliEMCALClusterizerv2.h"
57 #include "AliEMCALRecPoint.h"
58 #include "AliEMCALDigit.h"
59 #include "AliCaloCalibPedestal.h"
60 #include "AliEMCALCalibData.h"
61
62 #include "AliAnalysisTaskEMCALClusterize.h"
63
64 ClassImp(AliAnalysisTaskEMCALClusterize)
65
66 //______________________________________________________________________________
67 AliAnalysisTaskEMCALClusterize::AliAnalysisTaskEMCALClusterize(const char *name) 
68 : AliAnalysisTaskSE(name)
69 , fEvent(0)
70 , fGeom(0),               fGeomName("") 
71 , fGeomMatrixSet(kFALSE), fLoadGeomMatrices(kFALSE)
72 , fCalibData(0),          fPedestalData(0)
73 , fOCDBpath(""),          fAccessOCDB(kFALSE)
74 , fDigitsArr(0),          fClusterArr(0),             fCaloClusterArr(0)
75 , fRecParam(0),           fClusterizer(0)
76 , fUnfolder(0),           fJustUnfold(kFALSE) 
77 , fOutputAODBranch(0),    fOutputAODBranchName(""),   fOutputAODBranchSet(0)
78 , fFillAODFile(kFALSE),   fFillAODHeader(0)
79 , fFillAODCaloCells(0),   fRun(-1)
80 , fRecoUtils(0),          fConfigName("")
81 , fOrgClusterCellId()
82 , fCellLabels(),          fCellSecondLabels(),        fCellTime()
83 , fCellMatchdEta(),       fCellMatchdPhi()
84 , fRecalibrateWithClusterTime(0)
85 , fMaxEvent(0),           fDoTrackMatching(kFALSE)
86 , fSelectCell(kFALSE),    fSelectCellMinE(0),         fSelectCellMinFrac(0)
87 , fRejectBelowThreshold(kFALSE)
88 , fRemoveLEDEvents(kTRUE),fRemoveExoticEvents(kFALSE)
89 , fImportGeometryFromFile(kFALSE), fImportGeometryFilePath("") 
90 , fOADBSet(kFALSE),       fAccessOADB(kTRUE),         fOADBFilePath("")
91 , fCentralityClass(""),   fSelectEMCALEvent(0)
92 , fEMCALEnergyCut(0.),    fEMCALNcellsCut (0)
93 , fSetCellMCLabelFromCluster(0)
94 , fRemapMCLabelForAODs(0)
95 , fInputFromFilter(0)
96 {
97   // Constructor
98   
99   for(Int_t i = 0; i < 12;    i++)  fGeomMatrix[i] =  0;
100   
101   ResetArrays();
102   
103   fCentralityBin[0] = fCentralityBin[1]=-1;
104   
105 }
106
107 //______________________________________________________________
108 AliAnalysisTaskEMCALClusterize::AliAnalysisTaskEMCALClusterize() 
109 : AliAnalysisTaskSE("DefaultAnalysis_AliAnalysisTaskEMCALClusterize")
110 , fEvent(0)
111 , fGeom(0),                 fGeomName("") 
112 , fGeomMatrixSet(kFALSE),   fLoadGeomMatrices(kFALSE)
113 , fCalibData(0),            fPedestalData(0)
114 , fOCDBpath(""),            fAccessOCDB(kFALSE)
115 , fDigitsArr(0),            fClusterArr(0),             fCaloClusterArr(0)
116 , fRecParam(0),             fClusterizer(0)
117 , fUnfolder(0),             fJustUnfold(kFALSE) 
118 , fOutputAODBranch(0),      fOutputAODBranchName(""),   fOutputAODBranchSet(0)
119 , fFillAODFile(kFALSE),     fFillAODHeader(0)
120 , fFillAODCaloCells(0),     fRun(-1)
121 , fRecoUtils(0),            fConfigName("")
122 , fOrgClusterCellId()
123 , fCellLabels(),            fCellSecondLabels(),        fCellTime()
124 , fCellMatchdEta(),         fCellMatchdPhi()
125 , fRecalibrateWithClusterTime(0)
126 , fMaxEvent(0),             fDoTrackMatching(kFALSE)
127 , fSelectCell(kFALSE),      fSelectCellMinE(0),         fSelectCellMinFrac(0)
128 , fRejectBelowThreshold(kFALSE)
129 , fRemoveLEDEvents(kTRUE),  fRemoveExoticEvents(kFALSE)
130 , fImportGeometryFromFile(kFALSE), fImportGeometryFilePath("")
131 , fOADBSet(kFALSE),         fAccessOADB(kTRUE),        fOADBFilePath("")
132 , fCentralityClass(""),     fSelectEMCALEvent(0)
133 , fEMCALEnergyCut(0.),      fEMCALNcellsCut (0)
134 , fSetCellMCLabelFromCluster(0)
135 , fRemapMCLabelForAODs(0)
136 , fInputFromFilter(0)
137 {
138   // Constructor
139   
140   for(Int_t i = 0; i < 12;    i++)  fGeomMatrix[i] =  0;
141   
142   ResetArrays();
143   
144   fCentralityBin[0] = fCentralityBin[1]=-1;
145
146 }
147
148
149 //_______________________________________________________________
150 AliAnalysisTaskEMCALClusterize::~AliAnalysisTaskEMCALClusterize()
151 {
152   //dtor 
153   
154   if (fDigitsArr)
155   {
156     fDigitsArr->Clear("C");
157     delete fDigitsArr; 
158   }
159   
160   if (fClusterArr)
161   {
162     fClusterArr->Delete();
163     delete fClusterArr;
164   }
165   
166   if (fCaloClusterArr)
167   {
168     fCaloClusterArr->Delete();
169     delete fCaloClusterArr; 
170   }
171   
172   if(fClusterizer) delete fClusterizer;
173   if(fUnfolder)    delete fUnfolder;   
174   if(fRecoUtils)   delete fRecoUtils;
175   
176 }
177
178 //_______________________________________________________
179 Bool_t AliAnalysisTaskEMCALClusterize::AcceptEventEMCAL()
180 {
181   // Accept event given there is a EMCAL cluster with enough energy, and not noisy, exotic
182   
183   if(!fSelectEMCALEvent)   return kTRUE; // accept
184   
185   if(fEMCALEnergyCut <= 0) return kTRUE; // accept
186   
187   Int_t           nCluster = fEvent -> GetNumberOfCaloClusters();
188   AliVCaloCells * caloCell = fEvent -> GetEMCALCells();
189   Int_t           bc       = fEvent -> GetBunchCrossNumber();
190   
191   for(Int_t icalo = 0; icalo < nCluster; icalo++)
192   {
193     AliVCluster *clus = (AliVCluster*) (fEvent->GetCaloCluster(icalo));
194     
195     if( ( clus->IsEMCAL() ) && ( clus->GetNCells() > fEMCALNcellsCut ) && ( clus->E() > fEMCALEnergyCut ) &&
196        fRecoUtils->IsGoodCluster(clus,fGeom,caloCell,bc))
197     {
198       
199       if (fDebug > 0)
200         printf("AliAnalysisTaskEMCALClusterize::AcceptEventEMCAL() - Accept :  E %2.2f > %2.2f, nCells %d > %d \n",
201                clus->E(), fEMCALEnergyCut, clus->GetNCells(), fEMCALNcellsCut);
202       
203       return kTRUE;
204     }
205     
206   }// loop
207   
208   if (fDebug > 0)
209     printf("AliAnalysisTaskEMCALClusterize::AcceptEventEMCAL() - Reject \n");
210   
211   return kFALSE;
212   
213 }
214
215 //_______________________________________________
216 void AliAnalysisTaskEMCALClusterize::AccessOADB()
217 {
218   // Set the AODB calibration, bad channels etc. parameters at least once
219   // alignment matrices from OADB done in SetGeometryMatrices
220   
221   //Set it only once
222   if(fOADBSet) return ; 
223   
224   Int_t   runnumber = InputEvent()->GetRunNumber() ;
225   TString pass      = GetPass();
226   
227   printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Get AODB parameters from EMCAL in %s for run %d, and <%s> \n",fOADBFilePath.Data(),runnumber,pass.Data());
228   
229   Int_t nSM = fGeom->GetNumberOfSuperModules();
230   
231   // Bad map
232   if(fRecoUtils->IsBadChannelsRemovalSwitchedOn())
233   {
234     AliOADBContainer *contBC=new AliOADBContainer("");
235     contBC->InitFromFile(Form("%s/EMCALBadChannels.root",fOADBFilePath.Data()),"AliEMCALBadChannels"); 
236     
237     TObjArray *arrayBC=(TObjArray*)contBC->GetObject(runnumber);
238     
239     if(arrayBC)
240     {
241         fRecoUtils->SwitchOnDistToBadChannelRecalculation();
242         printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Remove EMCAL bad cells \n");
243         
244         for (Int_t i=0; i<nSM; ++i) 
245         {
246           TH2I *hbm = fRecoUtils->GetEMCALChannelStatusMap(i);
247           
248           if (hbm)
249             delete hbm;
250           
251           hbm=(TH2I*)arrayBC->FindObject(Form("EMCALBadChannelMap_Mod%d",i));
252           
253           if (!hbm) 
254           {
255             AliError(Form("Can not get EMCALBadChannelMap_Mod%d",i));
256             continue;
257           }
258           
259           hbm->SetDirectory(0);
260           fRecoUtils->SetEMCALChannelStatusMap(i,hbm);
261           
262         } // loop
263     } else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT remove EMCAL bad channels\n"); // run array
264   }  // Remove bad
265   
266   // Energy Recalibration
267   if(fRecoUtils->IsRecalibrationOn())
268   {
269     AliOADBContainer *contRF=new AliOADBContainer("");
270     
271     contRF->InitFromFile(Form("%s/EMCALRecalib.root",fOADBFilePath.Data()),"AliEMCALRecalib");
272     
273     TObjArray *recal=(TObjArray*)contRF->GetObject(runnumber); 
274     
275     if(recal)
276     {
277       TObjArray *recalpass=(TObjArray*)recal->FindObject(pass);
278       
279       if(recalpass)
280       {
281         TObjArray *recalib=(TObjArray*)recalpass->FindObject("Recalib");
282         
283         if(recalib)
284         {
285           printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Recalibrate EMCAL \n");
286           for (Int_t i=0; i<nSM; ++i) 
287           {
288             TH2F *h = fRecoUtils->GetEMCALChannelRecalibrationFactors(i);
289             
290             if (h)
291               delete h;
292             
293             h = (TH2F*)recalib->FindObject(Form("EMCALRecalFactors_SM%d",i));
294             
295             if (!h) 
296             {
297               AliError(Form("Could not load EMCALRecalFactors_SM%d",i));
298               continue;
299             }
300             
301             h->SetDirectory(0);
302             
303             fRecoUtils->SetEMCALChannelRecalibrationFactors(i,h);
304           } // SM loop
305         }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL, no params object array \n"); // array ok
306       }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL, no params for pass\n"); // array pass ok
307     }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL, no params for run\n");  // run number array ok
308         
309   } // Recalibration on
310   
311   // Energy Recalibration, apply on top of previous calibration factors
312   if(fRecoUtils->IsRunDepRecalibrationOn())
313   {
314     AliOADBContainer *contRFTD=new AliOADBContainer("");
315     
316     contRFTD->InitFromFile(Form("%s/EMCALTemperatureCorrCalib.root",fOADBFilePath.Data()),"AliEMCALRunDepTempCalibCorrections");
317     
318     TH1S *htd=(TH1S*)contRFTD->GetObject(runnumber); 
319     
320     //If it did not exist for this run, get closes one
321     if (!htd)
322     {
323       AliWarning(Form("No TemperatureCorrCalib Objects for run: %d",runnumber));
324       // let's get the closest runnumber instead then..
325       Int_t lower = 0;
326       Int_t ic = 0;
327       Int_t maxEntry = contRFTD->GetNumberOfEntries();
328       
329       while ( (ic < maxEntry) && (contRFTD->UpperLimit(ic) < runnumber) ) {
330         lower = ic;
331         ic++;
332       }
333       
334       Int_t closest = lower;
335       if ( (ic<maxEntry) &&
336           (contRFTD->LowerLimit(ic)-runnumber) < (runnumber - contRFTD->UpperLimit(lower)) ) {
337         closest = ic;
338       }
339       
340       AliWarning(Form("TemperatureCorrCalib Objects found closest id %d from run: %d", closest, contRFTD->LowerLimit(closest)));
341       htd = (TH1S*) contRFTD->GetObjectByIndex(closest);
342     } 
343
344     if(htd)
345     {
346       printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Recalibrate (Temperature) EMCAL \n");
347       
348       for (Int_t ism=0; ism<nSM; ++ism) 
349       {        
350         for (Int_t icol=0; icol<48; ++icol) 
351         {        
352           for (Int_t irow=0; irow<24; ++irow) 
353           {
354             Float_t factor = fRecoUtils->GetEMCALChannelRecalibrationFactor(ism,icol,irow);
355             
356             Int_t absID = fGeom->GetAbsCellIdFromCellIndexes(ism, irow, icol); // original calibration factor
357             factor *= htd->GetBinContent(absID) / 10000. ; // correction dependent on T
358             //printf("\t ism %d, icol %d, irow %d,absID %d, corrA %2.3f, corrB %2.3f, corrAB %2.3f\n",ism, icol, irow, absID, 
359             //      GetEMCALChannelRecalibrationFactor(ism,icol,irow) , htd->GetBinContent(absID) / 10000., factor);
360             fRecoUtils->SetEMCALChannelRecalibrationFactor(ism,icol,irow,factor);
361           } // columns
362         } // rows 
363       } // SM loop
364     }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate EMCAL with T variations, no params TH1 \n"); 
365   } // Run by Run T calibration
366   
367   // Time Recalibration
368   if(fRecoUtils->IsTimeRecalibrationOn())
369   {
370     AliOADBContainer *contTRF=new AliOADBContainer("");
371     
372     contTRF->InitFromFile(Form("%s/EMCALTimeCalib.root",fOADBFilePath.Data()),"AliEMCALTimeCalib");
373     
374     TObjArray *trecal=(TObjArray*)contTRF->GetObject(runnumber); 
375     
376     if(trecal)
377     {
378       TObjArray *trecalpass=(TObjArray*)trecal->FindObject(pass);
379
380       if(trecalpass)
381       {
382         printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Time Recalibrate EMCAL \n");
383         for (Int_t ibc = 0; ibc < 4; ++ibc) 
384         {
385           TH1F *h = fRecoUtils->GetEMCALChannelTimeRecalibrationFactors(ibc);
386           
387           if (h)
388             delete h;
389           
390           h = (TH1F*)trecalpass->FindObject(Form("hAllTimeAvBC%d",ibc));
391           
392           if (!h) 
393           {
394             AliError(Form("Could not load hAllTimeAvBC%d",ibc));
395             continue;
396           }
397           
398           h->SetDirectory(0);
399           
400           fRecoUtils->SetEMCALChannelTimeRecalibrationFactors(ibc,h);
401         } // bunch crossing loop
402       }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate time EMCAL, no params for pass\n"); // array pass ok
403     }else printf("AliAnalysisTaskEMCALClusterize::SetOADBParameters() - Do NOT recalibrate time EMCAL, no params for run\n");  // run number array ok
404     
405   } // Time recalibration on    
406   
407   // Parameters already set once, so do not it again
408   fOADBSet = kTRUE;
409   
410 }  
411
412 //_________________________________________________
413 Bool_t AliAnalysisTaskEMCALClusterize::AccessOCDB()
414 {
415   //Access to OCDB stuff
416   
417   fEvent = InputEvent();
418   if (!fEvent)
419   {
420     Warning("AccessOCDB","Event not available!!!");
421     return kFALSE;
422   }
423   
424   if (fEvent->GetRunNumber()==fRun)
425     return kTRUE;
426   fRun = fEvent->GetRunNumber();
427   
428   if(DebugLevel() > 1 )
429     printf("AliAnalysisTaksEMCALClusterize::AccessOCDB() - Begin");
430   
431   AliCDBManager *cdb = AliCDBManager::Instance();
432   
433   
434   if (fOCDBpath.Length())
435   {
436     cdb->SetDefaultStorage(fOCDBpath.Data());
437     printf("AliAnalysisTaksEMCALClusterize::AccessOCDB() - Default storage %s",fOCDBpath.Data());
438   }
439   
440   cdb->SetRun(fEvent->GetRunNumber());
441   
442   //
443   // EMCAL from RAW OCDB
444   if (fOCDBpath.Contains("alien:"))
445   {
446     cdb->SetSpecificStorage("EMCAL/Calib/Data","alien://Folder=/alice/data/2010/OCDB");
447     cdb->SetSpecificStorage("EMCAL/Calib/Pedestals","alien://Folder=/alice/data/2010/OCDB");
448   }
449   
450   TString path = cdb->GetDefaultStorage()->GetBaseFolder();
451   
452   // init parameters:
453   
454   //Get calibration parameters  
455   if(!fCalibData)
456   {
457     AliCDBEntry *entry = (AliCDBEntry*) 
458     AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
459     
460     if (entry) fCalibData =  (AliEMCALCalibData*) entry->GetObject();
461   }
462   
463   if(!fCalibData)
464     AliFatal("Calibration parameters not found in CDB!");
465   
466   //Get calibration parameters  
467   if(!fPedestalData)
468   {
469     AliCDBEntry *entry = (AliCDBEntry*) 
470     AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
471     
472     if (entry) fPedestalData =  (AliCaloCalibPedestal*) entry->GetObject();
473   }
474   
475   if(!fPedestalData)
476     AliFatal("Dead map not found in CDB!");
477   
478   return kTRUE;
479 }
480
481 //_____________________________________________________
482 void AliAnalysisTaskEMCALClusterize::CheckAndGetEvent()
483 {
484   // Get the input event, it can depend in embedded events what you want to get
485   // Also check if the quality of the event is good if not reject it
486   
487   fEvent = 0x0;
488       
489   AliAODInputHandler* aodIH = dynamic_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
490   Int_t eventN = Entry();
491   if(aodIH) eventN = aodIH->GetReadEntry(); 
492   
493   if (eventN > fMaxEvent) 
494     return ;
495   
496   //printf("Clusterizer --- Event %d-- \n",eventN);
497   
498   //Check if input event are embedded events
499   //If so, take output event
500   if (aodIH && aodIH->GetMergeEvents()) 
501   {
502     fEvent  = AODEvent();
503     
504     if(!aodIH->GetMergeEMCALCells()) 
505       AliFatal("Events merged but not EMCAL cells, check analysis settings!");
506     
507     if(DebugLevel() > 1)
508     {
509       printf("AliAnalysisTaksEMCALClusterize::UserExec() - Use embedded events\n");
510       
511       printf("\t InputEvent  N Clusters %d, N Cells %d\n",InputEvent()->GetNumberOfCaloClusters(),
512              InputEvent()->GetEMCALCells()->GetNumberOfCells());
513       
514       printf("\t MergedEvent  N Clusters %d, N Cells %d\n",aodIH->GetEventToMerge()->GetNumberOfCaloClusters(),
515              aodIH->GetEventToMerge()->GetEMCALCells()->GetNumberOfCells());
516       
517       for (Int_t icl=0; icl < aodIH->GetEventToMerge()->GetNumberOfCaloClusters(); icl++) 
518       {
519         AliAODCaloCluster *sigCluster = aodIH->GetEventToMerge()->GetCaloCluster(icl);
520         if(sigCluster->IsEMCAL()) printf("\t \t Signal cluster: i %d, E  %f\n",icl,sigCluster->E());
521       }
522       
523       printf("\t OutputEvent N Clusters %d, N Cells %d\n", AODEvent()->GetNumberOfCaloClusters(),
524              AODEvent()->GetEMCALCells()->GetNumberOfCells());
525     }
526   }
527   else if(fInputFromFilter)
528   {
529     //printf("Get Input From Filtered AOD\n");
530     fEvent =  AODEvent();
531   }
532   else 
533   {
534     fEvent =  InputEvent();
535     if(fFillAODCaloCells) FillAODCaloCells();   
536     if(fFillAODHeader)    FillAODHeader();
537   }
538   
539   if (!fEvent) 
540   {
541     Error("UserExec","Event not available");
542     return ;
543   }
544   
545   //Process events if there is a high energy cluster
546   if(!AcceptEventEMCAL())  { fEvent = 0x0 ; return ; }
547   
548   //-------------------------------------------------------------------------------------
549   // Reject events if LED was firing, use only for LHC11a data 
550   // Reject event if triggered by exotic cell and remove exotic cells if not triggered
551   //-------------------------------------------------------------------------------------
552   
553   if( IsLEDEvent( InputEvent()->GetRunNumber() ) ) { fEvent = 0x0 ; return ; }
554   
555   if( IsExoticEvent() )                            { fEvent = 0x0 ; return ; }
556   
557   //-------------------------------------------------------------------------------------
558   // Set the cluster array in the event (output or input)
559   //-------------------------------------------------------------------------------------
560   
561   if     ( fFillAODFile ) 
562   {
563     //Magic line to write events to AOD filem put after event rejection
564     AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kTRUE);
565   }
566   else if( !fOutputAODBranchSet )
567   {
568     // Create array and put it in the input event, if output AOD not selected, only once
569     InputEvent()->AddObject(fOutputAODBranch);
570     fOutputAODBranchSet = kTRUE;
571     printf("AliAnalysisTaskEMCALClusterize::UserExec() - Add AOD branch <%s> to input event\n",fOutputAODBranchName.Data());
572   }
573   
574 }
575
576 //____________________________________________________
577 void AliAnalysisTaskEMCALClusterize::ClusterizeCells()
578
579   // Recluster calocells, transform them into digits, 
580   // feed the clusterizer with them and get new list of clusters
581   
582   //In case of MC, first loop on the clusters and fill MC label to array  
583   Int_t nClusters     = fEvent->GetNumberOfCaloClusters();
584   Int_t nClustersOrg  = 0;
585   
586   AliAODInputHandler* aodIH = dynamic_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
587   if(aodIH && aodIH->GetEventToMerge())  //Embedding
588     nClusters = aodIH->GetEventToMerge()->GetNumberOfCaloClusters(); //Get clusters directly from embedded signal
589
590   ResetArrays();
591    
592   for (Int_t i = 0; i < nClusters; i++)
593   {
594     AliVCluster *clus = 0;
595     if(aodIH && aodIH->GetEventToMerge()) //Embedding
596       clus = aodIH->GetEventToMerge()->GetCaloCluster(i); //Get clusters directly from embedded signal
597     else      
598       clus = fEvent->GetCaloCluster(i);
599     
600     if(!clus) return;
601     
602     if(clus->IsEMCAL())
603     {
604       Int_t label = clus->GetLabel();
605       Int_t label2 = -1 ;
606       //printf("Org cluster E %f, Time  %e, Index = %d, ID %d, MC label %d\n", clus->E(), clus->GetTOF(),i, clus->GetID(),label );
607       //printf("Original list of labels from old cluster : \n");
608       //for(Int_t imc = 0; imc < clus->GetNLabels(); imc++) printf("\t Label %d\n",clus->GetLabelAt(imc));
609       
610       if (clus->GetNLabels()>=2) label2 = clus->GetLabelAt(1) ;
611       UShort_t * index    = clus->GetCellsAbsId() ;
612       for(Int_t icell=0; icell < clus->GetNCells(); icell++ )
613       {
614         //printf("\t cell %d, MC label %d\n",index[icell],fEvent->GetEMCALCells()->GetCellMCLabel(index[icell]));
615         fOrgClusterCellId[index[icell]] = i;
616         fCellLabels[index[icell]]       = label;
617         fCellSecondLabels[index[icell]] = label2;
618         fCellTime[index[icell]]         = clus->GetTOF();
619         fCellMatchdEta[index[icell]]    = clus->GetTrackDz();
620         fCellMatchdPhi[index[icell]]    = clus->GetTrackDx();
621       }
622       nClustersOrg++;
623     }
624     // printf("\n");
625   } 
626   
627   // Transform CaloCells into Digits
628   
629   Int_t    idigit =  0;
630   Int_t    id     = -1;
631   Float_t  amp    = -1; 
632   Double_t time   = -1; 
633   
634   AliVCaloCells *cells = fEvent->GetEMCALCells();
635   
636   Int_t bc = InputEvent()->GetBunchCrossNumber();
637
638   for (Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
639   {
640     // Get cell values, recalibrate and not include bad channels found in analysis, nor cells with too low energy, nor exotic cell
641     id = cells->GetCellNumber(icell);
642     Bool_t accept = fRecoUtils->AcceptCalibrateCell(id,bc,amp,time,cells);
643     
644     // Do not include cells with too low energy, nor exotic cell
645     if( amp  < fRecParam->GetMinECut() ||
646         time > fRecParam->GetTimeMax() ||
647         time < fRecParam->GetTimeMin()    ) accept = kFALSE;
648     
649     // In case of old AOD analysis cell time is -1 s, approximate replacing by time of the cluster the digit belongs.
650     if (fRecalibrateWithClusterTime)
651     { 
652       time = fCellTime[id];
653       //printf("cell %d time org %f - ",id, time*1.e9);
654       fRecoUtils->RecalibrateCellTime(id,bc,time);
655       //printf("recal %f\n",time*1.e9);
656     }
657     
658     //Exotic?
659     if (accept && fRecoUtils->IsExoticCell(id,cells,bc))
660         accept = kFALSE;
661     
662     if( !accept )
663     {
664       if( DebugLevel() > 2 )
665         printf("AliAnalysisTaksEMCALClusterize::ClusterizeCells() - Remove channel absId %d, index %d of %d, amp %f, time %f\n",
666                id,icell, cells->GetNumberOfCells(), amp, time*1.e9);
667       continue;
668     }
669     
670     Int_t mcLabel = cells->GetMCLabel(icell);
671     //printf("AliAnalysisTaksEMCALClusterize::ClusterizeCells() - cell %d, mc label %d\n",id,mcLabel);
672
673     //if(fCellLabels[id]!=mcLabel)printf("mcLabel %d - %d\n",mcLabel,fCellLabels[id]);
674     if     ( fSetCellMCLabelFromCluster == 1 ) mcLabel = fCellLabels[id]; // Older aliroot MC productions
675     else if( fSetCellMCLabelFromCluster == 0 && fRemapMCLabelForAODs) RemapMCLabelForAODs(mcLabel);
676     else mcLabel = -1; // found later
677     
678     //printf("\t new label %d\n",mcLabel);
679     
680     // Create the digit, put a fake primary deposited energy to trick the clusterizer
681     // when checking the most likely primary
682     
683     Float_t efrac = cells->GetEFraction(icell);
684     
685     //When checking the MC of digits, give weight to cells with embedded signal
686     if (mcLabel > 0 && efrac < 1.e-6) efrac = 1;
687     
688     //printf("******* Cell %d, id %d, e %f,  fraction %f, MC label %d, used MC label %d\n",icell,id,amp,cells->GetEFraction(icell),cells->GetMCLabel(icell),mcLabel);
689     
690     new((*fDigitsArr)[idigit]) AliEMCALDigit( mcLabel, mcLabel, id, amp, time,AliEMCALDigit::kHG,idigit, 0, 0, amp*efrac);
691     // Last parameter should be MC deposited energy, since it is not available, add just the cell amplitude so that
692     // we give more weight to the MC label of the cell with highest energy in the cluster
693         
694     idigit++;
695   }
696   
697   fDigitsArr->Sort();
698   
699   //-------------------------------------------------------------------------------------
700   //Do the clusterization
701   //-------------------------------------------------------------------------------------        
702   
703   fClusterizer->Digits2Clusters("");
704   
705   //-------------------------------------------------------------------------------------
706   //Transform the recpoints into AliVClusters
707   //-------------------------------------------------------------------------------------
708   
709   RecPoints2Clusters();
710   
711   if(!fCaloClusterArr)
712   { 
713     printf("AliAnalysisTaksEMCALClusterize::UserExec() - No array with CaloClusters, input RecPoints entries %d\n",fClusterArr->GetEntriesFast());
714     return;    
715   }
716   
717   if( DebugLevel() > 0 )
718   {
719     printf("AliAnalysisTaksEMCALClusterize::ClusterizeCells() - N clusters: before recluster %d, after recluster %d\n",nClustersOrg, fCaloClusterArr->GetEntriesFast());
720     
721     if(fCaloClusterArr->GetEntriesFast() != fClusterArr->GetEntriesFast())
722     {
723       printf("\t Some RecRoints not transformed into CaloClusters (clusterizer %d, unfold %d): Input entries %d - Output entries %d - %d (not fast)\n",
724              fRecParam->GetClusterizerFlag(),fRecParam->GetUnfold(),
725              fClusterArr->GetEntriesFast(), fCaloClusterArr->GetEntriesFast(), fCaloClusterArr->GetEntries());
726     }
727   }
728 }
729
730 //_____________________________________________________
731 void AliAnalysisTaskEMCALClusterize::ClusterUnfolding()
732 {
733   // Take the event clusters and unfold them
734   
735   AliVCaloCells *cells   = fEvent->GetEMCALCells();
736   Double_t cellAmplitude = 0;
737   Double_t cellTime      = 0;
738   Short_t  cellNumber    = 0;
739   Int_t    cellMCLabel   = 0;
740   Double_t cellEFrac     = 0;
741   Int_t    nClustersOrg  = 0;
742   
743   // Fill the array with the EMCAL clusters, copy them
744   for (Int_t i = 0; i < fEvent->GetNumberOfCaloClusters(); i++)
745   {
746     AliVCluster *clus = fEvent->GetCaloCluster(i);
747     if(clus->IsEMCAL())
748     {        
749       //recalibrate/remove bad channels/etc if requested
750       if(fRecoUtils->ClusterContainsBadChannel(fGeom,clus->GetCellsAbsId(), clus->GetNCells()))
751       {
752         continue;
753       } 
754       
755       if(fRecoUtils->IsRecalibrationOn())
756       {
757         //Calibrate cluster
758         fRecoUtils->RecalibrateClusterEnergy(fGeom, clus, cells);
759         
760         //CalibrateCells
761         for (Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
762         {
763           if (cells->GetCell(icell, cellNumber, cellAmplitude, cellTime, cellMCLabel, cellEFrac) != kTRUE)
764             break;
765           
766           Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1; 
767           fGeom->GetCellIndex(cellNumber,imod,iTower,iIphi,iIeta); 
768           fGeom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);       
769           
770           //Do not include bad channels found in analysis?
771           if( fRecoUtils->IsBadChannelsRemovalSwitchedOn() && 
772               fRecoUtils->GetEMCALChannelStatus(imod, ieta, iphi))
773             continue;
774           
775           cells->SetCell(icell, cellNumber, cellAmplitude*fRecoUtils->GetEMCALChannelRecalibrationFactor(imod,ieta,iphi),cellTime);
776           
777         }// cells loop            
778       }// recalibrate
779       
780       //Cast to ESD or AOD, needed to create the cluster array
781       AliESDCaloCluster * esdCluster = dynamic_cast<AliESDCaloCluster*> (clus);
782       AliAODCaloCluster * aodCluster = dynamic_cast<AliAODCaloCluster*> (clus);
783       
784       if     (esdCluster)
785       {
786         fCaloClusterArr->Add( new AliESDCaloCluster(*esdCluster) );   
787       }//ESD
788       else if(aodCluster)
789       {
790         fCaloClusterArr->Add( new AliAODCaloCluster(*aodCluster) );   
791       }//AOD
792       else 
793         Warning("UserExec()"," - Wrong CaloCluster type?");
794       
795       nClustersOrg++;
796     }
797   }
798   
799   //Do the unfolding
800   fUnfolder->UnfoldClusters(fCaloClusterArr, cells);
801   
802   //CLEAN-UP
803   fUnfolder->Clear();
804   
805 }
806
807 //_____________________________________________________
808 void AliAnalysisTaskEMCALClusterize::FillAODCaloCells() 
809 {
810   // Put calo cells in standard branch  
811   AliVCaloCells &eventEMcells = *(fEvent->GetEMCALCells());
812   Int_t nEMcell = eventEMcells.GetNumberOfCells() ;
813   
814   AliAODCaloCells &aodEMcells = *(AODEvent()->GetEMCALCells());
815   aodEMcells.CreateContainer(nEMcell);
816   aodEMcells.SetType(AliVCaloCells::kEMCALCell);
817   Double_t calibFactor = 1.;   
818   for (Int_t iCell = 0; iCell < nEMcell; iCell++) 
819   { 
820     Int_t imod = -1, iphi =-1, ieta=-1,iTower = -1, iIphi = -1, iIeta = -1; 
821     fGeom->GetCellIndex(eventEMcells.GetCellNumber(iCell),imod,iTower,iIphi,iIeta); 
822     fGeom->GetCellPhiEtaIndexInSModule(imod,iTower,iIphi, iIeta,iphi,ieta);     
823     
824     if(fRecoUtils->IsRecalibrationOn())
825     { 
826       calibFactor = fRecoUtils->GetEMCALChannelRecalibrationFactor(imod,ieta,iphi);
827     }
828     
829     if(!fRecoUtils->GetEMCALChannelStatus(imod, ieta, iphi))
830     { //Channel is not declared as bad
831       aodEMcells.SetCell(iCell,eventEMcells.GetCellNumber(iCell),eventEMcells.GetAmplitude(iCell)*calibFactor,
832                          eventEMcells.GetTime(iCell),eventEMcells.GetMCLabel(iCell),eventEMcells.GetEFraction(iCell));
833     }
834     else 
835     {
836       aodEMcells.SetCell(iCell,eventEMcells.GetCellNumber(iCell),0,-1,-1,0);
837     }
838   }
839   aodEMcells.Sort();
840   
841 }
842
843 //__________________________________________________
844 void AliAnalysisTaskEMCALClusterize::FillAODHeader() 
845 {
846   //Put event header information in standard AOD branch
847   
848   AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (fEvent);
849   AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (fEvent);
850   
851   Double_t pos[3]   ;
852   Double_t covVtx[6];
853   for (Int_t i = 0; i < 6; i++)  covVtx[i] = 0.;
854   
855   AliAODHeader* header = AODEvent()->GetHeader();
856   header->SetRunNumber(fEvent->GetRunNumber());
857   
858   if(esdevent)
859   {
860     TTree* tree = fInputHandler->GetTree();
861     if (tree) 
862     {
863       TFile* file = tree->GetCurrentFile();
864       if (file) header->SetESDFileName(file->GetName());
865     }
866   }
867   else if (aodevent) header->SetESDFileName(aodevent->GetHeader()->GetESDFileName());
868   
869   header->SetBunchCrossNumber(fEvent->GetBunchCrossNumber());
870   header->SetOrbitNumber(fEvent->GetOrbitNumber());
871   header->SetPeriodNumber(fEvent->GetPeriodNumber());
872   header->SetEventType(fEvent->GetEventType());
873   
874   //Centrality
875   if(fEvent->GetCentrality())
876     header->SetCentrality(new AliCentrality(*(fEvent->GetCentrality())));
877   else
878     header->SetCentrality(0);
879   
880   //Trigger  
881   header->SetOfflineTrigger(fInputHandler->IsEventSelected()); // propagate the decision of the physics selection
882
883   header->SetFiredTriggerClasses(fEvent->GetFiredTriggerClasses());
884   
885   header->SetTriggerMask(fEvent->GetTriggerMask()); 
886   header->SetTriggerCluster(fEvent->GetTriggerCluster());
887   
888   if      (esdevent)
889   {
890     header->SetL0TriggerInputs(esdevent->GetHeader()->GetL0TriggerInputs());    
891     header->SetL1TriggerInputs(esdevent->GetHeader()->GetL1TriggerInputs());    
892     header->SetL2TriggerInputs(esdevent->GetHeader()->GetL2TriggerInputs());    
893   }
894   else if (aodevent)
895   {
896     header->SetL0TriggerInputs(aodevent->GetHeader()->GetL0TriggerInputs());    
897     header->SetL1TriggerInputs(aodevent->GetHeader()->GetL1TriggerInputs());    
898     header->SetL2TriggerInputs(aodevent->GetHeader()->GetL2TriggerInputs());    
899   }
900   
901   header->SetMagneticField(fEvent->GetMagneticField());
902   //header->SetMuonMagFieldScale(esdevent->GetCurrentDip()/6000.); 
903   
904   header->SetZDCN1Energy(fEvent->GetZDCN1Energy());
905   header->SetZDCP1Energy(fEvent->GetZDCP1Energy());
906   header->SetZDCN2Energy(fEvent->GetZDCN2Energy());
907   header->SetZDCP2Energy(fEvent->GetZDCP2Energy());
908   header->SetZDCEMEnergy(fEvent->GetZDCEMEnergy(0),fEvent->GetZDCEMEnergy(1));
909   
910   Float_t diamxy[2]={(Float_t)fEvent->GetDiamondX(),(Float_t)fEvent->GetDiamondY()};
911   Float_t diamcov[3];
912   fEvent->GetDiamondCovXY(diamcov);
913   header->SetDiamond(diamxy,diamcov);
914   if      (esdevent) header->SetDiamondZ(esdevent->GetDiamondZ(),esdevent->GetSigma2DiamondZ());
915   else if (aodevent) header->SetDiamondZ(aodevent->GetDiamondZ(),aodevent->GetSigma2DiamondZ());
916   
917   //
918   Int_t nVertices = 1 ;/* = prim. vtx*/;
919   Int_t nCaloClus = fEvent->GetNumberOfCaloClusters();
920   
921   AODEvent()->ResetStd(0, nVertices, 0, 0, 0, nCaloClus, 0, 0);
922   
923   // Access to the AOD container of vertices
924   TClonesArray &vertices = *(AODEvent()->GetVertices());
925   Int_t jVertices=0;
926   
927   // Add primary vertex. The primary tracks will be defined
928   // after the loops on the composite objects (V0, cascades, kinks)
929   fEvent->GetPrimaryVertex()->GetXYZ(pos);
930   Float_t chi = 0;
931   if      (esdevent)
932   {
933     esdevent->GetPrimaryVertex()->GetCovMatrix(covVtx);
934     chi = esdevent->GetPrimaryVertex()->GetChi2toNDF();
935   }
936   else if (aodevent)
937   {
938     aodevent->GetPrimaryVertex()->GetCovMatrix(covVtx);
939     chi = aodevent->GetPrimaryVertex()->GetChi2perNDF();//Different from ESD?
940   }
941   
942   AliAODVertex * primary = new(vertices[jVertices++])
943   AliAODVertex(pos, covVtx, chi, NULL, -1, AliAODVertex::kPrimary);
944   primary->SetName(fEvent->GetPrimaryVertex()->GetName());
945   primary->SetTitle(fEvent->GetPrimaryVertex()->GetTitle());
946   
947 }
948
949 //___________________________________________________________
950 void AliAnalysisTaskEMCALClusterize::FillCaloClusterInEvent()
951 {
952   // Get the CaloClusters array, do some final calculations 
953   // and put the clusters in the output or input event 
954   // as a separate branch 
955   
956   //First recalculate track-matching for the new clusters
957   if(fDoTrackMatching) 
958   {
959     fRecoUtils->FindMatches(fEvent,fCaloClusterArr,fGeom);
960   }
961   //Put the new clusters in the AOD list
962   
963   Int_t kNumberOfCaloClusters   = fCaloClusterArr->GetEntriesFast();
964
965   for(Int_t i = 0; i < kNumberOfCaloClusters; i++)
966   {
967     AliAODCaloCluster *newCluster = (AliAODCaloCluster *) fCaloClusterArr->At(i);
968     
969     newCluster->SetID(i);
970     
971     // Correct cluster energy non linearity    
972     newCluster->SetE(fRecoUtils->CorrectClusterEnergyLinearity(newCluster));
973     
974     //Add matched track
975     if(fDoTrackMatching)
976     {
977       Int_t trackIndex = fRecoUtils->GetMatchedTrackIndex(i);
978       if(trackIndex >= 0)
979       {
980         newCluster->AddTrackMatched(fEvent->GetTrack(trackIndex));
981         if(DebugLevel() > 1) 
982           printf("AliAnalysisTaksEMCALClusterize::UserExec() - Matched Track index %d to new cluster %d \n",trackIndex,i);
983       }
984       
985       Float_t dR = 999., dZ = 999.;
986       fRecoUtils->GetMatchedResiduals(newCluster->GetID(),dR,dZ);
987       newCluster->SetTrackDistance(dR,dZ);
988       
989     }
990     else 
991     {// Assign previously assigned matched track in reco, very very rough
992       Int_t absId0 = newCluster->GetCellsAbsId()[0]; // Assign match of first cell in cluster
993       newCluster->SetTrackDistance(fCellMatchdPhi[absId0],fCellMatchdEta[absId0]);
994     }
995     
996     //printf("New cluster E %f, Time  %e, Id = ", newCluster->E(), newCluster->GetTOF() );
997     //for(Int_t icell=0; icell < newCluster->GetNCells(); icell++ ) printf(" %d,", newCluster->GetCellsAbsId() [icell] );
998     //printf("\n");
999     
1000     // Calculate distance to bad channel for new cluster. Make sure you give the list of bad channels.
1001     fRecoUtils->RecalculateClusterDistanceToBadChannel(fGeom, fEvent->GetEMCALCells(), newCluster);
1002     
1003     new((*fOutputAODBranch)[i])  AliAODCaloCluster(*newCluster);
1004     
1005     if(DebugLevel() > 1 )
1006       printf("AliAnalysisTaksEMCALClusterize::UserExec() - New cluster %d of %d, energy %f, mc label %d \n",newCluster->GetID(), kNumberOfCaloClusters, newCluster->E(), newCluster->GetLabel());
1007     
1008   } // cluster loop
1009   
1010   fOutputAODBranch->Expand(kNumberOfCaloClusters); // resize TObjArray to 'remove' slots
1011   
1012   
1013 }
1014
1015 //_______________________________________________
1016 TString AliAnalysisTaskEMCALClusterize::GetPass()
1017 {
1018   // Get passx from filename.
1019   
1020   if (!AliAnalysisManager::GetAnalysisManager()->GetTree()) 
1021   {
1022     AliError("AliAnalysisTaskEMCALClusterize::GetPass() - Pointer to tree = 0, returning null\n");
1023     return TString("");
1024   }
1025   
1026   if (!AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile()) 
1027   {
1028     AliError("AliAnalysisTaskEMCALClusterize::GetPass() - Null pointer input file, returning null\n");
1029     return TString("");
1030   }
1031   
1032   TString pass(AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile()->GetName());
1033   if      (pass.Contains("ass1")) return TString("pass1");
1034   else if (pass.Contains("ass2")) return TString("pass2");
1035   else if (pass.Contains("ass3")) return TString("pass3");
1036   else if (pass.Contains("ass4")) return TString("pass4");
1037   else if (pass.Contains("ass5")) return TString("pass5");
1038   else if (pass.Contains("calo") || pass.Contains("high_lumi"))
1039   {
1040     printf("AliAnalysisTaskEMCALClusterize::GetPass() - Path contains <calo> or <high-lumi>, set as <pass1>\n");
1041     return TString("pass1");
1042   }
1043   
1044   // No condition fullfilled, give a default value
1045   printf("AliAnalysisTaskEMCALClusterize::GetPass() - Pass number string not found \n");
1046   return TString("");            
1047   
1048 }
1049
1050 //_________________________________________
1051 void AliAnalysisTaskEMCALClusterize::Init()
1052 {
1053   //Init analysis with configuration macro if available
1054   
1055   fOADBSet           = kFALSE;
1056   if(fOADBFilePath == "") fOADBFilePath = "$ALICE_ROOT/OADB/EMCAL" ;          
1057   
1058   fBranchNames       = "ESD:AliESDHeader.,EMCALCells.";
1059   
1060   if(!fRecParam)     fRecParam  = new AliEMCALRecParam;
1061   if(!fRecoUtils)    fRecoUtils = new AliEMCALRecoUtils();  
1062   
1063   if(fMaxEvent          <= 0) fMaxEvent          = 1000000000;
1064   if(fSelectCellMinE    <= 0) fSelectCellMinE    = 0.005;     
1065   if(fSelectCellMinFrac <= 0) fSelectCellMinFrac = 0.001;
1066   fRejectBelowThreshold = kFALSE;
1067
1068   //Centrality
1069   if(fCentralityClass  == "") fCentralityClass  = "V0M";
1070   
1071   if (fOCDBpath            == "") fOCDBpath            = "raw://" ;
1072   if (fOutputAODBranchName == "") fOutputAODBranchName = "newEMCALClusters" ;
1073   
1074   if(gROOT->LoadMacro(fConfigName) >=0)
1075   {
1076     printf("AliAnalysisTaksEMCALClusterize::Init() - Configure analysis with %s\n",fConfigName.Data());
1077     AliAnalysisTaskEMCALClusterize *clus = (AliAnalysisTaskEMCALClusterize*)gInterpreter->ProcessLine("ConfigEMCALClusterize()");
1078     fGeomName         = clus->fGeomName; 
1079     fLoadGeomMatrices = clus->fLoadGeomMatrices;
1080     fOCDBpath         = clus->fOCDBpath;   
1081     fAccessOCDB       = clus->fAccessOCDB;
1082     fRecParam         = clus->fRecParam;
1083     fJustUnfold       = clus->fJustUnfold;
1084     fFillAODFile      = clus->fFillAODFile;
1085     fRecoUtils        = clus->fRecoUtils; 
1086     fConfigName       = clus->fConfigName;
1087     fMaxEvent         = clus->fMaxEvent;
1088     fDoTrackMatching  = clus->fDoTrackMatching;
1089     fOutputAODBranchName = clus->fOutputAODBranchName;
1090     for(Int_t i = 0; i < 12; i++) fGeomMatrix[i] = clus->fGeomMatrix[i] ;
1091     fCentralityClass  = clus->fCentralityClass;
1092     fCentralityBin[0] = clus->fCentralityBin[0];
1093     fCentralityBin[1] = clus->fCentralityBin[1];
1094   }
1095   
1096   // Init geometry, I do not like much to do it like this ...
1097   if(fImportGeometryFromFile && !gGeoManager) 
1098   {
1099     if (fImportGeometryFilePath == "") fImportGeometryFilePath = "$ALICE_ROOT/OADB/EMCAL/geometry_2011.root" ; // "$ALICE_ROOT/EVE/alice-data/default_geo.root"
1100     printf("AliAnalysisTaskEMCALClusterize::Init() - Import %s\n",fImportGeometryFilePath.Data());
1101     TGeoManager::Import(fImportGeometryFilePath) ; 
1102   }
1103   
1104 }  
1105
1106 //_______________________________________________________
1107 void AliAnalysisTaskEMCALClusterize::InitClusterization()
1108 {
1109   //Select clusterization/unfolding algorithm and set all the needed parameters
1110   
1111   if (fJustUnfold)
1112   {
1113     // init the unfolding afterburner 
1114     delete fUnfolder;
1115     fUnfolder =  new AliEMCALAfterBurnerUF(fRecParam->GetW0(),fRecParam->GetLocMaxCut(),fRecParam->GetMinECut());
1116     return;
1117   }
1118   
1119   //First init the clusterizer
1120   delete fClusterizer;
1121   if     (fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerv1)
1122     fClusterizer = new AliEMCALClusterizerv1 (fGeom, fCalibData, fPedestalData);
1123   else if(fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerv2) 
1124     fClusterizer = new AliEMCALClusterizerv2(fGeom, fCalibData, fPedestalData);
1125   else if(fRecParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerNxN)
1126   { 
1127     fClusterizer = new AliEMCALClusterizerNxN(fGeom, fCalibData, fPedestalData);
1128     fClusterizer->SetNRowDiff(fRecParam->GetNRowDiff());
1129     fClusterizer->SetNColDiff(fRecParam->GetNColDiff());
1130   } 
1131   else 
1132   {
1133     AliFatal(Form("Clusterizer < %d > not available", fRecParam->GetClusterizerFlag()));
1134   }
1135   
1136   //Now set the parameters
1137   fClusterizer->SetECAClusteringThreshold( fRecParam->GetClusteringThreshold() );
1138   fClusterizer->SetECALogWeight          ( fRecParam->GetW0()                  );
1139   fClusterizer->SetMinECut               ( fRecParam->GetMinECut()             );    
1140   fClusterizer->SetUnfolding             ( fRecParam->GetUnfold()              );
1141   fClusterizer->SetECALocalMaxCut        ( fRecParam->GetLocMaxCut()           );
1142   fClusterizer->SetTimeCut               ( fRecParam->GetTimeCut()             );
1143   fClusterizer->SetTimeMin               ( fRecParam->GetTimeMin()             );
1144   fClusterizer->SetTimeMax               ( fRecParam->GetTimeMax()             );
1145   fClusterizer->SetInputCalibrated       ( kTRUE                               );
1146   fClusterizer->SetJustClusters          ( kTRUE                               );  
1147   
1148   // Initialize the cluster rec points and digits arrays and get them.
1149   fClusterizer->SetOutput(0);
1150   fClusterArr = const_cast<TObjArray *>(fClusterizer->GetRecPoints());
1151   fDigitsArr  = fClusterizer->GetDigits();
1152   
1153   //In case of unfolding after clusterization is requested, set the corresponding parameters
1154   if(fRecParam->GetUnfold())
1155   {
1156     Int_t i=0;
1157     for (i = 0; i < 8; i++) 
1158     {
1159       fClusterizer->SetSSPars(i, fRecParam->GetSSPars(i));
1160     }//end of loop over parameters
1161     
1162     for (i = 0; i < 3; i++) 
1163     {
1164       fClusterizer->SetPar5  (i, fRecParam->GetPar5(i));
1165       fClusterizer->SetPar6  (i, fRecParam->GetPar6(i));
1166     }//end of loop over parameters
1167     fClusterizer->SetRejectBelowThreshold(fRejectBelowThreshold);//here we set option of unfolding: split or reject energy
1168     fClusterizer->InitClusterUnfolding();
1169     
1170   }// to unfold
1171 }
1172
1173 //_________________________________________________
1174 void AliAnalysisTaskEMCALClusterize::InitGeometry()
1175 {
1176   // Init geometry and set the geometry matrix, for the first event, skip the rest
1177   // Also set once the run dependent calibrations
1178   
1179   if(fGeomMatrixSet) return;
1180   
1181   Int_t runnumber = InputEvent()->GetRunNumber() ;
1182   if (!fGeom)
1183   {
1184     if(fGeomName=="")
1185     {
1186       if     (runnumber < 140000) fGeomName = "EMCAL_FIRSTYEARV1";
1187       else if(runnumber < 171000) fGeomName = "EMCAL_COMPLETEV1";
1188       else                        fGeomName = "EMCAL_COMPLETE12SMV1";  
1189       printf("AliAnalysisTaskEMCALClusterize::InitGeometry() - Set EMCAL geometry name to <%s> for run %d\n",fGeomName.Data(),runnumber);
1190     }
1191     
1192                 fGeom = AliEMCALGeometry::GetInstance(fGeomName);
1193     
1194                 if(fDebug > 0)
1195     {
1196                         printf("AliAnalysisTaskEMCALClusterize::InitGeometry(run=%d)",runnumber);
1197                         if (!gGeoManager) printf(" - Careful!, gGeoManager not loaded, load misalign matrices");
1198                         printf("\n");
1199                 }
1200         } // geometry pointer did not exist before
1201   
1202   if(fLoadGeomMatrices)
1203   {
1204     printf("AliAnalysisTaskEMCALClusterize::InitGeometry() - Load user defined EMCAL geometry matrices\n");
1205     
1206     // OADB if available
1207     AliOADBContainer emcGeoMat("AliEMCALgeo");
1208     emcGeoMat.InitFromFile(Form("%s/EMCALlocal2master.root",fOADBFilePath.Data()),"AliEMCALgeo");
1209     TObjArray *matEMCAL=(TObjArray*)emcGeoMat.GetObject(runnumber,"EmcalMatrices");
1210     
1211     
1212     for(Int_t mod=0; mod < (fGeom->GetEMCGeometry())->GetNumberOfSuperModules(); mod++)
1213     {
1214       
1215       if (!fGeomMatrix[mod]) // Get it from OADB
1216       {
1217         if(fDebug > 1 ) 
1218           printf("AliAnalysisTaskEMCALClusterize::InitGeometry() - EMCAL matrices SM %d, %p\n",
1219                  mod,((TGeoHMatrix*) matEMCAL->At(mod)));
1220         //((TGeoHMatrix*) matEMCAL->At(mod))->Print();
1221         
1222         fGeomMatrix[mod] = (TGeoHMatrix*) matEMCAL->At(mod) ;
1223       }        
1224       
1225       if(fGeomMatrix[mod])
1226       {
1227         if(DebugLevel() > 1) 
1228           fGeomMatrix[mod]->Print();
1229         
1230         fGeom->SetMisalMatrix(fGeomMatrix[mod],mod) ;  
1231       }
1232       
1233       fGeomMatrixSet=kTRUE;
1234       
1235     }//SM loop
1236   }//Load matrices
1237   else if(!gGeoManager)
1238   {
1239     printf("AliAnalysisTaksEMCALClusterize::InitGeometry() - Get geo matrices from data");
1240     //Still not implemented in AOD, just a workaround to be able to work at least with ESDs     
1241     if(!strcmp(fEvent->GetName(),"AliAODEvent")) 
1242     {
1243       if(DebugLevel() > 1) 
1244         Warning("UserExec","Use ideal geometry, values geometry matrix not kept in AODs.");
1245     }//AOD
1246     else 
1247     {   
1248       if(DebugLevel() > 1) 
1249         printf("AliAnalysisTaksEMCALClusterize::InitGeometry() - AliAnalysisTaskEMCALClusterize Load Misaligned matrices.");
1250       
1251       AliESDEvent* esd = dynamic_cast<AliESDEvent*>(fEvent) ;
1252       
1253       if(!esd)
1254       {
1255         Error("InitGeometry"," - This event does not contain ESDs?");
1256         if(fFillAODFile) AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kFALSE);
1257         return;
1258       }
1259       
1260       for(Int_t mod=0; mod < (fGeom->GetEMCGeometry())->GetNumberOfSuperModules(); mod++)
1261       {
1262         if(DebugLevel() > 1) 
1263           esd->GetEMCALMatrix(mod)->Print();
1264         
1265         if(esd->GetEMCALMatrix(mod)) fGeom->SetMisalMatrix(esd->GetEMCALMatrix(mod),mod) ;
1266         
1267       } 
1268       
1269       fGeomMatrixSet=kTRUE;
1270       
1271     }//ESD
1272   }//Load matrices from Data 
1273   
1274 }
1275
1276 //____________________________________________________
1277 Bool_t AliAnalysisTaskEMCALClusterize::IsExoticEvent()
1278 {
1279   
1280   // Check if event is exotic, get an exotic cell and compare with triggered patch
1281   // If there is a match remove event ... to be completed, filled with something provisional
1282   
1283   if(!fRemoveExoticEvents) return kFALSE;
1284   
1285   // Loop on cells
1286   AliVCaloCells * cells = fEvent->GetEMCALCells();
1287   Float_t totCellE = 0;
1288   Int_t bc = InputEvent()->GetBunchCrossNumber();
1289   for(Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
1290   {
1291     
1292     Float_t  ecell = 0 ;
1293     Double_t tcell = 0 ;
1294     
1295     Int_t absID   = cells->GetCellNumber(icell);
1296     Bool_t accept = fRecoUtils->AcceptCalibrateCell(absID,bc,ecell,tcell,cells);
1297     if(accept && !fRecoUtils->IsExoticCell(absID,cells,bc)) totCellE += ecell;
1298   }
1299   
1300   //  TString triggerclasses = event->GetFiredTriggerClasses();
1301   //    printf("AliAnalysisTaskEMCALClusterize - reject event %d with cluster  - reject event with ncells in SM3 %d and SM4 %d\n",(Int_t)Entry(),ncellsSM3, ncellsSM4);
1302   //    if(fFillAODFile) AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kFALSE);;
1303   //    return;
1304   //  
1305   
1306   //printf("TotE cell %f\n",totCellE);
1307   if(totCellE < 1) return kTRUE;
1308   
1309   return kFALSE;
1310   
1311
1312
1313 //________________________________________________________________
1314 Bool_t AliAnalysisTaskEMCALClusterize::IsLEDEvent(const Int_t run)
1315 {
1316   //Check if event is LED
1317   
1318   if(!fRemoveLEDEvents) return kFALSE;
1319   
1320   //check events of LHC11a period
1321   if(run < 146858 || run > 146860) return kFALSE ;
1322   
1323   // Count number of cells with energy larger than 0.1 in SM3, cut on this number
1324   Int_t ncellsSM3 = 0;
1325   AliVCaloCells * cells = fEvent->GetEMCALCells();
1326   for(Int_t icell = 0; icell < cells->GetNumberOfCells(); icell++)
1327   {
1328     if(cells->GetAmplitude(icell) > 0.1 && cells->GetCellNumber(icell)/(24*48)==3) ncellsSM3++;
1329   }
1330   
1331   TString triggerclasses = fEvent->GetFiredTriggerClasses();
1332   
1333   Int_t ncellcut = 21;
1334   if(triggerclasses.Contains("EMC")) ncellcut = 35;
1335   
1336   if( ncellsSM3 >= ncellcut)
1337   {
1338     printf("AliAnalysisTaksEMCALClusterize::IsLEDEvent() - reject event %d with ncells in SM3 %d\n",(Int_t)Entry(),ncellsSM3);
1339     if(fFillAODFile) AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillAOD(kFALSE);;
1340     return kTRUE;
1341   }
1342   
1343   return kFALSE;
1344   
1345
1346
1347 //_______________________________________________________
1348 void AliAnalysisTaskEMCALClusterize::RecPoints2Clusters()
1349 {
1350   // Restore clusters from recPoints
1351   // Cluster energy, global position, cells and their amplitude 
1352   // fractions are restored
1353   
1354   Int_t j = 0;
1355   for(Int_t i = 0; i < fClusterArr->GetEntriesFast(); i++)
1356   {
1357     AliEMCALRecPoint *recPoint = (AliEMCALRecPoint*) fClusterArr->At(i);
1358     
1359     const Int_t ncells = recPoint->GetMultiplicity();
1360     Int_t ncellsTrue = 0;
1361     
1362     if(recPoint->GetEnergy() < fRecParam->GetClusteringThreshold()) continue;
1363     
1364     // cells and their amplitude fractions
1365     UShort_t   absIds[ncells];  
1366     Double32_t ratios[ncells];
1367     
1368     //For later check embedding
1369     AliAODInputHandler* aodIH = dynamic_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
1370     
1371     Float_t clusterE = 0; 
1372     for (Int_t c = 0; c < ncells; c++) 
1373     {
1374       AliEMCALDigit *digit = (AliEMCALDigit*) fDigitsArr->At(recPoint->GetDigitsList()[c]);
1375       
1376       absIds[ncellsTrue] = digit->GetId();
1377       ratios[ncellsTrue] = recPoint->GetEnergiesList()[c]/digit->GetAmplitude();
1378             
1379       // In case of unfolding, remove digits with unfolded energy too low      
1380       if(fSelectCell)
1381       {
1382         if     (recPoint->GetEnergiesList()[c] < fSelectCellMinE || ratios[ncellsTrue] < fSelectCellMinFrac)  
1383         {
1384           if(DebugLevel() > 1)
1385           {
1386             printf("AliAnalysisTaksEMCALClusterize::RecPoints2Clusters() - Too small energy in cell of cluster: cluster cell %f, digit %f\n",
1387                    recPoint->GetEnergiesList()[c],digit->GetAmplitude());
1388           }
1389           
1390           continue;
1391           
1392         } // if cuts
1393       }// Select cells
1394       
1395       //Recalculate cluster energy and number of cluster cells in case any of the cells was rejected
1396       clusterE  +=recPoint->GetEnergiesList()[c];
1397       
1398       // In case of embedding, fill ratio with amount of signal, 
1399       if (aodIH && aodIH->GetMergeEvents()) 
1400       {
1401         //AliVCaloCells* inEMCALCells = InputEvent()->GetEMCALCells();
1402         AliVCaloCells* meEMCALCells = aodIH->GetEventToMerge()->GetEMCALCells();
1403         AliVCaloCells* ouEMCALCells = AODEvent()->GetEMCALCells();
1404         
1405         Float_t sigAmplitude = meEMCALCells->GetCellAmplitude(absIds[ncellsTrue]);
1406         //Float_t bkgAmplitude = inEMCALCells->GetCellAmplitude(absIds[ncellsTrue]);
1407         Float_t sumAmplitude = ouEMCALCells->GetCellAmplitude(absIds[ncellsTrue]);
1408         //printf("\t AbsID %d, amplitude : bkg %f, sigAmplitude %f, summed %f - %f\n",absIds[ncellsTrue], bkgAmplitude, sigAmplitude, sumAmplitude, digit->GetAmplitude());
1409         
1410         if(sumAmplitude > 0) ratios[ncellsTrue] = sigAmplitude/sumAmplitude;
1411         //printf("\t \t ratio %f\n",ratios[ncellsTrue]);
1412         
1413       }//Embedding
1414       
1415       ncellsTrue++;
1416       
1417     }// cluster cell loop
1418     
1419     if (ncellsTrue < 1) 
1420     {
1421       if (DebugLevel() > 1) 
1422         printf("AliAnalysisTaskEMCALClusterize::RecPoints2Clusters() - Skipping cluster with no cells avobe threshold E = %f, ncells %d\n",
1423                recPoint->GetEnergy(), ncells);
1424       continue;
1425     }
1426     
1427     //if(ncellsTrue != ncells) printf("Old E %f, ncells %d; New E %f, ncells %d\n",recPoint->GetEnergy(),ncells,clusterE,ncellsTrue);
1428     
1429     if(clusterE <  fRecParam->GetClusteringThreshold()) 
1430     {
1431       if (DebugLevel()>1)
1432         printf("AliAnalysisTaskEMCALClusterize::RecPoints2Clusters() - Remove cluster with energy below seed threshold %f\n",clusterE);
1433       continue;
1434     }
1435     
1436     TVector3 gpos;
1437     Float_t g[3];
1438     
1439     // calculate new cluster position
1440     
1441     recPoint->EvalGlobalPosition(fRecParam->GetW0(), fDigitsArr);
1442     recPoint->GetGlobalPosition(gpos);
1443     gpos.GetXYZ(g);
1444     
1445     // create a new cluster
1446     
1447     (*fCaloClusterArr)[j] = new AliAODCaloCluster() ;
1448     AliAODCaloCluster *clus = dynamic_cast<AliAODCaloCluster *>( fCaloClusterArr->At(j) ) ;
1449     j++;
1450     clus->SetType(AliVCluster::kEMCALClusterv1);
1451     clus->SetE(clusterE);
1452     clus->SetPosition(g);
1453     clus->SetNCells(ncellsTrue);
1454     clus->SetCellsAbsId(absIds);
1455     clus->SetCellsAmplitudeFraction(ratios);
1456     clus->SetChi2(-1); //not yet implemented
1457     clus->SetTOF(recPoint->GetTime()) ; //time-of-flight
1458     clus->SetNExMax(recPoint->GetNExMax()); //number of local maxima
1459     clus->SetDistanceToBadChannel(recPoint->GetDistanceToBadTower()); 
1460     
1461     if(ncells == ncellsTrue)
1462     {
1463       Float_t elipAxis[2];
1464       recPoint->GetElipsAxis(elipAxis);
1465       clus->SetM02(elipAxis[0]*elipAxis[0]) ;
1466       clus->SetM20(elipAxis[1]*elipAxis[1]) ;
1467       clus->SetDispersion(recPoint->GetDispersion());
1468     }
1469     else if(fSelectCell)
1470     {
1471       // In case some cells rejected, in unfolding case, recalculate
1472       // shower shape parameters and position
1473       if(DebugLevel() > 1) 
1474         printf("AliAnalysisTaskEMCALClusterize::RecPoints2Clusters() - Cells removed from cluster (ncells %d, ncellsTrue %d), recalculate Shower Shape\n",ncells,ncellsTrue);
1475       
1476       AliVCaloCells* cells = 0x0; 
1477       if (aodIH && aodIH->GetMergeEvents()) cells = AODEvent()  ->GetEMCALCells();
1478       else                                  cells = InputEvent()->GetEMCALCells();
1479       
1480       fRecoUtils->RecalculateClusterShowerShapeParameters(fGeom,cells,clus);
1481       fRecoUtils->RecalculateClusterPID(clus);
1482       fRecoUtils->RecalculateClusterPosition(fGeom,cells,clus); 
1483       
1484     }
1485     
1486     // MC
1487
1488     if     ( fSetCellMCLabelFromCluster == 1 ) SetClustersMCLabelFrom2SelectedLabels(recPoint,clus) ;
1489     else if( fSetCellMCLabelFromCluster == 2 ) SetClustersMCLabelFromOriginalClusters(clus) ;
1490     else
1491     {
1492       // Normal case, trust what the clusterizer has found
1493       Int_t  parentMult = 0;
1494       Int_t *parentList = recPoint->GetParents(parentMult);
1495       clus->SetLabel(parentList, parentMult);
1496 //      printf("Label list : ");
1497 //      for(Int_t ilabel = 0; ilabel < parentMult; ilabel++ ) printf(" %d ",parentList[ilabel]);
1498 //      printf("\n");
1499     }
1500     
1501   } // recPoints loop
1502   
1503 }
1504
1505 //___________________________________________________________________________
1506 void AliAnalysisTaskEMCALClusterize::RemapMCLabelForAODs(Int_t & label)
1507 {
1508   // MC label for Cells not remapped after ESD filtering, do it here.
1509
1510   if(label < 0) return ;
1511   
1512   AliAODEvent  * evt = dynamic_cast<AliAODEvent*> (fEvent) ;
1513   if(!evt) return ;
1514   
1515   TClonesArray * arr = dynamic_cast<TClonesArray*>(evt->FindListObject("mcparticles")) ;
1516   if(!arr) return ;
1517     
1518   if(label < arr->GetEntriesFast())
1519   {
1520     AliAODMCParticle * particle = dynamic_cast<AliAODMCParticle *>(arr->At(label));
1521     if(!particle) return ;
1522         
1523     if(label == particle->Label()) return ; // label already OK
1524     //else printf("AliAnalysisTaskEMCALClusterize::RemapMCLabelForAODs() - Label  %d - AOD stack %d \n",label, particle->Label());
1525   }
1526   //else printf("AliAnalysisTaskEMCALClusterize::RemapMCLabelForAODs() - Label  %d > AOD labels %d \n",label, arr->GetEntriesFast());
1527   
1528   // loop on the particles list and check if there is one with the same label
1529   for(Int_t ind = 0; ind < arr->GetEntriesFast(); ind++ )
1530   {
1531     AliAODMCParticle * particle = dynamic_cast<AliAODMCParticle *>(arr->At(ind));
1532     if(!particle) continue ;
1533
1534     if(label == particle->Label())
1535     {
1536       label = ind;
1537       //printf("AliAnalysisTaskEMCALClusterize::RemapMCLabelForAODs() - New Label Index  %d \n",label);
1538       return;
1539     }
1540   }
1541   
1542   label = -1;
1543   
1544   //printf("AliAnalysisTaskEMCALClusterize::RemapMCLabelForAODs() - Label not found set to -1 \n");
1545  
1546 }
1547
1548 //________________________________________________
1549 void AliAnalysisTaskEMCALClusterize::ResetArrays()
1550 {
1551   // Reset arrays containing information for all possible cells
1552   for(Int_t j = 0; j < 12672; j++)
1553   {
1554     fOrgClusterCellId[j] = -1;
1555     fCellLabels[j]       = -1;
1556     fCellSecondLabels[j] = -1;
1557     fCellTime[j]         =  0.;
1558     fCellMatchdEta[j]    = -999;
1559     fCellMatchdPhi[j]    = -999;
1560   }
1561 }
1562
1563 //_____________________________________________________________________________________________________
1564 void AliAnalysisTaskEMCALClusterize::SetClustersMCLabelFrom2SelectedLabels(AliEMCALRecPoint  * recPoint,
1565                                                                            AliAODCaloCluster * clus)
1566 {
1567   // Set the cluster MC label, the digizer was filled with most likely MC label for all cells in original cluster
1568   // Now check the second most likely MC label and add it to the new cluster
1569   
1570   Int_t  parentMult = 0;
1571   Int_t *parentList = recPoint->GetParents(parentMult);
1572   clus->SetLabel(parentList, parentMult);
1573   
1574   //Write the second major contributor to each MC cluster.
1575   Int_t iNewLabel ;
1576   for ( Int_t iLoopCell = 0 ; iLoopCell < clus->GetNCells() ; iLoopCell++ )
1577   {
1578     
1579     Int_t idCell = clus->GetCellAbsId(iLoopCell) ;
1580     if(idCell>=0)
1581     {
1582       iNewLabel = 1 ; //iNewLabel makes sure we  don't write twice the same label.
1583       for ( UInt_t iLoopLabels = 0 ; iLoopLabels < clus->GetNLabels() ; iLoopLabels++ )
1584       {
1585         if ( fCellSecondLabels[idCell] == -1 )  iNewLabel = 0;  // -1 is never a good second label.
1586           if ( fCellSecondLabels[idCell] == clus->GetLabelAt(iLoopLabels) )  iNewLabel = 0;
1587             }
1588       if (iNewLabel == 1)
1589       {
1590         Int_t * newLabelArray = new Int_t[clus->GetNLabels()+1] ;
1591         for ( UInt_t iLoopNewLabels = 0 ; iLoopNewLabels < clus->GetNLabels() ; iLoopNewLabels++ )
1592         {
1593           newLabelArray[iLoopNewLabels] = clus->GetLabelAt(iLoopNewLabels) ;
1594         }
1595         
1596         newLabelArray[clus->GetNLabels()] = fCellSecondLabels[idCell] ;
1597         clus->SetLabel(newLabelArray,clus->GetNLabels()+1) ;
1598         delete [] newLabelArray;
1599       }
1600     }//positive cell number
1601   }
1602 }
1603
1604 //___________________________________________________________________________________________________
1605 void AliAnalysisTaskEMCALClusterize::SetClustersMCLabelFromOriginalClusters(AliAODCaloCluster * clus)
1606 {
1607   // Get the original clusters that contribute to the new cluster, assign the labels of such clusters
1608   // to the new cluster.
1609   // Only approximatedly valid  when output are V1 clusters, handle with care
1610     
1611   TArrayI clArray(300) ; //Weird if more than a few clusters are in the origin ...
1612   clArray.Reset();
1613   Int_t nClu = 0;
1614   Int_t nLabTotOrg = 0;
1615   Float_t emax = -1;
1616   Int_t idMax = -1;
1617   
1618   AliVEvent * event = fEvent;
1619   
1620   //In case of embedding the MC signal is in the added event
1621   AliAODInputHandler* aodIH = dynamic_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
1622   if(aodIH && aodIH->GetEventToMerge())  //Embedding
1623     event = aodIH->GetEventToMerge(); //Get clusters directly from embedded signal
1624
1625   
1626   //Find the clusters that originally had the cells
1627   for ( Int_t iLoopCell = 0 ; iLoopCell < clus->GetNCells() ; iLoopCell++ )
1628   {
1629     Int_t idCell = clus->GetCellAbsId(iLoopCell) ;
1630     
1631     if(idCell>=0)
1632     {
1633       Int_t idCluster = fOrgClusterCellId[idCell];
1634       
1635       Bool_t set = kTRUE;
1636       for(Int_t icl =0; icl < nClu; icl++)
1637       {
1638         if(((Int_t)clArray.GetAt(icl))==-1 ) continue;
1639         if( idCluster == ((Int_t)clArray.GetAt(icl)) ) set = kFALSE;
1640         //   printf("\t \t icell %d  Cluster in array %d, IdCluster %d, in array %d, set %d\n",
1641         //          iLoopCell,                 icl,  idCluster,((Int_t)clArray.GetAt(icl)),set);
1642       }
1643       if( set && idCluster >= 0)
1644       {
1645         clArray.SetAt(idCluster,nClu++);
1646         //printf("******** idCluster %d \n",idCluster);
1647         nLabTotOrg+=(event->GetCaloCluster(idCluster))->GetNLabels();
1648
1649         //printf("Cluster in array %d, IdCluster %d\n",nClu-1,  idCluster);
1650
1651         //Search highest E cluster
1652         AliVCluster * clOrg = event->GetCaloCluster(idCluster);
1653         //printf("\t E %f\n",clOrg->E());
1654         if(emax < clOrg->E())
1655         {
1656           emax  = clOrg->E();
1657           idMax = idCluster;
1658         }        
1659       }
1660     }
1661   }// cell loop
1662
1663   if(nClu==0 || nLabTotOrg == 0)
1664   {
1665     //if(clus->E() > 0.25) printf("AliAnalysisTaskEMCALClusterize::SetClustersMCLabelFromOriginalClusters() - Check: N org clusters %d, n tot labels %d, cluster E %f,  n cells %d\n",nClu,nLabTotOrg,clus->E(), clus->GetNCells());
1666     //for(Int_t icell = 0; icell < clus->GetNCells(); icell++) printf("\t cell %d",clus->GetCellsAbsId()[icell]);
1667     //printf("\n");
1668   }
1669   
1670   // Put the first in the list the cluster with highest energy
1671   if(idMax != ((Int_t)clArray.GetAt(0))) // Max not at first position
1672   {
1673     Int_t maxIndex = -1;
1674     Int_t firstCluster = ((Int_t)clArray.GetAt(0));
1675     for ( Int_t iLoopCluster = 0 ; iLoopCluster < nClu ; iLoopCluster++ )
1676     {
1677       if(idMax == ((Int_t)clArray.GetAt(iLoopCluster))) maxIndex = iLoopCluster;
1678     }
1679     
1680     if(firstCluster >=0 && idMax >=0)
1681     {
1682       clArray.SetAt(idMax,0);
1683       clArray.SetAt(firstCluster,maxIndex);
1684     }
1685   }
1686   
1687   // Get the labels list in the original clusters, assign all to the new cluster
1688   TArrayI clMCArray(nLabTotOrg) ;
1689   clMCArray.Reset();
1690   
1691   Int_t nLabTot = 0;
1692   for ( Int_t iLoopCluster = 0 ; iLoopCluster < nClu ; iLoopCluster++ )
1693   {
1694     Int_t idCluster = (Int_t) clArray.GetAt(iLoopCluster);
1695     //printf("New Cluster in Array %d,  idCluster %d \n",iLoopCluster,idCluster);
1696     AliVCluster * clOrg = event->GetCaloCluster(idCluster);
1697     Int_t nLab = clOrg->GetNLabels();
1698     
1699     for ( Int_t iLab = 0 ; iLab < nLab ; iLab++ )
1700     {
1701       Int_t lab = clOrg->GetLabelAt(iLab) ;
1702       if(lab>=0)
1703       {
1704         Bool_t set = kTRUE;
1705         //printf("\t \t Set Label %d \n", lab);
1706         for(Int_t iLabTot =0; iLabTot < nLabTot; iLabTot++)
1707         {
1708           if( lab == ((Int_t)clMCArray.GetAt(iLabTot)) ) set = kFALSE;
1709           //printf("iLoopCluster %d, Label ID in Org Cluster %d,label %d  Label ID in array %d, label in array %d, set %d\n",
1710           //       iLoopCluster,                           iLab,     lab,           iLabTot,  ((Int_t)clMCArray.GetAt(iLabTot)),set);
1711         }
1712         if( set ) clMCArray.SetAt(lab,nLabTot++);
1713       }
1714     }
1715   }// cluster loop
1716   
1717   // Set the final list of labels
1718   
1719   clus->SetLabel(clMCArray.GetArray(), nLabTot);
1720   
1721 //  printf("Final list of labels for new cluster : \n");
1722 //  for(Int_t ice = 0; ice < clus->GetNCells() ; ice++)
1723 //  {
1724 //    printf("\t Cell %d ",clus->GetCellsAbsId()[ice]);
1725 //    Int_t label = InputEvent()->GetEMCALCells()->GetCellMCLabel(clus->GetCellsAbsId()[ice]);
1726 //    printf(" org %d ",label);
1727 //    RemapMCLabelForAODs(label);
1728 //    printf(" new %d \n",label);
1729 //  }
1730 //  for(Int_t imc = 0; imc < clus->GetNLabels(); imc++) printf("\t Label %d\n",clus->GetLabelAt(imc));
1731 }
1732
1733
1734 //____________________________________________________________
1735 void AliAnalysisTaskEMCALClusterize::UserCreateOutputObjects()
1736 {
1737   // Init geometry, create list of output clusters
1738   
1739
1740   fOutputAODBranch = new TClonesArray("AliAODCaloCluster", 0);
1741
1742   if(fOutputAODBranchName.Length()==0)
1743   {
1744     fOutputAODBranchName = "newEMCALClustersArray";
1745     printf("Cluster branch name not set, set it to newEMCALClustersArray \n");
1746   }
1747   
1748   fOutputAODBranch->SetName(fOutputAODBranchName);
1749   
1750   if( fFillAODFile )
1751   {
1752     //fOutputAODBranch = new TClonesArray("AliAODCaloCluster", 0);
1753
1754     
1755     //fOutputAODBranch->SetOwner(kFALSE);
1756     
1757     AddAODBranch("TClonesArray", &fOutputAODBranch);
1758   }
1759   
1760 }
1761
1762 //_______________________________________________________
1763 void AliAnalysisTaskEMCALClusterize::UserExec(Option_t *) 
1764 {
1765   // Do clusterization event by event, execute different steps
1766   // 1) Do some checks on the kind of events (ESD, AOD) or if some filtering is needed, initializations
1767   //    load things and clear arrays
1768   // 2) Clusterize a) just unfolding existing clusters (fJustUnfold)
1769   //               b) recluster cells
1770   //                   +  convert cells into digits (calibrating them)
1771   //                   +  recluster digits into recPoints with the chosen clusterizer (V1, V1+Unfold,V2, NxN) 
1772   //                      with methods in AliEMCALClusterizer
1773   //                   +  transform recPoints into CaloClusters
1774   // 3) Do some final calculations in the found clusters (track-matching) and put them in an AOD branch
1775   
1776   //-------
1777   // Step 1
1778   
1779
1780   //Remove the contents of AOD branch output list set in the previous event 
1781   fOutputAODBranch->Clear("C");
1782
1783   LoadBranches();
1784   
1785   //Check if there is a centrality value, PbPb analysis, and if a centrality bin selection is requested
1786   //If we need a centrality bin, we select only those events in the corresponding bin.
1787   if( GetCentrality() && fCentralityBin[0] >= 0 && fCentralityBin[1] >= 0 )
1788   {
1789     Float_t cen = GetEventCentrality();
1790     if(cen > fCentralityBin[1] || cen < fCentralityBin[0]) return ; //reject events out of bin.
1791   }  
1792     
1793   // intermediate array with new clusters : init the array only once or clear from previous event
1794   if(!fCaloClusterArr) fCaloClusterArr    = new TObjArray(10000);
1795   else                 fCaloClusterArr->Delete();//Clear("C"); it leaks?
1796
1797   InitGeometry(); // only once, must be done before OADB, geo OADB accessed here
1798   
1799   //Get the event, do some checks and settings
1800   CheckAndGetEvent() ;
1801   
1802   if (!fEvent) 
1803   {
1804     if(DebugLevel() > 0 ) printf("AliAnalysisTaksEMCALClusterize::UserExec() - Skip Event %d", (Int_t) Entry());
1805     return ;
1806   }
1807
1808   //Init pointers, geometry, clusterizer, ocdb, aodb
1809   
1810   
1811   if(fAccessOCDB) AccessOCDB();
1812   if(fAccessOADB) AccessOADB(); // only once
1813   
1814   InitClusterization();
1815   
1816   //-------
1817   // Step 2
1818   
1819   // Make clusters
1820   if (fJustUnfold) ClusterUnfolding();
1821   else             ClusterizeCells() ;
1822   
1823   //-------
1824   // Step 3
1825   
1826   FillCaloClusterInEvent();
1827   
1828 }      
1829