]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALClusterizerv1.cxx
Updates (Eulogio Serradilla <eulogio.serradilla@ciemat.es>)
[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 }
126
127 //____________________________________________________________________________
128 AliEMCALClusterizerv1::AliEMCALClusterizerv1(AliEMCALGeometry* geometry, AliEMCALCalibData * calib, AliCaloCalibPedestal * caloped)
129 : AliEMCALClusterizer(),
130 fGeom(geometry),
131 fDefaultInit(kFALSE),
132 fToUnfold(kFALSE),
133 fNumberOfECAClusters(0),fCalibData(calib), fCaloPed(caloped),
134 fADCchannelECA(0.),fADCpedestalECA(0.),fECAClusteringThreshold(0.),fECALocMaxCut(0.),
135 fECAW0(0.),fTimeCut(1.),fTimeMin(-1.),fTimeMax(1.),fMinECut(0.)
136 {
137         // ctor, geometry and calibration are initialized elsewhere.
138         
139         if (!fGeom)
140                 AliFatal("Geometry not initialized.");
141                         
142 }
143
144
145 //____________________________________________________________________________
146   AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
147 {
148   // dtor
149 }
150
151 //____________________________________________________________________________
152 Float_t  AliEMCALClusterizerv1::Calibrate(const Float_t amp, const Float_t time, const Int_t absId) 
153 {
154  
155   // Convert digitized amplitude into energy.
156   // Calibration parameters are taken from calibration data base for raw data,
157   // or from digitizer parameters for simulated data.
158
159   if(fCalibData){
160     
161     if (fGeom==0)
162       AliFatal("Did not get geometry from EMCALLoader") ;
163     
164     Int_t iSupMod = -1;
165     Int_t nModule = -1;
166     Int_t nIphi   = -1;
167     Int_t nIeta   = -1;
168     Int_t iphi    = -1;
169     Int_t ieta    = -1;
170     
171     Bool_t bCell = fGeom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
172     if(!bCell) {
173       fGeom->PrintGeometry();
174       Error("Calibrate()"," Wrong cell id number : %i", absId);
175       assert(0);
176     }
177
178     fGeom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
179           
180         // Check if channel is bad (dead or hot), in this case return 0.        
181         // Gustavo: 15-12-09 In case of RAW data this selection is already done, but not in simulation.
182         // for the moment keep it here but remember to do the selection at the sdigitizer level 
183         // and remove it from here
184         Int_t channelStatus = (Int_t)(fCaloPed->GetDeadMap(iSupMod))->GetBinContent(ieta,iphi);
185         if(channelStatus == AliCaloCalibPedestal::kHot || channelStatus == AliCaloCalibPedestal::kDead) {
186                   AliDebug(2,Form("Tower from SM %d, ieta %d, iphi %d is BAD : status %d !!!",iSupMod,ieta,iphi, channelStatus));
187                   return 0;
188         }
189         //Check if time is too large or too small, indication of a noisy channel, remove in this case
190         if(time > fTimeMax || time < fTimeMin) return 0;
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   fRecPoints->Delete();
252
253   if(strstr(option,"tim")){
254     gBenchmark->Stop("EMCALClusterizer");
255     printf("Exec took %f seconds for Clusterizing", 
256            gBenchmark->GetCpuTime("EMCALClusterizer"));
257   }    
258 }
259
260 //____________________________________________________________________________
261 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALRecPoint * recPoint, AliEMCALDigit ** maxAt, 
262                                       const Float_t* maxAtEnergy,
263                                       Int_t nPar, Float_t * fitparameters) const
264 {
265   // Calls TMinuit to fit the energy distribution of a cluster with several maxima
266   // The initial values for fitting procedure are set equal to the
267   // positions of local maxima.       
268   // Cluster will be fitted as a superposition of nPar/3
269   // electromagnetic showers
270
271   if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");
272         
273   if(!gMinuit)
274      gMinuit = new TMinuit(100) ;
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
429   fCalibData               = 0 ;
430   fCaloPed                 = 0 ;
431         
432   const AliEMCALRecParam* recParam = AliEMCALReconstructor::GetRecParam();
433   if(!recParam) {
434     AliFatal("Reconstruction parameters for EMCAL not set!");
435   } else {
436     fECAClusteringThreshold = recParam->GetClusteringThreshold();
437     fECAW0                  = recParam->GetW0();
438     fMinECut                = recParam->GetMinECut();    
439     fToUnfold               = recParam->GetUnfold();
440     if(fToUnfold) AliWarning("Cluster Unfolding ON. Implementing only for eta=0 case!!!"); 
441     fECALocMaxCut           = recParam->GetLocMaxCut();
442     fTimeCut                = recParam->GetTimeCut();
443     fTimeMin                = recParam->GetTimeMin();
444     fTimeMax                = recParam->GetTimeMax();
445
446     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",
447                     fECAClusteringThreshold,fECAW0,fMinECut,fToUnfold,fECALocMaxCut,fTimeCut, fTimeMin, fTimeMax));
448   }
449
450 }
451
452 //____________________________________________________________________________
453 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2, Bool_t & shared) const
454 {
455         // Gives the neighbourness of two digits = 0 are not neighbour ; continue searching 
456         //                                       = 1 are neighbour
457         //                                       = 2 is in different SM; continue searching 
458         // In case it is in different SM, but same phi rack, check if neigbours at eta=0
459         // neighbours are defined as digits having at least a common side 
460         // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster 
461         //                                      which is compared to a digit (d2)  not yet in a cluster  
462         
463         static Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0, iphi1=0, ieta1=0;
464         static Int_t nSupMod2=0, nModule2=0, nIphi2=0, nIeta2=0, iphi2=0, ieta2=0;
465         static Int_t rowdiff, coldiff;
466
467         shared = kFALSE;
468         
469         fGeom->GetCellIndex(d1->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
470         fGeom->GetCellIndex(d2->GetId(), nSupMod2,nModule2,nIphi2,nIeta2);
471         fGeom->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, iphi1,ieta1);
472         fGeom->GetCellPhiEtaIndexInSModule(nSupMod2,nModule2,nIphi2,nIeta2, iphi2,ieta2);
473         
474         //If different SM, check if they are in the same phi, then consider cells close to eta=0 as neighbours; May 2010
475         if(nSupMod1 != nSupMod2 ) {
476                 //Check if the 2 SM are in the same PHI position (0,1), (2,3), ...
477                 Float_t smPhi1 = fGeom->GetEMCGeometry()->GetPhiCenterOfSM(nSupMod1);
478                 Float_t smPhi2 = fGeom->GetEMCGeometry()->GetPhiCenterOfSM(nSupMod2);
479                 
480                 if(!TMath::AreEqualAbs(smPhi1, smPhi2, 1e-3)) return 2; //Not phi rack equal, not neighbours
481                                 
482                 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
483                 // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0
484                 if(nSupMod1%2) ieta1+=AliEMCALGeoParams::fgkEMCALCols;
485                 else           ieta2+=AliEMCALGeoParams::fgkEMCALCols;
486                 
487                 shared = kTRUE; // maybe a shared cluster, we know this later, set it for the moment.
488                 
489         }//Different SM, same phi
490         
491         rowdiff = TMath::Abs(iphi1 - iphi2);  
492         coldiff = TMath::Abs(ieta1 - ieta2) ;  
493
494         // neighbours with at least common side; May 11, 2007
495         if ((coldiff==0 && TMath::Abs(rowdiff)==1) || (rowdiff==0 && TMath::Abs(coldiff)==1)) {  
496         //Diagonal?
497         //if ((coldiff==0 && TMath::Abs(rowdiff==1)) || (rowdiff==0 && TMath::Abs(coldiff==1)) || (TMath::Abs(rowdiff)==1 && TMath::Abs(coldiff==1))) rv = 1;
498         
499         if (gDebug == 2) 
500                 printf("AliEMCALClusterizerv1::AreNeighbours(): id1=%d, (row %d, col %d) ; id2=%d, (row %d, col %d), shared %d \n",
501                                 d1->GetId(), iphi1,ieta1, d2->GetId(), iphi2,ieta2, shared);   
502         
503                 return 1;
504         }//Neighbours
505         else {
506                 shared = kFALSE;
507                 return 2 ; 
508         }//Not neighbours
509 }
510
511 //____________________________________________________________________________
512 void AliEMCALClusterizerv1::MakeClusters()
513 {
514   // Steering method to construct the clusters stored in a list of Reconstructed Points
515   // A cluster is defined as a list of neighbour digits
516   // Mar 03, 2007 by PAI
517
518   if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");
519
520   fRecPoints->Clear();
521
522   // Set up TObjArray with pointers to digits to work on 
523   TObjArray *digitsC = new TObjArray();
524   TIter nextdigit(fDigitsArr);
525   AliEMCALDigit *digit;
526   while ( (digit = dynamic_cast<AliEMCALDigit*>(nextdigit())) ) {
527     digitsC->AddLast(digit);
528   }
529
530   double e = 0.0, ehs = 0.0;
531   TIter nextdigitC(digitsC);
532   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // clean up digits
533     e = Calibrate(digit->GetAmplitude(), digit->GetTime(),digit->GetId());//Time or TimeR?
534     if ( e < fMinECut) //|| digit->GetTimeR() > fTimeCut ) // time window of cell checked in calibrate
535       digitsC->Remove(digit);
536     else    
537       ehs += e;
538   } 
539   AliDebug(1,Form("MakeClusters: Number of digits %d  -> (e %f), ehs %d\n",
540                   fDigitsArr->GetEntries(),fMinECut,ehs));
541
542   nextdigitC.Reset();
543
544   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // scan over the list of digitsC
545     TArrayI clusterECAdigitslist(fDigitsArr->GetEntries());
546
547     if(fGeom->CheckAbsCellId(digit->GetId()) && (Calibrate(digit->GetAmplitude(), digit->GetTime(),digit->GetId()) > fECAClusteringThreshold  ) ){
548       // start a new Tower RecPoint
549       if(fNumberOfECAClusters >= fRecPoints->GetSize()) fRecPoints->Expand(2*fNumberOfECAClusters+1) ;
550
551       AliEMCALRecPoint *recPoint = new  AliEMCALRecPoint("") ; 
552       fRecPoints->AddAt(recPoint, fNumberOfECAClusters) ;
553       recPoint = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(fNumberOfECAClusters)) ; 
554       fNumberOfECAClusters++ ; 
555
556       recPoint->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
557
558       recPoint->AddDigit(*digit, Calibrate(digit->GetAmplitude(), digit->GetTime(),digit->GetId()),kFALSE) ; //Time or TimeR?
559       TObjArray clusterDigits;
560       clusterDigits.AddLast(digit);     
561       digitsC->Remove(digit) ; 
562
563       AliDebug(1,Form("MakeClusters: OK id = %d, ene = %f , cell.th. = %f \n", digit->GetId(),
564       Calibrate(digit->GetAmplitude(),digit->GetTime(),digit->GetId()), fECAClusteringThreshold));  //Time or TimeR?
565           Float_t time = digit->GetTime();//Time or TimeR?
566       // Grow cluster by finding neighbours
567       TIter nextClusterDigit(&clusterDigits);
568       while ( (digit = dynamic_cast<AliEMCALDigit*>(nextClusterDigit())) ) { // scan over digits in cluster 
569         TIter nextdigitN(digitsC); 
570         AliEMCALDigit *digitN = 0; // digi neighbor
571         while ( (digitN = (AliEMCALDigit *)nextdigitN()) ) { // scan over all digits to look for neighbours
572                         
573                         //Do not add digits with too different time 
574                         Bool_t shared = kFALSE;//cluster shared by 2 SuperModules?
575                         if(TMath::Abs(time - digitN->GetTime()) > fTimeCut ) continue; //Time or TimeR?
576                         if (AreNeighbours(digit, digitN, shared)==1) {      // call (digit,digitN) in THAT order !!!!! 
577                                 recPoint->AddDigit(*digitN, Calibrate(digitN->GetAmplitude(), digitN->GetTime(), digitN->GetId()),shared) ;//Time or TimeR?
578                                 clusterDigits.AddLast(digitN) ; 
579                                 digitsC->Remove(digitN) ; 
580                         } // if(ineb==1)
581                 } // scan over digits
582       } // scan over digits already in cluster
583         
584       if(recPoint)
585         AliDebug(2,Form("MakeClusters: %d digitd, energy %f \n", clusterDigits.GetEntries(), recPoint->GetEnergy())); 
586     } // If seed found
587   } // while digit 
588
589   delete digitsC ;
590   
591   AliDebug(1,Form("total no of clusters %d from %d digits",fNumberOfECAClusters,fDigitsArr->GetEntriesFast())); 
592 }
593
594 //____________________________________________________________________________
595 void AliEMCALClusterizerv1::MakeUnfolding()
596 {
597   // Unfolds clusters using the shape of an ElectroMagnetic shower
598   // Performs unfolding of all clusters
599                 
600   if(fNumberOfECAClusters > 0){
601     if (fGeom==0)
602       AliFatal("Did not get geometry from EMCALLoader") ;
603     Int_t nModulesToUnfold = fGeom->GetNCells();
604
605     Int_t numberofNotUnfolded = fNumberOfECAClusters ;
606     Int_t index ;
607     for(index = 0 ; index < numberofNotUnfolded ; index++){
608
609       AliEMCALRecPoint * recPoint = dynamic_cast<AliEMCALRecPoint *>( fRecPoints->At(index) ) ;
610
611       TVector3 gpos;
612       Int_t absId;
613       recPoint->GetGlobalPosition(gpos);
614       fGeom->GetAbsCellIdFromEtaPhi(gpos.Eta(),gpos.Phi(),absId);
615       if(absId > nModulesToUnfold)
616         break ;
617
618       Int_t nMultipl = recPoint->GetMultiplicity() ;
619       AliEMCALDigit ** maxAt = new AliEMCALDigit*[nMultipl] ;
620       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
621       Int_t nMax = recPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fECALocMaxCut,fDigitsArr) ;
622
623       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0
624         UnfoldCluster(recPoint, nMax, maxAt, maxAtEnergy) ;
625         fRecPoints->Remove(recPoint);
626         fRecPoints->Compress() ;
627         index-- ;
628         fNumberOfECAClusters-- ;
629         numberofNotUnfolded-- ;
630       }
631       else{
632         recPoint->SetNExMax(1) ; //Only one local maximum
633       }
634
635       delete[] maxAt ;
636       delete[] maxAtEnergy ;
637     }
638   }
639   // End of Unfolding of clusters
640 }
641
642 //____________________________________________________________________________
643 Double_t  AliEMCALClusterizerv1::ShowerShape(Double_t x, Double_t y)
644
645   // Shape of the shower
646   // If you change this function, change also the gradient evaluation in ChiSquare()
647
648   Double_t r = sqrt(x*x+y*y);
649   Double_t r133  = TMath::Power(r, 1.33) ;
650   Double_t r669  = TMath::Power(r, 6.69) ;
651   Double_t shape = TMath::Exp( -r133 * (1. / (1.57 + 0.0860 * r133) - 0.55 / (1 + 0.000563 * r669) ) ) ;
652   return shape ;
653 }
654
655 //____________________________________________________________________________
656 void  AliEMCALClusterizerv1::UnfoldCluster(AliEMCALRecPoint * iniTower, 
657                                            Int_t nMax, 
658                                            AliEMCALDigit ** maxAt, 
659                                            Float_t * maxAtEnergy)
660 {
661   // Performs the unfolding of a cluster with nMax overlapping showers 
662   Int_t nPar = 3 * nMax ;
663   Float_t * fitparameters = new Float_t[nPar] ;
664
665   if (fGeom==0)
666     AliFatal("Did not get geometry from EMCALLoader") ;
667
668   Bool_t rv = FindFit(iniTower, maxAt, maxAtEnergy, nPar, fitparameters) ;
669   if( !rv ) {
670     // Fit failed, return and remove cluster
671     iniTower->SetNExMax(-1) ;
672     delete[] fitparameters ;
673     return ;
674   }
675
676   // create unfolded rec points and fill them with new energy lists
677   // First calculate energy deposited in each sell in accordance with
678   // fit (without fluctuations): efit[]
679   // and later correct this number in acordance with actual energy
680   // deposition
681
682   Int_t nDigits = iniTower->GetMultiplicity() ;
683   Float_t * efit = new Float_t[nDigits] ;
684   Double_t xDigit=0.,yDigit=0.,zDigit=0. ;
685   Float_t xpar=0.,zpar=0.,epar=0.  ;
686
687   AliEMCALDigit * digit = 0 ;
688   Int_t * digitsList = iniTower->GetDigitsList() ;
689
690   Int_t iparam ;
691   Int_t iDigit ;
692   for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
693     digit = dynamic_cast<AliEMCALDigit*>( fDigitsArr->At(digitsList[iDigit] ) ) ;
694     fGeom->RelPosCellInSModule(digit->GetId(), yDigit, xDigit, zDigit);
695     efit[iDigit] = 0;
696
697     iparam = 0 ;
698     while(iparam < nPar ){
699       xpar = fitparameters[iparam] ;
700       zpar = fitparameters[iparam+1] ;
701       epar = fitparameters[iparam+2] ;
702       iparam += 3 ;
703       efit[iDigit] += epar * ShowerShape(xDigit - xpar,zDigit - zpar) ;
704     }
705   }
706
707
708   // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
709   // so that energy deposited in each cell is distributed between new clusters proportionally
710   // to its contribution to efit
711
712   Float_t * energiesList = iniTower->GetEnergiesList() ;
713   Float_t ratio ;
714
715   iparam = 0 ;
716   while(iparam < nPar ){
717     xpar = fitparameters[iparam] ;
718     zpar = fitparameters[iparam+1] ;
719     epar = fitparameters[iparam+2] ;
720     iparam += 3 ;
721
722     AliEMCALRecPoint * recPoint = 0 ;
723
724     if(fNumberOfECAClusters >= fRecPoints->GetSize())
725       fRecPoints->Expand(2*fNumberOfECAClusters) ;
726
727     (*fRecPoints)[fNumberOfECAClusters] = new AliEMCALRecPoint("") ;
728     recPoint = dynamic_cast<AliEMCALRecPoint *>( fRecPoints->At(fNumberOfECAClusters) ) ;
729     fNumberOfECAClusters++ ;
730     recPoint->SetNExMax((Int_t)nPar/3) ;
731
732     Float_t eDigit ;
733     for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
734       digit = dynamic_cast<AliEMCALDigit*>( fDigitsArr->At( digitsList[iDigit] ) ) ;
735       fGeom->RelPosCellInSModule(digit->GetId(), yDigit, xDigit, zDigit);
736
737       ratio = epar * ShowerShape(xDigit - xpar,zDigit - zpar) / efit[iDigit] ;
738       eDigit = energiesList[iDigit] * ratio ;
739       recPoint->AddDigit( *digit, eDigit, kFALSE ) ; //FIXME, need to study the shared case
740     }
741   }
742
743   delete[] fitparameters ;
744   delete[] efit ;
745
746 }
747
748 //_____________________________________________________________________________
749 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad,
750                                                Double_t & fret,
751                                                Double_t * x, Int_t iflag)
752 {
753   // Calculates the Chi square for the cluster unfolding minimization
754   // Number of parameters, Gradient, Chi squared, parameters, what to do
755
756   TList * toMinuit = dynamic_cast<TList*>( gMinuit->GetObjectFit() ) ;
757
758   AliEMCALRecPoint * recPoint = dynamic_cast<AliEMCALRecPoint*>( toMinuit->At(0) )  ;
759   TClonesArray * digits = dynamic_cast<TClonesArray*>( toMinuit->At(1) )  ;
760   // A bit buggy way to get an access to the geometry
761   // To be revised!
762   AliEMCALGeometry *geom = dynamic_cast<AliEMCALGeometry *>(toMinuit->At(2));
763
764   Int_t * digitsList     = recPoint->GetDigitsList() ;
765
766   Int_t nOdigits = recPoint->GetDigitsMultiplicity() ;
767
768   Float_t * energiesList = recPoint->GetEnergiesList() ;
769
770   fret = 0. ;
771   Int_t iparam ;
772
773   if(iflag == 2)
774     for(iparam = 0 ; iparam < nPar ; iparam++)
775       Grad[iparam] = 0 ; // Will evaluate gradient
776
777   Double_t efit ;
778
779   AliEMCALDigit * digit ;
780   Int_t iDigit ;
781
782   for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
783
784     digit = dynamic_cast<AliEMCALDigit*>( digits->At( digitsList[iDigit] ) );
785
786     Double_t xDigit=0 ;
787     Double_t zDigit=0 ;
788     Double_t yDigit=0 ;//not used yet, assumed to be 0
789
790     geom->RelPosCellInSModule(digit->GetId(), yDigit, xDigit, zDigit);
791
792     if(iflag == 2){  // calculate gradient
793       Int_t iParam = 0 ;
794       efit = 0 ;
795       while(iParam < nPar ){
796         Double_t dx = (xDigit - x[iParam]) ;
797         iParam++ ;
798         Double_t dz = (zDigit - x[iParam]) ;
799         iParam++ ;
800         efit += x[iParam] * ShowerShape(dx,dz) ;
801         iParam++ ;
802       }
803       Double_t sum = 2. * (efit - energiesList[iDigit]) / energiesList[iDigit] ; // Here we assume, that sigma = sqrt(E)
804       iParam = 0 ;
805       while(iParam < nPar ){
806         Double_t xpar = x[iParam] ;
807         Double_t zpar = x[iParam+1] ;
808         Double_t epar = x[iParam+2] ;
809         Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
810         Double_t shape = sum * ShowerShape(xDigit - xpar,zDigit - zpar) ;
811         Double_t r133 =  TMath::Power(dr, 1.33);
812         Double_t r669 = TMath::Power(dr,6.69);
813         Double_t deriv =-1.33 * TMath::Power(dr,0.33)*dr * ( 1.57 / ( (1.57 + 0.0860 * r133) * (1.57 + 0.0860 * r133) )
814                                                              - 0.55 / (1 + 0.000563 * r669) / ( (1 + 0.000563 * r669) * (1 + 0.000563 * r669) ) ) ;
815
816         Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ;  // Derivative over x
817         iParam++ ;
818         Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ;  // Derivative over z
819         iParam++ ;
820         Grad[iParam] += shape ;                                  // Derivative over energy
821         iParam++ ;
822       }
823     }
824     efit = 0;
825     iparam = 0 ;
826
827
828     while(iparam < nPar ){
829       Double_t xpar = x[iparam] ;
830       Double_t zpar = x[iparam+1] ;
831       Double_t epar = x[iparam+2] ;
832       iparam += 3 ;
833       efit += epar * ShowerShape(xDigit - xpar,zDigit - zpar) ;
834     }
835
836     fret += (efit-energiesList[iDigit])*(efit-energiesList[iDigit])/energiesList[iDigit] ;
837     // Here we assume, that sigma = sqrt(E) 
838   }
839 }
840 //____________________________________________________________________________
841 void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
842 {
843   // Print clusterizer parameters
844
845   TString message("\n") ; 
846   
847   if( strcmp(GetName(), "") !=0 ){
848     
849     // Print parameters
850  
851     TString taskName(Version()) ;
852     
853     printf("--------------- "); 
854     printf("%s",taskName.Data()) ; 
855     printf(" "); 
856     printf("Clusterizing digits: "); 
857     printf("\n                       ECA Local Maximum cut    = %f", fECALocMaxCut); 
858     printf("\n                       ECA Logarithmic weight   = %f", fECAW0); 
859     if(fToUnfold)
860       printf("\nUnfolding on\n");
861     else
862       printf("\nUnfolding off\n");
863     
864     printf("------------------------------------------------------------------"); 
865   }
866   else
867     printf("AliEMCALClusterizerv1 not initialized ") ;
868 }
869
870 //____________________________________________________________________________
871 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
872 {
873   // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
874   if(strstr(option,"deb")) {
875     printf("PrintRecPoints: Clusterization result:") ; 
876   
877     printf("           Found %d ECA Rec Points\n ", 
878          fRecPoints->GetEntriesFast()) ; 
879   }
880
881   if(strstr(option,"all")) {
882     if(strstr(option,"deb")) {
883       printf("\n-----------------------------------------------------------------------\n") ;
884       printf("Clusters in ECAL section\n") ;
885       printf("Index    Ene(GeV) Multi Module     GX    GY   GZ  lX    lY   lZ   Dispersion Lambda 1   Lambda 2  # of prim  Primaries list\n") ;
886     }
887    Int_t index =0;
888
889     for (index = 0 ; index < fRecPoints->GetEntries() ; index++) {
890       AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint * >(fRecPoints->At(index)) ; 
891       TVector3  globalpos;  
892       //rp->GetGlobalPosition(globalpos);
893       TVector3  localpos;  
894       rp->GetLocalPosition(localpos);
895       Float_t lambda[2]; 
896       rp->GetElipsAxis(lambda);
897       Int_t * primaries; 
898       Int_t nprimaries;
899       primaries = rp->GetPrimaries(nprimaries);
900       if(strstr(option,"deb")) 
901       printf("\n%6d  %8.4f  %3d     %4.1f    %4.1f %4.1f  %4.1f %4.1f %4.1f    %4.1f   %4f  %4f    %2d     : ", 
902              rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(),
903              globalpos.X(), globalpos.Y(), globalpos.Z(), localpos.X(), localpos.Y(), localpos.Z(), 
904              rp->GetDispersion(), lambda[0], lambda[1], nprimaries) ; 
905       if(strstr(option,"deb")){ 
906         for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
907           printf("%d ", primaries[iprimary] ) ; 
908         }
909       }
910     }
911
912     if(strstr(option,"deb"))
913     printf("\n-----------------------------------------------------------------------\n");
914   }
915 }
916
917 //___________________________________________________________________
918 void  AliEMCALClusterizerv1::PrintRecoInfo()
919 {
920   printf(" AliEMCALClusterizerv1::PrintRecoInfo() : version %s \n", Version() );
921
922 }