21a52a8a237fb50c0b415c81c7af982235e34886
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPIDv2.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 v2 of the PHOS particle identifier 
20 // Particle identification based on the 
21 //     - RCPV: distance from CPV recpoint to EMCA recpoint.
22 //     - TOF 
23 //     - PCA: Principal Components Analisis..
24 // The identified particle has an identification number corresponding 
25 // to a 3 bits number:
26 //     -Bit 0: bit set if RCPV > fCpvEmcDistance 
27 //     -Bit 1: bit set if TOF  < fTimeGate
28 //     -Bit 2: bit set if Principal Components are 
29 //      inside an ellipse defined by fX_center, fY_center, fA, fB, fAngle
30 //
31 //
32 // use case:
33 //  root [0] AliPHOSPIDv2 * p1 = new AliPHOSPIDv2("galice.root")
34 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
35 //  root [1] p1->SetIdentificationMethod("disp ellipse")
36 //  root [2] p1->ExecuteTask()
37 //  root [3] AliPHOSPIDv2 * p2 = new AliPHOSPIDv2("galice1.root","ts1")
38 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
39 //                // reading headers from file galice1.root and TrackSegments 
40 //                // with title "ts1"
41 //  root [4] p2->SetRecParticlesBranch("rp1")
42 //                // set file name for the branch RecParticles
43 //  root [5] p2->ExecuteTask("deb all time")
44 //                // available options
45 //                // "deb" - prints # of reconstructed particles
46 //                // "deb all" -  prints # and list of RecParticles
47 //                // "time" - prints benchmarking results
48 //                  
49 //*-- Author: Yves Schutz (SUBATECH)  & Gines Martinez (SUBATECH) & 
50 //            Gustavo Conesa April 2002
51
52 // --- ROOT system ---
53 #include "TROOT.h"
54 #include "TTree.h"
55 #include "TFile.h"
56 #include "TF2.h"
57 #include "TFormula.h"
58 #include "TCanvas.h"
59 #include "TFolder.h"
60 #include "TSystem.h"
61 #include "TBenchmark.h"
62 #include "TEllipse.h"
63 #include "TPrincipal.h"
64 // --- Standard library ---
65
66 #include <iostream.h>
67 #include <iomanip.h>
68
69 // --- AliRoot header files ---
70
71 #include "AliRun.h"
72 #include "AliGenerator.h"
73 #include "AliPHOS.h"
74 #include "AliPHOSPIDv2.h"
75 #include "AliPHOSClusterizerv1.h"
76 #include "AliPHOSTrackSegment.h"
77 #include "AliPHOSTrackSegmentMakerv1.h"
78 #include "AliPHOSRecParticle.h"
79 #include "AliPHOSGeometry.h"
80 #include "AliPHOSGetter.h"
81
82 ClassImp( AliPHOSPIDv2) 
83
84 //____________________________________________________________________________
85 AliPHOSPIDv2::AliPHOSPIDv2():AliPHOSPID()
86
87   // default ctor
88  
89   fHeaderFileName    = "" ; 
90   fTrackSegmentsTitle= "" ; 
91   fRecPointsTitle    = "" ; 
92   fRecParticlesTitle = "" ; 
93   fRecParticlesInRun = 0 ;
94   fClusterizer       = 0 ; 
95   fTSMaker           = 0 ;
96   
97 }
98
99 //____________________________________________________________________________
100 AliPHOSPIDv2::AliPHOSPIDv2(const char * headerFile,const char * name) : AliPHOSPID(headerFile, name)
101
102   //ctor with the indication on where to look for the track segments
103  
104   fHeaderFileName     = GetTitle() ; 
105   fTrackSegmentsTitle = GetName() ; 
106   fRecPointsTitle     = GetName() ; 
107   fRecParticlesTitle  = GetName() ;
108   TString tempo(GetName()) ; 
109   tempo.Append(":") ;
110   tempo.Append(Version()) ; 
111   SetName(tempo) ; 
112   fRecParticlesInRun = 0 ; 
113
114   Init() ;
115
116 }
117
118 //____________________________________________________________________________
119 AliPHOSPIDv2::~AliPHOSPIDv2()
120
121   delete [] fX;
122   delete [] fP; 
123 }
124
125
126 //____________________________________________________________________________
127 Float_t  AliPHOSPIDv2::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t *  Axis)const
128 {
129   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
130  
131   const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ; 
132   TVector3 vecEmc ;
133   TVector3 vecCpv ;
134   
135   emc->GetLocalPosition(vecEmc) ;
136   cpv->GetLocalPosition(vecCpv) ; 
137   if(emc->GetPHOSMod() == cpv->GetPHOSMod()){ 
138     
139     // Correct to difference in CPV and EMC position due to different distance to center.
140     // we assume, that particle moves from center
141     Float_t dCPV = geom->GetIPtoOuterCoverDistance();
142     Float_t dEMC = geom->GetIPtoCrystalSurface() ;
143     dEMC         = dEMC / dCPV ;
144     vecCpv = dEMC * vecCpv  - vecEmc ; 
145     if (Axis == "X") return vecCpv.X();
146     if (Axis == "Y") return vecCpv.Y();
147     if (Axis == "Z") return vecCpv.Z();
148     if (Axis == "R") return vecCpv.Mag();
149   } 
150  
151   return 100000000 ;
152 }
153 //____________________________________________________________________________
154  void  AliPHOSPIDv2::SetEllipseParameters(Float_t x, Float_t y,Float_t a, Float_t b,Float_t angle)
155 {
156   fX_center = x ;
157   fY_center = y ;
158   fA = a  ;
159   fB = b ;
160   fAngle = angle ;
161 }
162
163 //____________________________________________________________________________
164 Int_t  AliPHOSPIDv2::GetPrincipalSign(Double_t* P )const
165 {
166   //This method gives if the PCA of the particle are inside a defined ellipse
167  
168   Int_t      prinsign;
169   Double_t   fDx        = 0. ; 
170   Double_t   fDelta     = 0. ; 
171   Double_t   fY         = 0. ; 
172   Double_t   fY_1       = 0. ; 
173   Double_t   fY_2       = 0. ;
174   Double_t   fPi        = TMath::Pi() ;
175   Double_t   fCos_Theta = TMath::Cos(fPi*fAngle/180.) ;
176   Double_t   fSin_Theta = TMath::Sin(fPi*fAngle/180.) ;   
177
178   fDx = P[0] - fX_center ; 
179   fDelta = 4.*fA*fA*fA*fB* (fA*fA*fCos_Theta*fCos_Theta + fB*fB*fSin_Theta*fSin_Theta - fDx*fDx) ; 
180   if (fDelta < 0.) 
181     {prinsign=0;} 
182   
183   else if (fDelta == 0.) 
184     { 
185       fY = fCos_Theta*fSin_Theta*(fA*fA - fB*fB)*fDx / (fA*fA*fCos_Theta*fCos_Theta +
186                                                         fB*fB*fSin_Theta*fSin_Theta) ; 
187       fY += fY_center ; 
188       if(P[1]==fY ) 
189         {prinsign=1;} 
190       else 
191         {prinsign=0;} 
192     } 
193   else 
194     { 
195       fY_1 = (fCos_Theta*fSin_Theta*(fA*fA - fB*fB) *fDx +
196              TMath::Sqrt(fDelta)/2.)/(fA*fA*fCos_Theta*fCos_Theta + fB*fB*fSin_Theta*fSin_Theta) ; 
197       fY_2 = (fCos_Theta*fSin_Theta*(fA*fA - fB*fB) *fDx -
198              TMath::Sqrt(fDelta)/2.)/(fA*fA*fCos_Theta*fCos_Theta + fB*fB*fSin_Theta*fSin_Theta) ; 
199       fY_1 += fY_center ; 
200       fY_2 += fY_center ; 
201       if ((P[1]<=fY_1) && (P[1]>=fY_2)) 
202         {prinsign=1;} 
203       else 
204         {prinsign=0;}  
205     } 
206   return prinsign;
207 }
208 //____________________________________________________________________________
209 void  AliPHOSPIDv2::Exec(Option_t * option) 
210 {
211   //Steering method
212   
213   if( strcmp(GetName(), "")== 0 ) 
214     Init() ;
215   
216   if(strstr(option,"tim"))
217     gBenchmark->Start("PHOSPID");
218   
219   if(strstr(option,"print")) {
220     Print("") ; 
221     return ; 
222   }
223
224   gAlice->GetEvent(0) ;
225   //check, if the branch with name of this" already exits?
226   TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
227   TIter next(lob) ; 
228   TBranch * branch = 0 ;  
229   Bool_t phospidfound = kFALSE, pidfound = kFALSE ; 
230   
231   TString taskName(GetName()) ; 
232   taskName.Remove(taskName.Index(Version())-1) ;
233
234   while ( (branch = (TBranch*)next()) && (!phospidfound || !pidfound) ) {
235     if ( (strcmp(branch->GetName(), "PHOSPID")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) 
236       phospidfound = kTRUE ;
237     
238     else if ( (strcmp(branch->GetName(), "AliPHOSPID")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) 
239       pidfound = kTRUE ; 
240   }
241
242   if ( phospidfound || pidfound ) {
243     cerr << "WARNING: AliPHOSPIDv2::Exec -> RecParticles and/or PIDtMaker branch with name " 
244          << taskName.Data() << " already exits" << endl ;
245     return ; 
246   }       
247   
248   Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
249   Int_t ievent ;
250   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
251   
252   for(ievent = 0; ievent < nevents; ievent++){
253     gime->Event(ievent,"R") ;
254     
255     MakeRecParticles() ;
256     
257     WriteRecParticles(ievent);
258     
259     if(strstr(option,"deb"))
260       PrintRecParticles(option) ;
261
262     //increment the total number of rec particles per run 
263     fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ; 
264
265   }
266   
267   if(strstr(option,"tim")){
268     gBenchmark->Stop("PHOSPID");
269     cout << "AliPHOSPID:" << endl ;
270     cout << "  took " << gBenchmark->GetCpuTime("PHOSPID") << " seconds for PID " 
271          <<  gBenchmark->GetCpuTime("PHOSPID")/nevents << " seconds per event " << endl ;
272     cout << endl ;
273   }
274   
275 }
276 //____________________________________________________________________________
277 void AliPHOSPIDv2::Init()
278 {
279   // Make all memory allocations that are not possible in default constructor
280   // Add the PID task to the list of PHOS tasks
281   
282   if ( strcmp(GetTitle(), "") == 0 )
283     SetTitle("galice.root") ;
284   
285   TString taskName(GetName()) ; 
286   taskName.Remove(taskName.Index(Version())-1) ;
287
288   fCpvEmcDistance = 4.0 ; 
289   fTimeGate       = 0.162e-7 ;
290
291   //PCA 
292   fX              = new double[7]; // Data for the PCA 
293   fP              = new double[7]; // Eigenvalues of the PCA  
294   fFileName       = "$ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root" ; 
295   TFile f( fFileName.Data(), "read" ) ;
296   fPrincipal      = dynamic_cast<TPrincipal*> (f.Get("principal")) ; 
297  
298   // Ellipse parameters 
299   fX_center          = 2.0 ; 
300   fY_center          = -0.35 ; 
301
302   fA                 = 1.9 ; 
303   fB                 = 1.0 ; 
304   fAngle             = -60. ;
305  
306   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), taskName.Data()) ; 
307   if ( gime == 0 ) {
308     cerr << "ERROR: AliPHOSPIDv2::Init -> Could not obtain the Getter object !" << endl ; 
309     return ;
310   } 
311    
312   gime->PostPID(this) ;
313   // create a folder on the white board //YSAlice/WhiteBoard/RecParticles/PHOS/recparticlesName
314   gime->PostRecParticles(taskName.Data() ) ; 
315   
316 }
317
318 //____________________________________________________________________________
319 void  AliPHOSPIDv2::MakeRecParticles(){
320
321   // Makes a RecParticle out of a TrackSegment
322   TString taskName(GetName()) ; 
323   taskName.Remove(taskName.Index(Version())-1) ;
324
325   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
326   TObjArray * emcRecPoints = gime->EmcRecPoints(taskName) ; 
327   TObjArray * cpvRecPoints = gime->CpvRecPoints(taskName) ; 
328   TClonesArray * trackSegments = gime->TrackSegments(taskName) ; 
329   TClonesArray * recParticles  = gime->RecParticles(taskName) ; 
330   recParticles->Clear();
331  
332
333   TIter next(trackSegments) ; 
334   AliPHOSTrackSegment * ts ; 
335   Int_t index = 0 ; 
336   AliPHOSRecParticle * rp ; 
337   
338   while ( (ts = (AliPHOSTrackSegment *)next()) ) {
339     
340     new( (*recParticles)[index] ) AliPHOSRecParticle() ;
341     rp = (AliPHOSRecParticle *)recParticles->At(index) ; 
342     rp->SetTraskSegment(index) ;
343     rp->SetIndexInList(index) ;
344     
345     AliPHOSEmcRecPoint * emc = 0 ;
346     if(ts->GetEmcIndex()>=0)
347       emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
348     
349     AliPHOSRecPoint    * cpv = 0 ;
350     if(ts->GetCpvIndex()>=0)
351       cpv = (AliPHOSRecPoint *)   cpvRecPoints->At(ts->GetCpvIndex()) ;
352     
353     //set momentum and energy first
354     Float_t    e = emc->GetEnergy() ;
355     TVector3 dir = GetMomentumDirection(emc,cpv) ; 
356     dir.SetMag(e) ;
357
358     rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),e) ;
359     rp->SetCalcMass(0);
360     
361     //now set type (reconstructed) of the particle    
362     
363     // Looking at the CPV detector. If RCPV grater than fCpvEmcDistance, first bit 1.
364     if(cpv)
365       if(GetDistance(emc, cpv,  "R") > fCpvEmcDistance )  
366         rp->SetPIDBit(0);
367     
368     //Looking the TOF. If TOF smaller than gate, second bit 1.
369     if(emc->GetTime()< fTimeGate) 
370       rp->SetPIDBit(1);                             
371     
372     
373     //Loking PCA. Define and calculate the data (X) introduce in the function 
374     //X2P that gives the components (P).  
375     Float_t    fSpher = 0. ;
376     Float_t    fEmaxdtotal = 0. ; 
377     Float_t lambda[2] ;
378     
379     emc->GetElipsAxis(lambda) ;
380     if((lambda[0]+lambda[1])!=0) fSpher=fabs(lambda[0]-lambda[1])/(lambda[0]+lambda[1]); 
381     
382     fEmaxdtotal=emc->GetMaximalEnergy()/emc->GetEnergy(); 
383     
384     fX[0] = lambda[0]; //cout <<" lambda0 "<<fX[0]; 
385     fX[1] = lambda[1]; //cout <<" lambda1 "<<fX[1]; 
386     fX[2] = emc->GetDispersion(); // cout <<" disper "<<fX[2]; 
387     fX[3] = fSpher; // cout <<" spher "<<fX[3]; 
388     fX[4] = emc->GetMultiplicity(); // cout <<" mult "<<fX[4]; 
389     fX[5] = fEmaxdtotal; // cout <<" emax "<<fX[5]; 
390     fX[6] = emc->GetCoreEnergy(); //  cout <<" core "<<fX[5]<< endl ; 
391     
392     // cout<<"Principal "<<fPrincipal<<endl;
393     fPrincipal->X2P(fX,fP);
394     
395     //If we are inside the ellipse, third bit 1
396     if(GetPrincipalSign(fP)== 1) 
397       rp->SetPIDBit(2) ;
398     
399     
400     rp->SetProductionVertex(0,0,0,0);
401     rp->SetFirstMother(-1);
402     rp->SetLastMother(-1);
403     rp->SetFirstDaughter(-1);
404     rp->SetLastDaughter(-1);
405     rp->SetPolarisation(0,0,0);
406     index++ ; 
407   }
408   
409 }
410
411 //____________________________________________________________________________
412 void  AliPHOSPIDv2:: Print(Option_t * option) const
413 {
414   // Print the parameters used for the particle type identification
415     cout <<  "=============== AliPHOSPID1 ================" << endl ;
416     cout <<  "Making PID "<< endl ;
417     cout <<  "    Headers file:               " << fHeaderFileName.Data() << endl ;
418     cout <<  "    RecPoints branch title:     " << fRecPointsTitle.Data() << endl ;
419     cout <<  "    TrackSegments Branch title: " << fTrackSegmentsTitle.Data() << endl ;
420     cout <<  "    RecParticles Branch title   " << fRecParticlesTitle.Data() << endl;
421     cout <<  "with parameters: " << endl ;
422     cout <<  "    Maximal EMC - CPV  distance (cm) " << fCpvEmcDistance << endl ;
423     cout <<  "    Time Gate used:                  " << fTimeGate <<  endl ;
424     cout <<  "    Principal Ellipse Parameters     " << endl ;
425     cout <<  "       Ellipse center   (x,y)           (" << fX_center<<","<<fY_center<<")"<< endl;
426     cout <<  "       Ellipse focus    (a,b)           (" << fA<<","<<fB<<")"<< endl;
427     cout <<  "       Ellipse angle                     " << fAngle<< endl;        
428     cout <<  "============================================" << endl ;
429 }
430
431 //____________________________________________________________________________
432 void  AliPHOSPIDv2::WriteRecParticles(Int_t event)
433 {
434  
435   AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ; 
436   TString taskName(GetName()) ; 
437   taskName.Remove(taskName.Index(Version())-1) ;
438   TClonesArray * recParticles = gime->RecParticles(taskName) ; 
439   recParticles->Expand(recParticles->GetEntriesFast() ) ;
440
441   //Make branch in TreeR for RecParticles 
442   char * filename = 0;
443   if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){   //generating file name
444     filename = new char[strlen(gAlice->GetBaseFile())+20] ;
445     sprintf(filename,"%s/PHOS.Reco.root",gAlice->GetBaseFile()) ; 
446   }
447   
448   TDirectory *cwd = gDirectory;
449   
450   //First rp
451   Int_t bufferSize = 32000 ;    
452   TBranch * rpBranch = gAlice->TreeR()->Branch("PHOSRP",&recParticles,bufferSize);
453   rpBranch->SetTitle(fRecParticlesTitle);
454   if (filename) {
455     rpBranch->SetFile(filename);
456     TIter next( rpBranch->GetListOfBranches());
457     TBranch * sb ;
458     while ((sb=(TBranch*)next())) {
459       sb->SetFile(filename);
460     }   
461     cwd->cd();
462   }
463   
464   //second, pid
465   Int_t splitlevel = 0 ; 
466   AliPHOSPIDv2 * pid = this ;
467   TBranch * pidBranch = gAlice->TreeR()->Branch("AliPHOSPID","AliPHOSPIDv2",&pid,bufferSize,splitlevel);
468   pidBranch->SetTitle(fRecParticlesTitle.Data());
469   if (filename) {
470     pidBranch->SetFile(filename);
471     TIter next( pidBranch->GetListOfBranches());
472     TBranch * sb ;
473     while ((sb=(TBranch*)next())) {
474       sb->SetFile(filename);
475     }   
476     cwd->cd();
477   }    
478   
479   rpBranch->Fill() ;
480   pidBranch->Fill() ;
481   
482   gAlice->TreeR()->Write(0,kOverwrite) ;  
483   //pidBranch->Write(0,kOverwrite) ;  
484
485   delete [] filename ; 
486 }
487
488 //____________________________________________________________________________
489 TVector3 AliPHOSPIDv2::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * cpv)const 
490
491   // Calculates the momentum direction:
492   //   1. if only a EMC RecPoint, direction is given by IP and this RecPoint
493   //   2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints 
494   //  However because of the poor position resolution of PPSD the direction is always taken as if we were 
495   //  in case 1.
496
497   TVector3 dir(0,0,0) ; 
498   
499   TVector3 emcglobalpos ;
500   TMatrix  dummy ;
501   
502   emc->GetGlobalPosition(emcglobalpos, dummy) ;
503   
504
505   dir = emcglobalpos ;  
506   dir.SetZ( -dir.Z() ) ;   // why ?  
507   dir.SetMag(1.) ;
508
509   //account correction to the position of IP
510   Float_t xo,yo,zo ; //Coordinates of the origin
511   gAlice->Generator()->GetOrigin(xo,yo,zo) ;
512   TVector3 origin(xo,yo,zo);
513   dir = dir - origin ;
514
515   return dir ;  
516 }
517 //____________________________________________________________________________
518 void AliPHOSPIDv2::PrintRecParticles(Option_t * option)
519 {
520   // Print table of reconstructed particles
521
522   AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ; 
523
524   TString taskName(GetName()) ; 
525   taskName.Remove(taskName.Index(Version())-1) ;
526   TClonesArray * recParticles = gime->RecParticles(taskName) ; 
527   
528   cout << "AliPHOSPIDv2: event "<<gAlice->GetEvNumber()  << endl ;
529   cout << "       found " << recParticles->GetEntriesFast() << " RecParticles " << endl ;
530   
531   if(strstr(option,"all")) {  // printing found TS
532     
533     cout << "  PARTICLE "   
534          << "  Index    "  << endl ;
535     
536     Int_t index ;
537     for (index = 0 ; index < recParticles->GetEntries() ; index++) {
538        AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;       
539       
540        Text_t particle[11];
541        switch(rp->GetType()) {
542        case  AliPHOSFastRecParticle::kCHARGEDHASLOW:
543          strcpy(particle, "CHARGED HA SLOW") ;
544          break ; 
545        case  AliPHOSFastRecParticle::kNEUTRALHASLOW: 
546          strcpy(particle, "NEUTRAL HA SLOW");
547          break ;    
548        case  AliPHOSFastRecParticle::kCHARGEDHAFAST:
549          strcpy(particle, "CHARGED HA FAST") ;
550          break ;        
551        case  AliPHOSFastRecParticle::kNEUTRALHAFAST:
552          strcpy(particle, "NEUTRAL HA FAST");
553          break;
554        case  AliPHOSFastRecParticle::kCHARGEDEMSLOW:
555          strcpy(particle, "CHARGED EM SLOW") ;
556          break ;
557        case  AliPHOSFastRecParticle::kNEUTRALEMSLOW:
558          strcpy(particle, "NEUTRAL EM SLOW");
559          break ;
560        case  AliPHOSFastRecParticle::kCHARGEDEMFAST:
561          strcpy(particle, "CHARGED EM FAST") ;
562          break ;
563        case  AliPHOSFastRecParticle::kNEUTRALEMFAST:
564          strcpy( particle, "NEUTRAL EM FAST");
565          break;
566       }
567
568       cout << setw(10) << particle << "  "
569            << setw(5) <<  rp->GetIndexInList() << " " <<endl;
570    
571     }
572     cout << "-------------------------------------------" << endl ;
573   }
574   
575 }
576
577
578