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