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