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