]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALClusterizerv1.cxx
583f522dc2941d9c7efe81108246b41f5169ba91
[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
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(0.),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(0.),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(0.),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(Int_t amp, 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, 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         if(fCaloPed->IsBadChannel(iSupMod,ieta,iphi)) {
191                   AliDebug(2,Form("Tower from SM %d, ieta %d, iphi %d is BAD!!!",iSupMod,ieta,iphi));
192                   return 0;
193         }
194           
195     fADCchannelECA  = fCalibData->GetADCchannel (iSupMod,ieta,iphi);
196     fADCpedestalECA = fCalibData->GetADCpedestal(iSupMod,ieta,iphi);
197   
198    return -fADCpedestalECA + amp * fADCchannelECA ;        
199  
200   }
201   else //Return energy with default parameters if calibration is not available
202     return -fADCpedestalECA + amp * fADCchannelECA ; 
203   
204 }
205
206 //____________________________________________________________________________
207 void AliEMCALClusterizerv1::Digits2Clusters(Option_t * option)
208 {
209   // Steering method to perform clusterization for the current event 
210   // in AliEMCALLoader
211
212   if(strstr(option,"tim"))
213     gBenchmark->Start("EMCALClusterizer"); 
214   
215   if(strstr(option,"print"))
216     Print("") ; 
217  
218   //Get calibration parameters from file or digitizer default values.
219   GetCalibrationParameters() ;
220
221   //Get dead channel map from file or digitizer default values.
222   GetCaloCalibPedestal() ;
223         
224   fNumberOfECAClusters = 0;
225
226   MakeClusters() ;  //only the real clusters
227
228   if(fToUnfold)
229     MakeUnfolding() ;
230
231   Int_t index ;
232
233   //Evaluate position, dispersion and other RecPoint properties for EC section                      
234   for(index = 0; index < fRecPoints->GetEntries(); index++) {
235       dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index))->EvalAll(fECAW0,fDigitsArr) ;
236           //For each rec.point set the distance to the nearest bad crystal
237           dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index))->EvalDistanceToBadChannels(fCaloPed);
238   }
239
240   fRecPoints->Sort() ;
241
242   for(index = 0; index < fRecPoints->GetEntries(); index++) {
243     (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index)))->SetIndexInList(index) ;
244     (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index)))->Print();
245   }
246
247   fTreeR->Fill();
248   
249   if(strstr(option,"deb") || strstr(option,"all"))  
250     PrintRecPoints(option) ;
251
252   AliDebug(1,Form("EMCAL Clusterizer found %d Rec Points",fRecPoints->GetEntriesFast()));
253
254   fRecPoints->Delete();
255
256   if(strstr(option,"tim")){
257     gBenchmark->Stop("EMCALClusterizer");
258     printf("Exec took %f seconds for Clusterizing", 
259            gBenchmark->GetCpuTime("EMCALClusterizer"));
260   }    
261 }
262
263 //____________________________________________________________________________
264 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALRecPoint * recPoint, AliEMCALDigit ** maxAt, 
265                                       const Float_t* maxAtEnergy,
266                                       Int_t nPar, Float_t * fitparameters) const
267 {
268   // Calls TMinuit to fit the energy distribution of a cluster with several maxima
269   // The initial values for fitting procedure are set equal to the
270   // positions of local maxima.       
271   // Cluster will be fitted as a superposition of nPar/3
272   // electromagnetic showers
273
274   if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");
275
276   gMinuit->mncler();                     // Reset Minuit's list of paramters
277   gMinuit->SetPrintLevel(-1) ;           // No Printout
278   gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;
279   // To set the address of the minimization function
280   TList * toMinuit = new TList();
281   toMinuit->AddAt(recPoint,0) ;
282   toMinuit->AddAt(fDigitsArr,1) ;
283   toMinuit->AddAt(fGeom,2) ;
284
285   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare
286
287   // filling initial values for fit parameters
288   AliEMCALDigit * digit ;
289
290   Int_t ierflg  = 0;
291   Int_t index   = 0 ;
292   Int_t nDigits = (Int_t) nPar / 3 ;
293
294   Int_t iDigit ;
295
296   for(iDigit = 0; iDigit < nDigits; iDigit++){
297     digit = maxAt[iDigit];
298     Double_t x = 0.;
299     Double_t y = 0.;
300     Double_t z = 0.;
301
302     fGeom->RelPosCellInSModule(digit->GetId(), y, x, z);
303
304     Float_t energy = maxAtEnergy[iDigit] ;
305
306     gMinuit->mnparm(index, "x",  x, 0.1, 0, 0, ierflg) ;
307     index++ ;
308     if(ierflg != 0){
309       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : x = %f", x ) ;
310       return kFALSE;
311     }
312     gMinuit->mnparm(index, "z",  z, 0.1, 0, 0, ierflg) ;
313     index++ ;
314     if(ierflg != 0){
315       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : z = %f", z) ;
316       return kFALSE;
317     }
318     gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;
319     index++ ;
320     if(ierflg != 0){
321       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : energy = %f", energy) ;
322       return kFALSE;
323     }
324   }
325
326   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; 
327                       // The number of function call slightly depends on it.
328   //Double_t p1 = 1.0 ;
329   Double_t p2 = 0.0 ;
330
331   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls
332   //  gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient
333   gMinuit->SetMaxIterations(5);
334   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
335   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize
336
337   if(ierflg == 4){  // Minimum not found
338     Error("FindFit", "EMCAL Unfolding  Fit not converged, cluster abandoned " ) ;
339     return kFALSE ;
340   }
341   for(index = 0; index < nPar; index++){
342     Double_t err ;
343     Double_t val ;
344     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index
345     fitparameters[index] = val ;
346   }
347
348   delete toMinuit ;
349   return kTRUE;
350
351 }
352
353 //____________________________________________________________________________
354 void AliEMCALClusterizerv1::GetCalibrationParameters() 
355 {
356   // Set calibration parameters:
357   // if calibration database exists, they are read from database,
358   // otherwise, they are taken from digitizer.
359   //
360   // It is a user responsilibity to open CDB before reconstruction, 
361   // for example: 
362   // AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage("local://CalibDB");
363
364   //Check if calibration is stored in data base
365
366   if(!fCalibData)
367     {
368       AliCDBEntry *entry = (AliCDBEntry*) 
369         AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
370       if (entry) fCalibData =  (AliEMCALCalibData*) entry->GetObject();
371     }
372   
373   if(!fCalibData)
374     AliFatal("Calibration parameters not found in CDB!");
375  
376 }
377
378 //____________________________________________________________________________
379 void AliEMCALClusterizerv1::GetCaloCalibPedestal() 
380 {
381         // Set calibration parameters:
382         // if calibration database exists, they are read from database,
383         // otherwise, they are taken from digitizer.
384         //
385         // It is a user responsilibity to open CDB before reconstruction, 
386         // for example: 
387         // AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage("local://CalibDB");
388         
389         //Check if calibration is stored in data base
390         
391         if(!fCaloPed)
392     {
393                 AliCDBEntry *entry = (AliCDBEntry*) 
394                 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
395                 if (entry) fCaloPed =  (AliCaloCalibPedestal*) entry->GetObject();
396     }
397         
398         if(!fCaloPed)
399                 AliFatal("Pedestal info not found in CDB!");
400         
401 }
402
403
404 //____________________________________________________________________________
405 void AliEMCALClusterizerv1::Init()
406 {
407   // Make all memory allocations which can not be done in default constructor.
408   // Attach the Clusterizer task to the list of EMCAL tasks
409   
410   AliRunLoader *rl = AliRunLoader::Instance();
411   if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL"))
412     fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
413   else 
414     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
415
416   AliDebug(1,Form("geom 0x%x",fGeom));
417
418   if(!gMinuit) 
419     gMinuit = new TMinuit(100) ;
420
421 }
422
423 //____________________________________________________________________________
424 void AliEMCALClusterizerv1::InitParameters()
425
426   // Initializes the parameters for the Clusterizer
427   fNumberOfECAClusters = 0;
428   fTimeCut = 1. ; // Originally 300 ns time cut, in data time found to be between 350 ns and 1500 ns, relax the cut for the moment.
429   //Gustavo, 17-12-09
430
431   fCalibData               = 0 ;
432   fCaloPed                 = 0 ;
433         
434   const AliEMCALRecParam* recParam = AliEMCALReconstructor::GetRecParam();
435   if(!recParam) {
436     AliFatal("Reconstruction parameters for EMCAL not set!");
437   } else {
438     fECAClusteringThreshold = recParam->GetClusteringThreshold();
439     fECAW0                  = recParam->GetW0();
440     fMinECut                = recParam->GetMinECut();    
441     fToUnfold               = recParam->GetUnfold();
442     if(fToUnfold) AliWarning("Cluster Unfolding ON. Implementing only for eta=0 case!!!"); 
443     fECALocMaxCut           = recParam->GetLocMaxCut();
444     
445     AliDebug(1,Form("Reconstruction parameters: fECAClusteringThreshold=%.3f, fECAW=%.3f, fMinECut=%.3f, fToUnfold=%d, fECALocMaxCut=%.3f",
446                     fECAClusteringThreshold,fECAW0,fMinECut,fToUnfold,fECALocMaxCut));
447   }
448
449 }
450
451 //____________________________________________________________________________
452 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2) const
453 {
454   // Gives the neighbourness of two digits = 0 are not neighbour ; continue searching 
455   //                                       = 1 are neighbour
456   //                                       = 2 is in different SM; continue searching 
457   // neighbours are defined as digits having at least a common vertex 
458   // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster 
459   //                                      which is compared to a digit (d2)  not yet in a cluster  
460
461   static Int_t rv; 
462   static Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0, iphi1=0, ieta1=0;
463   static Int_t nSupMod2=0, nModule2=0, nIphi2=0, nIeta2=0, iphi2=0, ieta2=0;
464   static Int_t rowdiff, coldiff;
465   rv = 0 ; 
466
467   fGeom->GetCellIndex(d1->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
468   fGeom->GetCellIndex(d2->GetId(), nSupMod2,nModule2,nIphi2,nIeta2);
469   if(nSupMod1 != nSupMod2) return 2; // different SM
470
471   fGeom->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, iphi1,ieta1);
472   fGeom->GetCellPhiEtaIndexInSModule(nSupMod2,nModule2,nIphi2,nIeta2, iphi2,ieta2);
473
474   rowdiff = TMath::Abs(iphi1 - iphi2);  
475   coldiff = TMath::Abs(ieta1 - ieta2) ;  
476   
477   // neighbours with at least commom side; May 11, 2007
478   if ((coldiff==0 && abs(rowdiff)==1) || (rowdiff==0 && abs(coldiff)==1)) rv = 1;  
479  
480   if (gDebug == 2 && rv==1) 
481   printf("AreNeighbours: neighbours=%d, id1=%d, relid1=%d,%d \n id2=%d, relid2=%d,%d \n", 
482          rv, d1->GetId(), iphi1,ieta1, d2->GetId(), iphi2,ieta2);   
483   
484   return rv ; 
485 }
486
487 //____________________________________________________________________________
488 void AliEMCALClusterizerv1::MakeClusters()
489 {
490   // Steering method to construct the clusters stored in a list of Reconstructed Points
491   // A cluster is defined as a list of neighbour digits
492   // Mar 03, 2007 by PAI
493
494   if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");
495
496   fRecPoints->Clear();
497
498   // Set up TObjArray with pointers to digits to work on 
499   TObjArray *digitsC = new TObjArray();
500   TIter nextdigit(fDigitsArr);
501   AliEMCALDigit *digit;
502   while ( (digit = dynamic_cast<AliEMCALDigit*>(nextdigit())) ) {
503     digitsC->AddLast(digit);
504   }
505
506   double e = 0.0, ehs = 0.0;
507   TIter nextdigitC(digitsC);
508
509   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // clean up digits
510     e = Calibrate(digit->GetAmp(), digit->GetId());
511     if ( e < fMinECut || digit->GetTimeR() > fTimeCut ) 
512       digitsC->Remove(digit);
513     else    
514       ehs += e;
515   } 
516   AliDebug(1,Form("MakeClusters: Number of digits %d  -> (e %f), ehs %d\n",
517                   fDigitsArr->GetEntries(),fMinECut,ehs));
518
519   nextdigitC.Reset();
520
521   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // scan over the list of digitsC
522     TArrayI clusterECAdigitslist(fDigitsArr->GetEntries());
523
524     if(fGeom->CheckAbsCellId(digit->GetId()) && (Calibrate(digit->GetAmp(), digit->GetId()) > fECAClusteringThreshold  ) ){
525       // start a new Tower RecPoint
526       if(fNumberOfECAClusters >= fRecPoints->GetSize()) fRecPoints->Expand(2*fNumberOfECAClusters+1) ;
527
528       AliEMCALRecPoint *recPoint = new  AliEMCALRecPoint("") ; 
529       fRecPoints->AddAt(recPoint, fNumberOfECAClusters) ;
530       recPoint = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(fNumberOfECAClusters)) ; 
531       fNumberOfECAClusters++ ; 
532
533       recPoint->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
534
535       recPoint->AddDigit(*digit, Calibrate(digit->GetAmp(), digit->GetId())) ; 
536       TObjArray clusterDigits;
537       clusterDigits.AddLast(digit);     
538       digitsC->Remove(digit) ; 
539
540       AliDebug(1,Form("MakeClusters: OK id = %d, ene = %f , cell.th. = %f \n", digit->GetId(),
541       Calibrate(digit->GetAmp(),digit->GetId()), fECAClusteringThreshold));  
542       
543       // Grow cluster by finding neighbours
544       TIter nextClusterDigit(&clusterDigits);
545       while ( (digit = dynamic_cast<AliEMCALDigit*>(nextClusterDigit())) ) { // scan over digits in cluster 
546         TIter nextdigitN(digitsC); 
547         AliEMCALDigit *digitN = 0; // digi neighbor
548         while ( (digitN = (AliEMCALDigit *)nextdigitN()) ) { // scan over all digits to look for neighbours
549           if (AreNeighbours(digit, digitN)==1) {      // call (digit,digitN) in THAT oder !!!!! 
550             recPoint->AddDigit(*digitN, Calibrate(digitN->GetAmp(),digitN->GetId()) ) ;
551             clusterDigits.AddLast(digitN) ; 
552             digitsC->Remove(digitN) ; 
553           } // if(ineb==1)
554         } // scan over digits
555       } // scan over digits already in cluster
556       if(recPoint)
557         AliDebug(2,Form("MakeClusters: %d digitd, energy %f \n", clusterDigits.GetEntries(), recPoint->GetEnergy())); 
558     } // If seed found
559   } // while digit 
560
561   delete digitsC ;
562   
563   AliDebug(1,Form("total no of clusters %d from %d digits",fNumberOfECAClusters,fDigitsArr->GetEntriesFast())); 
564 }
565
566 //____________________________________________________________________________
567 void AliEMCALClusterizerv1::MakeUnfolding()
568 {
569   // Unfolds clusters using the shape of an ElectroMagnetic shower
570   // Performs unfolding of all clusters
571
572   if(fNumberOfECAClusters > 0){
573     if (fGeom==0)
574       AliFatal("Did not get geometry from EMCALLoader") ;
575     Int_t nModulesToUnfold = fGeom->GetNCells();
576
577     Int_t numberofNotUnfolded = fNumberOfECAClusters ;
578     Int_t index ;
579     for(index = 0 ; index < numberofNotUnfolded ; index++){
580
581       AliEMCALRecPoint * recPoint = dynamic_cast<AliEMCALRecPoint *>( fRecPoints->At(index) ) ;
582
583       TVector3 gpos;
584       Int_t absId;
585       recPoint->GetGlobalPosition(gpos);
586       fGeom->GetAbsCellIdFromEtaPhi(gpos.Eta(),gpos.Phi(),absId);
587       if(absId > nModulesToUnfold)
588         break ;
589
590       Int_t nMultipl = recPoint->GetMultiplicity() ;
591       AliEMCALDigit ** maxAt = new AliEMCALDigit*[nMultipl] ;
592       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
593       Int_t nMax = recPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fECALocMaxCut,fDigitsArr) ;
594
595       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0
596         UnfoldCluster(recPoint, nMax, maxAt, maxAtEnergy) ;
597         fRecPoints->Remove(recPoint);
598         fRecPoints->Compress() ;
599         index-- ;
600         fNumberOfECAClusters-- ;
601         numberofNotUnfolded-- ;
602       }
603       else{
604         recPoint->SetNExMax(1) ; //Only one local maximum
605       }
606
607       delete[] maxAt ;
608       delete[] maxAtEnergy ;
609     }
610   }
611   // End of Unfolding of clusters
612 }
613
614 //____________________________________________________________________________
615 Double_t  AliEMCALClusterizerv1::ShowerShape(Double_t x, Double_t y)
616
617   // Shape of the shower
618   // If you change this function, change also the gradient evaluation in ChiSquare()
619
620   Double_t r = sqrt(x*x+y*y);
621   Double_t r133  = TMath::Power(r, 1.33) ;
622   Double_t r669  = TMath::Power(r, 6.69) ;
623   Double_t shape = TMath::Exp( -r133 * (1. / (1.57 + 0.0860 * r133) - 0.55 / (1 + 0.000563 * r669) ) ) ;
624   return shape ;
625 }
626
627 //____________________________________________________________________________
628 void  AliEMCALClusterizerv1::UnfoldCluster(AliEMCALRecPoint * iniTower, 
629                                            Int_t nMax, 
630                                            AliEMCALDigit ** maxAt, 
631                                            Float_t * maxAtEnergy)
632 {
633   // Performs the unfolding of a cluster with nMax overlapping showers 
634   Int_t nPar = 3 * nMax ;
635   Float_t * fitparameters = new Float_t[nPar] ;
636
637   if (fGeom==0)
638     AliFatal("Did not get geometry from EMCALLoader") ;
639
640   Bool_t rv = FindFit(iniTower, maxAt, maxAtEnergy, nPar, fitparameters) ;
641   if( !rv ) {
642     // Fit failed, return and remove cluster
643     iniTower->SetNExMax(-1) ;
644     delete[] fitparameters ;
645     return ;
646   }
647
648   // create unfolded rec points and fill them with new energy lists
649   // First calculate energy deposited in each sell in accordance with
650   // fit (without fluctuations): efit[]
651   // and later correct this number in acordance with actual energy
652   // deposition
653
654   Int_t nDigits = iniTower->GetMultiplicity() ;
655   Float_t * efit = new Float_t[nDigits] ;
656   Double_t xDigit=0.,yDigit=0.,zDigit=0. ;
657   Float_t xpar=0.,zpar=0.,epar=0.  ;
658
659   AliEMCALDigit * digit = 0 ;
660   Int_t * digitsList = iniTower->GetDigitsList() ;
661
662   Int_t iparam ;
663   Int_t iDigit ;
664   for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
665     digit = dynamic_cast<AliEMCALDigit*>( fDigitsArr->At(digitsList[iDigit] ) ) ;
666     fGeom->RelPosCellInSModule(digit->GetId(), yDigit, xDigit, zDigit);
667     efit[iDigit] = 0;
668
669     iparam = 0 ;
670     while(iparam < nPar ){
671       xpar = fitparameters[iparam] ;
672       zpar = fitparameters[iparam+1] ;
673       epar = fitparameters[iparam+2] ;
674       iparam += 3 ;
675       efit[iDigit] += epar * ShowerShape(xDigit - xpar,zDigit - zpar) ;
676     }
677   }
678
679
680   // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
681   // so that energy deposited in each cell is distributed between new clusters proportionally
682   // to its contribution to efit
683
684   Float_t * energiesList = iniTower->GetEnergiesList() ;
685   Float_t ratio ;
686
687   iparam = 0 ;
688   while(iparam < nPar ){
689     xpar = fitparameters[iparam] ;
690     zpar = fitparameters[iparam+1] ;
691     epar = fitparameters[iparam+2] ;
692     iparam += 3 ;
693
694     AliEMCALRecPoint * recPoint = 0 ;
695
696     if(fNumberOfECAClusters >= fRecPoints->GetSize())
697       fRecPoints->Expand(2*fNumberOfECAClusters) ;
698
699     (*fRecPoints)[fNumberOfECAClusters] = new AliEMCALRecPoint("") ;
700     recPoint = dynamic_cast<AliEMCALRecPoint *>( fRecPoints->At(fNumberOfECAClusters) ) ;
701     fNumberOfECAClusters++ ;
702     recPoint->SetNExMax((Int_t)nPar/3) ;
703
704     Float_t eDigit ;
705     for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
706       digit = dynamic_cast<AliEMCALDigit*>( fDigitsArr->At( digitsList[iDigit] ) ) ;
707       fGeom->RelPosCellInSModule(digit->GetId(), yDigit, xDigit, zDigit);
708
709       ratio = epar * ShowerShape(xDigit - xpar,zDigit - zpar) / efit[iDigit] ;
710       eDigit = energiesList[iDigit] * ratio ;
711       recPoint->AddDigit( *digit, eDigit ) ;
712     }
713   }
714
715   delete[] fitparameters ;
716   delete[] efit ;
717
718 }
719
720 //_____________________________________________________________________________
721 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad,
722                                                Double_t & fret,
723                                                Double_t * x, Int_t iflag)
724 {
725   // Calculates the Chi square for the cluster unfolding minimization
726   // Number of parameters, Gradient, Chi squared, parameters, what to do
727
728   TList * toMinuit = dynamic_cast<TList*>( gMinuit->GetObjectFit() ) ;
729
730   AliEMCALRecPoint * recPoint = dynamic_cast<AliEMCALRecPoint*>( toMinuit->At(0) )  ;
731   TClonesArray * digits = dynamic_cast<TClonesArray*>( toMinuit->At(1) )  ;
732   // A bit buggy way to get an access to the geometry
733   // To be revised!
734   AliEMCALGeometry *geom = dynamic_cast<AliEMCALGeometry *>(toMinuit->At(2));
735
736   Int_t * digitsList     = recPoint->GetDigitsList() ;
737
738   Int_t nOdigits = recPoint->GetDigitsMultiplicity() ;
739
740   Float_t * energiesList = recPoint->GetEnergiesList() ;
741
742   fret = 0. ;
743   Int_t iparam ;
744
745   if(iflag == 2)
746     for(iparam = 0 ; iparam < nPar ; iparam++)
747       Grad[iparam] = 0 ; // Will evaluate gradient
748
749   Double_t efit ;
750
751   AliEMCALDigit * digit ;
752   Int_t iDigit ;
753
754   for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
755
756     digit = dynamic_cast<AliEMCALDigit*>( digits->At( digitsList[iDigit] ) );
757
758     Double_t xDigit=0 ;
759     Double_t zDigit=0 ;
760     Double_t yDigit=0 ;//not used yet, assumed to be 0
761
762     geom->RelPosCellInSModule(digit->GetId(), yDigit, xDigit, zDigit);
763
764     if(iflag == 2){  // calculate gradient
765       Int_t iParam = 0 ;
766       efit = 0 ;
767       while(iParam < nPar ){
768         Double_t dx = (xDigit - x[iParam]) ;
769         iParam++ ;
770         Double_t dz = (zDigit - x[iParam]) ;
771         iParam++ ;
772         efit += x[iParam] * ShowerShape(dx,dz) ;
773         iParam++ ;
774       }
775       Double_t sum = 2. * (efit - energiesList[iDigit]) / energiesList[iDigit] ; // Here we assume, that sigma = sqrt(E)
776       iParam = 0 ;
777       while(iParam < nPar ){
778         Double_t xpar = x[iParam] ;
779         Double_t zpar = x[iParam+1] ;
780         Double_t epar = x[iParam+2] ;
781         Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
782         Double_t shape = sum * ShowerShape(xDigit - xpar,zDigit - zpar) ;
783         Double_t r133 =  TMath::Power(dr, 1.33);
784         Double_t r669 = TMath::Power(dr,6.69);
785         Double_t deriv =-1.33 * TMath::Power(dr,0.33)*dr * ( 1.57 / ( (1.57 + 0.0860 * r133) * (1.57 + 0.0860 * r133) )
786                                                              - 0.55 / (1 + 0.000563 * r669) / ( (1 + 0.000563 * r669) * (1 + 0.000563 * r669) ) ) ;
787
788         Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ;  // Derivative over x
789         iParam++ ;
790         Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ;  // Derivative over z
791         iParam++ ;
792         Grad[iParam] += shape ;                                  // Derivative over energy
793         iParam++ ;
794       }
795     }
796     efit = 0;
797     iparam = 0 ;
798
799
800     while(iparam < nPar ){
801       Double_t xpar = x[iparam] ;
802       Double_t zpar = x[iparam+1] ;
803       Double_t epar = x[iparam+2] ;
804       iparam += 3 ;
805       efit += epar * ShowerShape(xDigit - xpar,zDigit - zpar) ;
806     }
807
808     fret += (efit-energiesList[iDigit])*(efit-energiesList[iDigit])/energiesList[iDigit] ;
809     // Here we assume, that sigma = sqrt(E) 
810   }
811 }
812 //____________________________________________________________________________
813 void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
814 {
815   // Print clusterizer parameters
816
817   TString message("\n") ; 
818   
819   if( strcmp(GetName(), "") !=0 ){
820     
821     // Print parameters
822  
823     TString taskName(Version()) ;
824     
825     printf("--------------- "); 
826     printf("%s",taskName.Data()) ; 
827     printf(" "); 
828     printf("Clusterizing digits: "); 
829     printf("\n                       ECA Local Maximum cut    = %f", fECALocMaxCut); 
830     printf("\n                       ECA Logarithmic weight   = %f", fECAW0); 
831     if(fToUnfold)
832       printf("\nUnfolding on\n");
833     else
834       printf("\nUnfolding off\n");
835     
836     printf("------------------------------------------------------------------"); 
837   }
838   else
839     printf("AliEMCALClusterizerv1 not initialized ") ;
840 }
841
842 //____________________________________________________________________________
843 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
844 {
845   // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
846   if(strstr(option,"deb")) {
847     printf("PrintRecPoints: Clusterization result:") ; 
848   
849     printf("           Found %d ECA Rec Points\n ", 
850          fRecPoints->GetEntriesFast()) ; 
851   }
852
853   if(strstr(option,"all")) {
854     if(strstr(option,"deb")) {
855       printf("\n-----------------------------------------------------------------------\n") ;
856       printf("Clusters in ECAL section\n") ;
857       printf("Index    Ene(GeV) Multi Module     GX    GY   GZ  lX    lY   lZ   Dispersion Lambda 1   Lambda 2  # of prim  Primaries list\n") ;
858     }
859    Int_t index =0;
860
861     for (index = 0 ; index < fRecPoints->GetEntries() ; index++) {
862       AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint * >(fRecPoints->At(index)) ; 
863       TVector3  globalpos;  
864       //rp->GetGlobalPosition(globalpos);
865       TVector3  localpos;  
866       rp->GetLocalPosition(localpos);
867       Float_t lambda[2]; 
868       rp->GetElipsAxis(lambda);
869       Int_t * primaries; 
870       Int_t nprimaries;
871       primaries = rp->GetPrimaries(nprimaries);
872       if(strstr(option,"deb")) 
873       printf("\n%6d  %8.4f  %3d     %4.1f    %4.1f %4.1f  %4.1f %4.1f %4.1f    %4.1f   %4f  %4f    %2d     : ", 
874              rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(),
875              globalpos.X(), globalpos.Y(), globalpos.Z(), localpos.X(), localpos.Y(), localpos.Z(), 
876              rp->GetDispersion(), lambda[0], lambda[1], nprimaries) ; 
877       if(strstr(option,"deb")){ 
878         for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
879           printf("%d ", primaries[iprimary] ) ; 
880         }
881       }
882     }
883
884     if(strstr(option,"deb"))
885     printf("\n-----------------------------------------------------------------------\n");
886   }
887 }
888
889 //___________________________________________________________________
890 void  AliEMCALClusterizerv1::PrintRecoInfo()
891 {
892   printf(" AliEMCALClusterizerv1::PrintRecoInfo() : version %s \n", Version() );
893
894 }