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