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