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