merge conflicts solved
[u/mrichter/AliRoot.git] / PHOS / AliPHOSClusterizerv1.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 /* $Log:
19    1 October 2000. Yuri Kharlov:
20      AreNeighbours()
21      PPSD upper layer is considered if number of layers>1
22
23    18 October 2000. Yuri Kharlov:
24      AliPHOSClusterizerv1()
25      CPV clusterizing parameters added
26
27      MakeClusters()
28      After first PPSD digit remove EMC digits only once
29 */
30 //*-- Author: Yves Schutz (SUBATECH)  & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
31 //////////////////////////////////////////////////////////////////////////////
32 //  Clusterization class. Performs clusterization (collects neighbouring active cells) and 
33 //  unfolds the clusters having several local maxima.  
34 //  Results are stored in TreeR#, branches PHOSEmcRP (EMC recPoints),
35 //  PHOSCpvRP (CPV RecPoints) and AliPHOSClusterizer (Clusterizer with all 
36 //  parameters including input digits branch title, thresholds etc.)
37 //  This TTask is normally called from Reconstructioner, but can as well be used in 
38 //  standalone mode.
39 // Use Case:
40 //  root [0] AliPHOSClusterizerv1 * cl = new AliPHOSClusterizerv1("galice.root")  
41 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
42 //               //reads gAlice from header file "..."                      
43 //  root [1] cl->ExecuteTask()  
44 //               //finds RecPoints in all events stored in galice.root
45 //  root [2] cl->SetDigitsBranch("digits2") 
46 //               //sets another title for Digitis (input) branch
47 //  root [3] cl->SetRecPointsBranch("recp2")  
48 //               //sets another title four output branches
49 //  root [4] cl->SetEmcLocalMaxCut(0.03)  
50 //               //set clusterization parameters
51 //  root [5] cl->ExecuteTask("deb all time")  
52 //               //once more finds RecPoints options are 
53 //               // deb - print number of found rec points
54 //               // deb all - print number of found RecPoints and some their characteristics 
55 //               // time - print benchmarking results
56
57 // --- ROOT system ---
58
59 #include "TROOT.h" 
60 #include "TFile.h" 
61 #include "TFolder.h" 
62 #include "TMath.h" 
63 #include "TMinuit.h"
64 #include "TTree.h" 
65 #include "TSystem.h" 
66 #include "TBenchmark.h"
67
68 // --- Standard library ---
69
70 #include <iostream.h>
71 #include <iomanip.h>
72
73 // --- AliRoot header files ---
74
75 #include "AliPHOSClusterizerv1.h"
76 #include "AliPHOSCpvRecPoint.h"
77 #include "AliPHOSDigit.h"
78 #include "AliPHOSDigitizer.h"
79 #include "AliPHOSEmcRecPoint.h"
80 #include "AliPHOS.h"
81 #include "AliPHOSPpsdRecPoint.h"
82 #include "AliRun.h"
83
84 ClassImp(AliPHOSClusterizerv1)
85   
86 //____________________________________________________________________________
87   AliPHOSClusterizerv1::AliPHOSClusterizerv1():AliPHOSClusterizer()
88 {
89   // default ctor (to be used mainly by Streamer)
90   SetName("AliPHOSClusterizer"); 
91   SetTitle("Version 1") ;
92   
93   fNumberOfCpvClusters     = 0 ; 
94   fNumberOfEmcClusters     = 0 ; 
95   
96   fCpvClusteringThreshold  = 0.0;
97   fEmcClusteringThreshold  = 0.2;   
98   fPpsdClusteringThreshold = 0.0000002 ;
99   
100   fEmcLocMaxCut            = 0.03 ;
101   fCpvLocMaxCut            = 0.03 ;
102   
103   fW0                      = 4.5 ;
104   fW0CPV                   = 4.0 ;
105   
106   fGeom  = 0 ;
107   
108   fDigits = 0 ;
109   fDigitizer = 0 ;
110   fEmcRecPoints = 0 ;
111   fCpvRecPoints = 0 ;
112
113   fIsInitialized = kFALSE ;
114   
115 }
116 //____________________________________________________________________________
117 AliPHOSClusterizerv1::AliPHOSClusterizerv1(const char* headerFile,const char* digitsFile):AliPHOSClusterizer()
118 {
119   // ctor with the indication of the file where header Tree and digits Tree are stored
120   SetName("AliPHOSClusterizer");
121   SetTitle("Version 1") ;
122   
123   fNumberOfCpvClusters     = 0 ; 
124   fNumberOfEmcClusters     = 0 ; 
125   
126   fCpvClusteringThreshold  = 0.0;
127   fEmcClusteringThreshold  = 0.2;   
128   fPpsdClusteringThreshold = 0.0000002 ;
129   
130   fEmcLocMaxCut            = 0.03 ;
131   fCpvLocMaxCut            = 0.03 ;
132   
133   fW0                      = 4.5 ;
134   fW0CPV                   = 4.0 ;
135   
136   fToUnfold = kTRUE ;
137   
138   fHeaderFileName = headerFile  ;
139   fDigitsBranchTitle = digitsFile ;
140   
141   TFile * file = (TFile*) gROOT->GetFile(fHeaderFileName.Data() ) ;
142   
143   if(file == 0){
144     if(fHeaderFileName.Contains("rfio")) // if we read file using HPSS
145       file =    TFile::Open(fHeaderFileName.Data(),"update") ;
146     else
147       file = new TFile(fHeaderFileName.Data(),"update") ;
148     gAlice = (AliRun *) file->Get("gAlice") ;
149   }
150   
151   AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;    
152   fGeom  = AliPHOSGeometry::GetInstance(phos->GetGeometry()->GetName(),phos->GetGeometry()->GetTitle() );
153   
154   fDigits = new TClonesArray("AliPHOSDigit",10) ;
155   fDigitizer = new AliPHOSDigitizer() ;
156   fEmcRecPoints = new TObjArray(200) ;
157   fCpvRecPoints = new TObjArray(200) ;
158   fEmcRecPoints->SetOwner();          // This lets Clear() really detete rec.points in array
159   fCpvRecPoints->SetOwner();
160   
161   if(!gMinuit) gMinuit = new TMinuit(100) ;
162   
163   // add Task to //root/Tasks folder
164   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
165   TTask * phostasks = (TTask*)(roottasks->GetListOfTasks()->FindObject("PHOS"));
166   phostasks->Add(this) ; 
167
168
169   fIsInitialized = kTRUE ;
170
171 }
172 //____________________________________________________________________________
173 void AliPHOSClusterizerv1::Exec(Option_t * option)
174 {
175   // Steering function
176
177   if(!fIsInitialized) Init() ;
178
179   if(strstr(option,"tim"))
180     gBenchmark->Start("PHOSClusterizer");  
181   
182   Int_t nEvents = (Int_t) gAlice->TreeE()->GetEntries() ;
183   
184   for(fEvent = 0;fEvent< nEvents; fEvent++){
185     if(!ReadDigits())  //reads digits for event fEvent
186       return;
187     MakeClusters() ;
188     
189     if(fToUnfold) MakeUnfolding() ;
190     WriteRecPoints() ;
191     if(strstr(option,"deb"))
192       PrintRecPoints(option) ;
193   }
194   
195   if(strstr(option,"tim")){
196     gBenchmark->Stop("PHOSClusterizer");
197     cout << "AliPHOSClusterizer:" << endl ;
198     cout << "  took " << gBenchmark->GetCpuTime("PHOSClusterizer") << " seconds for Clusterizing " 
199          <<  gBenchmark->GetCpuTime("PHOSClusterizer")/nEvents << " seconds per event " << endl ;
200     cout << endl ;
201   }
202   
203   
204 }
205
206 //____________________________________________________________________________
207 Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, int * maxAt, Float_t * maxAtEnergy,
208                                     Int_t nPar, Float_t * fitparameters) const
209
210   // Calls TMinuit to fit the energy distribution of a cluster with several maxima 
211   // The initial values for fitting procedure are set equal to the positions of local maxima.
212   // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
213
214   gMinuit->mncler();                     // Reset Minuit's list of paramters
215   gMinuit->SetPrintLevel(-1) ;           // No Printout
216   gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;  
217                                          // To set the address of the minimization function 
218
219   TList * toMinuit = new TList();
220   toMinuit->AddAt(emcRP,0) ;
221   toMinuit->AddAt(fDigits,1) ;
222   
223   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare
224
225   // filling initial values for fit parameters
226   AliPHOSDigit * digit ;
227
228   Int_t ierflg  = 0; 
229   Int_t index   = 0 ;
230   Int_t nDigits = (Int_t) nPar / 3 ;
231
232   Int_t iDigit ;
233
234
235   for(iDigit = 0; iDigit < nDigits; iDigit++){
236     digit = (AliPHOSDigit *) maxAt[iDigit]; 
237
238     Int_t relid[4] ;
239     Float_t x ;
240     Float_t z ;
241     fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
242     fGeom->RelPosInModule(relid, x, z) ;
243
244     Float_t energy = maxAtEnergy[iDigit] ;
245
246     gMinuit->mnparm(index, "x",  x, 0.1, 0, 0, ierflg) ;
247     index++ ;   
248     if(ierflg != 0){ 
249       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : x = " << x << endl ;
250       return kFALSE;
251     }
252     gMinuit->mnparm(index, "z",  z, 0.1, 0, 0, ierflg) ;
253     index++ ;   
254     if(ierflg != 0){
255       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : z = " << z << endl ;
256       return kFALSE;
257     }
258     gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;
259     index++ ;   
260     if(ierflg != 0){
261       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : energy = " << energy << endl ;      
262       return kFALSE;
263     }
264   }
265
266   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
267                       //  depends on it. 
268   Double_t p1 = 1.0 ;
269   Double_t p2 = 0.0 ;
270
271   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
272   gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
273   gMinuit->SetMaxIterations(5);
274   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
275
276   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
277
278   if(ierflg == 4){  // Minimum not found   
279     cout << "PHOS Unfolding>  Fit not converged, cluster abandoned "<< endl ;      
280     return kFALSE ;
281   }            
282   for(index = 0; index < nPar; index++){
283     Double_t err ;
284     Double_t val ;
285     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index
286     fitparameters[index] = val ;
287    }
288
289   delete toMinuit ;
290   return kTRUE;
291
292 }
293
294 //____________________________________________________________________________
295 void AliPHOSClusterizerv1::Init()
296 {
297   //Make all memory allocations which can not be done in default constructor.
298   if(!fIsInitialized){
299     if(fHeaderFileName.IsNull())
300       fHeaderFileName = "galice.root" ;
301     
302     
303     TFile * file = (TFile*) gROOT->GetFile(fHeaderFileName.Data() ) ;
304     
305     if(file == 0){
306       file = new TFile(fHeaderFileName.Data(),"update") ;
307       gAlice = (AliRun *) file->Get("gAlice") ;
308     }
309      
310     AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;    
311     fGeom  = AliPHOSGeometry::GetInstance(phos->GetGeometry()->GetName(),phos->GetGeometry()->GetTitle() );
312     
313     fDigits = new TClonesArray("AliPHOSDigit",10) ;
314     fDigitizer = new AliPHOSDigitizer() ;
315     fEmcRecPoints = new TObjArray(200) ;
316     fCpvRecPoints = new TObjArray(200) ;
317     fEmcRecPoints->SetOwner();          // This lets Clear() really detete rec.points in array
318     fCpvRecPoints->SetOwner();
319     
320     if(!gMinuit) gMinuit = new TMinuit(100) ;
321
322     // add Task to //root/Tasks folder
323     TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
324     TTask * phostasks = (TTask*)(roottasks->GetListOfTasks()->FindObject("PHOS"));
325     phostasks->Add(this) ; 
326
327     fIsInitialized = kTRUE ;
328   }
329 }
330 //____________________________________________________________________________
331 Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
332 {
333   // Gives the neighbourness of two digits = 0 are not neighbour but continue searching 
334   //                                       = 1 are neighbour
335   //                                       = 2 are not neighbour but do not continue searching
336   // neighbours are defined as digits having at least a common vertex 
337   // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster 
338   //                                      which is compared to a digit (d2)  not yet in a cluster  
339
340   Int_t rv = 0 ; 
341
342   Int_t relid1[4] ; 
343   fGeom->AbsToRelNumbering(d1->GetId(), relid1) ; 
344
345   Int_t relid2[4] ; 
346   fGeom->AbsToRelNumbering(d2->GetId(), relid2) ; 
347  
348   if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module and the same PPSD Module 
349     Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
350     Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
351     
352     if (( coldiff <= 1 )  && ( rowdiff <= 1 )){
353       rv = 1 ; 
354     }
355     else {
356       if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1)) 
357         rv = 2; //  Difference in row numbers is too large to look further 
358     }
359
360   } 
361   else {
362     
363     if( (relid1[0] < relid2[0]) || (relid1[1] < relid2[1]) )  
364       rv=2 ;
365
366   }
367
368   //Do NOT clusterize upper PPSD  
369   if( IsInPpsd(d1) && IsInPpsd(d2) &&
370      relid1[1] > 0                 &&
371      relid1[1] < fGeom->GetNumberOfPadsPhi()*fGeom->GetNumberOfPadsPhi() ) rv = 2 ;
372
373   return rv ; 
374 }
375
376
377 //____________________________________________________________________________
378 Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
379 {
380   // Tells if (true) or not (false) the digit is in a PHOS-EMC module
381  
382   Bool_t rv = kFALSE ; 
383
384   Int_t relid[4] ; 
385   fGeom->AbsToRelNumbering(digit->GetId(), relid) ; 
386
387   if ( relid[1] == 0  ) rv = kTRUE; 
388
389   return rv ; 
390 }
391
392 //____________________________________________________________________________
393 Bool_t AliPHOSClusterizerv1::IsInPpsd(AliPHOSDigit * digit) const
394 {
395   // Tells if (true) or not (false) the digit is in a PHOS-PPSD module
396  
397   Bool_t rv = kFALSE ; 
398
399   Int_t relid[4] ; 
400   fGeom->AbsToRelNumbering(digit->GetId(), relid) ; 
401
402   if ( relid[1] > 0 && relid[0] > fGeom->GetNCPVModules() ) rv = kTRUE; 
403
404   return rv ; 
405 }
406
407 //____________________________________________________________________________
408 Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
409 {
410   // Tells if (true) or not (false) the digit is in a PHOS-CPV module
411  
412   Bool_t rv = kFALSE ; 
413
414   Int_t relid[4] ; 
415   fGeom->AbsToRelNumbering(digit->GetId(), relid) ; 
416
417   if ( relid[1] > 0 && relid[0] <= fGeom->GetNCPVModules() ) rv = kTRUE; 
418
419   return rv ; 
420 }
421 //____________________________________________________________________________
422 Bool_t AliPHOSClusterizerv1::ReadDigits()
423  {
424    // reads digitis with specified title from TreeD
425
426   fNumberOfEmcClusters  = 0 ;
427   fNumberOfCpvClusters  = 0 ;
428
429   // Get Digits Tree header from file
430   gAlice->GetEvent(fEvent) ;
431   gAlice->SetEvent(fEvent) ;
432
433   TTree * treeD = gAlice->TreeD()  ;
434
435   if(treeD==0){
436     char treeName[20]; 
437     sprintf(treeName,"TreeD%d",fEvent);
438     cout << "Error in AliPHOSClusterizerv1 : no "<<treeName << endl  ;
439     cout << "    Do nothing " << endl ;
440     return kFALSE ;
441   }
442
443   TBranch * digitsBranch = 0;
444   TBranch * digitizerBranch = 0;
445
446   TObjArray * branches = treeD->GetListOfBranches() ;
447   Int_t ibranch;
448   Bool_t phosNotFound = kTRUE ;
449   Bool_t digitizerNotFound = kTRUE ;
450   
451   for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
452
453     if(phosNotFound){
454       digitsBranch=(TBranch *) branches->At(ibranch) ;
455       if( fDigitsBranchTitle.CompareTo(digitsBranch->GetTitle())==0 )
456         if( strcmp(digitsBranch->GetName(),"PHOS") == 0)
457           phosNotFound = kFALSE ;
458     }
459     
460     if(digitizerNotFound){
461       digitizerBranch = (TBranch *) branches->At(ibranch) ;
462       if( fDigitsBranchTitle.CompareTo(digitizerBranch->GetTitle()) == 0)
463         if( strcmp(digitizerBranch->GetName(),"AliPHOSDigitizer") == 0) 
464           digitizerNotFound = kFALSE ;
465     }
466    
467   }
468   
469   if(digitizerNotFound || phosNotFound){
470     cout << "ERROR in AliPHOSClusterizerv1: " << endl ;
471     cout << "      Can't find Branch with digits or Digitizer "<< endl ; ;
472     cout << "      Do nothing" <<endl  ;
473     return kFALSE ;
474   }
475   
476   digitsBranch->SetAddress(&fDigits) ;
477   digitizerBranch->SetAddress(&fDigitizer) ;
478
479   digitsBranch->GetEntry(0) ;
480   digitizerBranch->GetEntry(0) ;
481   
482   fPedestal = fDigitizer->GetPedestal() ;
483   fSlope    = fDigitizer->GetSlope() ;
484   return kTRUE ;
485 }
486
487 //____________________________________________________________________________
488 void AliPHOSClusterizerv1::WriteRecPoints()
489 {
490   // Checks, if PHOSEmcRP etc. branches with given title already exist, 
491   // if yes exits without writing, 
492   // else creates new branches with given title
493   // fills and writes into TreeR.
494   
495   Int_t index ;
496   //Evaluate poisition, dispersion and other RecPoint properties...
497   for(index = 0; index < fEmcRecPoints->GetEntries(); index++)
498     ((AliPHOSEmcRecPoint *)fEmcRecPoints->At(index))->EvalAll(fW0,fDigits) ;
499
500   fEmcRecPoints->Sort() ;
501
502   for(index = 0; index < fEmcRecPoints->GetEntries(); index++)
503     ((AliPHOSEmcRecPoint *)fEmcRecPoints->At(index))->SetIndexInList(index) ;
504
505   fEmcRecPoints->Expand(fEmcRecPoints->GetEntriesFast()) ; 
506
507   //Now the same for CPV
508   for(index = 0; index < fCpvRecPoints->GetEntries(); index++)
509     ((AliPHOSRecPoint *)fCpvRecPoints->At(index))->EvalAll(fW0CPV,fDigits)  ;
510
511   fCpvRecPoints->Sort() ;
512
513   for(index = 0; index < fCpvRecPoints->GetEntries(); index++)
514     ((AliPHOSRecPoint *)fCpvRecPoints->At(index))->SetIndexInList(index) ;
515
516   fCpvRecPoints->Expand(fCpvRecPoints->GetEntriesFast()) ;
517
518   if(gAlice->TreeR()==0)
519     gAlice->MakeTree("R") ;
520   
521
522   //Check, if branches already exist
523   TBranch * emcBranch = 0;
524   TBranch * cpvBranch = 0;
525   TBranch * clusterizerBranch = 0;
526   
527   TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
528   Int_t ibranch;
529   Bool_t emcNotFound = kTRUE ;
530   Bool_t cpvNotFound = kTRUE ;  
531   Bool_t clusterizerNotFound = kTRUE ;
532   
533   for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
534     
535     if(emcNotFound){
536       emcBranch=(TBranch *) branches->At(ibranch) ;
537       if( fRecPointsBranchTitle.CompareTo(emcBranch->GetTitle())==0 ){
538         if( strcmp(emcBranch->GetName(),"PHOSEmcRP") == 0) {
539           emcNotFound = kFALSE ;
540         }
541       }
542     }
543     if(cpvNotFound){
544       cpvBranch=(TBranch *) branches->At(ibranch) ;
545       if( fRecPointsBranchTitle.CompareTo(cpvBranch->GetTitle())==0 ){
546         if( strcmp(cpvBranch->GetName(),"PHOSCpvRP") == 0) {
547           cpvNotFound = kFALSE ;
548         }
549       }
550     }
551
552     if(clusterizerNotFound){
553       clusterizerBranch = (TBranch *) branches->At(ibranch) ;
554       if( fRecPointsBranchTitle.CompareTo(clusterizerBranch->GetTitle()) == 0){
555         if( strcmp(clusterizerBranch->GetName(),"AliPHOSClusterizer") == 0) {
556           clusterizerNotFound = kFALSE ;
557         }
558       }
559     }
560     
561   }
562   
563   if(!(clusterizerNotFound && emcNotFound && cpvNotFound)){
564     cout << "AliPHOSClusterizer error" << endl;
565     cout << "       Branches PHOSEmcRP, PHOSCpvRP and AliPHOSClusterizer " << endl ;
566     cout << "       with title '" << fRecPointsBranchTitle.Data() <<"' already exist" << endl ;
567     cout << "       can not overwrite " << endl ;
568     return ;
569   }
570
571   //Make branches in TreeR for RecPoints and Clusterizer
572   char * filename = 0;
573   if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){   //generating file name
574     filename = new char[strlen(gAlice->GetBaseFile())+20] ;
575     sprintf(filename,"%s/PHOS.Reco.root",gAlice->GetBaseFile()) ;
576   }
577   
578   //Make new branches
579   TDirectory *cwd = gDirectory;
580   
581   //First EMC
582   Int_t bufferSize = 32000 ;    
583   Int_t splitlevel = 0 ;
584   emcBranch = gAlice->TreeR()->Branch("PHOSEmcRP","TObjArray",&fEmcRecPoints,bufferSize,splitlevel);
585   emcBranch->SetTitle(fRecPointsBranchTitle.Data());
586   if (filename) {
587     emcBranch->SetFile(filename);
588     TIter next( emcBranch->GetListOfBranches());
589     TBranch * sb ;
590     while ((sb=(TBranch*)next())) {
591       sb->SetFile(filename);
592     }   
593
594     cwd->cd();
595   }
596     
597   //Now CPV branch
598   cpvBranch = gAlice->TreeR()->Branch("PHOSCpvRP","TObjArray",&fCpvRecPoints,bufferSize,splitlevel);
599   cpvBranch->SetTitle(fRecPointsBranchTitle.Data());
600   if (filename) {
601     cpvBranch->SetFile(filename);
602     TIter next( cpvBranch->GetListOfBranches());
603     TBranch * sb;
604     while ((sb=(TBranch*)next())) {
605       sb->SetFile(filename);
606     }   
607     cwd->cd();
608   } 
609     
610   //And Finally  clusterizer branch
611   AliPHOSClusterizerv1 * cl = this ;
612   clusterizerBranch = gAlice->TreeR()->Branch("AliPHOSClusterizer","AliPHOSClusterizerv1",
613                                               &cl,bufferSize,splitlevel);
614   clusterizerBranch->SetTitle(fRecPointsBranchTitle.Data());
615   if (filename) {
616     clusterizerBranch->SetFile(filename);
617     TIter next( clusterizerBranch->GetListOfBranches());
618     TBranch * sb ;
619     while ((sb=(TBranch*)next())) {
620       sb->SetFile(filename);
621     }   
622     cwd->cd();
623   }
624   
625   emcBranch->Fill() ;
626   cpvBranch->Fill() ;
627   clusterizerBranch->Fill() ;
628
629   gAlice->TreeR()->Write(0,kOverwrite) ;  
630   
631 }
632
633 //____________________________________________________________________________
634 void AliPHOSClusterizerv1::MakeClusters()
635 {
636   // Steering method to construct the clusters stored in a list of Reconstructed Points
637   // A cluster is defined as a list of neighbour digits
638   fEmcRecPoints->Clear() ;
639   fCpvRecPoints->Clear() ;
640   
641
642   // Clusterization starts  
643   TClonesArray * digits =  (TClonesArray*)fDigits->Clone() ;
644
645   TIter nextdigit(digits) ; 
646   AliPHOSDigit * digit ; 
647   Bool_t notremoved = kTRUE ;
648
649   while ( (digit = (AliPHOSDigit *)nextdigit()) ) { // scan over the list of digits
650     AliPHOSRecPoint * clu = 0 ; 
651
652     TArrayI clusterdigitslist(1500) ;   
653     Int_t index ;
654
655     if (( IsInEmc (digit) && Calibrate(digit->GetAmp()) > fEmcClusteringThreshold  ) || 
656         ( IsInPpsd(digit) && Calibrate(digit->GetAmp()) > fPpsdClusteringThreshold ) ||
657         ( IsInCpv (digit) && Calibrate(digit->GetAmp()) > fCpvClusteringThreshold  ) ) {
658       
659       Int_t iDigitInCluster = 0 ; 
660
661       if  ( IsInEmc(digit) ) {   
662         // start a new EMC RecPoint
663         if(fNumberOfEmcClusters >= fEmcRecPoints->GetSize()) fEmcRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
664         fEmcRecPoints->AddAt(new  AliPHOSEmcRecPoint(), fNumberOfEmcClusters) ;
665         clu = (AliPHOSEmcRecPoint *) fEmcRecPoints->At(fNumberOfEmcClusters) ; 
666         fNumberOfEmcClusters++ ; 
667         clu->AddDigit(*digit, Calibrate(digit->GetAmp())) ; 
668         clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;  
669         iDigitInCluster++ ; 
670         digits->Remove(digit) ; 
671
672       } else { 
673         
674         // start a new PPSD/CPV cluster
675         if(fNumberOfCpvClusters >= fCpvRecPoints->GetSize()) fCpvRecPoints->Expand(2*fNumberOfCpvClusters+1);
676         if(IsInPpsd(digit)) 
677           fCpvRecPoints->AddAt(new AliPHOSPpsdRecPoint(),fNumberOfCpvClusters) ;
678         else
679           fCpvRecPoints->AddAt(new AliPHOSCpvRecPoint(), fNumberOfCpvClusters) ;
680         clu =  (AliPHOSPpsdRecPoint *) fCpvRecPoints->At(fNumberOfCpvClusters)  ;  
681         fNumberOfCpvClusters++ ; 
682
683         clu->AddDigit(*digit, Calibrate(digit->GetAmp()) ) ;    
684         clusterdigitslist[iDigitInCluster] = digit->GetIndexInList()  ; 
685         iDigitInCluster++ ; 
686         digits->Remove(digit) ; 
687         nextdigit.Reset() ;
688         
689         // Here we remove resting EMC digits, which cannot make cluster
690         
691         if( notremoved ) { 
692           while( ( digit = (AliPHOSDigit *)nextdigit() ) ) {
693             if( IsInEmc(digit) ) 
694               digits->Remove(digit) ;
695             else 
696               break ;
697           }
698           notremoved = kFALSE ;
699         }
700         
701       } // else        
702       
703       nextdigit.Reset() ;
704       
705       AliPHOSDigit * digitN ; 
706       index = 0 ;
707       while (index < iDigitInCluster){ // scan over digits already in cluster 
708         digit =  (AliPHOSDigit*)fDigits->At(clusterdigitslist[index])  ;      
709         index++ ; 
710         while ( (digitN = (AliPHOSDigit *)nextdigit()) ) { // scan over the reduced list of digits 
711           Int_t ineb = AreNeighbours(digit, digitN);       // call (digit,digitN) in THAT oder !!!!!
712           switch (ineb ) {
713           case 0 :   // not a neighbour
714             break ;
715           case 1 :   // are neighbours 
716             clu->AddDigit(*digitN, Calibrate( digitN->GetAmp() ) ) ;
717             clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ; 
718             iDigitInCluster++ ; 
719             digits->Remove(digitN) ;
720             break ;
721           case 2 :   // too far from each other
722             goto endofloop;   
723           } // switch
724           
725         } // while digitN
726         
727       endofloop: ;
728         nextdigit.Reset() ; 
729         
730       } // loop over cluster     
731
732     } // energy theshold  
733
734     
735   } // while digit
736
737   delete digits ;
738
739 }
740
741 //____________________________________________________________________________
742 void AliPHOSClusterizerv1::MakeUnfolding()
743 {
744   // Unfolds clusters using the shape of an ElectroMagnetic shower
745   // Performs unfolding of all EMC/CPV but NOT ppsd clusters
746
747   // Unfold first EMC clusters 
748   if(fNumberOfEmcClusters > 0){
749
750     Int_t nModulesToUnfold = fGeom->GetNModules() ; 
751
752     Int_t numberofNotUnfolded = fNumberOfEmcClusters ; 
753     Int_t index ;   
754     for(index = 0 ; index < numberofNotUnfolded ; index++){
755       
756       AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint *) fEmcRecPoints->At(index) ;
757       if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
758         break ;
759       
760       Int_t nMultipl = emcRecPoint->GetMultiplicity() ; 
761       Int_t * maxAt = new Int_t[nMultipl] ;
762       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
763       Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,fDigits) ;
764       
765       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0       
766         UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
767         fEmcRecPoints->Remove(emcRecPoint); 
768         fEmcRecPoints->Compress() ;
769         index-- ;
770         fNumberOfEmcClusters -- ;
771         numberofNotUnfolded-- ;
772       }
773       
774       delete[] maxAt ; 
775       delete[] maxAtEnergy ; 
776     }
777   } 
778   // Unfolding of EMC clusters finished
779
780
781   // Unfold now CPV clusters
782   if(fNumberOfCpvClusters > 0){
783     
784     Int_t nModulesToUnfold = fGeom->GetNCPVModules() ;
785
786     Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;     
787     Int_t index ;   
788     for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
789       
790       AliPHOSRecPoint * recPoint = (AliPHOSRecPoint *) fCpvRecPoints->At(index) ;
791
792       if(recPoint->GetPHOSMod()> nModulesToUnfold)
793         break ;
794       
795       AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint*) recPoint ; 
796       
797       Int_t nMultipl = emcRecPoint->GetMultiplicity() ; 
798       Int_t * maxAt = new Int_t[nMultipl] ;
799       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
800       Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fCpvLocMaxCut,fDigits) ;
801       
802       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0       
803         UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
804         fCpvRecPoints->Remove(emcRecPoint); 
805         fCpvRecPoints->Compress() ;
806         index-- ;
807         numberofCpvNotUnfolded-- ;
808         fNumberOfCpvClusters-- ;
809       }
810       
811       delete[] maxAt ; 
812       delete[] maxAtEnergy ; 
813     } 
814   }
815   //Unfolding of Cpv clusters finished
816   
817 }
818
819 //____________________________________________________________________________
820 Double_t  AliPHOSClusterizerv1::ShowerShape(Double_t r)
821
822   // Shape of the shower (see PHOS TDR)
823   // If you change this function, change also the gradient evaluation in ChiSquare()
824
825   Double_t r4    = r*r*r*r ;
826   Double_t r295  = TMath::Power(r, 2.95) ;
827   Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
828   return shape ;
829 }
830
831 //____________________________________________________________________________
832 void  AliPHOSClusterizerv1::UnfoldCluster(AliPHOSEmcRecPoint * iniEmc, 
833                                                  Int_t nMax, 
834                                                  int * maxAt, 
835                                                  Float_t * maxAtEnergy)
836
837   // Performs the unfolding of a cluster with nMax overlapping showers 
838
839   Int_t nPar = 3 * nMax ;
840   Float_t * fitparameters = new Float_t[nPar] ;
841
842   Bool_t rv = FindFit(iniEmc, maxAt, maxAtEnergy, nPar, fitparameters) ;
843   if( !rv ) {
844     // Fit failed, return and remove cluster
845     delete[] fitparameters ; 
846     return ;
847   }
848
849   // create ufolded rec points and fill them with new energy lists
850   // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
851   // and later correct this number in acordance with actual energy deposition
852
853   Int_t nDigits = iniEmc->GetMultiplicity() ;  
854   Float_t * efit = new Float_t[nDigits] ;
855   Float_t xDigit,zDigit,distance ;
856   Float_t xpar,zpar,epar  ;
857   Int_t relid[4] ;
858   AliPHOSDigit * digit ;
859   Int_t * emcDigits = iniEmc->GetDigitsList() ;
860
861   Int_t iparam ;
862   Int_t iDigit ;
863   for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
864     digit = (AliPHOSDigit*) fDigits->At(emcDigits[iDigit] ) ;   
865     fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
866     fGeom->RelPosInModule(relid, xDigit, zDigit) ;
867     efit[iDigit] = 0;
868
869     iparam = 0 ;    
870     while(iparam < nPar ){
871       xpar = fitparameters[iparam] ;
872       zpar = fitparameters[iparam+1] ;
873       epar = fitparameters[iparam+2] ;
874       iparam += 3 ;
875       distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
876       distance =  TMath::Sqrt(distance) ;
877       efit[iDigit] += epar * ShowerShape(distance) ;
878     }
879   }
880   
881
882   // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
883   // so that energy deposited in each cell is distributed betwin new clusters proportionally
884   // to its contribution to efit
885
886   Float_t * emcEnergies = iniEmc->GetEnergiesList() ;
887   Float_t ratio ;
888
889   iparam = 0 ;
890   while(iparam < nPar ){
891     xpar = fitparameters[iparam] ;
892     zpar = fitparameters[iparam+1] ;
893     epar = fitparameters[iparam+2] ;
894     iparam += 3 ;    
895     
896     AliPHOSEmcRecPoint * emcRP ;  
897
898     if(iniEmc->IsEmc()){ //create new entries in fEmcRecPoints...
899       
900       if(fNumberOfEmcClusters >= fEmcRecPoints->GetSize())
901         fEmcRecPoints->Expand(2*fNumberOfEmcClusters) ;
902       
903       (*fEmcRecPoints)[fNumberOfEmcClusters] = new AliPHOSEmcRecPoint() ;
904       emcRP = (AliPHOSEmcRecPoint *) fEmcRecPoints->At(fNumberOfEmcClusters);
905       fNumberOfEmcClusters++ ;
906     }
907     else{//create new entries in fCpvRecPoints
908       if(fNumberOfCpvClusters >= fCpvRecPoints->GetSize())
909         fCpvRecPoints->Expand(2*fNumberOfCpvClusters) ;
910       
911       (*fCpvRecPoints)[fNumberOfCpvClusters] = new AliPHOSCpvRecPoint() ;
912       emcRP = (AliPHOSEmcRecPoint *) fCpvRecPoints->At(fNumberOfCpvClusters);
913       fNumberOfCpvClusters++ ;
914     }
915     
916     Float_t eDigit ;
917     for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
918       digit = (AliPHOSDigit*) fDigits->At( emcDigits[iDigit] ) ; 
919       fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
920       fGeom->RelPosInModule(relid, xDigit, zDigit) ;
921       distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
922       distance =  TMath::Sqrt(distance) ;
923       ratio = epar * ShowerShape(distance) / efit[iDigit] ; 
924       eDigit = emcEnergies[iDigit] * ratio ;
925       emcRP->AddDigit( *digit, eDigit ) ;
926     }   
927   }
928  
929   delete[] fitparameters ; 
930   delete[] efit ; 
931   
932 }
933
934 //_____________________________________________________________________________
935 void AliPHOSClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
936 {
937   // Calculates the Chi square for the cluster unfolding minimization
938   // Number of parameters, Gradient, Chi squared, parameters, what to do
939
940
941   TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
942
943   AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint*) toMinuit->At(0)  ;
944   TClonesArray * digits = (TClonesArray*)toMinuit->At(1)  ;
945
946
947   
948   //  AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint *) gMinuit->GetObjectFit() ; // EmcRecPoint to fit
949
950   Int_t * emcDigits     = emcRP->GetDigitsList() ;
951
952   Int_t nOfDigits = emcRP->GetDigitsMultiplicity() ; 
953
954   Float_t * emcEnergies = emcRP->GetEnergiesList() ;
955
956   AliPHOSGeometry * geom = AliPHOSGeometry::GetInstance() ;
957
958   fret = 0. ;     
959   Int_t iparam ;
960
961   if(iflag == 2)
962     for(iparam = 0 ; iparam < nPar ; iparam++)    
963       Grad[iparam] = 0 ; // Will evaluate gradient
964   
965   Double_t efit ;    
966
967   AliPHOSDigit * digit ;
968   Int_t iDigit ;
969
970   for( iDigit = 0 ; iDigit < nOfDigits ; iDigit++) {
971
972     digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ; 
973
974     Int_t relid[4] ;
975     Float_t xDigit ;
976     Float_t zDigit ;
977
978     geom->AbsToRelNumbering(digit->GetId(), relid) ;
979
980     geom->RelPosInModule(relid, xDigit, zDigit) ;
981
982      if(iflag == 2){  // calculate gradient
983        Int_t iParam = 0 ;
984        efit = 0 ;
985        while(iParam < nPar ){
986          Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
987          iParam++ ; 
988          distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ; 
989          distance = TMath::Sqrt( distance ) ; 
990          iParam++ ;      
991          efit += x[iParam] * ShowerShape(distance) ;
992          iParam++ ;
993        }
994        Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E) 
995        iParam = 0 ;
996        while(iParam < nPar ){
997          Double_t xpar = x[iParam] ;
998          Double_t zpar = x[iParam+1] ;
999          Double_t epar = x[iParam+2] ;
1000          Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
1001          Double_t shape = sum * ShowerShape(dr) ;
1002          Double_t r4 = dr*dr*dr*dr ;
1003          Double_t r295 = TMath::Power(dr,2.95) ;
1004          Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
1005                                          0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
1006          
1007          Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ;  // Derivative over x    
1008          iParam++ ; 
1009          Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ;  // Derivative over z         
1010          iParam++ ; 
1011          Grad[iParam] += shape ;                                  // Derivative over energy             
1012          iParam++ ; 
1013        }
1014      }
1015      efit = 0;
1016      iparam = 0 ;
1017
1018      while(iparam < nPar ){
1019        Double_t xpar = x[iparam] ;
1020        Double_t zpar = x[iparam+1] ;
1021        Double_t epar = x[iparam+2] ;
1022        iparam += 3 ;
1023        Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
1024        distance =  TMath::Sqrt(distance) ;
1025        efit += epar * ShowerShape(distance) ;
1026      }
1027
1028      fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ; 
1029      // Here we assume, that sigma = sqrt(E)
1030   }
1031
1032 }
1033
1034 //____________________________________________________________________________
1035 void AliPHOSClusterizerv1::Print(Option_t * option)const
1036 {
1037   // Print clusterizer parameters
1038
1039   if(fIsInitialized){
1040     
1041     // Print parameters
1042     
1043     cout << "---------------"<< GetName() << " " << GetTitle()<< "-----------" << endl 
1044          << "Clusterizing digits from the file: " << fHeaderFileName.Data() << endl 
1045          << "                           Branch: " << fDigitsBranchTitle.Data() << endl 
1046          << endl 
1047          << "                       EMC Clustering threshold = " << fEmcClusteringThreshold << endl
1048          << "                       EMC Local Maximum cut    = " << fEmcLocMaxCut << endl
1049          << "                       EMC Logarothmic weight   = " << fW0 << endl
1050          << endl
1051          << "                       CPV Clustering threshold = " << fCpvClusteringThreshold << endl
1052          << "                       CPV Local Maximum cut    = " << fCpvLocMaxCut << endl
1053        << "                       CPV Logarothmic weight   = " << fW0CPV << endl
1054          << endl
1055          << "                      PPSD  Clustering threshold = " << fPpsdClusteringThreshold << endl;
1056     if(fToUnfold)
1057       cout << " Unfolding on " << endl ;
1058     else
1059       cout << " Unfolding off " << endl ;
1060     
1061     cout << "------------------------------------------------------------------" <<endl ;
1062   }
1063   else
1064     cout << " AliPHOSClusterizerv1 not initialized " << endl ;
1065 }
1066 //____________________________________________________________________________
1067 void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option)
1068 {
1069   // Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
1070
1071   cout << "AliPHOSClusterizerv1: " << endl ;
1072   cout << "       Found "<< fEmcRecPoints->GetEntriesFast() << " EMC Rec Points and " 
1073            << fCpvRecPoints->GetEntriesFast() << " CPV RecPoints" << endl ;
1074
1075   if(strstr(option,"all")) {
1076     cout << "EMC clusters " << endl ;
1077     cout << "  Index    " 
1078          << "  Ene(MeV) " 
1079          << "  Multi    "
1080          << "  Module   "  
1081          << "    X      "
1082          << "    Y      "
1083          << "    Z      "
1084          << " Lambda 1  "
1085          << " Lambda 2  "
1086          << " MaxEnergy "
1087          << " # of prim "
1088          << " Primaries list "      <<  endl;      
1089     
1090     Int_t index ;
1091     for (index = 0 ; index < fEmcRecPoints->GetEntries() ; index++) {
1092       AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )fEmcRecPoints->At(index) ; 
1093
1094       cout << setw(6) << rp->GetIndexInList() << "     ";
1095       cout << setw(6) << rp->GetEnergy()      << "     ";
1096       cout << setw(6) << rp->GetMultiplicity()<< "     ";
1097       cout << setw(6) << rp->GetPHOSMod()     << "     ";
1098
1099       TVector3  locpos;  
1100       rp->GetLocalPosition(locpos);
1101       cout << setw(8) <<  locpos.X()          << "     ";
1102       cout << setw(8) <<  locpos.Y()          << "     ";
1103       cout << setw(8) <<  locpos.Z()          << "     ";
1104
1105       Float_t lambda[2]; 
1106       rp->GetElipsAxis(lambda);
1107       cout << setw(10)<<  lambda[0]           << "     ";
1108       cout << setw(10)<<  lambda[1]           << "     ";
1109       
1110       
1111       Int_t * primaries; 
1112       Int_t nprimaries;
1113       primaries = rp->GetPrimaries(nprimaries);
1114       cout << setw(8) <<    primaries         << "     ";
1115
1116       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1117         cout << setw(4)  <<  primaries[iprimary] << " ";
1118       cout << endl;      
1119     }
1120     cout << endl ;
1121
1122     //Now plot PCV/PPSD recPoints
1123     cout << "EMC clusters " << endl ;
1124     cout << "  Index    " 
1125          << "  Multi    "
1126          << "  Module   "  
1127          << "  Layer    "  
1128          << "    X      "
1129          << "    Y      "
1130          << "    Z      "
1131          << " # of prim "
1132          << " Primaries list "      <<  endl;      
1133     
1134     for (index = 0 ; index < fEmcRecPoints->GetEntries() ; index++) {
1135       AliPHOSRecPoint * rp = (AliPHOSRecPoint * )fCpvRecPoints->At(index) ; 
1136       cout << setw(6) << rp->GetIndexInList() << "     ";
1137       cout << setw(6) << rp->GetPHOSMod()     << "     ";
1138
1139       if( (strcmp(rp->ClassName() , "AliPHOSPpsdRecPoint" )) == 0){
1140         AliPHOSPpsdRecPoint * ppsd = (AliPHOSPpsdRecPoint*) rp ;
1141         if(ppsd->GetUp())
1142           cout <<"        CPV     ";
1143         else
1144           cout <<"       PPSD     ";
1145       }
1146       else
1147         cout <<"        CPV     ";
1148       
1149       TVector3  locpos;  
1150       rp->GetLocalPosition(locpos);
1151       cout << setw(8) <<  locpos.X()          << "     ";
1152       cout << setw(8) <<  locpos.Y()          << "     ";
1153       cout << setw(8) <<  locpos.Z()          << "     ";
1154       
1155       Int_t * primaries; 
1156       Int_t nprimaries;
1157       primaries = rp->GetPrimaries(nprimaries);
1158       cout << setw(8) <<    primaries         << "     ";
1159
1160       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1161         cout << setw(4)  <<  primaries[iprimary] << " ";
1162       cout << endl;      
1163     }
1164
1165
1166     cout << "-------------------------------------------------"<<endl ;
1167   }
1168 }
1169