]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSTrackSegmentMakerv1.cxx
No cast to AliMagFCM.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSTrackSegmentMakerv1.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 /* $Id$ */
16 //_________________________________________________________________________
17 // Implementation version 1 of algorithm class to construct PHOS track segments
18 // Track segment for PHOS is list of 
19 //        EMC RecPoint + (possibly) CPV RecPoint + (possibly) PPSD RecPoint
20 // To find TrackSegments we do the following: 
21 //  for each EMC RecPoints we look at
22 //   CPV/PPSD RecPoints in the radious fRcpv. 
23 //  If there is such a CPV RecPoint, 
24 //   we make "Link" it is just indexes of EMC and CPV/PPSD RecPoint and distance
25 //   between them in the PHOS plane. 
26 //  Then we sort "Links" and starting from the 
27 //   least "Link" pointing to the unassined EMC and CPV RecPoints assing them to 
28 //   new TrackSegment. 
29 // If there is no CPV/PPSD RecPoint we make TrackSegment 
30 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
31 //
32 // In principle this class should be called from AliPHOSReconstructioner, but 
33 // one can use it as well in standalone mode.
34 // Use  case:
35 //  root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
36 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
37 //               // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
38 //               // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")                       
39 //  root [1] t->ExecuteTask()
40 //  root [2] t->SetMaxEmcPpsdDistance(5)
41 //  root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
42 //  root [4] t->ExecuteTask("deb all time") 
43 //                 
44 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH) 
45 //
46
47 // --- ROOT system ---
48 #include "TROOT.h"
49 #include "TFile.h"
50 #include "TFolder.h"
51 #include "TTree.h"
52 #include "TSystem.h"
53 #include "TBenchmark.h"
54 // --- Standard library ---
55
56 #include <iostream.h>
57 #include <iomanip.h>
58
59 // --- AliRoot header files ---
60
61 #include "AliPHOSTrackSegmentMakerv1.h"
62 #include "AliPHOSClusterizerv1.h"
63 #include "AliPHOSTrackSegment.h"
64 #include "AliPHOSCpvRecPoint.h"
65 #include "AliPHOSLink.h"
66 #include "AliPHOSGetter.h"
67 #include "AliPHOS.h"
68 #include "AliRun.h"
69
70 ClassImp( AliPHOSTrackSegmentMakerv1) 
71
72
73 //____________________________________________________________________________
74   AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() : AliPHOSTrackSegmentMaker()
75 {
76   // default ctor (to be used mainly by Streamer)
77
78   InitParameters() ; 
79 //   fHeaderFileName           = "" ;
80 //   fRecPointsBranchTitle     = "" ;
81 //   fTrackSegmentsBranchTitle = "" ; 
82 //   fFrom                     = "" ; 
83
84   fTrackSegmentsInRun       = 0 ; 
85
86   fDefaultInit = kTRUE ; 
87 }
88
89 //____________________________________________________________________________
90  AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const char * headerFile, const char * name, const Bool_t toSplit) : AliPHOSTrackSegmentMaker(headerFile, name, toSplit)
91 {
92   // ctor
93
94   InitParameters() ; 
95 //   fHeaderFileName           = GetTitle() ;
96 //   fRecPointsBranchTitle     = GetName() ;
97 //   fTrackSegmentsBranchTitle = GetName() ; 
98   fTrackSegmentsInRun       = 0 ; 
99
100 //   if ( from == 0 ) 
101 //     fFrom = name ; 
102 //   else
103 //     fFrom = from ; 
104   Init() ;
105
106   fDefaultInit = kFALSE ; 
107
108 }
109
110 //____________________________________________________________________________
111  AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
112
113   // dtor
114   // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
115   
116   if (!fDefaultInit) {
117     delete fLinkUpArray  ;
118     
119 //     AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
120     
121 //     // remove the task from the folder list
122 //     gime->RemoveTask("T",GetName()) ;
123 //     TString name(GetName()) ; 
124 //     name.ReplaceAll("tsm", "clu") ; 
125 //     gime->RemoveTask("C",name) ;
126     
127 //     // remove the data from the folder list
128 //     name = GetName() ; 
129 //     name.Remove(name.Index(":")) ; 
130 //     gime->RemoveObjects("RE", name) ; // EMCARecPoints
131 //     gime->RemoveObjects("RC", name) ; // CPVRecPoints
132 //     gime->RemoveObjects("T", name) ;  // TrackSegments
133     
134 //     // Delete gAlice
135 //     gime->CloseFile() ; 
136     
137     fSplitFile = 0 ; 
138   }
139 }
140
141
142 //____________________________________________________________________________
143 const TString AliPHOSTrackSegmentMakerv1::BranchName() const 
144 {  
145   TString branchName(GetName() ) ;
146   branchName.Remove(branchName.Index(Version())-1) ;
147   return branchName ;
148 }
149
150 //____________________________________________________________________________
151 void  AliPHOSTrackSegmentMakerv1::FillOneModule()
152 {
153   // Finds first and last indexes between which 
154   // clusters from one PHOS module are
155   
156   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
157   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
158   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
159  
160   //First EMC clusters
161   Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
162   for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&  
163         ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule ); 
164       fEmcLast ++)  ;
165   
166   //Now CPV clusters
167   Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
168
169     for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) && 
170           ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule ); 
171         fCpvLast ++) ;
172       
173 }
174
175 //____________________________________________________________________________
176 Float_t  AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSRecPoint * cpvClu, Bool_t &toofar)const
177 {
178   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
179   // Clusters are sorted in "rows" and "columns" of width 1 cm
180
181   Float_t delta = 1 ;  // Width of the rows in sorting of RecPoints (in cm)
182                        // if you change this value, change it as well in xxxRecPoint::Compare()
183   Float_t r = fRcpv ;
184  
185   TVector3 vecEmc ;
186   TVector3 vecCpv ;
187   
188   emcClu->GetLocalPosition(vecEmc) ;
189   cpvClu->GetLocalPosition(vecCpv)  ; 
190
191   if(emcClu->GetPHOSMod() == cpvClu->GetPHOSMod()){ 
192     if(vecCpv.X() <= vecEmc.X() + fRcpv + 2*delta ){ 
193
194       vecCpv = vecCpv  - vecEmc ; 
195       r = vecCpv.Mag() ;
196       toofar = kFALSE ;
197
198     } // if  xPpsd >= xEmc + ...
199     else 
200       toofar = kTRUE ;
201   } 
202   else 
203     toofar = kTRUE ;
204
205   //toofar = kFALSE ;
206  
207   
208   return r ;
209 }
210
211 //____________________________________________________________________________
212 void  AliPHOSTrackSegmentMakerv1::Init()
213 {
214   // Make all memory allocations that are not possible in default constructor
215   
216   if ( strcmp(GetTitle(), "") == 0 )
217     SetTitle("galice.root") ;
218     
219   TString branchname = GetName() ;
220   branchname.Remove(branchname.Index(Version())-1) ;
221   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(),branchname.Data(), fToSplit ) ; 
222   if ( gime == 0 ) {
223     cerr << "ERROR: AliPHOSTrackSegmentMakerv1::Init -> Could not obtain the Getter object !" << endl ; 
224     return ;
225   } 
226   
227   fSplitFile = 0 ;
228   if(fToSplit){
229     //First - extract full path if necessary
230     TString fileName(GetTitle()) ;
231     Ssiz_t islash = fileName.Last('/') ;
232     if(islash<fileName.Length())
233       fileName.Remove(islash+1,fileName.Length()) ;
234     else
235       fileName="" ;
236     fileName+="PHOS.RecData." ;
237     if((strcmp(branchname.Data(),"Default")!=0)&&(strcmp(branchname.Data(),"")!=0)){
238       fileName+=branchname ;
239       fileName+="." ;
240     }
241     fileName+="root" ;
242     fSplitFile = static_cast<TFile*>(gROOT->GetFile(fileName.Data()));   
243     if(!fSplitFile)
244       fSplitFile =  TFile::Open(fileName.Data(),"update") ;
245   }
246   
247   fLinkUpArray  = new TClonesArray("AliPHOSLink", 1000); 
248   
249
250   gime->PostTrackSegmentMaker(this) ;
251   gime->PostTrackSegments(BranchName()) ; 
252
253 }
254
255 //____________________________________________________________________________
256 void  AliPHOSTrackSegmentMakerv1::InitParameters()
257 {
258   fRcpv      = 10. ;   
259   fEmcFirst  = 0 ;    
260   fEmcLast   = 0 ;   
261   fCpvFirst  = 0 ;   
262   fCpvLast   = 0 ;   
263   fLinkUpArray = 0 ;
264   TString tsmName( GetName()) ; 
265   if (tsmName.IsNull() ) 
266     tsmName = "Default" ; 
267   tsmName.Append(":") ; 
268   tsmName.Append(Version()) ; 
269   SetName(tsmName) ;
270 }
271
272
273 //____________________________________________________________________________
274 void  AliPHOSTrackSegmentMakerv1::MakeLinks()const
275
276   // Finds distances (links) between all EMC and PPSD clusters, 
277   // which are not further apart from each other than fRcpv 
278   // and sort them in accordance with this distance
279   
280   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
281   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
282   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
283
284   fLinkUpArray->Clear() ;    
285
286   AliPHOSRecPoint * cpv ;
287   AliPHOSEmcRecPoint * emcclu ;
288
289   Int_t iLinkUp  = 0 ;
290   
291   Int_t iEmcRP;
292   for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
293     emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
294
295     Bool_t toofar ;        
296     Int_t iCpv = 0 ;    
297     for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) { 
298       
299       cpv = dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(iCpv)) ;
300       Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, toofar) ;
301       
302       if(toofar)
303         break ;  
304       if(r < fRcpv) { 
305         new ((*fLinkUpArray)[iLinkUp++])  AliPHOSLink(r, iEmcRP, iCpv) ;
306       }      
307     }
308   } 
309   
310   fLinkUpArray->Sort() ;  //first links with smallest distances
311 }
312
313 //____________________________________________________________________________
314 void  AliPHOSTrackSegmentMakerv1::MakePairs()
315
316   // Using the previously made list of "links", we found the smallest link - i.e. 
317   // link with the least distance between EMC and CPV and pointing to still 
318   // unassigned RecParticles. We assign these RecPoints to TrackSegment and 
319   // remove them from the list of "unassigned". 
320   
321   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
322   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
323   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
324   TClonesArray * trackSegments = gime->TrackSegments(BranchName()) ;   
325     
326   //Make arrays to mark clusters already chosen
327   Int_t * emcExist = 0;
328   if(fEmcLast > fEmcFirst)
329     emcExist = new Int_t[fEmcLast-fEmcFirst] ;
330   
331   Int_t index;
332   for(index = 0; index <fEmcLast-fEmcFirst; index ++)
333     emcExist[index] = 1 ;
334   
335   Bool_t * cpvExist = 0;
336   if(fCpvLast > fCpvFirst)
337     cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
338   for(index = 0; index <fCpvLast-fCpvFirst; index ++)
339     cpvExist[index] = kTRUE ;
340   
341   
342   // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance 
343   TIter nextUp(fLinkUpArray) ;
344   
345   AliPHOSLink * linkUp ;
346   
347   AliPHOSRecPoint * nullpointer = 0 ;
348   
349   while ( (linkUp =  static_cast<AliPHOSLink *>(nextUp()) ) ){  
350
351     if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){ //without ppsd Up yet 
352
353       if(cpvExist[linkUp->GetPpsd()-fCpvFirst]){ //CPV still exist
354         
355         new ((* trackSegments)[fNTrackSegments]) 
356           AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) , 
357                               dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(linkUp->GetPpsd()))) ;
358         (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
359         fNTrackSegments++ ;
360         
361         emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc  that Cpv was found 
362         //mark CPV recpoint as already used 
363         cpvExist[linkUp->GetPpsd()-fCpvFirst] = kFALSE ;
364       } //if ppsdUp still exist
365     } 
366   }      
367
368   //look through emc recPoints left without CPV/PPSD
369   if(emcExist){ //if there is emc rec point
370     Int_t iEmcRP ;
371     for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst  ; iEmcRP++ ){
372       if(emcExist[iEmcRP] > 0 ){
373         new ((*trackSegments)[fNTrackSegments])  
374           AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)), 
375                               nullpointer) ;
376         (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
377         fNTrackSegments++;    
378       } 
379     }
380   }
381   delete [] emcExist ; 
382   delete [] cpvExist ; 
383 }
384
385 //____________________________________________________________________________
386 void  AliPHOSTrackSegmentMakerv1::Exec(Option_t * option)
387 {
388   // STEERing method
389
390   if( strcmp(GetName(), "")== 0 ) 
391     Init() ;
392
393   if(strstr(option,"tim"))
394     gBenchmark->Start("PHOSTSMaker");
395  
396   if(strstr(option,"print")) {
397     Print("") ; 
398     return ; 
399   }
400
401 //   gAlice->GetEvent(0) ;
402 //   //check, if the branch with name of this" already exits?
403 //   if (gAlice->TreeR()) { 
404 //     TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeR()->GetListOfBranches()) ;
405 //     TIter next(lob) ; 
406 //     TBranch * branch = 0 ;  
407 //     Bool_t phostsfound = kFALSE, tracksegmentmakerfound = kFALSE ; 
408     
409 //     TString branchname = GetName() ;
410 //     branchname.Remove(branchname.Index(Version())-1) ;
411     
412 //     while ( (branch = static_cast<TBranch*>(next())) && (!phostsfound || !tracksegmentmakerfound) ) {
413 //       if ( (strcmp(branch->GetName(), "PHOSTS")==0) && (strcmp(branch->GetTitle(), branchname.Data())==0) ) 
414 //      phostsfound = kTRUE ;
415       
416 //       else if ( (strcmp(branch->GetName(), "AliPHOSTrackSegmentMaker")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
417 //      tracksegmentmakerfound = kTRUE ; 
418 //     }
419     
420 //     if ( phostsfound || tracksegmentmakerfound ) {
421 //       cerr << "WARNING: AliPHOSTrackSegmentMakerv1::Exec -> TrackSegments and/or TrackSegmentMaker branch with name " 
422 //         << branchname.Data() << " already exits" << endl ;
423 //       return ; 
424 //     }       
425 //   }
426
427   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
428   if(gime->BranchExists("TrackSegments") )
429     return ;
430   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
431   Int_t nevents = gime->MaxEvent() ;       //(Int_t) gAlice->TreeE()->GetEntries() ;
432   Int_t ievent ;
433
434   for(ievent = 0; ievent < nevents; ievent++){
435
436     gime->Event(ievent,"R") ;
437     //Make some initializations 
438     fNTrackSegments = 0 ;
439     fEmcFirst = 0 ;    
440     fEmcLast  = 0 ;   
441     fCpvFirst = 0 ;   
442     fCpvLast  = 0 ;   
443     gime->TrackSegments(BranchName())->Clear() ; 
444
445     //    if(!ReadRecPoints(ievent))   continue; //reads RecPoints for event ievent
446     
447     for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ){
448       
449       FillOneModule() ; 
450       
451       MakeLinks() ;
452       
453       MakePairs() ;
454       
455     }
456
457     WriteTrackSegments(ievent) ;
458
459     if(strstr(option,"deb"))
460       PrintTrackSegments(option) ;
461     
462     //increment the total number of track segments per run 
463     fTrackSegmentsInRun += gime->TrackSegments(BranchName())->GetEntriesFast() ; 
464
465   }
466
467   if(strstr(option,"tim")){
468     gBenchmark->Stop("PHOSTSMaker");
469     cout << "AliPHOSTSMaker:" << endl ;
470     cout << "  took " << gBenchmark->GetCpuTime("PHOSTSMaker") << " seconds for making TS " 
471          <<  gBenchmark->GetCpuTime("PHOSTSMaker")/nevents << " seconds per event " << endl ;
472     cout << endl ;
473   }
474     
475 }
476
477 //____________________________________________________________________________
478 void AliPHOSTrackSegmentMakerv1::Print(Option_t * option)const
479 {
480   //  Print TrackSegmentMaker parameters
481
482   if( strcmp(GetName(), "") != 0 ) {
483     cout <<  "======== AliPHOSTrackSegmentMakerv1 ========" << endl ;
484     cout <<  "Making Track segments "<< endl ;
485     //    cout <<  "    Headers file:                   " << fHeaderFileName.Data() << endl ;
486     //    cout <<  "    RecPoints branch file name:     " << fRecPointsBranchTitle.Data() << endl ;
487     //    cout <<  "    TrackSegments Branch file name: " << fTrackSegmentsBranchTitle.Data() << endl ;
488     cout <<  "with parameters: " << endl ;
489     cout <<  "    Maximal EMC - CPV (PPSD) distance (cm)" << fRcpv << endl ;
490     cout <<  "============================================" << endl ;
491   }
492   else
493     cout << "AliPHOSTrackSegmentMakerv1 not initialized " << endl ;
494 }
495
496 //____________________________________________________________________________
497 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments(Int_t event)
498 {
499   // Writes found TrackSegments to TreeR. Creates branches 
500   // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
501   // In the former branch found TrackSegments are stored, while 
502   // in the latter all parameters, with which TS were made. 
503   // ROOT does not allow overwriting existing branches, therefore
504   // first we check, if branches with the same title already exist.
505   // If yes - exits without writing.
506   
507   AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ; 
508
509   TClonesArray * trackSegments = gime->TrackSegments() ; 
510   trackSegments->Expand(trackSegments->GetEntriesFast()) ;
511   TTree * treeR ;
512   
513   if(fToSplit){
514     if(!fSplitFile)
515       return ;
516     fSplitFile->cd() ;
517     char name[10] ;
518     sprintf(name,"%s%d", "TreeR",event) ;
519     treeR = dynamic_cast<TTree*>(fSplitFile->Get(name)); 
520   }
521   else{
522     treeR = gAlice->TreeR();
523   }
524   
525   if(!treeR){
526     gAlice->MakeTree("R", fSplitFile);
527     treeR = gAlice->TreeR() ;
528   }
529   
530   //First TS
531   Int_t bufferSize = 32000 ;    
532   TBranch * tsBranch = treeR->Branch("PHOSTS",&trackSegments,bufferSize);
533   tsBranch->SetTitle(BranchName());
534
535   //Second -TSMaker
536   Int_t splitlevel = 0 ;
537   AliPHOSTrackSegmentMakerv1 * ts = this ;
538   TBranch * tsMakerBranch = treeR->Branch("AliPHOSTrackSegmentMaker","AliPHOSTrackSegmentMakerv1",
539                                           &ts,bufferSize,splitlevel);
540   tsMakerBranch->SetTitle(BranchName());
541
542   tsBranch->Fill() ;  
543   tsMakerBranch->Fill() ;
544
545   treeR->AutoSave() ; //Write(0,kOverwrite) ;  
546   if(gAlice->TreeR()!=treeR)
547     treeR->Delete();
548 }
549
550
551 //____________________________________________________________________________
552 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
553 {
554   // option deb - prints # of found TrackSegments
555   // option deb all - prints as well indexed of found RecParticles assigned to the TS
556   TString taskName(GetName()) ; 
557   taskName.Remove(taskName.Index(Version())-1) ;
558   
559   TClonesArray * trackSegments = AliPHOSGetter::GetInstance()->TrackSegments(taskName) ; 
560
561   
562   cout << "AliPHOSTrackSegmentMakerv1: event "<<gAlice->GetEvNumber()  << endl ;
563   cout << "       Found " << trackSegments->GetEntriesFast() << "  trackSegments " << endl ;
564   
565   if(strstr(option,"all")) {  // printing found TS
566     cout << "TrackSegment # " << "    EMC RP#    " << "    CPV RP#    " << endl ; 
567     
568     Int_t index;
569     for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
570       AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ; 
571       cout<<"   "<< setw(4) << ts->GetIndexInList() << "            " 
572           <<setw(4) << ts->GetEmcIndex()<< "            " 
573           <<setw(4) << ts->GetCpvIndex()<< "            " << endl ;
574     }   
575     
576     cout << "-------------------------------------------------------"<< endl ;
577   }
578 }