]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSClusterizerv1.cxx
add a data member to count the number of produced objects by run
[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   fRecPointsInRun          = 0 ; 
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   fRecPointsInRun          = 0 ; 
139
140   Init() ;
141
142 }
143
144 //____________________________________________________________________________
145 void AliPHOSClusterizerv1::Exec(Option_t * option)
146 {
147   // Steering method
148
149   if( strcmp(GetName(), "")== 0 ) 
150     Init() ;
151
152   if(strstr(option,"tim"))
153     gBenchmark->Start("PHOSClusterizer"); 
154   
155   if(strstr(option,"print"))
156     Print("") ; 
157
158  //check, if the branch with name of this" already exits?
159   TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
160   TIter next(lob) ; 
161   TBranch * branch = 0 ;  
162   Bool_t phosemcfound = kFALSE, phoscpvfound = kFALSE, clusterizerfound = kFALSE ; 
163   
164   TString taskName(GetName()) ; 
165   taskName.ReplaceAll(Version(), "") ;
166
167   while ( (branch = (TBranch*)next()) && (!phosemcfound || !phoscpvfound || !clusterizerfound) ) {
168     if ( (strcmp(branch->GetName(), "PHOSEmcRP")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) 
169       phosemcfound = kTRUE ;
170     
171     else if ( (strcmp(branch->GetName(), "PHOSCpvRP")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) 
172       phoscpvfound = kTRUE ;
173    
174     else if ( (strcmp(branch->GetName(), "AliPHOSClusterizer")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) 
175       clusterizerfound = kTRUE ; 
176   }
177
178   if ( phoscpvfound || phosemcfound || clusterizerfound ) {
179     cerr << "WARNING: AliPHOSClusterizer::Exec -> Emc(Cpv)RecPoints and/or Clusterizer branch with name " 
180          << taskName.Data() << " already exits" << endl ;
181     return ; 
182   }       
183
184   Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
185   Int_t ievent ;
186
187   for(ievent = 0; ievent < nevents; ievent++){
188     
189     gAlice->GetEvent(ievent) ;
190     gAlice->SetEvent(ievent) ;
191
192     if(!ReadDigits(ievent))  //reads digits for event fEvent
193       continue;
194     
195     MakeClusters() ;
196     
197     if(fToUnfold) 
198       MakeUnfolding() ;
199
200     WriteRecPoints(ievent) ;
201
202     if(strstr(option,"deb"))
203       PrintRecPoints(option) ;
204   }
205   
206   if(strstr(option,"tim")){
207     gBenchmark->Stop("PHOSClusterizer");
208     cout << "AliPHOSClusterizer:" << endl ;
209     cout << "  took " << gBenchmark->GetCpuTime("PHOSClusterizer") << " seconds for Clusterizing " 
210          <<  gBenchmark->GetCpuTime("PHOSClusterizer")/nevents << " seconds per event " << endl ;
211     cout << endl ;
212   }
213   
214 }
215
216 //____________________________________________________________________________
217 Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, int * maxAt, Float_t * maxAtEnergy,
218                                     Int_t nPar, Float_t * fitparameters) const
219
220   // Calls TMinuit to fit the energy distribution of a cluster with several maxima 
221   // The initial values for fitting procedure are set equal to the positions of local maxima.
222   // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
223
224   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
225   TClonesArray * digits = gime->Digits() ; 
226   
227
228   gMinuit->mncler();                     // Reset Minuit's list of paramters
229   gMinuit->SetPrintLevel(-1) ;           // No Printout
230   gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;  
231                                          // To set the address of the minimization function 
232
233   TList * toMinuit = new TList();
234   toMinuit->AddAt(emcRP,0) ;
235   toMinuit->AddAt(digits,1) ;
236   
237   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare
238
239   // filling initial values for fit parameters
240   AliPHOSDigit * digit ;
241
242   Int_t ierflg  = 0; 
243   Int_t index   = 0 ;
244   Int_t nDigits = (Int_t) nPar / 3 ;
245
246   Int_t iDigit ;
247
248   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
249
250   for(iDigit = 0; iDigit < nDigits; iDigit++){
251     digit = (AliPHOSDigit *) maxAt[iDigit]; 
252
253     Int_t relid[4] ;
254     Float_t x = 0.;
255     Float_t z = 0.;
256     geom->AbsToRelNumbering(digit->GetId(), relid) ;
257     geom->RelPosInModule(relid, x, z) ;
258
259     Float_t energy = maxAtEnergy[iDigit] ;
260
261     gMinuit->mnparm(index, "x",  x, 0.1, 0, 0, ierflg) ;
262     index++ ;   
263     if(ierflg != 0){ 
264       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : x = " << x << endl ;
265       return kFALSE;
266     }
267     gMinuit->mnparm(index, "z",  z, 0.1, 0, 0, ierflg) ;
268     index++ ;   
269     if(ierflg != 0){
270       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : z = " << z << endl ;
271       return kFALSE;
272     }
273     gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;
274     index++ ;   
275     if(ierflg != 0){
276       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : energy = " << energy << endl ;      
277       return kFALSE;
278     }
279   }
280
281   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
282                       //  depends on it. 
283   Double_t p1 = 1.0 ;
284   Double_t p2 = 0.0 ;
285
286   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
287   gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
288   gMinuit->SetMaxIterations(5);
289   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
290
291   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
292
293   if(ierflg == 4){  // Minimum not found   
294     cout << "PHOS Unfolding>  Fit not converged, cluster abandoned "<< endl ;      
295     return kFALSE ;
296   }            
297   for(index = 0; index < nPar; index++){
298     Double_t err ;
299     Double_t val ;
300     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index
301     fitparameters[index] = val ;
302    }
303
304   delete toMinuit ;
305   return kTRUE;
306
307 }
308
309 //____________________________________________________________________________
310 void AliPHOSClusterizerv1::Init()
311 {
312   // Make all memory allocations which can not be done in default constructor.
313   // Attach the Clusterizer task to the list of PHOS tasks
314   
315   if ( strcmp(GetTitle(), "") == 0 )
316     SetTitle("galice.root") ;
317
318   TString taskName(GetName()) ; 
319   taskName.ReplaceAll(Version(), "") ;
320
321   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), taskName.Data()) ; 
322   if ( gime == 0 ) {
323     cerr << "ERROR: AliPHOSClusterizerv1::Init -> Could not obtain the Getter object !" << endl ; 
324     return ;
325   } 
326     
327   if(!gMinuit) 
328     gMinuit = new TMinuit(100) ;
329
330   //add Task to //YSAlice/tasks/Reconstructioner/PHOS
331   TTask * aliceRe  = (TTask*)gROOT->FindObjectAny("YSAlice/tasks/Reconstructioner") ; 
332   TTask * phosRe   = (TTask*)aliceRe->GetListOfTasks()->FindObject("PHOS") ;
333   phosRe->Add(this) ; 
334   // create a folder on the white board //YSAlice/WhiteBoard/RecPoints/PHOS/recpointsName
335   gime->Post(GetTitle(), "R", taskName.Data() ) ; 
336   
337 }
338
339 //____________________________________________________________________________
340 Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
341 {
342   // Gives the neighbourness of two digits = 0 are not neighbour but continue searching 
343   //                                       = 1 are neighbour
344   //                                       = 2 are not neighbour but do not continue searching
345   // neighbours are defined as digits having at least a common vertex 
346   // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster 
347   //                                      which is compared to a digit (d2)  not yet in a cluster  
348
349   const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
350
351   Int_t rv = 0 ; 
352
353   Int_t relid1[4] ; 
354   geom->AbsToRelNumbering(d1->GetId(), relid1) ; 
355
356   Int_t relid2[4] ; 
357   geom->AbsToRelNumbering(d2->GetId(), relid2) ; 
358  
359   if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module and the same PPSD Module 
360     Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
361     Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
362     
363     if (( coldiff <= 1 )  && ( rowdiff <= 1 )){
364       rv = 1 ; 
365     }
366     else {
367       if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1)) 
368         rv = 2; //  Difference in row numbers is too large to look further 
369     }
370
371   } 
372   else {
373     
374     if( (relid1[0] < relid2[0]) || (relid1[1] < relid2[1]) )  
375       rv=2 ;
376
377   }
378
379   //Do NOT clusterize upper PPSD  
380   if( IsInPpsd(d1) && IsInPpsd(d2) &&
381      relid1[1] > 0                 &&
382      relid1[1] < geom->GetNumberOfPadsPhi()*geom->GetNumberOfPadsPhi() ) rv = 2 ;
383
384   return rv ; 
385 }
386
387
388 //____________________________________________________________________________
389 Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
390 {
391   // Tells if (true) or not (false) the digit is in a PHOS-EMC module
392  
393   Bool_t rv = kFALSE ; 
394   const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
395
396   Int_t relid[4] ; 
397   geom->AbsToRelNumbering(digit->GetId(), relid) ; 
398
399   if ( relid[1] == 0  ) rv = kTRUE; 
400
401   return rv ; 
402 }
403
404 //____________________________________________________________________________
405 Bool_t AliPHOSClusterizerv1::IsInPpsd(AliPHOSDigit * digit) const
406 {
407   // Tells if (true) or not (false) the digit is in a PHOS-PPSD module
408  
409   Bool_t rv = kFALSE ; 
410   const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
411
412   Int_t relid[4] ; 
413   geom->AbsToRelNumbering(digit->GetId(), relid) ; 
414
415   if ( relid[1] > 0 && relid[0] > geom->GetNCPVModules() ) rv = kTRUE; 
416
417   return rv ; 
418 }
419
420 //____________________________________________________________________________
421 Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
422 {
423   // Tells if (true) or not (false) the digit is in a PHOS-CPV module
424  
425   Bool_t rv = kFALSE ; 
426   const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
427
428   Int_t relid[4] ; 
429   geom->AbsToRelNumbering(digit->GetId(), relid) ; 
430
431   if ( relid[1] > 0 && relid[0] <= geom->GetNCPVModules() ) rv = kTRUE; 
432
433   return rv ; 
434 }
435 //____________________________________________________________________________
436 Bool_t AliPHOSClusterizerv1::ReadDigits(Int_t event)
437  {
438    // reads digitis with specified title from TreeD
439
440   fNumberOfEmcClusters  = 0 ;
441   fNumberOfCpvClusters  = 0 ;
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 position, 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   Int_t bufferSize = 32000 ;    
547   Int_t splitlevel = 0 ;
548
549   //First EMC
550   TBranch * emcBranch = gAlice->TreeR()->Branch("PHOSEmcRP","TObjArray",&emcRecPoints,bufferSize,splitlevel);
551   emcBranch->SetTitle(branchName);
552   if (filename) {
553     emcBranch->SetFile(filename);
554     TIter next( emcBranch->GetListOfBranches());
555     TBranch * sb ;
556     while ((sb=(TBranch*)next())) {
557       sb->SetFile(filename);
558     }   
559     
560     cwd->cd();
561   }
562     
563   //Now CPV branch
564   TBranch * cpvBranch = gAlice->TreeR()->Branch("PHOSCpvRP","TObjArray",&cpvRecPoints,bufferSize,splitlevel);
565   cpvBranch->SetTitle(branchName);
566   if (filename) {
567     cpvBranch->SetFile(filename);
568     TIter next( cpvBranch->GetListOfBranches());
569     TBranch * sb;
570     while ((sb=(TBranch*)next())) {
571       sb->SetFile(filename);
572     }   
573     cwd->cd();
574   } 
575     
576   //And Finally  clusterizer branch
577   AliPHOSClusterizerv1 * cl = (AliPHOSClusterizerv1*)gime->Clusterizer(GetName()) ;
578   TBranch * clusterizerBranch = gAlice->TreeR()->Branch("AliPHOSClusterizer","AliPHOSClusterizerv1",
579                                               &cl,bufferSize,splitlevel);
580   clusterizerBranch->SetTitle(branchName);
581   if (filename) {
582     clusterizerBranch->SetFile(filename);
583     TIter next( clusterizerBranch->GetListOfBranches());
584     TBranch * sb ;
585     while ((sb=(TBranch*)next())) {
586       sb->SetFile(filename);
587     }   
588     cwd->cd();
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->Delete() ;
609   cpvRecPoints->Delete() ;
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   const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ; 
946   fret = 0. ;     
947   Int_t iparam ;
948
949   if(iflag == 2)
950     for(iparam = 0 ; iparam < nPar ; iparam++)    
951       Grad[iparam] = 0 ; // Will evaluate gradient
952   
953   Double_t efit ;    
954
955   AliPHOSDigit * digit ;
956   Int_t iDigit ;
957
958   for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
959
960     digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ; 
961
962     Int_t relid[4] ;
963     Float_t xDigit ;
964     Float_t zDigit ;
965
966     geom->AbsToRelNumbering(digit->GetId(), relid) ;
967
968     geom->RelPosInModule(relid, xDigit, zDigit) ;
969
970      if(iflag == 2){  // calculate gradient
971        Int_t iParam = 0 ;
972        efit = 0 ;
973        while(iParam < nPar ){
974          Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
975          iParam++ ; 
976          distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ; 
977          distance = TMath::Sqrt( distance ) ; 
978          iParam++ ;      
979          efit += x[iParam] * ShowerShape(distance) ;
980          iParam++ ;
981        }
982        Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E) 
983        iParam = 0 ;
984        while(iParam < nPar ){
985          Double_t xpar = x[iParam] ;
986          Double_t zpar = x[iParam+1] ;
987          Double_t epar = x[iParam+2] ;
988          Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
989          Double_t shape = sum * ShowerShape(dr) ;
990          Double_t r4 = dr*dr*dr*dr ;
991          Double_t r295 = TMath::Power(dr,2.95) ;
992          Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
993                                          0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
994          
995          Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ;  // Derivative over x    
996          iParam++ ; 
997          Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ;  // Derivative over z         
998          iParam++ ; 
999          Grad[iParam] += shape ;                                  // Derivative over energy             
1000          iParam++ ; 
1001        }
1002      }
1003      efit = 0;
1004      iparam = 0 ;
1005
1006      while(iparam < nPar ){
1007        Double_t xpar = x[iparam] ;
1008        Double_t zpar = x[iparam+1] ;
1009        Double_t epar = x[iparam+2] ;
1010        iparam += 3 ;
1011        Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
1012        distance =  TMath::Sqrt(distance) ;
1013        efit += epar * ShowerShape(distance) ;
1014      }
1015
1016      fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ; 
1017      // Here we assume, that sigma = sqrt(E)
1018   }
1019
1020 }
1021
1022 //____________________________________________________________________________
1023 void AliPHOSClusterizerv1::Print(Option_t * option)const
1024 {
1025   // Print clusterizer parameters
1026
1027   if( strcmp(GetName(), "") !=0 ){
1028     
1029     // Print parameters
1030  
1031     TString taskName(GetName()) ; 
1032     taskName.ReplaceAll(Version(), "") ;
1033
1034     cout << "---------------"<< taskName.Data() << " " << GetTitle()<< "-----------" << endl 
1035          << "Clusterizing digits from the file: " << fHeaderFileName.Data() << endl 
1036          << "                           Branch: " << fDigitsBranchTitle.Data() << endl 
1037          << endl 
1038          << "                       EMC Clustering threshold = " << fEmcClusteringThreshold << endl
1039          << "                       EMC Local Maximum cut    = " << fEmcLocMaxCut << endl
1040          << "                       EMC Logarothmic weight   = " << fW0 << endl
1041          << endl
1042          << "                       CPV Clustering threshold = " << fCpvClusteringThreshold << endl
1043          << "                       CPV Local Maximum cut    = " << fCpvLocMaxCut << endl
1044        << "                       CPV Logarothmic weight   = " << fW0CPV << endl
1045          << endl
1046          << "                      PPSD  Clustering threshold = " << fPpsdClusteringThreshold << endl;
1047     if(fToUnfold)
1048       cout << " Unfolding on " << endl ;
1049     else
1050       cout << " Unfolding off " << endl ;
1051     
1052     cout << "------------------------------------------------------------------" <<endl ;
1053   }
1054   else
1055     cout << " AliPHOSClusterizerv1 not initialized " << endl ;
1056 }
1057 //____________________________________________________________________________
1058 void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option)
1059 {
1060   // Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
1061
1062   TObjArray * emcRecPoints = AliPHOSGetter::GetInstance()->EmcRecPoints() ; 
1063   TObjArray * cpvRecPoints = AliPHOSGetter::GetInstance()->CpvRecPoints() ; 
1064
1065   cout << "AliPHOSClusterizerv1: : event "<<gAlice->GetEvNumber() << endl ;
1066   cout << "       Found "<< emcRecPoints->GetEntriesFast() << " EMC Rec Points and " 
1067            << cpvRecPoints->GetEntriesFast() << " CPV RecPoints" << endl ;
1068
1069   fRecPointsInRun +=  emcRecPoints->GetEntriesFast() ; 
1070   fRecPointsInRun +=  cpvRecPoints->GetEntriesFast() ; 
1071
1072   if(strstr(option,"all")) {
1073     cout << "EMC clusters " << endl ;
1074     cout << "  Index    " 
1075          << "  Ene(MeV) " 
1076          << "  Multi    "
1077          << "  Module   "  
1078          << "    X      "
1079          << "    Y      "
1080          << "    Z      "
1081          << " Lambda 1  "
1082          << " Lambda 2  "
1083          << " MaxEnergy "
1084          << " # of prim "
1085          << " Primaries list "      <<  endl;      
1086     
1087     Int_t index ;
1088     for (index = 0 ; index < emcRecPoints->GetEntries() ; index++) {
1089       AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )emcRecPoints->At(index) ; 
1090
1091       cout << setw(6) << rp->GetIndexInList() << "     ";
1092       cout << setw(6) << rp->GetEnergy()      << "     ";
1093       cout << setw(6) << rp->GetMultiplicity()<< "     ";
1094       cout << setw(6) << rp->GetPHOSMod()     << "     ";
1095
1096       TVector3  locpos;  
1097       rp->GetLocalPosition(locpos);
1098       cout << setw(8) <<  locpos.X()          << "     ";
1099       cout << setw(8) <<  locpos.Y()          << "     ";
1100       cout << setw(8) <<  locpos.Z()          << "     ";
1101
1102       Float_t lambda[2]; 
1103       rp->GetElipsAxis(lambda);
1104       cout << setw(10)<<  lambda[0]           << "     ";
1105       cout << setw(10)<<  lambda[1]           << "     ";
1106       
1107       
1108       Int_t * primaries; 
1109       Int_t nprimaries;
1110       primaries = rp->GetPrimaries(nprimaries);
1111       cout << setw(8) <<    primaries         << "     ";
1112
1113       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1114         cout << setw(4)  <<  primaries[iprimary] << " ";
1115       cout << endl;      
1116     }
1117     cout << endl ;
1118
1119     //Now plot CPV/PPSD recPoints
1120     cout << "EMC clusters " << endl ;
1121     cout << "  Index    " 
1122          << "  Multi    "
1123          << "  Module   "  
1124          << "  Layer    "  
1125          << "    X      "
1126          << "    Y      "
1127          << "    Z      "
1128          << " # of prim "
1129          << " Primaries list "      <<  endl;      
1130     
1131     for (index = 0 ; index < cpvRecPoints->GetEntries() ; index++) {
1132       AliPHOSRecPoint * rp = (AliPHOSRecPoint * )cpvRecPoints->At(index) ; 
1133       cout << setw(6) << rp->GetIndexInList() << "     ";
1134       cout << setw(6) << rp->GetPHOSMod()     << "     ";
1135
1136       if( (strcmp(rp->ClassName() , "AliPHOSPpsdRecPoint" )) == 0){
1137         AliPHOSPpsdRecPoint * ppsd = (AliPHOSPpsdRecPoint*) rp ;
1138         if(ppsd->GetUp())
1139           cout <<"        CPV     ";
1140         else
1141           cout <<"       PPSD     ";
1142       }
1143       else
1144         cout <<"        CPV     ";
1145       
1146       TVector3  locpos;  
1147       rp->GetLocalPosition(locpos);
1148       cout << setw(8) <<  locpos.X()          << "     ";
1149       cout << setw(8) <<  locpos.Y()          << "     ";
1150       cout << setw(8) <<  locpos.Z()          << "     ";
1151       
1152       Int_t * primaries; 
1153       Int_t nprimaries;
1154       primaries = rp->GetPrimaries(nprimaries);
1155       cout << setw(8) <<    primaries         << "     ";
1156
1157       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1158         cout << setw(4)  <<  primaries[iprimary] << " ";
1159       cout << endl;      
1160     }
1161
1162
1163     cout << "-------------------------------------------------"<<endl ;
1164   }
1165 }
1166