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