]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALClusterizerv1.cxx
include call to EvalAll() before writing RecPoints to tree
[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 "AliEMCALHistoUtilities.h"
77 #include "AliEMCALRecParam.h"
78 #include "AliEMCALReconstructor.h"
79 #include "AliCDBManager.h"
80
81 class AliCDBStorage;
82 #include "AliCDBEntry.h"
83
84 ClassImp(AliEMCALClusterizerv1)
85
86 //____________________________________________________________________________
87 AliEMCALClusterizerv1::AliEMCALClusterizerv1()
88   : AliEMCALClusterizer(),
89     fHists(0),fPointE(0),fPointL1(0),fPointL2(0),
90     fPointDis(0),fPointMult(0),fDigitAmp(0),fMaxE(0),
91     fMaxL1(0),fMaxL2(0),fMaxDis(0),fGeom(0),
92     fDefaultInit(kFALSE),
93     fToUnfold(kFALSE),
94     fNumberOfECAClusters(0),fNTowerInGroup(0),fCalibData(0),
95     fADCchannelECA(0.),fADCpedestalECA(0.),fECAClusteringThreshold(0.),fECALocMaxCut(0.),
96     fECAW0(0.),fTimeCut(0.),fMinECut(0.)
97 {
98   // ctor with the indication of the file where header Tree and digits Tree are stored
99   
100   InitParameters() ; 
101   Init() ;
102 }
103
104 //____________________________________________________________________________
105   AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
106 {
107   // dtor
108 }
109
110 //____________________________________________________________________________
111 Float_t  AliEMCALClusterizerv1::Calibrate(Int_t amp, Int_t AbsId) 
112 {
113  
114   // Convert digitized amplitude into energy.
115   // Calibration parameters are taken from calibration data base for raw data,
116   // or from digitizer parameters for simulated data.
117
118   if(fCalibData){
119     
120     if (fGeom==0)
121       AliFatal("Did not get geometry from EMCALLoader") ;
122     
123     Int_t iSupMod = -1;
124     Int_t nModule  = -1;
125     Int_t nIphi   = -1;
126     Int_t nIeta   = -1;
127     Int_t iphi    = -1;
128     Int_t ieta    = -1;
129     
130     Bool_t bCell = fGeom->GetCellIndex(AbsId, iSupMod, nModule, nIphi, nIeta) ;
131     if(!bCell) {
132       fGeom->PrintGeometry();
133       Error("Calibrate()"," Wrong cell id number : %i", AbsId);
134       assert(0);
135     }
136
137     fGeom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
138
139     fADCchannelECA  = fCalibData->GetADCchannel (iSupMod,ieta,iphi);
140     fADCpedestalECA = fCalibData->GetADCpedestal(iSupMod,ieta,iphi);
141   
142    return -fADCpedestalECA + amp * fADCchannelECA ;        
143  
144   }
145   else //Return energy with default parameters if calibration is not available
146     return -fADCpedestalECA + amp * fADCchannelECA ; 
147   
148 }
149
150 //____________________________________________________________________________
151 void AliEMCALClusterizerv1::Digits2Clusters(Option_t * option)
152 {
153   // Steering method to perform clusterization for the current event 
154   // in AliEMCALLoader
155
156   if(strstr(option,"tim"))
157     gBenchmark->Start("EMCALClusterizer"); 
158   
159   if(strstr(option,"print"))
160     Print("") ; 
161  
162   //Get calibration parameters from file or digitizer default values.
163   GetCalibrationParameters() ;
164
165
166   fNumberOfECAClusters = 0;
167
168   if(strstr(option,"pseudo"))
169     MakeClusters("pseudo") ;  //both types
170   else
171     MakeClusters("") ;  //only the real clusters
172
173   if(fToUnfold)
174     MakeUnfolding() ;
175
176   Int_t index ;
177
178   //Evaluate position, dispersion and other RecPoint properties for EC section                      
179   for(index = 0; index < fRecPoints->GetEntries(); index++) {
180     if (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index))->GetClusterType() != AliESDCaloCluster::kPseudoCluster)
181       dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index))->EvalAll(fECAW0,fDigitsArr) ;
182   }
183
184   fRecPoints->Sort() ;
185
186   for(index = 0; index < fRecPoints->GetEntries(); index++) {
187     (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index)))->SetIndexInList(index) ;
188     (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index)))->Print();
189   }
190
191   fTreeR->Fill();
192   
193   if(strstr(option,"deb") || strstr(option,"all"))  
194     PrintRecPoints(option) ;
195
196   AliDebug(1,Form("EMCAL Clusterizer found %d Rec Points",fRecPoints->GetEntriesFast()));
197
198   if(strstr(option,"tim")){
199     gBenchmark->Stop("EMCALClusterizer");
200     printf("Exec took %f seconds for Clusterizing", 
201            gBenchmark->GetCpuTime("EMCALClusterizer"));
202   }    
203 }
204
205 //____________________________________________________________________________
206 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
207                                     Int_t nPar, Float_t * fitparameters) const
208
209   // Calls TMinuit to fit the energy distribution of a cluster with several maxima 
210   // The initial values for fitting procedure are set equal to the positions of local maxima.
211   // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
212
213   gMinuit->mncler();                     // Reset Minuit's list of paramters
214   gMinuit->SetPrintLevel(-1) ;           // No Printout
215   gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;  
216                                          // To set the address of the minimization function 
217   TList * toMinuit = new TList();
218   toMinuit->AddAt(emcRP,0) ;
219   toMinuit->AddAt(fDigitsArr,1) ;
220   
221   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare
222
223   // filling initial values for fit parameters
224   AliEMCALDigit * digit ;
225
226   Int_t ierflg  = 0; 
227   Int_t index   = 0 ;
228   Int_t nDigits = (Int_t) nPar / 3 ;
229
230   Int_t iDigit ;
231
232   for(iDigit = 0; iDigit < nDigits; iDigit++){
233     digit = maxAt[iDigit]; 
234
235     Float_t x = 0.;
236     Float_t z = 0.;
237     //   have to be tune for TRD1; May 31,06
238     //   Int_t relid[2] ;
239     //   fGeom->AbsToRelNumbering(digit->GetId(), relid) ; // obsolete method
240     //   fGeom->PosInAlice(relid, x, z) ;                  // obsolete method
241
242     Float_t energy = maxAtEnergy[iDigit] ;
243
244     gMinuit->mnparm(index, "x",  x, 0.1, 0, 0, ierflg) ;
245     index++ ;   
246     if(ierflg != 0){ 
247       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : x = %f",  x ) ;
248       return kFALSE;
249     }
250     gMinuit->mnparm(index, "z",  z, 0.1, 0, 0, ierflg) ;
251     index++ ;   
252     if(ierflg != 0){
253        Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : z = %f", z) ;
254       return kFALSE;
255     }
256     gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;
257     index++ ;   
258     if(ierflg != 0){
259      Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : energy = %f", energy) ;      
260       return kFALSE;
261     }
262   }
263
264   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
265                       //  depends on it. 
266   Double_t p1 = 1.0 ;
267   Double_t p2 = 0.0 ;
268
269   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
270   gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
271   gMinuit->SetMaxIterations(5);
272   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
273
274   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
275
276   if(ierflg == 4){  // Minimum not found   
277     Error("FindFit", "EMCAL Unfolding  Fit not converged, cluster abandoned " ) ;      
278     return kFALSE ;
279   }            
280   for(index = 0; index < nPar; index++){
281     Double_t err ;
282     Double_t val ;
283     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index
284     fitparameters[index] = val ;
285    }
286
287   delete toMinuit ;
288   return kTRUE;
289
290 }
291
292 //____________________________________________________________________________
293 void AliEMCALClusterizerv1::GetCalibrationParameters() 
294 {
295   // Set calibration parameters:
296   // if calibration database exists, they are read from database,
297   // otherwise, they are taken from digitizer.
298   //
299   // It is a user responsilibity to open CDB before reconstruction, 
300   // for example: 
301   // AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage("local://CalibDB");
302
303   //Check if calibration is stored in data base
304
305   if(!fCalibData && (AliCDBManager::Instance()->IsDefaultStorageSet()))
306     {
307       AliCDBEntry *entry = (AliCDBEntry*) 
308         AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
309       if (entry) fCalibData =  (AliEMCALCalibData*) entry->GetObject();
310     }
311   
312   if(!fCalibData)
313     AliFatal("Calibration parameters not found in CDB!");
314  
315   //  Please fix it!! Or better just remove it...
316 //    if(!fCalibData)
317 //      {
318 //        //If calibration is not available use default parameters
319 //        //Loader
320 //        if ( !emcalLoader->Digitizer() ) 
321 //       emcalLoader->LoadDigitizer();
322 //        AliEMCALDigitizer * dig = dynamic_cast<AliEMCALDigitizer*>(emcalLoader->Digitizer());
323        
324 //        fADCchannelECA   = dig->GetECAchannel() ;
325 //        fADCpedestalECA  = dig->GetECApedestal();
326 //      }
327 }
328
329 //____________________________________________________________________________
330 void AliEMCALClusterizerv1::Init()
331 {
332   // Make all memory allocations which can not be done in default constructor.
333   // Attach the Clusterizer task to the list of EMCAL tasks
334   
335   AliRunLoader *rl = AliRunLoader::GetRunLoader();
336   if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL"))
337     fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
338   else 
339     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaulGeometryName());
340
341   fGeom->GetTransformationForSM(); // Global <-> Local
342   AliInfo(Form("geom 0x%x",fGeom));
343
344   if(!gMinuit) 
345     gMinuit = new TMinuit(100) ;
346
347   fHists = BookHists();
348 }
349
350 //____________________________________________________________________________
351 void AliEMCALClusterizerv1::InitParameters()
352
353   // Initializes the parameters for the Clusterizer
354   fNumberOfECAClusters = 0;
355
356   fNTowerInGroup = 36;  //Produces maximum of 80 pseudoclusters per event
357
358   fECALocMaxCut = 0.03; // ??
359
360   fTimeCut = 300e-9 ; // 300 ns time cut (to be tuned) 
361   fToUnfold = kFALSE ;
362
363   fCalibData               = 0 ;
364
365   const AliEMCALRecParam* recParam = AliEMCALReconstructor::GetRecParam();
366   if(!recParam) {
367     AliFatal("Reconstruction parameters for EMCAL not set!");
368   }
369   else {
370     fECAClusteringThreshold = recParam->GetClusteringThreshold();
371     fECAW0                  = recParam->GetW0();
372     fMinECut                = recParam->GetMinECut();
373     AliDebug(1,Form("Reconstruction parameters: fECAClusteringThreshold=%.3f, fECAW=%.3f, fMinECut=%.3f",
374                  fECAClusteringThreshold,fECAW0,fMinECut));
375   }
376
377 }
378
379 //____________________________________________________________________________
380 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2) const
381 {
382   // Gives the neighbourness of two digits = 0 are not neighbour ; continue searching 
383   //                                       = 1 are neighbour
384   //                                       = 2 is in different SM; continue searching 
385   // neighbours are defined as digits having at least a common vertex 
386   // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster 
387   //                                      which is compared to a digit (d2)  not yet in a cluster  
388
389   static Int_t rv; 
390   static Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0, iphi1=0, ieta1=0;
391   static Int_t nSupMod2=0, nModule2=0, nIphi2=0, nIeta2=0, iphi2=0, ieta2=0;
392   static Int_t rowdiff, coldiff;
393   rv = 0 ; 
394
395   fGeom->GetCellIndex(d1->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
396   fGeom->GetCellIndex(d2->GetId(), nSupMod2,nModule2,nIphi2,nIeta2);
397   if(nSupMod1 != nSupMod2) return 2; // different SM
398
399   fGeom->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, iphi1,ieta1);
400   fGeom->GetCellPhiEtaIndexInSModule(nSupMod2,nModule2,nIphi2,nIeta2, iphi2,ieta2);
401
402   rowdiff = TMath::Abs(iphi1 - iphi2);  
403   coldiff = TMath::Abs(ieta1 - ieta2) ;  
404   
405   // neighbours with at least commom side; May 11, 2007
406   if ((coldiff==0 && abs(rowdiff)==1) || (rowdiff==0 && abs(coldiff)==1)) rv = 1;  
407  
408   if (gDebug == 2 && rv==1) 
409   printf("AreNeighbours: neighbours=%d, id1=%d, relid1=%d,%d \n id2=%d, relid2=%d,%d \n", 
410          rv, d1->GetId(), iphi1,ieta1, d2->GetId(), iphi2,ieta2);   
411   
412   return rv ; 
413 }
414
415 //____________________________________________________________________________
416 Int_t AliEMCALClusterizerv1::AreInGroup(AliEMCALDigit * d1, AliEMCALDigit * d2) const
417 {
418   // Tells whether two digits fall within the same supermodule and
419   // tower grouping.  The number of towers in a group is controlled by
420   // the parameter nTowersInGroup
421   //    = 0 are not in same group but continue searching 
422   //    = 1 same group
423   //    = 2 is in different SM, quit from searching
424   //    = 3 different tower group, quit from searching
425   //
426   // The order of d1 and d2 is important: first (d1) should be a digit 
427   // already in a cluster which is compared to a digit (d2)  not yet in a cluster  
428
429   //JLK Question: does the quit from searching assume that the digits
430   //are ordered, so that once you are in a different SM, you'll not
431   //find another in the list that will match?  How about my TowerGroup search?
432
433   static Int_t rv;
434   static Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0, iphi1=0, ieta1=0;
435   static Int_t nSupMod2=0, nModule2=0, nIphi2=0, nIeta2=0, iphi2=0, ieta2=0;
436   static Int_t towerGroup1 = -1, towerGroup2 = -1;
437   rv = 0 ;
438
439   fGeom->GetCellIndex(d1->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
440   fGeom->GetCellIndex(d2->GetId(), nSupMod2,nModule2,nIphi2,nIeta2);
441   if(nSupMod1 != nSupMod2) return 2; // different SM
442
443   static Int_t nTowerInSM = fGeom->GetNCellsInSupMod()/fGeom->GetNCellsInModule();
444
445   //figure out which tower grouping each digit belongs to
446   for(int it = 0; it < nTowerInSM/fNTowerInGroup; it++) {
447     if(nModule1 <= nTowerInSM - it*fNTowerInGroup) towerGroup1 = it;
448     if(nModule2 <= nTowerInSM - it*fNTowerInGroup) towerGroup2 = it;
449   }
450   if(towerGroup1 != towerGroup2) return 3; //different Towergroup
451
452   //same SM, same towergroup, we're happy
453   if(towerGroup1 == towerGroup2 && towerGroup2 >= 0)
454     rv = 1;
455
456   if (gDebug == 2 && rv==1)
457   printf("AreInGroup: neighbours=%d, id1=%d, relid1=%d,%d \n id2=%d, relid2=%d,%d \n",
458          rv, d1->GetId(), iphi1,ieta1, d2->GetId(), iphi2,ieta2);
459
460   return rv ;
461 }
462  
463 //____________________________________________________________________________
464 void AliEMCALClusterizerv1::MakeClusters(char* option)
465 {
466   // Steering method to construct the clusters stored in a list of Reconstructed Points
467   // A cluster is defined as a list of neighbour digits
468   // Mar 03, 2007 by PAI
469
470   if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");
471
472   fRecPoints->Clear();
473
474   // Set up TObjArray with pointers to digits to work on 
475   TObjArray *digitsC = new TObjArray();
476   TIter nextdigit(fDigitsArr);
477   AliEMCALDigit *digit;
478   while ( (digit = dynamic_cast<AliEMCALDigit*>(nextdigit())) ) {
479     digitsC->AddLast(digit);
480   }
481
482   //Start with pseudoclusters, if option
483   if(strstr(option,"pseudo")) { 
484     //
485     // New algorithm : will be created one pseudo cluster per module  
486     //
487     AliDebug(1,Form("Pseudo clustering #digits : %i ",fDigitsArr->GetEntries()));
488
489     AliEMCALRecPoint *recPoints[12]; // max size is 12 : see fGeom->GetNumberOfSuperModules();
490     for(int i=0; i<12; i++) recPoints[i] = 0;
491     TIter nextdigitC(digitsC) ;
492
493     // PseudoClusterization starts  
494     int nSM = 0; // # of SM
495     while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // scan over the list of digitsC
496       if(fGeom->CheckAbsCellId(digit->GetId()) ) { //Is this an EMCAL digit? Just maing sure...
497         nSM = fGeom->GetSuperModuleNumber(digit->GetId());
498         if(recPoints[nSM] == 0) {
499           recPoints[nSM] = new AliEMCALRecPoint(Form("PC%2.2i", nSM));
500           recPoints[nSM]->SetClusterType(AliESDCaloCluster::kPseudoCluster);
501         }
502         recPoints[nSM]->AddDigit(*digit, Calibrate(digit->GetAmp(), digit->GetId()));
503       }
504     }
505     fNumberOfECAClusters = 0;
506     for(int i=0; i<fGeom->GetNumberOfSuperModules(); i++) { // put non empty rec.points to container
507       if(recPoints[i]) fRecPoints->AddAt(recPoints[i], fNumberOfECAClusters++);
508     }
509     AliDebug(1,Form(" Number of PC %d ", fNumberOfECAClusters));
510   }
511
512   //
513   // Now do real clusters
514   //
515
516   double e = 0.0, ehs = 0.0;
517   TIter nextdigitC(digitsC);
518
519   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // clean up digits
520     e = Calibrate(digit->GetAmp(), digit->GetId());
521     AliEMCALHistoUtilities::FillH1(fHists, 10, digit->GetAmp());
522     AliEMCALHistoUtilities::FillH1(fHists, 11, e);
523     if ( e < fMinECut || digit->GetTimeR() > fTimeCut ) 
524       digitsC->Remove(digit);
525     else    
526       ehs += e;
527   } 
528   AliDebug(1,Form("MakeClusters: Number of digits %d  -> (e %f), ehs %d\n",
529                   fDigitsArr->GetEntries(),fMinECut,ehs));
530
531   nextdigitC.Reset();
532
533   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // scan over the list of digitsC
534     TArrayI clusterECAdigitslist(fDigitsArr->GetEntries());
535
536     if(fGeom->CheckAbsCellId(digit->GetId()) && (Calibrate(digit->GetAmp(), digit->GetId()) > fECAClusteringThreshold  ) ){
537       // start a new Tower RecPoint
538       if(fNumberOfECAClusters >= fRecPoints->GetSize()) fRecPoints->Expand(2*fNumberOfECAClusters+1) ;
539       AliEMCALRecPoint *recPoint = new  AliEMCALRecPoint("") ; 
540       fRecPoints->AddAt(recPoint, fNumberOfECAClusters) ;
541       recPoint = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(fNumberOfECAClusters)) ; 
542       fNumberOfECAClusters++ ; 
543
544       recPoint->SetClusterType(AliESDCaloCluster::kClusterv1);
545
546       recPoint->AddDigit(*digit, Calibrate(digit->GetAmp(), digit->GetId())) ; 
547       TObjArray clusterDigits;
548       clusterDigits.AddLast(digit);     
549       digitsC->Remove(digit) ; 
550
551       AliDebug(1,Form("MakeClusters: OK id = %d, ene = %f , cell.th. = %f \n", digit->GetId(),
552       Calibrate(digit->GetAmp(),digit->GetId()), fECAClusteringThreshold));  
553       
554       // Grow cluster by finding neighbours
555       TIter nextClusterDigit(&clusterDigits);
556       while ( (digit = dynamic_cast<AliEMCALDigit*>(nextClusterDigit())) ) { // scan over digits in cluster 
557         TIter nextdigitN(digitsC); 
558         AliEMCALDigit *digitN = 0; // digi neighbor
559         while ( (digitN = (AliEMCALDigit *)nextdigitN()) ) { // scan over all digits to look for neighbours
560           if (AreNeighbours(digit, digitN)==1) {      // call (digit,digitN) in THAT oder !!!!! 
561             recPoint->AddDigit(*digitN, Calibrate(digitN->GetAmp(),digitN->GetId()) ) ;
562             clusterDigits.AddLast(digitN) ; 
563             digitsC->Remove(digitN) ; 
564           } // if(ineb==1)
565         } // scan over digits
566       } // scan over digits already in cluster
567       if(recPoint)
568         AliDebug(2,Form("MakeClusters: %d digitd, energy %f \n", clusterDigits.GetEntries(), recPoint->GetEnergy())); 
569     } // If seed found
570   } // while digit 
571
572   delete digitsC ;
573
574   AliDebug(1,Form("total no of clusters %d from %d digits",fNumberOfECAClusters,fDigitsArr->GetEntriesFast())); 
575 }
576
577 void AliEMCALClusterizerv1::MakeUnfolding() const
578 {
579   Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
580 }
581
582 //____________________________________________________________________________
583 Double_t  AliEMCALClusterizerv1::ShowerShape(Double_t r)
584
585   // Shape of the shower (see EMCAL TDR)
586   // If you change this function, change also the gradient evaluation in ChiSquare()
587
588   Double_t r4    = r*r*r*r ;
589   Double_t r295  = TMath::Power(r, 2.95) ;
590   Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
591   return shape ;
592 }
593
594 //____________________________________________________________________________
595 void  AliEMCALClusterizerv1::UnfoldCluster(AliEMCALRecPoint * /*iniTower*/, 
596                                            Int_t /*nMax*/, 
597                                            AliEMCALDigit ** /*maxAt*/, 
598                                            Float_t * /*maxAtEnergy*/) const
599 {
600   // Performs the unfolding of a cluster with nMax overlapping showers 
601   
602   Fatal("UnfoldCluster", "--> Unfolding not implemented") ;
603
604 }
605
606 //_____________________________________________________________________________
607 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * /*Grad*/,
608                                                Double_t & /*fret*/,
609                                                Double_t * /*x*/, Int_t /*iflag*/)
610 {
611   // Calculates the Chi square for the cluster unfolding minimization
612   // Number of parameters, Gradient, Chi squared, parameters, what to do
613   
614   ::Fatal("UnfoldingChiSquare","Unfolding not implemented") ;
615 }
616 //____________________________________________________________________________
617 void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
618 {
619   // Print clusterizer parameters
620
621   TString message("\n") ; 
622   
623   if( strcmp(GetName(), "") !=0 ){
624     
625     // Print parameters
626  
627     TString taskName(Version()) ;
628     
629     printf("--------------- "); 
630     printf(taskName.Data()) ; 
631     printf(" "); 
632     printf("Clusterizing digits: "); 
633     printf("\n                       ECA Local Maximum cut    = %f", fECALocMaxCut); 
634     printf("\n                       ECA Logarithmic weight   = %f", fECAW0); 
635     if(fToUnfold)
636       printf("\nUnfolding on\n");
637     else
638       printf("\nUnfolding off\n");
639     
640     printf("------------------------------------------------------------------"); 
641   }
642   else
643     printf("AliEMCALClusterizerv1 not initialized ") ;
644 }
645
646 //____________________________________________________________________________
647 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
648 {
649   // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
650   if(strstr(option,"deb")) {
651     printf("PrintRecPoints: Clusterization result:") ; 
652   
653     printf("           Found %d ECA Rec Points\n ", 
654          fRecPoints->GetEntriesFast()) ; 
655   }
656
657   if(strstr(option,"all")) {
658     if(strstr(option,"deb")) {
659       printf("\n-----------------------------------------------------------------------\n") ;
660       printf("Clusters in ECAL section\n") ;
661       printf("Index    Ene(GeV) Multi Module     GX    GY   GZ  lX    lY   lZ   Dispersion Lambda 1   Lambda 2  # of prim  Primaries list\n") ;
662     }
663    Int_t index =0;
664    Float_t maxE=0; 
665    Float_t maxL1=0; 
666    Float_t maxL2=0; 
667    Float_t maxDis=0; 
668
669     AliEMCALHistoUtilities::FillH1(fHists, 12, double(fRecPoints->GetEntries()));
670
671     for (index = 0 ; index < fRecPoints->GetEntries() ; index++) {
672       AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint * >(fRecPoints->At(index)) ; 
673       TVector3  globalpos;  
674       //rp->GetGlobalPosition(globalpos);
675       TVector3  localpos;  
676       rp->GetLocalPosition(localpos);
677       Float_t lambda[2]; 
678       rp->GetElipsAxis(lambda);
679       Int_t * primaries; 
680       Int_t nprimaries;
681       primaries = rp->GetPrimaries(nprimaries);
682       if(strstr(option,"deb")) 
683       printf("\n%6d  %8.4f  %3d     %4.1f    %4.1f %4.1f  %4.1f %4.1f %4.1f    %4.1f   %4f  %4f    %2d     : ", 
684              rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(),
685              globalpos.X(), globalpos.Y(), globalpos.Z(), localpos.X(), localpos.Y(), localpos.Z(), 
686              rp->GetDispersion(), lambda[0], lambda[1], nprimaries) ; 
687   /////////////
688       if(rp->GetEnergy()>maxE){
689               maxE=rp->GetEnergy();
690               maxL1=lambda[0];
691               maxL2=lambda[1];
692               maxDis=rp->GetDispersion();
693       }
694       fPointE->Fill(rp->GetEnergy());
695       fPointL1->Fill(lambda[0]);
696       fPointL2->Fill(lambda[1]);
697       fPointDis->Fill(rp->GetDispersion());
698       fPointMult->Fill(rp->GetMultiplicity());
699       ///////////// 
700       if(strstr(option,"deb")){ 
701         for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
702           printf("%d ", primaries[iprimary] ) ; 
703         }
704       }
705     }
706
707       fMaxE->Fill(maxE);
708       fMaxL1->Fill(maxL1);
709       fMaxL2->Fill(maxL2);
710       fMaxDis->Fill(maxDis);
711
712     if(strstr(option,"deb"))
713     printf("\n-----------------------------------------------------------------------\n");
714   }
715 }
716 TList* AliEMCALClusterizerv1::BookHists()
717 {
718   //set up histograms for monitoring clusterizer performance
719
720   gROOT->cd();
721
722         fPointE = new TH1F("00_pointE","point energy", 2000, 0.0, 150.);
723         fPointL1 = new TH1F("01_pointL1","point L1", 1000, 0.0, 3.);
724         fPointL2 = new TH1F("02_pointL2","point L2", 1000, 0.0, 3.);
725         fPointDis = new TH1F("03_pointDisp","point dispersion", 1000, 0.0, 10.);
726         fPointMult = new TH1F("04_pointMult","#cell in point(cluster)", 101, -0.5, 100.5);
727         fDigitAmp = new TH1F("05_digitAmp","Digit Amplitude", 2000, 0.0, 5000.);
728         fMaxE = new TH1F("06_maxE","Max point energy", 2000, 0.0, 150.);
729         fMaxL1 = new TH1F("07_maxL1","Largest (first) of eigenvalue of covariance matrix", 1000, 0.0, 3.);
730         fMaxL2 = new TH1F("08_maxL2","Smalest (second) of eigenvalue of covariace matrix", 1000, 0.0, 3.);
731         fMaxDis = new TH1F("09_maxDis","Point dispersion", 1000, 0.0, 10.); // 9
732         //
733         new TH1F("10_adcOfDigits","adc of digits(threshold control)", 1001, -0.5, 1000.5);   // 10
734         new TH1F("11_energyOfDigits","energy of digits(threshold control)", 1000, 0.0, 1.);  // 11
735         new TH1F("12_numberOfPoints","number of points(clusters)", 101, -0.5, 100.5);        // 12
736
737   return AliEMCALHistoUtilities::MoveHistsToList("EmcalClusterizerv1ControlHists", kFALSE);
738 }
739
740 void AliEMCALClusterizerv1::SaveHists(const char *fn)
741 {
742   AliEMCALHistoUtilities::SaveListOfHists(fHists, fn, kTRUE);
743 }
744
745 void  AliEMCALClusterizerv1::PrintRecoInfo()
746 {
747   printf(" AliEMCALClusterizerv1::PrintRecoInfo() : version %s \n", Version() );
748   TH1F *h = (TH1F*)fHists->At(12);
749   if(h) {
750     printf(" ## Multiplicity of RecPoints ## \n");
751     for(int i=1; i<=h->GetNbinsX(); i++) {
752       int nbin = int((*h)[i]);
753       int mult = int(h->GetBinCenter(i));
754       if(nbin > 0) printf(" %i : %5.5i %6.3f %% \n", mult, nbin, 100.*nbin/h->GetEntries()); 
755     }    
756   }
757 }
758
759 void AliEMCALClusterizerv1::DrawLambdasHists()
760 {
761   if(fMaxL1) {
762     fMaxL1->Draw();
763     if(fMaxL2) fMaxL2->Draw("same");
764     if(fMaxDis) {
765       fMaxDis->Draw("same");
766     }
767   }
768 }