]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALClusterizerv1.cxx
correction of trivial typo preventing compilation
[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),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   MakeClusters() ;  //only the real clusters
169
170   if(fToUnfold)
171     MakeUnfolding() ;
172
173   Int_t index ;
174
175   //Evaluate position, dispersion and other RecPoint properties for EC section                      
176   for(index = 0; index < fRecPoints->GetEntries(); index++) {
177       dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index))->EvalAll(fECAW0,fDigitsArr) ;
178   }
179
180   fRecPoints->Sort() ;
181
182   for(index = 0; index < fRecPoints->GetEntries(); index++) {
183     (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index)))->SetIndexInList(index) ;
184     (dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(index)))->Print();
185   }
186
187   fTreeR->Fill();
188   
189   if(strstr(option,"deb") || strstr(option,"all"))  
190     PrintRecPoints(option) ;
191
192   AliDebug(1,Form("EMCAL Clusterizer found %d Rec Points",fRecPoints->GetEntriesFast()));
193
194   if(strstr(option,"tim")){
195     gBenchmark->Stop("EMCALClusterizer");
196     printf("Exec took %f seconds for Clusterizing", 
197            gBenchmark->GetCpuTime("EMCALClusterizer"));
198   }    
199 }
200
201 //____________________________________________________________________________
202 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
203                                     Int_t nPar, Float_t * fitparameters) const
204
205   // Calls TMinuit to fit the energy distribution of a cluster with several maxima 
206   // The initial values for fitting procedure are set equal to the positions of local maxima.
207   // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
208
209   gMinuit->mncler();                     // Reset Minuit's list of paramters
210   gMinuit->SetPrintLevel(-1) ;           // No Printout
211   gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;  
212                                          // To set the address of the minimization function 
213   TList * toMinuit = new TList();
214   toMinuit->AddAt(emcRP,0) ;
215   toMinuit->AddAt(fDigitsArr,1) ;
216   
217   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare
218
219   // filling initial values for fit parameters
220   AliEMCALDigit * digit ;
221
222   Int_t ierflg  = 0; 
223   Int_t index   = 0 ;
224   Int_t nDigits = (Int_t) nPar / 3 ;
225
226   Int_t iDigit ;
227
228   for(iDigit = 0; iDigit < nDigits; iDigit++){
229     digit = maxAt[iDigit]; 
230
231     Float_t x = 0.;
232     Float_t z = 0.;
233     //   have to be tune for TRD1; May 31,06
234     //   Int_t relid[2] ;
235     //   fGeom->AbsToRelNumbering(digit->GetId(), relid) ; // obsolete method
236     //   fGeom->PosInAlice(relid, x, z) ;                  // obsolete method
237
238     Float_t energy = maxAtEnergy[iDigit] ;
239
240     gMinuit->mnparm(index, "x",  x, 0.1, 0, 0, ierflg) ;
241     index++ ;   
242     if(ierflg != 0){ 
243       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : x = %f",  x ) ;
244       return kFALSE;
245     }
246     gMinuit->mnparm(index, "z",  z, 0.1, 0, 0, ierflg) ;
247     index++ ;   
248     if(ierflg != 0){
249        Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : z = %f", z) ;
250       return kFALSE;
251     }
252     gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;
253     index++ ;   
254     if(ierflg != 0){
255      Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : energy = %f", energy) ;      
256       return kFALSE;
257     }
258   }
259
260   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
261                       //  depends on it. 
262   Double_t p1 = 1.0 ;
263   Double_t p2 = 0.0 ;
264
265   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
266   gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
267   gMinuit->SetMaxIterations(5);
268   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
269
270   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
271
272   if(ierflg == 4){  // Minimum not found   
273     Error("FindFit", "EMCAL Unfolding  Fit not converged, cluster abandoned " ) ;      
274     return kFALSE ;
275   }            
276   for(index = 0; index < nPar; index++){
277     Double_t err ;
278     Double_t val ;
279     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index
280     fitparameters[index] = val ;
281    }
282
283   delete toMinuit ;
284   return kTRUE;
285
286 }
287
288 //____________________________________________________________________________
289 void AliEMCALClusterizerv1::GetCalibrationParameters() 
290 {
291   // Set calibration parameters:
292   // if calibration database exists, they are read from database,
293   // otherwise, they are taken from digitizer.
294   //
295   // It is a user responsilibity to open CDB before reconstruction, 
296   // for example: 
297   // AliCDBStorage* storage = AliCDBManager::Instance()->GetStorage("local://CalibDB");
298
299   //Check if calibration is stored in data base
300
301   if(!fCalibData && (AliCDBManager::Instance()->IsDefaultStorageSet()))
302     {
303       AliCDBEntry *entry = (AliCDBEntry*) 
304         AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
305       if (entry) fCalibData =  (AliEMCALCalibData*) entry->GetObject();
306     }
307   
308   if(!fCalibData)
309     AliFatal("Calibration parameters not found in CDB!");
310  
311   //  Please fix it!! Or better just remove it...
312 //    if(!fCalibData)
313 //      {
314 //        //If calibration is not available use default parameters
315 //        //Loader
316 //        if ( !emcalLoader->Digitizer() ) 
317 //       emcalLoader->LoadDigitizer();
318 //        AliEMCALDigitizer * dig = dynamic_cast<AliEMCALDigitizer*>(emcalLoader->Digitizer());
319        
320 //        fADCchannelECA   = dig->GetECAchannel() ;
321 //        fADCpedestalECA  = dig->GetECApedestal();
322 //      }
323 }
324
325 //____________________________________________________________________________
326 void AliEMCALClusterizerv1::Init()
327 {
328   // Make all memory allocations which can not be done in default constructor.
329   // Attach the Clusterizer task to the list of EMCAL tasks
330   
331   AliRunLoader *rl = AliRunLoader::GetRunLoader();
332   if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL"))
333     fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
334   else 
335     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaulGeometryName());
336
337   AliDebug(1,Form("geom 0x%x",fGeom));
338
339   if(!gMinuit) 
340     gMinuit = new TMinuit(100) ;
341
342   fHists = BookHists();
343 }
344
345 //____________________________________________________________________________
346 void AliEMCALClusterizerv1::InitParameters()
347
348   // Initializes the parameters for the Clusterizer
349   fNumberOfECAClusters = 0;
350
351   fECALocMaxCut = 0.03; // ??
352
353   fTimeCut = 300e-9 ; // 300 ns time cut (to be tuned) 
354   fToUnfold = kFALSE ;
355
356   fCalibData               = 0 ;
357
358   const AliEMCALRecParam* recParam = AliEMCALReconstructor::GetRecParam();
359   if(!recParam) {
360     AliFatal("Reconstruction parameters for EMCAL not set!");
361   }
362   else {
363     fECAClusteringThreshold = recParam->GetClusteringThreshold();
364     fECAW0                  = recParam->GetW0();
365     fMinECut                = recParam->GetMinECut();
366     AliDebug(1,Form("Reconstruction parameters: fECAClusteringThreshold=%.3f, fECAW=%.3f, fMinECut=%.3f",
367                  fECAClusteringThreshold,fECAW0,fMinECut));
368   }
369
370 }
371
372 //____________________________________________________________________________
373 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2) const
374 {
375   // Gives the neighbourness of two digits = 0 are not neighbour ; continue searching 
376   //                                       = 1 are neighbour
377   //                                       = 2 is in different SM; continue searching 
378   // neighbours are defined as digits having at least a common vertex 
379   // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster 
380   //                                      which is compared to a digit (d2)  not yet in a cluster  
381
382   static Int_t rv; 
383   static Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0, iphi1=0, ieta1=0;
384   static Int_t nSupMod2=0, nModule2=0, nIphi2=0, nIeta2=0, iphi2=0, ieta2=0;
385   static Int_t rowdiff, coldiff;
386   rv = 0 ; 
387
388   fGeom->GetCellIndex(d1->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
389   fGeom->GetCellIndex(d2->GetId(), nSupMod2,nModule2,nIphi2,nIeta2);
390   if(nSupMod1 != nSupMod2) return 2; // different SM
391
392   fGeom->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, iphi1,ieta1);
393   fGeom->GetCellPhiEtaIndexInSModule(nSupMod2,nModule2,nIphi2,nIeta2, iphi2,ieta2);
394
395   rowdiff = TMath::Abs(iphi1 - iphi2);  
396   coldiff = TMath::Abs(ieta1 - ieta2) ;  
397   
398   // neighbours with at least commom side; May 11, 2007
399   if ((coldiff==0 && abs(rowdiff)==1) || (rowdiff==0 && abs(coldiff)==1)) rv = 1;  
400  
401   if (gDebug == 2 && rv==1) 
402   printf("AreNeighbours: neighbours=%d, id1=%d, relid1=%d,%d \n id2=%d, relid2=%d,%d \n", 
403          rv, d1->GetId(), iphi1,ieta1, d2->GetId(), iphi2,ieta2);   
404   
405   return rv ; 
406 }
407
408 //____________________________________________________________________________
409 void AliEMCALClusterizerv1::MakeClusters()
410 {
411   // Steering method to construct the clusters stored in a list of Reconstructed Points
412   // A cluster is defined as a list of neighbour digits
413   // Mar 03, 2007 by PAI
414
415   if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");
416
417   fRecPoints->Clear();
418
419   // Set up TObjArray with pointers to digits to work on 
420   TObjArray *digitsC = new TObjArray();
421   TIter nextdigit(fDigitsArr);
422   AliEMCALDigit *digit;
423   while ( (digit = dynamic_cast<AliEMCALDigit*>(nextdigit())) ) {
424     digitsC->AddLast(digit);
425   }
426
427   double e = 0.0, ehs = 0.0;
428   TIter nextdigitC(digitsC);
429
430   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // clean up digits
431     e = Calibrate(digit->GetAmp(), digit->GetId());
432     AliEMCALHistoUtilities::FillH1(fHists, 10, digit->GetAmp());
433     AliEMCALHistoUtilities::FillH1(fHists, 11, e);
434     if ( e < fMinECut || digit->GetTimeR() > fTimeCut ) 
435       digitsC->Remove(digit);
436     else    
437       ehs += e;
438   } 
439   AliDebug(1,Form("MakeClusters: Number of digits %d  -> (e %f), ehs %d\n",
440                   fDigitsArr->GetEntries(),fMinECut,ehs));
441
442   nextdigitC.Reset();
443
444   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigitC())) ) { // scan over the list of digitsC
445     TArrayI clusterECAdigitslist(fDigitsArr->GetEntries());
446
447     if(fGeom->CheckAbsCellId(digit->GetId()) && (Calibrate(digit->GetAmp(), digit->GetId()) > fECAClusteringThreshold  ) ){
448       // start a new Tower RecPoint
449       if(fNumberOfECAClusters >= fRecPoints->GetSize()) fRecPoints->Expand(2*fNumberOfECAClusters+1) ;
450       AliEMCALRecPoint *recPoint = new  AliEMCALRecPoint("") ; 
451       fRecPoints->AddAt(recPoint, fNumberOfECAClusters) ;
452       recPoint = dynamic_cast<AliEMCALRecPoint *>(fRecPoints->At(fNumberOfECAClusters)) ; 
453       fNumberOfECAClusters++ ; 
454
455       recPoint->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
456
457       recPoint->AddDigit(*digit, Calibrate(digit->GetAmp(), digit->GetId())) ; 
458       TObjArray clusterDigits;
459       clusterDigits.AddLast(digit);     
460       digitsC->Remove(digit) ; 
461
462       AliDebug(1,Form("MakeClusters: OK id = %d, ene = %f , cell.th. = %f \n", digit->GetId(),
463       Calibrate(digit->GetAmp(),digit->GetId()), fECAClusteringThreshold));  
464       
465       // Grow cluster by finding neighbours
466       TIter nextClusterDigit(&clusterDigits);
467       while ( (digit = dynamic_cast<AliEMCALDigit*>(nextClusterDigit())) ) { // scan over digits in cluster 
468         TIter nextdigitN(digitsC); 
469         AliEMCALDigit *digitN = 0; // digi neighbor
470         while ( (digitN = (AliEMCALDigit *)nextdigitN()) ) { // scan over all digits to look for neighbours
471           if (AreNeighbours(digit, digitN)==1) {      // call (digit,digitN) in THAT oder !!!!! 
472             recPoint->AddDigit(*digitN, Calibrate(digitN->GetAmp(),digitN->GetId()) ) ;
473             clusterDigits.AddLast(digitN) ; 
474             digitsC->Remove(digitN) ; 
475           } // if(ineb==1)
476         } // scan over digits
477       } // scan over digits already in cluster
478       if(recPoint)
479         AliDebug(2,Form("MakeClusters: %d digitd, energy %f \n", clusterDigits.GetEntries(), recPoint->GetEnergy())); 
480     } // If seed found
481   } // while digit 
482
483   delete digitsC ;
484
485   AliDebug(1,Form("total no of clusters %d from %d digits",fNumberOfECAClusters,fDigitsArr->GetEntriesFast())); 
486 }
487
488 void AliEMCALClusterizerv1::MakeUnfolding() const
489 {
490   Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
491 }
492
493 //____________________________________________________________________________
494 Double_t  AliEMCALClusterizerv1::ShowerShape(Double_t r)
495
496   // Shape of the shower (see EMCAL TDR)
497   // If you change this function, change also the gradient evaluation in ChiSquare()
498
499   Double_t r4    = r*r*r*r ;
500   Double_t r295  = TMath::Power(r, 2.95) ;
501   Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
502   return shape ;
503 }
504
505 //____________________________________________________________________________
506 void  AliEMCALClusterizerv1::UnfoldCluster(AliEMCALRecPoint * /*iniTower*/, 
507                                            Int_t /*nMax*/, 
508                                            AliEMCALDigit ** /*maxAt*/, 
509                                            Float_t * /*maxAtEnergy*/) const
510 {
511   // Performs the unfolding of a cluster with nMax overlapping showers 
512   
513   Fatal("UnfoldCluster", "--> Unfolding not implemented") ;
514
515 }
516
517 //_____________________________________________________________________________
518 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * /*Grad*/,
519                                                Double_t & /*fret*/,
520                                                Double_t * /*x*/, Int_t /*iflag*/)
521 {
522   // Calculates the Chi square for the cluster unfolding minimization
523   // Number of parameters, Gradient, Chi squared, parameters, what to do
524   
525   ::Fatal("UnfoldingChiSquare","Unfolding not implemented") ;
526 }
527 //____________________________________________________________________________
528 void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
529 {
530   // Print clusterizer parameters
531
532   TString message("\n") ; 
533   
534   if( strcmp(GetName(), "") !=0 ){
535     
536     // Print parameters
537  
538     TString taskName(Version()) ;
539     
540     printf("--------------- "); 
541     printf(taskName.Data()) ; 
542     printf(" "); 
543     printf("Clusterizing digits: "); 
544     printf("\n                       ECA Local Maximum cut    = %f", fECALocMaxCut); 
545     printf("\n                       ECA Logarithmic weight   = %f", fECAW0); 
546     if(fToUnfold)
547       printf("\nUnfolding on\n");
548     else
549       printf("\nUnfolding off\n");
550     
551     printf("------------------------------------------------------------------"); 
552   }
553   else
554     printf("AliEMCALClusterizerv1 not initialized ") ;
555 }
556
557 //____________________________________________________________________________
558 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
559 {
560   // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
561   if(strstr(option,"deb")) {
562     printf("PrintRecPoints: Clusterization result:") ; 
563   
564     printf("           Found %d ECA Rec Points\n ", 
565          fRecPoints->GetEntriesFast()) ; 
566   }
567
568   if(strstr(option,"all")) {
569     if(strstr(option,"deb")) {
570       printf("\n-----------------------------------------------------------------------\n") ;
571       printf("Clusters in ECAL section\n") ;
572       printf("Index    Ene(GeV) Multi Module     GX    GY   GZ  lX    lY   lZ   Dispersion Lambda 1   Lambda 2  # of prim  Primaries list\n") ;
573     }
574    Int_t index =0;
575    Float_t maxE=0; 
576    Float_t maxL1=0; 
577    Float_t maxL2=0; 
578    Float_t maxDis=0; 
579
580     AliEMCALHistoUtilities::FillH1(fHists, 12, double(fRecPoints->GetEntries()));
581
582     for (index = 0 ; index < fRecPoints->GetEntries() ; index++) {
583       AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint * >(fRecPoints->At(index)) ; 
584       TVector3  globalpos;  
585       //rp->GetGlobalPosition(globalpos);
586       TVector3  localpos;  
587       rp->GetLocalPosition(localpos);
588       Float_t lambda[2]; 
589       rp->GetElipsAxis(lambda);
590       Int_t * primaries; 
591       Int_t nprimaries;
592       primaries = rp->GetPrimaries(nprimaries);
593       if(strstr(option,"deb")) 
594       printf("\n%6d  %8.4f  %3d     %4.1f    %4.1f %4.1f  %4.1f %4.1f %4.1f    %4.1f   %4f  %4f    %2d     : ", 
595              rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(),
596              globalpos.X(), globalpos.Y(), globalpos.Z(), localpos.X(), localpos.Y(), localpos.Z(), 
597              rp->GetDispersion(), lambda[0], lambda[1], nprimaries) ; 
598   /////////////
599       if(rp->GetEnergy()>maxE){
600               maxE=rp->GetEnergy();
601               maxL1=lambda[0];
602               maxL2=lambda[1];
603               maxDis=rp->GetDispersion();
604       }
605       fPointE->Fill(rp->GetEnergy());
606       fPointL1->Fill(lambda[0]);
607       fPointL2->Fill(lambda[1]);
608       fPointDis->Fill(rp->GetDispersion());
609       fPointMult->Fill(rp->GetMultiplicity());
610       ///////////// 
611       if(strstr(option,"deb")){ 
612         for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
613           printf("%d ", primaries[iprimary] ) ; 
614         }
615       }
616     }
617
618       fMaxE->Fill(maxE);
619       fMaxL1->Fill(maxL1);
620       fMaxL2->Fill(maxL2);
621       fMaxDis->Fill(maxDis);
622
623     if(strstr(option,"deb"))
624     printf("\n-----------------------------------------------------------------------\n");
625   }
626 }
627 TList* AliEMCALClusterizerv1::BookHists()
628 {
629   //set up histograms for monitoring clusterizer performance
630
631   gROOT->cd();
632
633         fPointE = new TH1F("00_pointE","point energy", 2000, 0.0, 150.);
634         fPointL1 = new TH1F("01_pointL1","point L1", 1000, 0.0, 3.);
635         fPointL2 = new TH1F("02_pointL2","point L2", 1000, 0.0, 3.);
636         fPointDis = new TH1F("03_pointDisp","point dispersion", 1000, 0.0, 10.);
637         fPointMult = new TH1F("04_pointMult","#cell in point(cluster)", 101, -0.5, 100.5);
638         fDigitAmp = new TH1F("05_digitAmp","Digit Amplitude", 2000, 0.0, 5000.);
639         fMaxE = new TH1F("06_maxE","Max point energy", 2000, 0.0, 150.);
640         fMaxL1 = new TH1F("07_maxL1","Largest (first) of eigenvalue of covariance matrix", 1000, 0.0, 3.);
641         fMaxL2 = new TH1F("08_maxL2","Smalest (second) of eigenvalue of covariace matrix", 1000, 0.0, 3.);
642         fMaxDis = new TH1F("09_maxDis","Point dispersion", 1000, 0.0, 10.); // 9
643         //
644         new TH1F("10_adcOfDigits","adc of digits(threshold control)", 1001, -0.5, 1000.5);   // 10
645         new TH1F("11_energyOfDigits","energy of digits(threshold control)", 1000, 0.0, 1.);  // 11
646         new TH1F("12_numberOfPoints","number of points(clusters)", 101, -0.5, 100.5);        // 12
647
648   return AliEMCALHistoUtilities::MoveHistsToList("EmcalClusterizerv1ControlHists", kFALSE);
649 }
650
651 void AliEMCALClusterizerv1::SaveHists(const char *fn)
652 {
653   AliEMCALHistoUtilities::SaveListOfHists(fHists, fn, kTRUE);
654 }
655
656 void  AliEMCALClusterizerv1::PrintRecoInfo()
657 {
658   printf(" AliEMCALClusterizerv1::PrintRecoInfo() : version %s \n", Version() );
659   TH1F *h = (TH1F*)fHists->At(12);
660   if(h) {
661     printf(" ## Multiplicity of RecPoints ## \n");
662     for(int i=1; i<=h->GetNbinsX(); i++) {
663       int nbin = int((*h)[i]);
664       int mult = int(h->GetBinCenter(i));
665       if(nbin > 0) printf(" %i : %5.5i %6.3f %% \n", mult, nbin, 100.*nbin/h->GetEntries()); 
666     }    
667   }
668 }
669
670 void AliEMCALClusterizerv1::DrawLambdasHists()
671 {
672   if(fMaxL1) {
673     fMaxL1->Draw();
674     if(fMaxL2) fMaxL2->Draw("same");
675     if(fMaxDis) {
676       fMaxDis->Draw("same");
677     }
678   }
679 }