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