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