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