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