]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALClusterizerv1.cxx
1. Style modifications to clone PHOS
[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 /* $Id$ */
16
17 /* $Log:
18    1 October 2000. Yuri Kharlov:
19      AreNeighbours()
20      PPSD upper layer is considered if number of layers>1
21
22    18 October 2000. Yuri Kharlov:
23      AliEMCALClusterizerv1()
24      CPV clusterizing parameters added
25
26      MakeClusters()
27      After first PPSD digit remove EMC digits only once
28 */
29 //*-- Author: Yves Schutz (SUBATECH)  & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
30 //  August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
31 //                           of new  IO (à la PHOS)
32 //////////////////////////////////////////////////////////////////////////////
33 //  Clusterization class. Performs clusterization (collects neighbouring active cells) and 
34 //  unfolds the clusters having several local maxima.  
35 //  Results are stored in TreeR#, branches EMCALTowerRP (EMC recPoints),
36 //  EMCALPreShoRP (CPV RecPoints) and AliEMCALClusterizer (Clusterizer with all 
37 //  parameters including input digits branch title, thresholds etc.)
38 //  This TTask is normally called from Reconstructioner, but can as well be used in 
39 //  standalone mode.
40 // Use Case:
41 //  root [0] AliEMCALClusterizerv1 * cl = new AliEMCALClusterizerv1("galice.root")  
42 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
43 //               //reads gAlice from header file "..."                      
44 //  root [1] cl->ExecuteTask()  
45 //               //finds RecPoints in all events stored in galice.root
46 //  root [2] cl->SetDigitsBranch("digits2") 
47 //               //sets another title for Digitis (input) branch
48 //  root [3] cl->SetRecPointsBranch("recp2")  
49 //               //sets another title four output branches
50 //  root [4] cl->SetTowerLocalMaxCut(0.03)  
51 //               //set clusterization parameters
52 //  root [5] cl->ExecuteTask("deb all time")  
53 //               //once more finds RecPoints options are 
54 //               // deb - print number of found rec points
55 //               // deb all - print number of found RecPoints and some their characteristics 
56 //               // time - print benchmarking results
57
58 // --- ROOT system ---
59
60 #include "TROOT.h" 
61 #include "TFile.h" 
62 #include "TFolder.h" 
63 #include "TMath.h" 
64 #include "TMinuit.h"
65 #include "TTree.h" 
66 #include "TSystem.h" 
67 #include "TBenchmark.h"
68
69 // --- Standard library ---
70
71 #include <iostream.h>
72 #include <iomanip.h>
73
74 // --- AliRoot header files ---
75
76 #include "AliEMCALClusterizerv1.h"
77 #include "AliEMCALDigit.h"
78 #include "AliEMCALDigitizer.h"
79 #include "AliEMCALTowerRecPoint.h"
80 #include "AliEMCAL.h"
81 #include "AliEMCALGetter.h"
82 #include "AliEMCALGeometry.h"
83 #include "AliRun.h"
84
85 ClassImp(AliEMCALClusterizerv1)
86   
87 //____________________________________________________________________________
88   AliEMCALClusterizerv1::AliEMCALClusterizerv1() : AliEMCALClusterizer()
89 {
90   // default ctor (to be used mainly by Streamer)
91   
92   InitParameters() ; 
93   fDefaultInit = kTRUE ; 
94 }
95
96 //____________________________________________________________________________
97 AliEMCALClusterizerv1::AliEMCALClusterizerv1(const char* headerFile, const char* name, const Bool_t toSplit)
98 :AliEMCALClusterizer(headerFile, name, toSplit)
99 {
100   // ctor with the indication of the file where header Tree and digits Tree are stored
101   
102   InitParameters() ; 
103   Init() ;
104   fDefaultInit = kFALSE ; 
105
106 }
107
108 //____________________________________________________________________________
109   AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
110 {
111   // dtor
112   fSplitFile = 0 ; 
113   
114 }
115
116 //____________________________________________________________________________
117 const TString AliEMCALClusterizerv1::BranchName() const 
118 {  
119   TString branchName(GetName() ) ;
120   branchName.Remove(branchName.Index(Version())-1) ;
121   return branchName ;
122 }
123
124 //____________________________________________________________________________
125 Float_t  AliEMCALClusterizerv1::Calibrate(Int_t amp, Bool_t inpresho) const
126 {//To be replased later by the method, reading individual parameters from the database
127
128   if ( inpresho ) // calibrate as pre shower
129      return -fADCpedestalPreSho + amp * fADCchannelPreSho ; 
130
131   else //calibrate as tower 
132     return -fADCpedestalTower + amp * fADCchannelTower ;                
133 }
134
135 //____________________________________________________________________________
136 void AliEMCALClusterizerv1::Exec(Option_t * option)
137 {
138   // Steering method
139
140   if( strcmp(GetName(), "")== 0 ) 
141     Init() ;
142
143   if(strstr(option,"tim"))
144     gBenchmark->Start("EMCALClusterizer"); 
145   
146   if(strstr(option,"print"))
147     Print("") ; 
148
149   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
150   if(gime->BranchExists("RecPoints"))
151     return ;
152   Int_t nevents = gime->MaxEvent() ;
153   Int_t ievent ;
154
155   for(ievent = 0; ievent < nevents; ievent++){
156
157     gime->Event(ievent,"D") ;
158
159     if(ievent == 0)
160       GetCalibrationParameters() ;
161
162     fNumberOfTowerClusters = fNumberOfPreShoClusters = 0 ;
163            
164     MakeClusters() ;
165     
166     if(fToUnfold)
167       MakeUnfolding() ;
168
169     WriteRecPoints(ievent) ;
170
171     if(strstr(option,"deb"))  
172       PrintRecPoints(option) ;
173
174     //increment the total number of digits per run 
175     fRecPointsInRun += gime->TowerRecPoints()->GetEntriesFast() ;  
176     fRecPointsInRun += gime->PreShowerRecPoints()->GetEntriesFast() ;  
177  }
178   
179   if(strstr(option,"tim")){
180     gBenchmark->Stop("EMCALClusterizer");
181     cout << "AliEMCALClusterizer:" << endl ;
182     cout << "  took " << gBenchmark->GetCpuTime("EMCALClusterizer") << " seconds for Clusterizing " 
183          <<  gBenchmark->GetCpuTime("EMCALClusterizer")/nevents << " seconds per event " << endl ;
184     cout << endl ;
185   }
186   
187 }
188
189 //____________________________________________________________________________
190 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALTowerRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
191                                     Int_t nPar, Float_t * fitparameters) const
192
193   // Calls TMinuit to fit the energy distribution of a cluster with several maxima 
194   // The initial values for fitting procedure are set equal to the positions of local maxima.
195   // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
196
197   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
198   TClonesArray * digits = gime->Digits() ; 
199   
200
201   gMinuit->mncler();                     // Reset Minuit's list of paramters
202   gMinuit->SetPrintLevel(-1) ;           // No Printout
203   gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;  
204                                          // To set the address of the minimization function 
205
206   TList * toMinuit = new TList();
207   toMinuit->AddAt(emcRP,0) ;
208   toMinuit->AddAt(digits,1) ;
209   
210   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare
211
212   // filling initial values for fit parameters
213   AliEMCALDigit * digit ;
214
215   Int_t ierflg  = 0; 
216   Int_t index   = 0 ;
217   Int_t nDigits = (Int_t) nPar / 3 ;
218
219   Int_t iDigit ;
220
221   AliEMCALGeometry * geom = gime->EMCALGeometry() ; 
222
223   for(iDigit = 0; iDigit < nDigits; iDigit++){
224     digit = maxAt[iDigit]; 
225
226     Int_t relid[4] ;
227     Float_t x = 0.;
228     Float_t z = 0.;
229     geom->AbsToRelNumbering(digit->GetId(), relid) ;
230     geom->PosInAlice(relid, x, z) ;
231
232     Float_t energy = maxAtEnergy[iDigit] ;
233
234     gMinuit->mnparm(index, "x",  x, 0.1, 0, 0, ierflg) ;
235     index++ ;   
236     if(ierflg != 0){ 
237       cout << "EMCAL Unfolding>  Unable to set initial value for fit procedure : x = " << x << endl ;
238       return kFALSE;
239     }
240     gMinuit->mnparm(index, "z",  z, 0.1, 0, 0, ierflg) ;
241     index++ ;   
242     if(ierflg != 0){
243       cout << "EMCAL Unfolding>  Unable to set initial value for fit procedure : z = " << z << endl ;
244       return kFALSE;
245     }
246     gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;
247     index++ ;   
248     if(ierflg != 0){
249       cout << "EMCAL Unfolding>  Unable to set initial value for fit procedure : energy = " << energy << endl ;      
250       return kFALSE;
251     }
252   }
253
254   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
255                       //  depends on it. 
256   Double_t p1 = 1.0 ;
257   Double_t p2 = 0.0 ;
258
259   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
260   gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
261   gMinuit->SetMaxIterations(5);
262   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
263
264   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
265
266   if(ierflg == 4){  // Minimum not found   
267     cout << "EMCAL Unfolding>  Fit not converged, cluster abandoned "<< endl ;      
268     return kFALSE ;
269   }            
270   for(index = 0; index < nPar; index++){
271     Double_t err ;
272     Double_t val ;
273     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index
274     fitparameters[index] = val ;
275    }
276
277   delete toMinuit ;
278   return kTRUE;
279
280 }
281
282 //____________________________________________________________________________
283 void AliEMCALClusterizerv1::GetCalibrationParameters() 
284 {
285   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
286   const AliEMCALDigitizer * dig = gime->Digitizer(BranchName()) ;
287
288   fADCchannelTower   = dig->GetTowerchannel() ;
289   fADCpedestalTower  = dig->GetTowerpedestal();
290
291   fADCchannelPreSho  = dig->GetPreShochannel() ;
292   fADCpedestalPreSho = dig->GetPreShopedestal() ; 
293
294 }
295
296 //____________________________________________________________________________
297 void AliEMCALClusterizerv1::Init()
298 {
299   // Make all memory allocations which can not be done in default constructor.
300   // Attach the Clusterizer task to the list of EMCAL tasks
301   
302   if ( strcmp(GetTitle(), "") == 0 )
303     SetTitle("galice.root") ;
304
305   TString branchname = GetName() ;
306   branchname.Remove(branchname.Index(Version())-1) ;
307
308   AliEMCALGetter * gime = AliEMCALGetter::GetInstance(GetTitle(), branchname.Data(), fToSplit ) ; 
309   if ( gime == 0 ) {
310     cerr << "ERROR: AliEMCALClusterizerv1::Init -> Could not obtain the Getter object !" << endl ; 
311     return ;
312   } 
313
314   fSplitFile = 0 ;
315   if(fToSplit){
316     // construct the name of the file as /path/EMCAL.SDigits.root
317     //First - extract full path if necessary
318     TString fileName(GetTitle()) ;
319     Ssiz_t islash = fileName.Last('/') ;
320     if(islash<fileName.Length())
321       fileName.Remove(islash+1,fileName.Length()) ;
322     else
323       fileName="" ;
324     // Next - append the file name 
325     fileName+="EMCAL.RecData." ;
326     if((strcmp(branchname.Data(),"Default")!=0)&&(strcmp(branchname.Data(),"")!=0)){
327       fileName+=branchname ;
328       fileName+="." ;
329     }
330     fileName+="root" ;
331     // Finally - check if the file already opened or open the file
332     fSplitFile = static_cast<TFile*>(gROOT->GetFile(fileName.Data()));   
333     if(!fSplitFile)
334       fSplitFile =  TFile::Open(fileName.Data(),"update") ;
335   }
336
337     
338   const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
339   fNTowers = geom->GetNZ() *  geom->GetNPhi() ;
340
341   if(!gMinuit) 
342     gMinuit = new TMinuit(100) ;
343
344   gime->PostClusterizer(this) ;
345   gime->PostRecPoints(branchname ) ;
346  
347 }
348
349 //____________________________________________________________________________
350 void AliEMCALClusterizerv1::InitParameters()
351 {
352   fNumberOfPreShoClusters = fNumberOfTowerClusters = 0 ; 
353
354  
355   
356   fPreShoClusteringThreshold  = 0.0001;
357   fTowerClusteringThreshold   = 0.2;   
358   
359   fTowerLocMaxCut  = 0.03 ;
360   fPreShoLocMaxCut = 0.03 ;
361   
362   fW0     = 4.5 ;
363   fW0CPV  = 4.0 ;
364
365   fTimeGate = 1.e-8 ; 
366   
367   fToUnfold = kFALSE ;
368    
369   TString clusterizerName( GetName()) ; 
370   if (clusterizerName.IsNull() ) 
371     clusterizerName = "Default" ; 
372   clusterizerName.Append(":") ; 
373   clusterizerName.Append(Version()) ; 
374   SetName(clusterizerName) ;
375   fRecPointsInRun          = 0 ; 
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 but continue searching 
383   //                                       = 1 are neighbour
384   //                                       = 2 are not neighbour but do not 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    AliEMCALGeometry * geom = AliEMCALGetter::GetInstance()->EMCALGeometry() ;
390
391   Int_t rv = 0 ; 
392
393   Int_t relid1[4] ; 
394   geom->AbsToRelNumbering(d1->GetId(), relid1) ; 
395
396   Int_t relid2[4] ; 
397   geom->AbsToRelNumbering(d2->GetId(), relid2) ; 
398  
399   if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same EMCAL Arm 
400     Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
401     Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
402     
403     if (( coldiff <= 1 )  && ( rowdiff <= 1 )){
404       if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fTimeGate))
405       rv = 1 ; 
406     }
407     else {
408       if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1)) 
409         rv = 2; //  Difference in row numbers is too large to look further 
410     }
411
412   } 
413   else {
414     
415     if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )  
416       rv=2 ;
417
418   }
419
420   return rv ; 
421 }
422
423
424 //____________________________________________________________________________
425 Bool_t AliEMCALClusterizerv1::IsInTower(AliEMCALDigit * digit) const
426 {
427   // Tells if (true) or not (false) the digit is in a EMCAL-Tower 
428  
429   Bool_t rv = kFALSE ; 
430   if (!digit->IsInPreShower()) 
431     rv = kTRUE; 
432   return rv ; 
433 }
434
435 //____________________________________________________________________________
436 Bool_t AliEMCALClusterizerv1::IsInPreShower(AliEMCALDigit * digit) const
437 {
438   // Tells if (true) or not (false) the digit is in a EMCAL-PreShower
439  
440   Bool_t rv = kFALSE ; 
441   if (digit->IsInPreShower()) 
442     rv = kTRUE; 
443   return rv ; 
444 }
445
446 //____________________________________________________________________________
447 void AliEMCALClusterizerv1::WriteRecPoints(Int_t event)
448 {
449
450   // Creates new branches with given title
451   // fills and writes into TreeR.
452
453   AliEMCALGetter *gime = AliEMCALGetter::GetInstance() ; 
454   TObjArray * towerRecPoints = gime->TowerRecPoints() ; 
455   TObjArray * preshoRecPoints = gime->PreShowerRecPoints() ; 
456   TClonesArray * digits = gime->Digits() ; 
457   TTree * treeR ; 
458   
459   if(fToSplit){
460     if(!fSplitFile)
461       return ;
462     fSplitFile->cd() ;
463     TString name("TreeR") ;
464     name += event ; 
465     treeR = dynamic_cast<TTree*>(fSplitFile->Get(name)); 
466   }
467   else{
468     treeR = gAlice->TreeR();
469   }
470
471   if(!treeR){
472     gAlice->MakeTree("R", fSplitFile);
473     treeR = gAlice->TreeR() ;
474   }
475  
476   Int_t index ;
477   //Evaluate position, dispersion and other RecPoint properties...
478   for(index = 0; index < towerRecPoints->GetEntries(); index++)
479     (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->EvalAll(fW0,digits) ;
480   
481   towerRecPoints->Sort() ;
482
483   for(index = 0; index < towerRecPoints->GetEntries(); index++)
484     (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->SetIndexInList(index) ;
485
486   towerRecPoints->Expand(towerRecPoints->GetEntriesFast()) ; 
487
488   //Now the same for pre shower
489   for(index = 0; index < preshoRecPoints->GetEntries(); index++)
490     (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->EvalAll(fW0CPV,digits)  ;
491
492   preshoRecPoints->Sort() ;
493
494   for(index = 0; index < preshoRecPoints->GetEntries(); index++)
495     (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->SetIndexInList(index) ;
496
497   preshoRecPoints->Expand(preshoRecPoints->GetEntriesFast()) ;
498   
499   Int_t bufferSize = 32000 ;    
500   Int_t splitlevel = 0 ;
501
502   //First Tower branch
503   TBranch * towerBranch = treeR->Branch("EMCALTowerRP","TObjArray",&towerRecPoints,bufferSize,splitlevel);
504   towerBranch->SetTitle(BranchName());
505   
506   //Now Pre Shower branch 
507   TBranch * preshoBranch = treeR->Branch("EMCALPreShoRP","TObjArray",&preshoRecPoints,bufferSize,splitlevel);
508   preshoBranch->SetTitle(BranchName());
509     
510   //And Finally  clusterizer branch
511   AliEMCALClusterizerv1 * cl = (AliEMCALClusterizerv1*)gime->Clusterizer(BranchName()) ;
512   TBranch * clusterizerBranch = treeR->Branch("AliEMCALClusterizer","AliEMCALClusterizerv1",
513                                               &cl,bufferSize,splitlevel);
514   clusterizerBranch->SetTitle(BranchName());
515
516   towerBranch        ->Fill() ;
517   preshoBranch        ->Fill() ;
518   clusterizerBranch->Fill() ;
519
520   treeR->AutoSave() ; //Write(0,kOverwrite) ;  
521   if(gAlice->TreeR()!=treeR)
522     treeR->Delete(); 
523 }
524
525 //____________________________________________________________________________
526 void AliEMCALClusterizerv1::MakeClusters()
527 {
528   // Steering method to construct the clusters stored in a list of Reconstructed Points
529   // A cluster is defined as a list of neighbour digits
530     
531   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
532   
533   TObjArray * towerRecPoints  = gime->TowerRecPoints(BranchName()) ; 
534   TObjArray * preshoRecPoints = gime->PreShowerRecPoints(BranchName()) ; 
535   towerRecPoints->Delete() ;
536   preshoRecPoints->Delete() ;
537   
538   TClonesArray * digits = gime->Digits() ; 
539   if ( !digits ) {
540     cerr << "ERROR:  AliEMCALClusterizerv1::MakeClusters -> Digits with name " 
541          << GetName() << " not found ! " << endl ; 
542     abort() ; 
543   } 
544   TClonesArray * digitsC =  dynamic_cast<TClonesArray*>(digits->Clone()) ;
545   
546   
547   // Clusterization starts  
548   
549   TIter nextdigit(digitsC) ; 
550   AliEMCALDigit * digit ; 
551   Bool_t notremoved = kTRUE ;
552   
553   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
554     AliEMCALRecPoint * clu = 0 ; 
555     
556     TArrayI clusterdigitslist(1500) ;   
557     Int_t index ;
558  
559     if (( IsInTower (digit)  && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fTowerClusteringThreshold  ) || 
560         ( IsInPreShower (digit) && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fPreShoClusteringThreshold  ) ) {
561       
562       Int_t iDigitInCluster = 0 ; 
563       
564       if  ( IsInTower(digit) ) {   
565         // start a new Tower RecPoint
566         if(fNumberOfTowerClusters >= towerRecPoints->GetSize()) 
567           towerRecPoints->Expand(2*fNumberOfTowerClusters+1) ;
568         
569         towerRecPoints->AddAt(new  AliEMCALTowerRecPoint(""), fNumberOfTowerClusters) ;
570         clu = dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(fNumberOfTowerClusters)) ; 
571         fNumberOfTowerClusters++ ; 
572         clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower())) ; 
573         clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;  
574         iDigitInCluster++ ; 
575         digitsC->Remove(digit) ; 
576         
577       } else { 
578         
579         // start a new Pre Shower cluster
580         if(fNumberOfPreShoClusters >= preshoRecPoints->GetSize()) 
581           preshoRecPoints->Expand(2*fNumberOfPreShoClusters+1);
582         
583         preshoRecPoints->AddAt(new AliEMCALTowerRecPoint(""), fNumberOfPreShoClusters) ;
584         
585         clu =  dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(fNumberOfPreShoClusters))  ;  
586         fNumberOfPreShoClusters++ ; 
587         clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower() ) );     
588         clusterdigitslist[iDigitInCluster] = digit->GetIndexInList()  ; 
589         iDigitInCluster++ ; 
590         digitsC->Remove(digit) ; 
591         nextdigit.Reset() ;
592         
593         // Here we remove remaining Tower digits, which cannot make a cluster
594         
595         if( notremoved ) { 
596           while( ( digit = dynamic_cast<AliEMCALDigit *>(nextdigit()) ) ) {
597             if( IsInTower(digit) )
598               digitsC->Remove(digit) ;
599             else 
600               break ; 
601           }
602           notremoved = kFALSE ;
603         }
604         
605       } // else        
606       
607       nextdigit.Reset() ;
608       
609       AliEMCALDigit * digitN ; 
610       index = 0 ;
611       while (index < iDigitInCluster){ // scan over digits already in cluster 
612         digit =  (AliEMCALDigit*)digits->At(clusterdigitslist[index])  ;      
613         index++ ; 
614         while ( (digitN = (AliEMCALDigit *)nextdigit()) ) { // scan over the reduced list of digits 
615           Int_t ineb = AreNeighbours(digit, digitN);       // call (digit,digitN) in THAT oder !!!!!
616          switch (ineb ) {
617           case 0 :   // not a neighbour
618             break ;
619           case 1 :   // are neighbours 
620             clu->AddDigit(*digitN, Calibrate( digitN->GetAmp(), digitN->IsInPreShower() ) ) ;
621             clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ; 
622             iDigitInCluster++ ; 
623             digitsC->Remove(digitN) ;
624             break ;
625           case 2 :   // too far from each other
626             goto endofloop;   
627           } // switch
628           
629         } // while digitN
630         
631       endofloop: ;
632         nextdigit.Reset() ; 
633         
634       } // loop over cluster     
635       
636     } // energy theshold  
637     
638     
639   } // while digit
640   
641   delete digitsC ;
642   
643 }
644
645 //____________________________________________________________________________
646 void AliEMCALClusterizerv1::MakeUnfolding()
647 {
648   Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
649   
650 //   // Unfolds clusters using the shape of an ElectroMagnetic shower
651 //   // Performs unfolding of all EMC/CPV clusters
652
653 //   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
654   
655 //   const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
656 //   TObjArray * emcRecPoints = gime->TowerRecPoints() ; 
657 //   TObjArray * cpvRecPoints = gime->PreShoRecPoints() ; 
658 //   TClonesArray * digits = gime->Digits() ; 
659   
660 //   // Unfold first EMC clusters 
661 //   if(fNumberOfTowerClusters > 0){
662
663 //     Int_t nModulesToUnfold = geom->GetNModules() ; 
664
665 //     Int_t numberofNotUnfolded = fNumberOfTowerClusters ; 
666 //     Int_t index ;   
667 //     for(index = 0 ; index < numberofNotUnfolded ; index++){
668       
669 //       AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint *) emcRecPoints->At(index) ;
670 //       if(emcRecPoint->GetEMCALMod()> nModulesToUnfold)
671 //      break ;
672       
673 //       Int_t nMultipl = emcRecPoint->GetMultiplicity() ; 
674 //       Int_t * maxAt = new Int_t[nMultipl] ;
675 //       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
676 //       Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fTowerLocMaxCut,digits) ;
677       
678 //       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0       
679 //      UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
680 //      emcRecPoints->Remove(emcRecPoint); 
681 //      emcRecPoints->Compress() ;
682 //      index-- ;
683 //      fNumberOfTowerClusters -- ;
684 //      numberofNotUnfolded-- ;
685 //       }
686       
687 //       delete[] maxAt ; 
688 //       delete[] maxAtEnergy ; 
689 //     }
690 //   } 
691 //   // Unfolding of EMC clusters finished
692
693
694 //   // Unfold now CPV clusters
695 //   if(fNumberOfPreShoClusters > 0){
696     
697 //     Int_t nModulesToUnfold = geom->GetNModules() ;
698
699 //     Int_t numberofPreShoNotUnfolded = fNumberOfPreShoClusters ;     
700 //     Int_t index ;   
701 //     for(index = 0 ; index < numberofPreShoNotUnfolded ; index++){
702       
703 //       AliEMCALRecPoint * recPoint = (AliEMCALRecPoint *) cpvRecPoints->At(index) ;
704
705 //       if(recPoint->GetEMCALMod()> nModulesToUnfold)
706 //      break ;
707       
708 //       AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint*) recPoint ; 
709       
710 //       Int_t nMultipl = emcRecPoint->GetMultiplicity() ; 
711 //       Int_t * maxAt = new Int_t[nMultipl] ;
712 //       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
713 //       Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fPreShoLocMaxCut,digits) ;
714       
715 //       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0       
716 //      UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
717 //      cpvRecPoints->Remove(emcRecPoint); 
718 //      cpvRecPoints->Compress() ;
719 //      index-- ;
720 //      numberofPreShoNotUnfolded-- ;
721 //      fNumberOfPreShoClusters-- ;
722 //       }
723       
724 //       delete[] maxAt ; 
725 //       delete[] maxAtEnergy ; 
726 //     } 
727 //   }
728 //   //Unfolding of PreSho clusters finished
729   
730 }
731
732 //____________________________________________________________________________
733 Double_t  AliEMCALClusterizerv1::ShowerShape(Double_t r)
734
735   // Shape of the shower (see EMCAL TDR)
736   // If you change this function, change also the gradient evaluation in ChiSquare()
737
738   Double_t r4    = r*r*r*r ;
739   Double_t r295  = TMath::Power(r, 2.95) ;
740   Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
741   return shape ;
742 }
743
744 //____________________________________________________________________________
745 void  AliEMCALClusterizerv1::UnfoldCluster(AliEMCALTowerRecPoint * iniTower, 
746                                                  Int_t nMax, 
747                                                  AliEMCALDigit ** maxAt, 
748                                                  Float_t * maxAtEnergy)
749 {
750   // Performs the unfolding of a cluster with nMax overlapping showers 
751   
752   Fatal("AliEMCALClusterizerv1::UnfoldCluster", "--> Unfolding not implemented") ;
753
754  //  AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
755 //   const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
756 //   const TClonesArray * digits = gime->Digits() ; 
757 //   TObjArray * emcRecPoints = gime->TowerRecPoints() ; 
758 //   TObjArray * cpvRecPoints = gime->PreShoRecPoints() ; 
759
760 //   Int_t nPar = 3 * nMax ;
761 //   Float_t * fitparameters = new Float_t[nPar] ;
762
763 //   Bool_t rv = FindFit(iniTower, maxAt, maxAtEnergy, nPar, fitparameters) ;
764 //   if( !rv ) {
765 //     // Fit failed, return and remove cluster
766 //     delete[] fitparameters ; 
767 //     return ;
768 //   }
769
770 //   // create ufolded rec points and fill them with new energy lists
771 //   // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
772 //   // and later correct this number in acordance with actual energy deposition
773
774 //   Int_t nDigits = iniTower->GetMultiplicity() ;  
775 //   Float_t * efit = new Float_t[nDigits] ;
776 //   Float_t xDigit=0.,zDigit=0.,distance=0. ;
777 //   Float_t xpar=0.,zpar=0.,epar=0.  ;
778 //   Int_t relid[4] ;
779 //   AliEMCALDigit * digit = 0 ;
780 //   Int_t * emcDigits = iniTower->GetDigitsList() ;
781
782 //   Int_t iparam ;
783 //   Int_t iDigit ;
784 //   for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
785 //     digit = (AliEMCALDigit*) digits->At(emcDigits[iDigit] ) ;   
786 //     geom->AbsToRelNumbering(digit->GetId(), relid) ;
787 //     geom->RelPosInModule(relid, xDigit, zDigit) ;
788 //     efit[iDigit] = 0;
789
790 //     iparam = 0 ;    
791 //     while(iparam < nPar ){
792 //       xpar = fitparameters[iparam] ;
793 //       zpar = fitparameters[iparam+1] ;
794 //       epar = fitparameters[iparam+2] ;
795 //       iparam += 3 ;
796 //       distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
797 //       distance =  TMath::Sqrt(distance) ;
798 //       efit[iDigit] += epar * ShowerShape(distance) ;
799 //     }
800 //   }
801   
802
803 //   // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
804 //   // so that energy deposited in each cell is distributed betwin new clusters proportionally
805 //   // to its contribution to efit
806
807 //   Float_t * emcEnergies = iniTower->GetEnergiesList() ;
808 //   Float_t ratio ;
809
810 //   iparam = 0 ;
811 //   while(iparam < nPar ){
812 //     xpar = fitparameters[iparam] ;
813 //     zpar = fitparameters[iparam+1] ;
814 //     epar = fitparameters[iparam+2] ;
815 //     iparam += 3 ;    
816     
817 //     AliEMCALTowerRecPoint * emcRP = 0 ;  
818
819 //     if(iniTower->IsTower()){ //create new entries in fTowerRecPoints...
820       
821 //       if(fNumberOfTowerClusters >= emcRecPoints->GetSize())
822 //      emcRecPoints->Expand(2*fNumberOfTowerClusters) ;
823       
824 //       (*emcRecPoints)[fNumberOfTowerClusters] = new AliEMCALTowerRecPoint("") ;
825 //       emcRP = (AliEMCALTowerRecPoint *) emcRecPoints->At(fNumberOfTowerClusters);
826 //       fNumberOfTowerClusters++ ;
827 //     }
828 //     else{//create new entries in fPreShoRecPoints
829 //       if(fNumberOfPreShoClusters >= cpvRecPoints->GetSize())
830 //      cpvRecPoints->Expand(2*fNumberOfPreShoClusters) ;
831       
832 //       (*cpvRecPoints)[fNumberOfPreShoClusters] = new AliEMCALPreShoRecPoint("") ;
833 //       emcRP = (AliEMCALTowerRecPoint *) cpvRecPoints->At(fNumberOfPreShoClusters);
834 //       fNumberOfPreShoClusters++ ;
835 //     }
836     
837 //     Float_t eDigit ;
838 //     for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
839 //       digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ; 
840 //       geom->AbsToRelNumbering(digit->GetId(), relid) ;
841 //       geom->RelPosInModule(relid, xDigit, zDigit) ;
842 //       distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
843 //       distance =  TMath::Sqrt(distance) ;
844 //       ratio = epar * ShowerShape(distance) / efit[iDigit] ; 
845 //       eDigit = emcEnergies[iDigit] * ratio ;
846 //       emcRP->AddDigit( *digit, eDigit ) ;
847 //     }        
848 //   }
849  
850 //   delete[] fitparameters ; 
851 //   delete[] efit ; 
852   
853 }
854
855 //_____________________________________________________________________________
856 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
857 {
858   // Calculates the Chi square for the cluster unfolding minimization
859   // Number of parameters, Gradient, Chi squared, parameters, what to do
860
861   abort() ; 
862  //  Fatal("AliEMCALClusterizerv1::UnfoldingChiSquare","-->Unfolding not implemented") ;
863
864 //   TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
865
866 //   AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint*) toMinuit->At(0)  ;
867 //   TClonesArray * digits = (TClonesArray*)toMinuit->At(1)  ;
868
869
870   
871 //   //  AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint *) gMinuit->GetObjectFit() ; // TowerRecPoint to fit
872
873 //   Int_t * emcDigits     = emcRP->GetDigitsList() ;
874
875 //   Int_t nOdigits = emcRP->GetDigitsMultiplicity() ; 
876
877 //   Float_t * emcEnergies = emcRP->GetEnergiesList() ;
878
879 //   const AliEMCALGeometry * geom = AliEMCALGetter::GetInstance()->EMCALGeometry() ; 
880 //   fret = 0. ;     
881 //   Int_t iparam ;
882
883 //   if(iflag == 2)
884 //     for(iparam = 0 ; iparam < nPar ; iparam++)    
885 //       Grad[iparam] = 0 ; // Will evaluate gradient
886   
887 //   Double_t efit ;    
888
889 //   AliEMCALDigit * digit ;
890 //   Int_t iDigit ;
891
892 //   for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
893
894 //     digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ; 
895
896 //     Int_t relid[4] ;
897 //     Float_t xDigit ;
898 //     Float_t zDigit ;
899
900 //     geom->AbsToRelNumbering(digit->GetId(), relid) ;
901
902 //     geom->RelPosInModule(relid, xDigit, zDigit) ;
903
904 //      if(iflag == 2){  // calculate gradient
905 //        Int_t iParam = 0 ;
906 //        efit = 0 ;
907 //        while(iParam < nPar ){
908 //       Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
909 //       iParam++ ; 
910 //       distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ; 
911 //       distance = TMath::Sqrt( distance ) ; 
912 //       iParam++ ;      
913 //       efit += x[iParam] * ShowerShape(distance) ;
914 //       iParam++ ;
915 //        }
916 //        Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E) 
917 //        iParam = 0 ;
918 //        while(iParam < nPar ){
919 //       Double_t xpar = x[iParam] ;
920 //       Double_t zpar = x[iParam+1] ;
921 //       Double_t epar = x[iParam+2] ;
922 //       Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
923 //       Double_t shape = sum * ShowerShape(dr) ;
924 //       Double_t r4 = dr*dr*dr*dr ;
925 //       Double_t r295 = TMath::Power(dr,2.95) ;
926 //       Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
927 //                                       0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
928          
929 //       Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ;  // Derivative over x    
930 //       iParam++ ; 
931 //       Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ;  // Derivative over z         
932 //       iParam++ ; 
933 //       Grad[iParam] += shape ;                                  // Derivative over energy             
934 //       iParam++ ; 
935 //        }
936 //      }
937 //      efit = 0;
938 //      iparam = 0 ;
939
940 //      while(iparam < nPar ){
941 //        Double_t xpar = x[iparam] ;
942 //        Double_t zpar = x[iparam+1] ;
943 //        Double_t epar = x[iparam+2] ;
944 //        iparam += 3 ;
945 //        Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
946 //        distance =  TMath::Sqrt(distance) ;
947 //        efit += epar * ShowerShape(distance) ;
948 //      }
949
950 //      fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ; 
951 //      // Here we assume, that sigma = sqrt(E)
952 //   }
953
954 }
955
956 //____________________________________________________________________________
957 void AliEMCALClusterizerv1::Print(Option_t * option)const
958 {
959   // Print clusterizer parameters
960
961   if( strcmp(GetName(), "") !=0 ){
962     
963     // Print parameters
964  
965     TString taskName(GetName()) ; 
966     taskName.ReplaceAll(Version(), "") ;
967
968     cout << "---------------"<< taskName.Data() << " " << GetTitle()<< "-----------" << endl 
969          << "Clusterizing digits from the file: " << taskName.Data() << endl 
970          << "                           Branch: " << GetName() << endl 
971          << endl 
972          << "                       EMC Clustering threshold = " << fTowerClusteringThreshold << endl
973          << "                       EMC Local Maximum cut    = " << fTowerLocMaxCut << endl
974          << "                       EMC Logarothmic weight   = " << fW0 << endl
975          << endl
976          << "                       CPV Clustering threshold = " << fPreShoClusteringThreshold << endl
977          << "                       CPV Local Maximum cut    = " << fPreShoLocMaxCut << endl
978        << "                       CPV Logarothmic weight   = " << fW0CPV << endl
979          << endl ;
980     if(fToUnfold)
981       cout << " Unfolding on " << endl ;
982     else
983       cout << " Unfolding off " << endl ;
984     
985     cout << "------------------------------------------------------------------" <<endl ;
986   }
987   else
988     cout << " AliEMCALClusterizerv1 not initialized " << endl ;
989 }
990 //____________________________________________________________________________
991 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
992 {
993   // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
994
995   TObjArray * towerRecPoints = AliEMCALGetter::GetInstance()->TowerRecPoints() ; 
996   TObjArray * preshoRecPoints = AliEMCALGetter::GetInstance()->PreShowerRecPoints() ; 
997
998   cout << "AliEMCALClusterizerv1: : event "<<gAlice->GetEvNumber() << endl ;
999   cout << "       Found "<< towerRecPoints->GetEntriesFast() << " TOWER Rec Points and " 
1000            << preshoRecPoints->GetEntriesFast() << " PRE SHOWER RecPoints" << endl ;
1001
1002   fRecPointsInRun +=  towerRecPoints->GetEntriesFast() ; 
1003   fRecPointsInRun +=  preshoRecPoints->GetEntriesFast() ; 
1004
1005   if(strstr(option,"all")) {
1006
1007     cout << "Tower clusters " << endl ;
1008     cout << " Index  Ene(MeV)   Multi  Module     phi     r  theta    Lambda 1   Lambda 2  # of prim  Primaries list "      <<  endl;      
1009     
1010     Int_t index ;
1011     for (index = 0 ; index < towerRecPoints->GetEntries() ; index++) {
1012       AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint * >(towerRecPoints->At(index)) ; 
1013       TVector3  globalpos;  
1014       rp->GetGlobalPosition(globalpos);
1015       Float_t lambda[2]; 
1016       rp->GetElipsAxis(lambda);
1017       Int_t * primaries; 
1018       Int_t nprimaries;
1019       primaries = rp->GetPrimaries(nprimaries);
1020
1021       cout << setw(4) << rp->GetIndexInList() << "   " 
1022            << setw(7) << setprecision(3) << rp->GetEnergy() << "      "
1023            << setw(3) << rp->GetMultiplicity() << "      " 
1024            << setw(1) << rp->GetEMCALArm() << "     " 
1025            << setw(5) << setprecision(4) << globalpos.X() << "  " 
1026            << setw(5) << setprecision(4) << globalpos.Y() << "  " 
1027            << setw(5) << setprecision(4) << globalpos.Z() << "     "
1028            << setw(4) << setprecision(2) << lambda[0]  << "  "
1029            << setw(4) << setprecision(2) << lambda[1]  << "  "
1030            << setw(2) << nprimaries << "  " ;
1031      
1032       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1033         cout << setw(4) <<   primaries[iprimary] << "  "  ;
1034       cout << endl ;     
1035     }
1036
1037     //Now plot Pre shower recPoints
1038
1039     cout << "-----------------------------------------------------------------------"<<endl ;
1040
1041     cout << "PreShower clusters " << endl ;
1042     cout << " Index  Ene(MeV)   Multi  Module     phi     r  theta    Lambda 1   Lambda 2  # of prim  Primaries list "      <<  endl;      
1043     
1044     for (index = 0 ; index < preshoRecPoints->GetEntries() ; index++) {
1045       AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(index)) ; 
1046       TVector3  globalpos;  
1047       rp->GetGlobalPosition(globalpos);
1048       Float_t lambda[2]; 
1049       rp->GetElipsAxis(lambda);
1050       Int_t * primaries; 
1051       Int_t nprimaries;
1052       primaries = rp->GetPrimaries(nprimaries);
1053
1054       cout << setw(4) << rp->GetIndexInList() << "   " 
1055            << setw(7) << setprecision(3) << rp->GetEnergy() << "      "
1056            << setw(3) << rp->GetMultiplicity() << "      " 
1057            << setw(1) << rp->GetEMCALArm() << "     " 
1058            << setw(5) << setprecision(4) << globalpos.X() << "  " 
1059            << setw(5) << setprecision(4) << globalpos.Y() << "  " 
1060            << setw(5) << setprecision(4) << globalpos.Z() << "     "
1061            << setw(4) << setprecision(2) << lambda[0]  << "  "
1062            << setw(4) << setprecision(2) << lambda[1]  << "  "
1063            << setw(2) << nprimaries << "  " ;
1064      
1065       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1066         cout << setw(4) <<   primaries[iprimary] << "  "  ;
1067       cout << endl ;     
1068     }
1069
1070     cout << "-----------------------------------------------------------------------"<<endl ;
1071   }
1072 }
1073