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