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