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