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