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