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