]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSPIDv1.cxx
e1f9c4adcf3517b046248ffd808016327c1a02e3
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPIDv1.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 //_________________________________________________________________________
19 // Implementation version v1 of the PHOS particle identifier 
20 // Particle identification based on the 
21 //     - CPV information, 
22 //     - Preshower information (in MIXT or GPS2 geometries)
23 //     - shower width.
24 //
25 // CPV or Preshower clusters should be closer in PHOS plane than fCpvEmcDistance (in cm).
26 // This parameter can be set by method SetCpvtoEmcDistanceCut(Float_t cut)  
27 //
28 // One can set desirable ID method by the function SetIdentificationMethod(option).
29 // Presently the following options can be used together or separately :
30 //     - "disp": use dispersion cut on shower width 
31 //               (width can be set by method SetDispersionCut(Float_t cut)
32 //     - "ell" : use cut on the axis of the ellipse, drawn around shower 
33 //       (this cut can be changed by SetShowerProfileCut(char* formula), 
34 //        where formula - any function of two variables f(lambda[0],lambda[1]).
35 //        Shower is considered as EM if f() > 0 )
36 // One can visualize current cuts calling method PlotDispersionCuts().    
37 //
38 // use case:
39 //  root [0] AliPHOSPIDv1 * p1 = new AliPHOSPIDv1("galice.root")
40 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
41 //  root [1] p1->SetIdentificationMethod("disp ellipse")
42 //  root [2] p1->ExecuteTask()
43 //  root [3] AliPHOSPIDv1 * p2 = new AliPHOSPIDv1("galice1.root","ts1")
44 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
45 //                // reading headers from file galice1.root and TrackSegments 
46 //                // with title "ts1"
47 //  root [4] p2->SetRecParticlesBranch("rp1")
48 //                // set file name for the branch RecParticles
49 //  root [5] p2->ExecuteTask("deb all time")
50 //                // available options
51 //                // "deb" - prints # of reconstructed particles
52 //                // "deb all" -  prints # and list of RecParticles
53 //                // "time" - prints benchmarking results
54 //                  
55 //*-- Author: Yves Schutz (SUBATECH)  & Gines Martinez (SUBATECH) & 
56 //            Dmitri Peressounko (SUBATECH & Kurchatov Institute)
57 //            Completely redesined by Dmitri Peressounko, March 2001
58
59 // --- ROOT system ---
60 #include "TROOT.h"
61 #include "TTree.h"
62 #include "TFile.h"
63 #include "TF2.h"
64 #include "TFormula.h"
65 #include "TCanvas.h"
66 #include "TFolder.h"
67 #include "TSystem.h"
68 #include "TBenchmark.h"
69 // --- Standard library ---
70
71 #include <iostream.h>
72 #include <iomanip.h>
73
74 // --- AliRoot header files ---
75
76 #include "AliRun.h"
77 #include "AliGenerator.h"
78 #include "AliPHOS.h"
79 #include "AliPHOSPIDv1.h"
80 #include "AliPHOSClusterizerv1.h"
81 #include "AliPHOSTrackSegment.h"
82 #include "AliPHOSTrackSegmentMakerv1.h"
83 #include "AliPHOSRecParticle.h"
84 #include "AliPHOSGeometry.h"
85 #include "AliPHOSGetter.h"
86
87 ClassImp( AliPHOSPIDv1) 
88
89 //____________________________________________________________________________
90 AliPHOSPIDv1::AliPHOSPIDv1():AliPHOSPID()
91
92   // default ctor
93   fFormula           = 0 ;
94   fDispersion        = 0. ; 
95   fCpvEmcDistance    = 0 ; 
96   fTimeGate          = 2.e-9 ;
97   fHeaderFileName    = "" ; 
98   fTrackSegmentsTitle= "" ; 
99   fRecPointsTitle    = "" ; 
100   fRecParticlesTitle = "" ; 
101   fIDOptions         = "dis time" ; 
102   fRecParticlesInRun = 0 ;
103   fEmcRecPoints = 0;
104   fCpvRecPoints = 0;
105   fTrackSegments = 0;
106   fRecParticles = 0;
107   fClusterizer = 0;
108   fTSMaker = 0;
109 }
110
111 //____________________________________________________________________________
112 AliPHOSPIDv1::AliPHOSPIDv1(const char * headerFile,const char * name) : AliPHOSPID(headerFile, name)
113
114   //ctor with the indication on where to look for the track segments
115
116   fFormula        = new TFormula("LambdaCuts","(x>1)*(x<2.5)*(y>0)*(y<x)") ;   
117   fDispersion     = 2.0 ; 
118   fCpvEmcDistance = 3.0 ;
119   fTimeGate          = 2.e-9 ;
120  
121   fHeaderFileName     = GetTitle() ; 
122   fTrackSegmentsTitle = GetName() ; 
123   fRecPointsTitle     = GetName() ; 
124   fRecParticlesTitle  = GetName() ; 
125   fIDOptions          = "dis time" ;
126     
127   TString tempo(GetName()) ; 
128   tempo.Append(":") ;
129   tempo.Append(Version()) ; 
130   SetName(tempo) ; 
131   fRecParticlesInRun = 0 ; 
132
133   Init() ;
134
135 }
136
137 //____________________________________________________________________________
138 AliPHOSPIDv1::~AliPHOSPIDv1()
139
140   if(fTrackSegments) fTrackSegments->Delete() ; 
141   if(fEmcRecPoints) fEmcRecPoints->Delete() ; 
142   if(fCpvRecPoints) fCpvRecPoints->Delete() ;
143   if(fRecParticles) fRecParticles->Delete() ;
144 }
145
146
147 //____________________________________________________________________________
148 Float_t  AliPHOSPIDv1::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t *  Axis)const
149 {
150   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
151  
152   const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ; 
153   TVector3 vecEmc ;
154   TVector3 vecCpv ;
155   
156   emc->GetLocalPosition(vecEmc) ;
157   cpv->GetLocalPosition(vecCpv) ; 
158   if(emc->GetPHOSMod() == cpv->GetPHOSMod()){ 
159     
160     // Correct to difference in CPV and EMC position due to different distance to center.
161     // we assume, that particle moves from center
162     Float_t dCPV = geom->GetIPtoOuterCoverDistance();
163     Float_t dEMC = geom->GetIPtoCrystalSurface() ;
164     dEMC         = dEMC / dCPV ;
165     vecCpv = dEMC * vecCpv  - vecEmc ; 
166     if (Axis == "X") return vecCpv.X();
167     if (Axis == "Y") return vecCpv.Y();
168     if (Axis == "Z") return vecCpv.Z();
169     if (Axis == "R") return vecCpv.Mag();
170   } 
171  
172   return 100000000 ;
173 }
174
175 //____________________________________________________________________________
176 void  AliPHOSPIDv1::Exec(Option_t * option) 
177 {
178   //Steering method
179   
180   if( strcmp(GetName(), "")== 0 ) 
181     Init() ;
182   
183   if(strstr(option,"tim"))
184     gBenchmark->Start("PHOSPID");
185   
186   if(strstr(option,"print")) {
187     Print("") ; 
188     return ; 
189   }
190
191   gAlice->GetEvent(0) ;
192   //check, if the branch with name of this" already exits?
193   TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
194   TIter next(lob) ; 
195   TBranch * branch = 0 ;  
196   Bool_t phospidfound = kFALSE, pidfound = kFALSE ; 
197   
198   TString taskName(GetName()) ; 
199   taskName.Remove(taskName.Index(Version())-1) ;
200
201   while ( (branch = (TBranch*)next()) && (!phospidfound || !pidfound) ) {
202     if ( (strcmp(branch->GetName(), "PHOSPID")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) 
203       phospidfound = kTRUE ;
204     
205     else if ( (strcmp(branch->GetName(), "AliPHOSPID")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) 
206       pidfound = kTRUE ; 
207   }
208
209   if ( phospidfound || pidfound ) {
210     cerr << "WARNING: AliPHOSPIDv1::Exec -> RecParticles and/or PIDtMaker branch with name " 
211          << taskName.Data() << " already exits" << endl ;
212     return ; 
213   }       
214   
215   Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
216   Int_t ievent ;
217   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
218   
219   for(ievent = 0; ievent < nevents; ievent++){
220     gime->Event(ievent,"R") ;
221     
222     MakeRecParticles() ;
223     
224     WriteRecParticles(ievent);
225     
226     if(strstr(option,"deb"))
227       PrintRecParticles(option) ;
228   }
229   
230   if(strstr(option,"tim")){
231     gBenchmark->Stop("PHOSPID");
232     cout << "AliPHOSPID:" << endl ;
233     cout << "  took " << gBenchmark->GetCpuTime("PHOSPID") << " seconds for PID " 
234          <<  gBenchmark->GetCpuTime("PHOSPID")/nevents << " seconds per event " << endl ;
235     cout << endl ;
236   }
237   
238 }
239 //____________________________________________________________________________
240 void AliPHOSPIDv1::Init()
241 {
242   // Make all memory allocations that are not possible in default constructor
243   // Add the PID task to the list of PHOS tasks
244   
245   if ( strcmp(GetTitle(), "") == 0 )
246     SetTitle("galice.root") ;
247   
248   TString taskName(GetName()) ; 
249   taskName.Remove(taskName.Index(Version())-1) ;
250   
251   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), taskName.Data()) ; 
252   if ( gime == 0 ) {
253     cerr << "ERROR: AliPHOSPIDv1::Init -> Could not obtain the Getter object !" << endl ; 
254     return ;
255   } 
256    
257   gime->PostPID(this) ;
258   // create a folder on the white board //YSAlice/WhiteBoard/RecParticles/PHOS/recparticlesName
259   gime->PostRecParticles(taskName.Data() ) ; 
260   
261 }
262
263 //____________________________________________________________________________
264 void  AliPHOSPIDv1::MakeRecParticles(){
265
266   // Makes a RecParticle out of a TrackSegment
267   TString taskName(GetName()) ; 
268   taskName.Remove(taskName.Index(Version())-1) ;
269
270   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
271   TObjArray * emcRecPoints = gime->EmcRecPoints(taskName) ; 
272   TObjArray * cpvRecPoints = gime->CpvRecPoints(taskName) ; 
273   TClonesArray * trackSegments = gime->TrackSegments(taskName) ; 
274   TClonesArray * recParticles  = gime->RecParticles(taskName) ; 
275   recParticles->Clear();
276   
277   TIter next(trackSegments) ; 
278   AliPHOSTrackSegment * ts ; 
279   Int_t index = 0 ; 
280   AliPHOSRecParticle * rp ; 
281   
282   Bool_t ellips = fIDOptions.Contains("ell",TString::kIgnoreCase ) ;
283   Bool_t disp   = fIDOptions.Contains("dis",TString::kIgnoreCase ) ;
284   Bool_t time   = fIDOptions.Contains("tim",TString::kIgnoreCase ) ;
285   
286   while ( (ts = (AliPHOSTrackSegment *)next()) ) {
287     
288     new( (*recParticles)[index] ) AliPHOSRecParticle() ;
289     rp = (AliPHOSRecParticle *)recParticles->At(index) ; 
290     rp->SetTraskSegment(index) ;
291     rp->SetIndexInList(index) ;
292     
293     AliPHOSEmcRecPoint * emc = 0 ;
294     if(ts->GetEmcIndex()>=0)
295       emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
296     
297     AliPHOSRecPoint    * cpv = 0 ;
298     if(ts->GetCpvIndex()>=0)
299       cpv = (AliPHOSRecPoint *)   cpvRecPoints->At(ts->GetCpvIndex()) ;
300     
301     //set momentum and energy first
302     Float_t    e = emc->GetEnergy() ;
303     TVector3 dir = GetMomentumDirection(emc,cpv) ; 
304     dir.SetMag(e) ;
305
306     rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),e) ;
307     rp->SetCalcMass(0);
308
309     //now set type (reconstructed) of the particle    
310     Int_t showerprofile = 0;  // 0 narrow and 1 wide
311     
312     if(ellips){
313       Float_t lambda[2] ;
314       emc->GetElipsAxis(lambda) ;
315       if(fFormula->Eval(lambda[0],lambda[1]) <= 0 )
316         showerprofile = 1 ;  // not narrow
317     }
318     
319     if(disp)
320       if(emc->GetDispersion() > fDispersion )
321         showerprofile = 1 ;  // not narrow
322     
323     Int_t slow = 0 ;
324     if(time)
325       if(emc->GetTime() > fTimeGate )
326         slow = 0 ; 
327         
328     // Looking at the CPV detector
329     Int_t cpvdetector= 0 ;  //1 hit and 0 no hit     
330     if(cpv)
331       if(GetDistance(emc, cpv,  "R") < fCpvEmcDistance) 
332         cpvdetector = 1 ;  
333     
334     Int_t type = showerprofile + 2 * slow  + 4 * cpvdetector ;
335     rp->SetType(type) ; 
336     index++ ; 
337   }
338   
339 }
340
341 //____________________________________________________________________________
342 void  AliPHOSPIDv1:: Print(Option_t * option) const
343 {
344   // Print the parameters used for the particle type identification
345     cout <<  "=============== AliPHOSPID1 ================" << endl ;
346     cout <<  "Making PID "<< endl ;
347     cout <<  "    Headers file:               " << fHeaderFileName.Data() << endl ;
348     cout <<  "    RecPoints branch title:     " << fRecPointsTitle.Data() << endl ;
349     cout <<  "    TrackSegments Branch title: " << fTrackSegmentsTitle.Data() << endl ;
350     cout <<  "    RecParticles Branch title   " << fRecParticlesTitle.Data() << endl;
351     cout <<  "with parameters: " << endl ;
352     cout <<  "    Maximal EMC - CPV  distance (cm) " << fCpvEmcDistance << endl ;
353     if(fIDOptions.Contains("dis",TString::kIgnoreCase ))
354       cout <<  "                    dispersion cut " << fDispersion << endl ;
355     if(fIDOptions.Contains("ell",TString::kIgnoreCase )){
356       cout << "             Eliptic cuts function: " << endl ;
357       cout << fFormula->GetTitle() << endl ;
358     }
359     if(fIDOptions.Contains("tim",TString::kIgnoreCase ))
360       cout << "             Time Gate uzed: " << fTimeGate <<  endl ;
361     cout <<  "============================================" << endl ;
362 }
363
364 //____________________________________________________________________________
365 void  AliPHOSPIDv1::SetShowerProfileCut(char * formula)
366 {
367   //set shape of the cut on the axis of ellipce, drown around shouer
368   //shower considered "narrow" if Formula(lambda[0],lambda[1]) > 0.
369   if(fFormula) 
370     delete fFormula; 
371   fFormula = new TFormula("Lambda Cut",formula) ;
372 }
373 //____________________________________________________________________________
374 void  AliPHOSPIDv1::WriteRecParticles(Int_t event)
375 {
376  
377   AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ; 
378   TString taskName(GetName()) ; 
379   taskName.Remove(taskName.Index(Version())-1) ;
380   TClonesArray * recParticles = gime->RecParticles(taskName) ; 
381   recParticles->Expand(recParticles->GetEntriesFast() ) ;
382
383   //Make branch in TreeR for RecParticles 
384   char * filename = 0;
385   if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){   //generating file name
386     filename = new char[strlen(gAlice->GetBaseFile())+20] ;
387     sprintf(filename,"%s/PHOS.Reco.root",gAlice->GetBaseFile()) ; 
388   }
389   
390   TDirectory *cwd = gDirectory;
391   
392   //First rp
393   Int_t bufferSize = 32000 ;    
394   TBranch * rpBranch = gAlice->TreeR()->Branch("PHOSRP",&recParticles,bufferSize);
395   rpBranch->SetTitle(fRecParticlesTitle);
396   if (filename) {
397     rpBranch->SetFile(filename);
398     TIter next( rpBranch->GetListOfBranches());
399     TBranch * sb ;
400     while ((sb=(TBranch*)next())) {
401       sb->SetFile(filename);
402     }   
403     cwd->cd();
404   }
405   
406   //second, pid
407   Int_t splitlevel = 0 ; 
408   AliPHOSPIDv1 * pid = this ;
409   TBranch * pidBranch = gAlice->TreeR()->Branch("AliPHOSPID","AliPHOSPIDv1",&pid,bufferSize,splitlevel);
410   pidBranch->SetTitle(fRecParticlesTitle.Data());
411   if (filename) {
412     pidBranch->SetFile(filename);
413     TIter next( pidBranch->GetListOfBranches());
414     TBranch * sb ;
415     while ((sb=(TBranch*)next())) {
416       sb->SetFile(filename);
417     }   
418     cwd->cd();
419   }    
420   
421   rpBranch->Fill() ;
422   pidBranch->Fill() ;
423   
424   gAlice->TreeR()->Write(0,kOverwrite) ;  
425   
426 }
427
428 //____________________________________________________________________________
429 void  AliPHOSPIDv1::PlotDispersionCuts()const
430 {
431   // produces a plot of the dispersion cut
432   TCanvas*  lambdas = new TCanvas("lambdas","Cuts on the ellipse axis",200,10,700,500);
433  
434   if(fIDOptions.Contains("ell",TString::kIgnoreCase ) ){
435     TF2 * ell = new TF2("Elliptic Cuts",fFormula->GetName(),0,3,0,3) ;
436     ell->SetMinimum(0.0000001) ;
437     ell->SetMaximum(0.001) ;
438     ell->SetLineStyle(1) ;
439     ell->SetLineWidth(2) ;
440     ell->Draw() ;
441   }
442   
443   if( fIDOptions.Contains("dis",TString::kIgnoreCase ) ){
444     TF2 * dsp = new TF2("dispersion","(y<x)*(x*x+y*y < [0]*[0])",0,3,0,3) ;
445     dsp->SetParameter(0,fDispersion) ;
446     dsp->SetMinimum(0.0000001) ;
447     dsp->SetMaximum(0.001) ;
448     dsp->SetLineStyle(1) ;
449     dsp->SetLineColor(2) ;
450     dsp->SetLineWidth(2) ;
451     dsp->SetNpx(200) ;
452     dsp->SetNpy(200) ;
453     if(fIDOptions.Contains("ell",TString::kIgnoreCase ) )
454       dsp->Draw("same") ;
455     else
456       dsp->Draw() ;
457   }
458   lambdas->Update();
459 }
460
461 //____________________________________________________________________________
462 TVector3 AliPHOSPIDv1::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * cpv)const 
463
464   // Calculates the momentum direction:
465   //   1. if only a EMC RecPoint, direction is given by IP and this RecPoint
466   //   2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints 
467   //  However because of the poor position resolution of PPSD the direction is always taken as if we were 
468   //  in case 1.
469
470   TVector3 dir(0,0,0) ; 
471   
472   TVector3 emcglobalpos ;
473   TMatrix  dummy ;
474   
475   emc->GetGlobalPosition(emcglobalpos, dummy) ;
476   
477  
478   // The following commented code becomes valid once the PPSD provides 
479   // a reasonable position resolution, at least as good as EMC ! 
480   //   TVector3 ppsdlglobalpos ;
481   //   TVector3 ppsduglobalpos ;
482   //   if( fPpsdLowRecPoint ){ // certainly a photon that has concerted
483   //     fPpsdLowRecPoint->GetGlobalPosition(ppsdlglobalpos, mdummy) ; 
484   //     dir = emcglobalpos -  ppsdlglobalpos ; 
485   //     if( fPpsdUpRecPoint ){ // not looks like a charged       
486   //        fPpsdUpRecPoint->GetGlobalPosition(ppsduglobalpos, mdummy) ; 
487   //        dir = ( dir +  emcglobalpos -  ppsduglobalpos ) * 0.5 ; 
488   //      }
489   //   }
490   //   else { // looks like a neutral
491   //    dir = emcglobalpos ;  
492   //  }
493   
494   dir = emcglobalpos ;  
495   dir.SetZ( -dir.Z() ) ;   // why ?  
496   dir.SetMag(1.) ;
497
498   //account correction to the position of IP
499   Float_t xo,yo,zo ; //Coordinates of the origin
500   gAlice->Generator()->GetOrigin(xo,yo,zo) ;
501   TVector3 origin(xo,yo,zo);
502   dir = dir - origin ;
503
504   return dir ;  
505 }
506 //____________________________________________________________________________
507 void AliPHOSPIDv1::PrintRecParticles(Option_t * option)
508 {
509   // Print table of reconstructed particles
510
511   AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ; 
512
513   TString taskName(GetName()) ; 
514   taskName.Remove(taskName.Index(Version())-1) ;
515   TClonesArray * recParticles = gime->RecParticles(taskName) ; 
516   
517   cout << "AliPHOSPIDv1: event "<<gAlice->GetEvNumber()  << endl ;
518   cout << "       found " << recParticles->GetEntriesFast() << " RecParticles " << endl ;
519
520   fRecParticlesInRun += recParticles->GetEntriesFast() ; 
521   
522   if(strstr(option,"all")) {  // printing found TS
523     
524     cout << "  PARTICLE "   
525          << "  Index    "  << endl ;
526       //         << "  X        "     
527       //         << "  Y        " 
528       //         << "  Z        "    
529       //         << " # of primaries "          
530       //         << " Primaries list "    <<  endl;      
531     
532     Int_t index ;
533     for (index = 0 ; index < recParticles->GetEntries() ; index++) {
534       AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;       
535       
536       Text_t particle[11];
537       switch(rp->GetType()) {
538       case  AliPHOSFastRecParticle::kNEUTRALEMFAST:
539         strcpy( particle, "NEUTRAL EM FAST");
540         break;
541       case  AliPHOSFastRecParticle::kNEUTRALHAFAST:
542         strcpy(particle, "NEUTRAL HA FAST");
543         break;
544       case  AliPHOSFastRecParticle::kNEUTRALEMSLOW:
545         strcpy(particle, "NEUTRAL EM SLOW");
546         break ;
547       case  AliPHOSFastRecParticle::kNEUTRALHASLOW: 
548         strcpy(particle, "NEUTRAL HA SLOW");
549         break ;
550       case  AliPHOSFastRecParticle::kCHARGEDEMFAST:
551         strcpy(particle, "CHARGED EM FAST") ;
552         break ;
553       case  AliPHOSFastRecParticle::kCHARGEDHAFAST:
554         strcpy(particle, "CHARGED HA FAST") ;
555         break ; 
556       case  AliPHOSFastRecParticle::kCHARGEDEMSLOW:
557         strcpy(particle, "CHARGEDEMSLOW") ;
558         break ;
559       case  AliPHOSFastRecParticle::kCHARGEDHASLOW:
560         strcpy(particle, "CHARGED HA SLOW") ;
561         break ; 
562       }
563       
564       //    Int_t * primaries; 
565       //    Int_t nprimaries;
566       //    primaries = rp->GetPrimaries(nprimaries);
567       
568       cout << setw(10) << particle << "  "
569            << setw(5) <<  rp->GetIndexInList() << " "  ;
570         //         << setw(4) <<  nprimaries << "  ";
571         //      for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
572         //      cout << setw(4)  <<  primaries[iprimary] << " ";
573       cout << endl;      
574     }
575     cout << "-------------------------------------------" << endl ;
576   }
577   
578 }
579
580
581