]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALClusterizerv1.cxx
Allow to have clusters in 2 SuperModules with common eta =0.
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALClusterizerv1.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 /* $Id$ */
17
18 //-- Author: Yves Schutz (SUBATECH)  & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
19 //  August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
20 //                           of new  IO (à la PHOS)
21 //  Mar 2007, Aleksei Pavlinov - new algoritmh of pseudo clusters
22 //////////////////////////////////////////////////////////////////////////////
23 //  Clusterization class. Performs clusterization (collects neighbouring active cells) and 
24 //  unfolds the clusters having several local maxima.  
25 //  Results are stored in TreeR#, branches EMCALTowerRP (EMC recPoints),
26 //  EMCALPreShoRP (CPV RecPoints) and AliEMCALClusterizer (Clusterizer with all 
27 //  parameters including input digits branch title, thresholds etc.)
28 //  This TTask is normally called from Reconstructioner, but can as well be used in 
29 //  standalone mode.
30 // Use Case:
31 //  root [0] AliEMCALClusterizerv1 * cl = new AliEMCALClusterizerv1("galice.root")  
32 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
33 //               //reads gAlice from header file "..."                      
34 //  root [1] cl->ExecuteTask()  
35 //               //finds RecPoints in all events stored in galice.root
36 //  root [2] cl->SetDigitsBranch("digits2") 
37 //               //sets another title for Digitis (input) branch
38 //  root [3] cl->SetRecPointsBranch("recp2")  
39 //               //sets another title four output branches
40 //  root [4] cl->SetTowerLocalMaxCut(0.03)  
41 //               //set clusterization parameters
42 //  root [5] cl->ExecuteTask("deb all time")  
43 //               //once more finds RecPoints options are 
44 //               // deb - print number of found rec points
45 //               // deb all - print number of found RecPoints and some their characteristics 
46 //               // time - print benchmarking results
47
48 // --- ROOT system ---
49 #include <cassert>
50
51 class TROOT;
52 #include <TH1.h>
53 #include <TFile.h> 
54 class TFolder;
55 #include <TMath.h> 
56 #include <TMinuit.h>
57 #include <TTree.h> 
58 class TSystem; 
59 #include <TBenchmark.h>
60 #include <TBrowser.h>
61 #include <TROOT.h>
62
63 // --- Standard library ---
64
65
66 // --- AliRoot header files ---
67 #include "AliRunLoader.h"
68 #include "AliRun.h"
69 #include "AliESD.h"
70 #include "AliEMCALClusterizerv1.h"
71 #include "AliEMCALRecPoint.h"
72 #include "AliEMCALDigit.h"
73 #include "AliEMCALDigitizer.h"
74 #include "AliEMCAL.h"
75 #include "AliEMCALGeometry.h"
76 #include "AliEMCALRecParam.h"
77 #include "AliEMCALReconstructor.h"
78 #include "AliCDBManager.h"
79 #include "AliCaloCalibPedestal.h"
80 #include "AliEMCALCalibData.h"
81 class AliCDBStorage;
82 #include "AliCDBEntry.h"
83
84 ClassImp(AliEMCALClusterizerv1)
85
86 //____________________________________________________________________________
87 AliEMCALClusterizerv1::AliEMCALClusterizerv1()
88   : AliEMCALClusterizer(),
89     fGeom(0),
90     fDefaultInit(kFALSE),
91     fToUnfold(kFALSE),
92     fNumberOfECAClusters(0),fCalibData(0),fCaloPed(0),
93     fADCchannelECA(0.),fADCpedestalECA(0.),fECAClusteringThreshold(0.),fECALocMaxCut(0.),
94     fECAW0(0.),fTimeCut(1.),fTimeMin(-1.),fTimeMax(1.),fMinECut(0.)
95 {
96   // ctor with the indication of the file where header Tree and digits Tree are stored
97   
98   Init() ;
99 }
100
101 //____________________________________________________________________________
102 AliEMCALClusterizerv1::AliEMCALClusterizerv1(AliEMCALGeometry* geometry)
103   : AliEMCALClusterizer(),
104     fGeom(geometry),
105     fDefaultInit(kFALSE),
106     fToUnfold(kFALSE),
107     fNumberOfECAClusters(0),fCalibData(0), fCaloPed(0),
108     fADCchannelECA(0.),fADCpedestalECA(0.),fECAClusteringThreshold(0.),fECALocMaxCut(0.),
109     fECAW0(0.),fTimeCut(1.),fTimeMin(-1.),fTimeMax(1.),fMinECut(0.)
110 {
111   // ctor with the indication of the file where header Tree and digits Tree are stored
112   // use this contructor to avoid usage of Init() which uses runloader
113   // change needed by HLT - MP
114
115   // Note for the future: the use on runloader should be avoided or optional at least
116   // another way is to make Init virtual and protected at least such that the deriving classes can overload
117   // Init() ;
118   //
119
120   if (!fGeom)
121     {
122       AliFatal("Geometry not initialized.");
123     }
124
125   if(!gMinuit)
126     gMinuit = new TMinuit(100) ;
127
128 }
129
130 //____________________________________________________________________________
131 AliEMCALClusterizerv1::AliEMCALClusterizerv1(AliEMCALGeometry* geometry, AliEMCALCalibData * calib, AliCaloCalibPedestal * caloped)
132 : AliEMCALClusterizer(),
133 fGeom(geometry),
134 fDefaultInit(kFALSE),
135 fToUnfold(kFALSE),
136 fNumberOfECAClusters(0),fCalibData(calib), fCaloPed(caloped),
137 fADCchannelECA(0.),fADCpedestalECA(0.),fECAClusteringThreshold(0.),fECALocMaxCut(0.),
138 fECAW0(0.),fTimeCut(1.),fTimeMin(-1.),fTimeMax(1.),fMinECut(0.)
139 {
140         // ctor, geometry and calibration are initialized elsewhere.
141         
142         if (!fGeom)
143                 AliFatal("Geometry not initialized.");
144                 
145         if(!gMinuit)
146                 gMinuit = new TMinuit(100) ;
147         
148 }
149
150
151 //____________________________________________________________________________
152   AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
153 {
154   // dtor
155 }
156
157 //____________________________________________________________________________
158 Float_t  AliEMCALClusterizerv1::Calibrate(const Float_t amp, const Float_t time, const Int_t absId) 
159 {
160  
161   // Convert digitized amplitude into energy.
162   // Calibration parameters are taken from calibration data base for raw data,
163   // or from digitizer parameters for simulated data.
164
165   if(fCalibData){
166     
167     if (fGeom==0)
168       AliFatal("Did not get geometry from EMCALLoader") ;
169     
170     Int_t iSupMod = -1;
171     Int_t nModule = -1;
172     Int_t nIphi   = -1;
173     Int_t nIeta   = -1;
174     Int_t iphi    = -1;
175     Int_t ieta    = -1;
176     
177     Bool_t bCell = fGeom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
178     if(!bCell) {
179       fGeom->PrintGeometry();
180       Error("Calibrate()"," Wrong cell id number : %i", absId);
181       assert(0);
182     }
183
184     fGeom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
185           
186         // Check if channel is bad (dead or hot), in this case return 0.        
187         // Gustavo: 15-12-09 In case of RAW data this selection is already done, but not in simulation.
188         // for the moment keep it here but remember to do the selection at the sdigitizer level 
189         // and remove it from here
190         Int_t channelStatus = (Int_t)(fCaloPed->GetDeadMap(iSupMod))->GetBinContent(ieta,iphi);
191         if(channelStatus == AliCaloCalibPedestal::kHot || channelStatus == AliCaloCalibPedestal::kDead) {
192                   AliDebug(2,Form("Tower from SM %d, ieta %d, iphi %d is BAD : status %d !!!",iSupMod,ieta,iphi, channelStatus));
193                   return 0;
194         }
195         //Check if time is too large or too small, indication of a noisy channel, remove in this case
196         if(time > fTimeMax || time < fTimeMin) return 0;
197           
198     fADCchannelECA  = fCalibData->GetADCchannel (iSupMod,ieta,iphi);
199     fADCpedestalECA = fCalibData->GetADCpedestal(iSupMod,ieta,iphi);
200   
201    return -fADCpedestalECA + amp * fADCchannelECA ;        
202  
203   }
204   else //Return energy with default parameters if calibration is not available
205     return -fADCpedestalECA + amp * fADCchannelECA ; 
206   
207 }
208
209 //____________________________________________________________________________
210 void AliEMCALClusterizerv1::Digits2Clusters(Option_t * option)
211 {
212   // Steering method to perform clusterization for the current event 
213   // in AliEMCALLoader
214
215   if(strstr(option,"tim"))
216     gBenchmark->Start("EMCALClusterizer"); 
217   
218   if(strstr(option,"print"))
219     Print("") ; 
220  
221   //Get calibration parameters from file or digitizer default values.
222   GetCalibrationParameters() ;
223
224   //Get dead channel map from file or digitizer default values.
225   GetCaloCalibPedestal() ;
226         
227   fNumberOfECAClusters = 0;
228
229   MakeClusters() ;  //only the real clusters
230
231   if(fToUnfold)
232     MakeUnfolding() ;
233
234   Int_t index ;
235
236   //Evaluate position, dispersion and other RecPoint properties for EC section                      
237   for(index = 0; index < fRecPoints->GetEntries(); index++) {
238       dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index))->EvalAll(fECAW0,fDigitsArr) ;
239           //For each rec.point set the distance to the nearest bad crystal
240           dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index))->EvalDistanceToBadChannels(fCaloPed);
241   }
242
243   fRecPoints->Sort() ;
244
245   for(index = 0; index < fRecPoints->GetEntries(); index++) {
246     (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index)))->SetIndexInList(index) ;
247     (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index)))->Print();
248   }
249
250   fTreeR->Fill();
251   
252   if(strstr(option,"deb") || strstr(option,"all"))  
253     PrintRecPoints(option) ;
254
255   AliDebug(1,Form("EMCAL Clusterizer found %d Rec Points",fRecPoints->GetEntriesFast()));
256
257   fRecPoints->Delete();
258
259   if(strstr(option,"tim")){
260     gBenchmark->Stop("EMCALClusterizer");
261     printf("Exec took %f seconds for Clusterizing", 
262            gBenchmark->GetCpuTime("EMCALClusterizer"));
263   }    
264 }
265
266 //____________________________________________________________________________
267 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALRecPoint * recPoint, AliEMCALDigit ** maxAt, 
268                                       const Float_t* maxAtEnergy,
269                                       Int_t nPar, Float_t * fitparameters) const
270 {
271   // Calls TMinuit to fit the energy distribution of a cluster with several maxima
272   // The initial values for fitting procedure are set equal to the
273   // positions of local maxima.       
274   // Cluster will be fitted as a superposition of nPar/3
275   // electromagnetic showers
276
277   if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");
278
279   gMinuit->mncler();                     // Reset Minuit's list of paramters
280   gMinuit->SetPrintLevel(-1) ;           // No Printout
281   gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;
282   // To set the address of the minimization function
283   TList * toMinuit = new TList();
284   toMinuit->AddAt(recPoint,0) ;
285   toMinuit->AddAt(fDigitsArr,1) ;
286   toMinuit->AddAt(fGeom,2) ;
287
288   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare
289
290   // filling initial values for fit parameters
291   AliEMCALDigit * digit ;
292
293   Int_t ierflg  = 0;
294   Int_t index   = 0 ;
295   Int_t nDigits = (Int_t) nPar / 3 ;
296
297   Int_t iDigit ;
298
299   for(iDigit = 0; iDigit < nDigits; iDigit++){
300     digit = maxAt[iDigit];
301     Double_t x = 0.;
302     Double_t y = 0.;
303     Double_t z = 0.;
304
305     fGeom->RelPosCellInSModule(digit->GetId(), y, x, z);
306
307     Float_t energy = maxAtEnergy[iDigit] ;
308
309     gMinuit->mnparm(index, "x",  x, 0.1, 0, 0, ierflg) ;
310     index++ ;
311     if(ierflg != 0){
312       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : x = %f", x ) ;
313       return kFALSE;
314     }
315     gMinuit->mnparm(index, "z",  z, 0.1, 0, 0, ierflg) ;
316     index++ ;
317     if(ierflg != 0){
318       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : z = %f", z) ;
319       return kFALSE;
320     }
321     gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;
322     index++ ;
323     if(ierflg != 0){
324       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : energy = %f", energy) ;
325       return kFALSE;
326     }
327   }
328
329   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; 
330                       // The number of function call slightly depends on it.
331   //Double_t p1 = 1.0 ;
332   Double_t p2 = 0.0 ;
333
334   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls
335   //  gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient
336   gMinuit->SetMaxIterations(5);
337   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
338   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize
339
340   if(ierflg == 4){  // Minimum not found
341     Error("FindFit", "EMCAL Unfolding  Fit not converged, cluster abandoned " ) ;
342     return kFALSE ;
343   }
344   for(index = 0; index < nPar; index++){
345     Double_t err ;
346     Double_t val ;
347     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index
348     fitparameters[index] = val ;
349   }
350
351   delete toMinuit ;
352   return kTRUE;
353
354 }
355
356 //____________________________________________________________________________
357 void AliEMCALClusterizerv1::GetCalibrationParameters() 
358 {
359   // Set calibration parameters:
360   // if calibration database exists, they are read from database,
361   // otherwise, they are taken from digitizer.
362   //
363   // It is a user responsilibity to open CDB before reconstruction, 
364   // for example: 
365   // AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage("local://CalibDB");
366
367   //Check if calibration is stored in data base
368
369   if(!fCalibData)
370     {
371       AliCDBEntry *entry = (AliCDBEntry*) 
372         AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
373       if (entry) fCalibData =  (AliEMCALCalibData*) entry->GetObject();
374     }
375   
376   if(!fCalibData)
377     AliFatal("Calibration parameters not found in CDB!");
378  
379 }
380
381 //____________________________________________________________________________
382 void AliEMCALClusterizerv1::GetCaloCalibPedestal() 
383 {
384         // Set calibration parameters:
385         // if calibration database exists, they are read from database,
386         // otherwise, they are taken from digitizer.
387         //
388         // It is a user responsilibity to open CDB before reconstruction, 
389         // for example: 
390         // AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage("local://CalibDB");
391         
392         //Check if calibration is stored in data base
393         
394         if(!fCaloPed)
395     {
396                 AliCDBEntry *entry = (AliCDBEntry*) 
397                 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
398                 if (entry) fCaloPed =  (AliCaloCalibPedestal*) entry->GetObject();
399     }
400         
401         if(!fCaloPed)
402                 AliFatal("Pedestal info not found in CDB!");
403         
404 }
405
406
407 //____________________________________________________________________________
408 void AliEMCALClusterizerv1::Init()
409 {
410   // Make all memory allocations which can not be done in default constructor.
411   // Attach the Clusterizer task to the list of EMCAL tasks
412   
413   AliRunLoader *rl = AliRunLoader::Instance();
414   if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL"))
415     fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
416   else 
417     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
418
419   AliDebug(1,Form("geom 0x%x",fGeom));
420
421   if(!gMinuit) 
422     gMinuit = new TMinuit(100) ;
423
424 }
425
426 //____________________________________________________________________________
427 void AliEMCALClusterizerv1::InitParameters()
428
429   // Initializes the parameters for the Clusterizer
430   fNumberOfECAClusters = 0;
431
432   fCalibData               = 0 ;
433   fCaloPed                 = 0 ;
434         
435   const AliEMCALRecParam* recParam = AliEMCALReconstructor::GetRecParam();
436   if(!recParam) {
437     AliFatal("Reconstruction parameters for EMCAL not set!");
438   } else {
439     fECAClusteringThreshold = recParam->GetClusteringThreshold();
440     fECAW0                  = recParam->GetW0();
441     fMinECut                = recParam->GetMinECut();    
442     fToUnfold               = recParam->GetUnfold();
443     if(fToUnfold) AliWarning("Cluster Unfolding ON. Implementing only for eta=0 case!!!"); 
444     fECALocMaxCut           = recParam->GetLocMaxCut();
445     fTimeCut                = recParam->GetTimeCut();
446     fTimeMin                = recParam->GetTimeMin();
447     fTimeMax                = recParam->GetTimeMax();
448
449     AliDebug(1,Form("Reconstruction parameters: fECAClusteringThreshold=%.3f GeV, fECAW=%.3f, fMinECut=%.3f GeV, fToUnfold=%d, fECALocMaxCut=%.3f GeV, fTimeCut=%e s,fTimeMin=%e s,fTimeMax=%e s",
450                     fECAClusteringThreshold,fECAW0,fMinECut,fToUnfold,fECALocMaxCut,fTimeCut, fTimeMin, fTimeMax));
451   }
452
453 }
454
455 //____________________________________________________________________________
456 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2, Bool_t & shared) const
457 {
458         // Gives the neighbourness of two digits = 0 are not neighbour ; continue searching 
459         //                                       = 1 are neighbour
460         //                                       = 2 is in different SM; continue searching 
461         // In case it is in different SM, but same phi rack, check if neigbours at eta=0
462         // neighbours are defined as digits having at least a common side 
463         // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster 
464         //                                      which is compared to a digit (d2)  not yet in a cluster  
465         
466         static Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0, iphi1=0, ieta1=0;
467         static Int_t nSupMod2=0, nModule2=0, nIphi2=0, nIeta2=0, iphi2=0, ieta2=0;
468         static Int_t rowdiff, coldiff;
469
470         shared = kFALSE;
471         
472         fGeom->GetCellIndex(d1->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
473         fGeom->GetCellIndex(d2->GetId(), nSupMod2,nModule2,nIphi2,nIeta2);
474         fGeom->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, iphi1,ieta1);
475         fGeom->GetCellPhiEtaIndexInSModule(nSupMod2,nModule2,nIphi2,nIeta2, iphi2,ieta2);
476         
477         //If different SM, check if they are in the same phi, then consider cells close to eta=0 as neighbours; May 2010
478         if(nSupMod1 != nSupMod2 ) {
479                 //Check if the 2 SM are in the same PHI position (0,1), (2,3), ...
480                 Float_t smPhi1 = fGeom->GetEMCGeometry()->GetPhiCenterOfSM(nSupMod1);
481                 Float_t smPhi2 = fGeom->GetEMCGeometry()->GetPhiCenterOfSM(nSupMod2);
482                 
483                 if(!TMath::AreEqualAbs(smPhi1, smPhi2, 1e-3)) return 2; //Not phi rack equal, not neighbours
484                                 
485                 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
486                 // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0
487                 if(nSupMod1%2) ieta1+=AliEMCALGeoParams::fgkEMCALCols;
488                 else           ieta2+=AliEMCALGeoParams::fgkEMCALCols;
489                 
490                 shared = kTRUE; // maybe a shared cluster, we know this later, set it for the moment.
491                 
492         }//Different SM, same phi
493         
494         rowdiff = TMath::Abs(iphi1 - iphi2);  
495         coldiff = TMath::Abs(ieta1 - ieta2) ;  
496
497         // neighbours with at least common side; May 11, 2007
498         if ((coldiff==0 && TMath::Abs(rowdiff)==1) || (rowdiff==0 && TMath::Abs(coldiff)==1)) {  
499         //Diagonal?
500         //if ((coldiff==0 && TMath::Abs(rowdiff==1)) || (rowdiff==0 && TMath::Abs(coldiff==1)) || (TMath::Abs(rowdiff)==1 && TMath::Abs(coldiff==1))) rv = 1;
501         
502         if (gDebug == 2) 
503                 printf("AliEMCALClusterizerv1::AreNeighbours(): id1=%d, (row %d, col %d) ; id2=%d, (row %d, col %d), shared %d \n",
504                                 d1->GetId(), iphi1,ieta1, d2->GetId(), iphi2,ieta2, shared);   
505         
506                 return 1;
507         }//Neighbours
508         else {
509                 shared = kFALSE;
510                 return 2 ; 
511         }//Not neighbours
512 }
513
514 //____________________________________________________________________________
515 void AliEMCALClusterizerv1::MakeClusters()
516 {
517   // Steering method to construct the clusters stored in a list of Reconstructed Points
518   // A cluster is defined as a list of neighbour digits
519   // Mar 03, 2007 by PAI
520
521   if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");
522
523   fRecPoints->Clear();
524
525   // Set up TObjArray with pointers to digits to work on 
526   TObjArray *digitsC = new TObjArray();
527   TIter nextdigit(fDigitsArr);
528   AliEMCALDigit *digit;
529   while ( (digit = dynamic_cast<AliEMCALDigit*>(nextdigit())) ) {
530     digitsC->AddLast(digit);
531   }
532
533   double e = 0.0, ehs = 0.0;
534   TIter nextdigitC(digitsC);
535   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // clean up digits
536     e = Calibrate(digit->GetAmplitude(), digit->GetTime(),digit->GetId());//Time or TimeR?
537     if ( e < fMinECut) //|| digit->GetTimeR() > fTimeCut ) // time window of cell checked in calibrate
538       digitsC->Remove(digit);
539     else    
540       ehs += e;
541   } 
542   AliDebug(1,Form("MakeClusters: Number of digits %d  -> (e %f), ehs %d\n",
543                   fDigitsArr->GetEntries(),fMinECut,ehs));
544
545   nextdigitC.Reset();
546
547   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // scan over the list of digitsC
548     TArrayI clusterECAdigitslist(fDigitsArr->GetEntries());
549
550     if(fGeom->CheckAbsCellId(digit->GetId()) && (Calibrate(digit->GetAmplitude(), digit->GetTime(),digit->GetId()) > fECAClusteringThreshold  ) ){
551       // start a new Tower RecPoint
552       if(fNumberOfECAClusters >= fRecPoints->GetSize()) fRecPoints->Expand(2*fNumberOfECAClusters+1) ;
553
554       AliEMCALRecPoint *recPoint = new  AliEMCALRecPoint("") ; 
555       fRecPoints->AddAt(recPoint, fNumberOfECAClusters) ;
556       recPoint = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(fNumberOfECAClusters)) ; 
557       fNumberOfECAClusters++ ; 
558
559       recPoint->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
560
561       recPoint->AddDigit(*digit, Calibrate(digit->GetAmplitude(), digit->GetTime(),digit->GetId()),kFALSE) ; //Time or TimeR?
562       TObjArray clusterDigits;
563       clusterDigits.AddLast(digit);     
564       digitsC->Remove(digit) ; 
565
566       AliDebug(1,Form("MakeClusters: OK id = %d, ene = %f , cell.th. = %f \n", digit->GetId(),
567       Calibrate(digit->GetAmplitude(),digit->GetTime(),digit->GetId()), fECAClusteringThreshold));  //Time or TimeR?
568           Float_t time = digit->GetTime();//Time or TimeR?
569       // Grow cluster by finding neighbours
570       TIter nextClusterDigit(&clusterDigits);
571       while ( (digit = dynamic_cast<AliEMCALDigit*>(nextClusterDigit())) ) { // scan over digits in cluster 
572         TIter nextdigitN(digitsC); 
573         AliEMCALDigit *digitN = 0; // digi neighbor
574         while ( (digitN = (AliEMCALDigit *)nextdigitN()) ) { // scan over all digits to look for neighbours
575                         
576                         //Do not add digits with too different time 
577                         Bool_t shared = kFALSE;//cluster shared by 2 SuperModules?
578                         if(TMath::Abs(time - digitN->GetTime()) > fTimeCut ) continue; //Time or TimeR?
579                         if (AreNeighbours(digit, digitN, shared)==1) {      // call (digit,digitN) in THAT order !!!!! 
580                                 recPoint->AddDigit(*digitN, Calibrate(digitN->GetAmplitude(), digitN->GetTime(), digitN->GetId()),shared) ;//Time or TimeR?
581                                 clusterDigits.AddLast(digitN) ; 
582                                 digitsC->Remove(digitN) ; 
583                         } // if(ineb==1)
584                 } // scan over digits
585       } // scan over digits already in cluster
586         
587       if(recPoint)
588         AliDebug(2,Form("MakeClusters: %d digitd, energy %f \n", clusterDigits.GetEntries(), recPoint->GetEnergy())); 
589     } // If seed found
590   } // while digit 
591
592   delete digitsC ;
593   
594   AliDebug(1,Form("total no of clusters %d from %d digits",fNumberOfECAClusters,fDigitsArr->GetEntriesFast())); 
595 }
596
597 //____________________________________________________________________________
598 void AliEMCALClusterizerv1::MakeUnfolding()
599 {
600   // Unfolds clusters using the shape of an ElectroMagnetic shower
601   // Performs unfolding of all clusters
602
603   if(fNumberOfECAClusters > 0){
604     if (fGeom==0)
605       AliFatal("Did not get geometry from EMCALLoader") ;
606     Int_t nModulesToUnfold = fGeom->GetNCells();
607
608     Int_t numberofNotUnfolded = fNumberOfECAClusters ;
609     Int_t index ;
610     for(index = 0 ; index < numberofNotUnfolded ; index++){
611
612       AliEMCALRecPoint * recPoint = dynamic_cast<AliEMCALRecPoint *>( fRecPoints->At(index) ) ;
613
614       TVector3 gpos;
615       Int_t absId;
616       recPoint->GetGlobalPosition(gpos);
617       fGeom->GetAbsCellIdFromEtaPhi(gpos.Eta(),gpos.Phi(),absId);
618       if(absId > nModulesToUnfold)
619         break ;
620
621       Int_t nMultipl = recPoint->GetMultiplicity() ;
622       AliEMCALDigit ** maxAt = new AliEMCALDigit*[nMultipl] ;
623       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
624       Int_t nMax = recPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fECALocMaxCut,fDigitsArr) ;
625
626       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0
627         UnfoldCluster(recPoint, nMax, maxAt, maxAtEnergy) ;
628         fRecPoints->Remove(recPoint);
629         fRecPoints->Compress() ;
630         index-- ;
631         fNumberOfECAClusters-- ;
632         numberofNotUnfolded-- ;
633       }
634       else{
635         recPoint->SetNExMax(1) ; //Only one local maximum
636       }
637
638       delete[] maxAt ;
639       delete[] maxAtEnergy ;
640     }
641   }
642   // End of Unfolding of clusters
643 }
644
645 //____________________________________________________________________________
646 Double_t  AliEMCALClusterizerv1::ShowerShape(Double_t x, Double_t y)
647
648   // Shape of the shower
649   // If you change this function, change also the gradient evaluation in ChiSquare()
650
651   Double_t r = sqrt(x*x+y*y);
652   Double_t r133  = TMath::Power(r, 1.33) ;
653   Double_t r669  = TMath::Power(r, 6.69) ;
654   Double_t shape = TMath::Exp( -r133 * (1. / (1.57 + 0.0860 * r133) - 0.55 / (1 + 0.000563 * r669) ) ) ;
655   return shape ;
656 }
657
658 //____________________________________________________________________________
659 void  AliEMCALClusterizerv1::UnfoldCluster(AliEMCALRecPoint * iniTower, 
660                                            Int_t nMax, 
661                                            AliEMCALDigit ** maxAt, 
662                                            Float_t * maxAtEnergy)
663 {
664   // Performs the unfolding of a cluster with nMax overlapping showers 
665   Int_t nPar = 3 * nMax ;
666   Float_t * fitparameters = new Float_t[nPar] ;
667
668   if (fGeom==0)
669     AliFatal("Did not get geometry from EMCALLoader") ;
670
671   Bool_t rv = FindFit(iniTower, maxAt, maxAtEnergy, nPar, fitparameters) ;
672   if( !rv ) {
673     // Fit failed, return and remove cluster
674     iniTower->SetNExMax(-1) ;
675     delete[] fitparameters ;
676     return ;
677   }
678
679   // create unfolded rec points and fill them with new energy lists
680   // First calculate energy deposited in each sell in accordance with
681   // fit (without fluctuations): efit[]
682   // and later correct this number in acordance with actual energy
683   // deposition
684
685   Int_t nDigits = iniTower->GetMultiplicity() ;
686   Float_t * efit = new Float_t[nDigits] ;
687   Double_t xDigit=0.,yDigit=0.,zDigit=0. ;
688   Float_t xpar=0.,zpar=0.,epar=0.  ;
689
690   AliEMCALDigit * digit = 0 ;
691   Int_t * digitsList = iniTower->GetDigitsList() ;
692
693   Int_t iparam ;
694   Int_t iDigit ;
695   for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
696     digit = dynamic_cast<AliEMCALDigit*>( fDigitsArr->At(digitsList[iDigit] ) ) ;
697     fGeom->RelPosCellInSModule(digit->GetId(), yDigit, xDigit, zDigit);
698     efit[iDigit] = 0;
699
700     iparam = 0 ;
701     while(iparam < nPar ){
702       xpar = fitparameters[iparam] ;
703       zpar = fitparameters[iparam+1] ;
704       epar = fitparameters[iparam+2] ;
705       iparam += 3 ;
706       efit[iDigit] += epar * ShowerShape(xDigit - xpar,zDigit - zpar) ;
707     }
708   }
709
710
711   // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
712   // so that energy deposited in each cell is distributed between new clusters proportionally
713   // to its contribution to efit
714
715   Float_t * energiesList = iniTower->GetEnergiesList() ;
716   Float_t ratio ;
717
718   iparam = 0 ;
719   while(iparam < nPar ){
720     xpar = fitparameters[iparam] ;
721     zpar = fitparameters[iparam+1] ;
722     epar = fitparameters[iparam+2] ;
723     iparam += 3 ;
724
725     AliEMCALRecPoint * recPoint = 0 ;
726
727     if(fNumberOfECAClusters >= fRecPoints->GetSize())
728       fRecPoints->Expand(2*fNumberOfECAClusters) ;
729
730     (*fRecPoints)[fNumberOfECAClusters] = new AliEMCALRecPoint("") ;
731     recPoint = dynamic_cast<AliEMCALRecPoint *>( fRecPoints->At(fNumberOfECAClusters) ) ;
732     fNumberOfECAClusters++ ;
733     recPoint->SetNExMax((Int_t)nPar/3) ;
734
735     Float_t eDigit ;
736     for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
737       digit = dynamic_cast<AliEMCALDigit*>( fDigitsArr->At( digitsList[iDigit] ) ) ;
738       fGeom->RelPosCellInSModule(digit->GetId(), yDigit, xDigit, zDigit);
739
740       ratio = epar * ShowerShape(xDigit - xpar,zDigit - zpar) / efit[iDigit] ;
741       eDigit = energiesList[iDigit] * ratio ;
742       recPoint->AddDigit( *digit, eDigit, kFALSE ) ; //FIXME, need to study the shared case
743     }
744   }
745
746   delete[] fitparameters ;
747   delete[] efit ;
748
749 }
750
751 //_____________________________________________________________________________
752 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad,
753                                                Double_t & fret,
754                                                Double_t * x, Int_t iflag)
755 {
756   // Calculates the Chi square for the cluster unfolding minimization
757   // Number of parameters, Gradient, Chi squared, parameters, what to do
758
759   TList * toMinuit = dynamic_cast<TList*>( gMinuit->GetObjectFit() ) ;
760
761   AliEMCALRecPoint * recPoint = dynamic_cast<AliEMCALRecPoint*>( toMinuit->At(0) )  ;
762   TClonesArray * digits = dynamic_cast<TClonesArray*>( toMinuit->At(1) )  ;
763   // A bit buggy way to get an access to the geometry
764   // To be revised!
765   AliEMCALGeometry *geom = dynamic_cast<AliEMCALGeometry *>(toMinuit->At(2));
766
767   Int_t * digitsList     = recPoint->GetDigitsList() ;
768
769   Int_t nOdigits = recPoint->GetDigitsMultiplicity() ;
770
771   Float_t * energiesList = recPoint->GetEnergiesList() ;
772
773   fret = 0. ;
774   Int_t iparam ;
775
776   if(iflag == 2)
777     for(iparam = 0 ; iparam < nPar ; iparam++)
778       Grad[iparam] = 0 ; // Will evaluate gradient
779
780   Double_t efit ;
781
782   AliEMCALDigit * digit ;
783   Int_t iDigit ;
784
785   for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
786
787     digit = dynamic_cast<AliEMCALDigit*>( digits->At( digitsList[iDigit] ) );
788
789     Double_t xDigit=0 ;
790     Double_t zDigit=0 ;
791     Double_t yDigit=0 ;//not used yet, assumed to be 0
792
793     geom->RelPosCellInSModule(digit->GetId(), yDigit, xDigit, zDigit);
794
795     if(iflag == 2){  // calculate gradient
796       Int_t iParam = 0 ;
797       efit = 0 ;
798       while(iParam < nPar ){
799         Double_t dx = (xDigit - x[iParam]) ;
800         iParam++ ;
801         Double_t dz = (zDigit - x[iParam]) ;
802         iParam++ ;
803         efit += x[iParam] * ShowerShape(dx,dz) ;
804         iParam++ ;
805       }
806       Double_t sum = 2. * (efit - energiesList[iDigit]) / energiesList[iDigit] ; // Here we assume, that sigma = sqrt(E)
807       iParam = 0 ;
808       while(iParam < nPar ){
809         Double_t xpar = x[iParam] ;
810         Double_t zpar = x[iParam+1] ;
811         Double_t epar = x[iParam+2] ;
812         Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
813         Double_t shape = sum * ShowerShape(xDigit - xpar,zDigit - zpar) ;
814         Double_t r133 =  TMath::Power(dr, 1.33);
815         Double_t r669 = TMath::Power(dr,6.69);
816         Double_t deriv =-1.33 * TMath::Power(dr,0.33)*dr * ( 1.57 / ( (1.57 + 0.0860 * r133) * (1.57 + 0.0860 * r133) )
817                                                              - 0.55 / (1 + 0.000563 * r669) / ( (1 + 0.000563 * r669) * (1 + 0.000563 * r669) ) ) ;
818
819         Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ;  // Derivative over x
820         iParam++ ;
821         Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ;  // Derivative over z
822         iParam++ ;
823         Grad[iParam] += shape ;                                  // Derivative over energy
824         iParam++ ;
825       }
826     }
827     efit = 0;
828     iparam = 0 ;
829
830
831     while(iparam < nPar ){
832       Double_t xpar = x[iparam] ;
833       Double_t zpar = x[iparam+1] ;
834       Double_t epar = x[iparam+2] ;
835       iparam += 3 ;
836       efit += epar * ShowerShape(xDigit - xpar,zDigit - zpar) ;
837     }
838
839     fret += (efit-energiesList[iDigit])*(efit-energiesList[iDigit])/energiesList[iDigit] ;
840     // Here we assume, that sigma = sqrt(E) 
841   }
842 }
843 //____________________________________________________________________________
844 void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
845 {
846   // Print clusterizer parameters
847
848   TString message("\n") ; 
849   
850   if( strcmp(GetName(), "") !=0 ){
851     
852     // Print parameters
853  
854     TString taskName(Version()) ;
855     
856     printf("--------------- "); 
857     printf("%s",taskName.Data()) ; 
858     printf(" "); 
859     printf("Clusterizing digits: "); 
860     printf("\n                       ECA Local Maximum cut    = %f", fECALocMaxCut); 
861     printf("\n                       ECA Logarithmic weight   = %f", fECAW0); 
862     if(fToUnfold)
863       printf("\nUnfolding on\n");
864     else
865       printf("\nUnfolding off\n");
866     
867     printf("------------------------------------------------------------------"); 
868   }
869   else
870     printf("AliEMCALClusterizerv1 not initialized ") ;
871 }
872
873 //____________________________________________________________________________
874 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
875 {
876   // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
877   if(strstr(option,"deb")) {
878     printf("PrintRecPoints: Clusterization result:") ; 
879   
880     printf("           Found %d ECA Rec Points\n ", 
881          fRecPoints->GetEntriesFast()) ; 
882   }
883
884   if(strstr(option,"all")) {
885     if(strstr(option,"deb")) {
886       printf("\n-----------------------------------------------------------------------\n") ;
887       printf("Clusters in ECAL section\n") ;
888       printf("Index    Ene(GeV) Multi Module     GX    GY   GZ  lX    lY   lZ   Dispersion Lambda 1   Lambda 2  # of prim  Primaries list\n") ;
889     }
890    Int_t index =0;
891
892     for (index = 0 ; index < fRecPoints->GetEntries() ; index++) {
893       AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint * >(fRecPoints->At(index)) ; 
894       TVector3  globalpos;  
895       //rp->GetGlobalPosition(globalpos);
896       TVector3  localpos;  
897       rp->GetLocalPosition(localpos);
898       Float_t lambda[2]; 
899       rp->GetElipsAxis(lambda);
900       Int_t * primaries; 
901       Int_t nprimaries;
902       primaries = rp->GetPrimaries(nprimaries);
903       if(strstr(option,"deb")) 
904       printf("\n%6d  %8.4f  %3d     %4.1f    %4.1f %4.1f  %4.1f %4.1f %4.1f    %4.1f   %4f  %4f    %2d     : ", 
905              rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(),
906              globalpos.X(), globalpos.Y(), globalpos.Z(), localpos.X(), localpos.Y(), localpos.Z(), 
907              rp->GetDispersion(), lambda[0], lambda[1], nprimaries) ; 
908       if(strstr(option,"deb")){ 
909         for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
910           printf("%d ", primaries[iprimary] ) ; 
911         }
912       }
913     }
914
915     if(strstr(option,"deb"))
916     printf("\n-----------------------------------------------------------------------\n");
917   }
918 }
919
920 //___________________________________________________________________
921 void  AliEMCALClusterizerv1::PrintRecoInfo()
922 {
923   printf(" AliEMCALClusterizerv1::PrintRecoInfo() : version %s \n", Version() );
924
925 }