c1a231b1e1b989aa5f8407041056d0fd356176f0
[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 /* History of cvs commits:
18  *
19  * $Log$
20  * Revision 1.84  2007/03/07 07:01:21  hristov
21  * Fixing copy/paste erro. Additional protections
22  *
23  * Revision 1.83  2007/03/06 21:07:37  kharlov
24  * DP: xz CPV-EMC distance filled to TS
25  *
26  * Revision 1.82  2007/03/06 06:54:48  kharlov
27  * DP:Calculation of cluster properties dep. on vertex added
28  *
29  * Revision 1.81  2007/02/05 10:02:40  kharlov
30  * Module numbering is corrected
31  *
32  * Revision 1.80  2006/08/28 10:01:56  kharlov
33  * Effective C++ warnings fixed (Timur Pocheptsov)
34  *
35  * Revision 1.79  2006/04/25 12:41:15  hristov
36  * Moving non-persistent data to AliESDfriend (Yu.Belikov)
37  *
38  * Revision 1.78  2005/11/18 13:04:51  hristov
39  * Bug fix
40  *
41  * Revision 1.77  2005/11/17 23:34:36  hristov
42  * Corrected logics
43  *
44  * Revision 1.76  2005/11/17 22:29:12  hristov
45  * Faster version, no attempt to match tracks outside the PHOS acceptance
46  *
47  * Revision 1.75  2005/11/17 12:35:27  hristov
48  * Use references instead of objects. Avoid to create objects when they are not really needed
49  *
50  * Revision 1.74  2005/07/08 14:01:36  hristov
51  * Tracking in non-uniform nmagnetic field (Yu.Belikov)
52  *
53  * Revision 1.73  2005/05/28 14:19:05  schutz
54  * Compilation warnings fixed by T.P.
55  *
56  */
57
58 //_________________________________________________________________________
59 // Implementation version 1 of algorithm class to construct PHOS track segments
60 // Track segment for PHOS is list of 
61 //        EMC RecPoint + (possibly) CPV RecPoint
62 // To find TrackSegments we do the following: 
63 //  for each EMC RecPoints we look at
64 //   CPV RecPoints in the radious fRcpv. 
65 //  If there is such a CPV RecPoint, 
66 //   we make "Link" it is just indexes of EMC and CPV RecPoint and distance
67 //   between them in the PHOS plane. 
68 //  Then we sort "Links" and starting from the 
69 //   least "Link" pointing to the unassined EMC and CPV RecPoints assing them to 
70 //   new TrackSegment. 
71 // If there is no CPV RecPoint we make TrackSegment 
72 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
73 //// In principle this class should be called from AliPHOSReconstructor, but 
74 // one can use it as well in standalone mode.
75 // Use  case:
76 //  root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
77 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
78 //               // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
79 //               // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")                       
80 //  root [1] t->ExecuteTask()
81 //  root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
82 //  root [4] t->ExecuteTask("deb all time") 
83 //                 
84 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH) 
85 //
86
87 // --- ROOT system ---
88 #include "TTree.h"
89 #include "TBenchmark.h"
90
91 // --- Standard library ---
92 #include "Riostream.h"
93 // --- AliRoot header files ---
94 #include "AliPHOSGeometry.h"
95 #include "AliPHOSTrackSegmentMakerv1.h"
96 #include "AliPHOSTrackSegment.h"
97 #include "AliPHOSLink.h"
98 #include "AliPHOSGetter.h"
99 #include "AliESD.h"
100 #include "AliESDtrack.h"
101
102 ClassImp( AliPHOSTrackSegmentMakerv1) 
103
104
105 //____________________________________________________________________________
106 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() :
107   AliPHOSTrackSegmentMaker(),
108   fDefaultInit(kTRUE),
109   fWrite(kFALSE),
110   fNTrackSegments(0),
111   fRcpv(0.f),
112   fRtpc(0.f),
113   fLinkUpArray(0),
114   fEmcFirst(0),
115   fEmcLast(0),
116   fCpvFirst(0),
117   fCpvLast(0),
118   fModule(0),
119   fTrackSegmentsInRun(0)
120   
121 {
122   // default ctor (to be used mainly by Streamer)
123   InitParameters() ; 
124 }
125
126 //____________________________________________________________________________
127 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const TString & alirunFileName, const TString & eventFolderName) :
128   AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName),
129   fDefaultInit(kFALSE),
130   fWrite(kFALSE),
131   fNTrackSegments(0),
132   fRcpv(0.f),
133   fRtpc(0.f),
134   fLinkUpArray(0),
135   fEmcFirst(0),
136   fEmcLast(0),
137   fCpvFirst(0),
138   fCpvLast(0),
139   fModule(0),
140   fTrackSegmentsInRun(0)
141 {
142   // ctor
143   InitParameters() ; 
144   Init() ;
145   fESD = 0;
146 }
147
148
149 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const AliPHOSTrackSegmentMakerv1 & tsm) :
150   AliPHOSTrackSegmentMaker(tsm),
151   fDefaultInit(kFALSE),
152   fWrite(kFALSE),
153   fNTrackSegments(0),
154   fRcpv(0.f),
155   fRtpc(0.f),
156   fLinkUpArray(0),
157   fEmcFirst(0),
158   fEmcLast(0),
159   fCpvFirst(0),
160   fCpvLast(0),
161   fModule(0),
162   fTrackSegmentsInRun(0)
163 {
164   // cpy ctor: no implementation yet
165   // requested by the Coding Convention
166   Fatal("cpy ctor", "not implemented") ;
167 }
168
169
170 //____________________________________________________________________________
171  AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
172
173   // dtor
174   // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
175   if (!fDefaultInit)  
176     delete fLinkUpArray ;
177 }
178
179
180 //____________________________________________________________________________
181 const TString AliPHOSTrackSegmentMakerv1::BranchName() const 
182 {  
183  
184   return GetName() ;
185 }
186
187 //____________________________________________________________________________
188 void  AliPHOSTrackSegmentMakerv1::FillOneModule()
189 {
190   // Finds first and last indexes between which 
191   // clusters from one PHOS module are
192
193   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
194   
195   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
196   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
197  
198   //First EMC clusters
199   Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
200   for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&  
201         ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule ); 
202       fEmcLast ++)  ;
203   
204   //Now CPV clusters
205   Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
206
207     for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) && 
208          ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule ); 
209        fCpvLast ++) ;
210       
211 }
212
213 //____________________________________________________________________________
214 void  AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,
215                                                          AliPHOSCpvRecPoint * cpvClu, 
216                                                          Int_t &trackindex, 
217                                                          Float_t &dx, Float_t &dz) const
218 {
219   // Calculates the distance between the EMC RecPoint and the CPV RecPoint
220   // Clusters are sorted in "rows" and "columns" of width 1 cm
221
222 //  Float_t delta = 1 ;  // Width of the rows in sorting of RecPoints (in cm)
223 //                       // if you change this value, change it as well in xxxRecPoint::Compare()
224   Float_t distance2Track = fRtpc ; 
225
226   trackindex = -1 ; // closest track within fRCpv 
227
228   TVector3 vecEmc ;   // Local position of EMC recpoint
229   TVector3 vecP ;     // Momentum direction at CPV plain
230   TVector3 vecPloc ;     // Momentum direction at CPV plain
231   
232   //toofar = kTRUE ;
233   if(emcClu->GetPHOSMod() != cpvClu->GetPHOSMod()){
234     dx=999. ;
235     dz=999. ;
236     return ;
237   }
238
239   emcClu->GetLocalPosition(vecEmc) ;
240
241   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
242   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
243     
244   Double_t xCPV,zCPV ; //EMC-projected coordinates of CPV cluster 
245   TVector3 cpvGlobal; // Global position of the CPV recpoint
246   geom->GetGlobal((AliRecPoint*)cpvClu,cpvGlobal);
247   Double_t vtxCPV[3]={cpvGlobal.X(),cpvGlobal.Y(),cpvGlobal.Z()} ;
248
249   if (fESD == 0x0) {
250      //if no track information available, assume straight line from IP to emcal
251      Int_t dummyMod ;
252      geom->ImpactOnEmc(vtxCPV,cpvGlobal.Theta(),cpvGlobal.Phi(),dummyMod,xCPV,zCPV) ;
253      dx=xCPV - vecEmc.X() ;
254      dz=zCPV - vecEmc.Z() ;
255      return ;
256   } 
257
258   //if there is ESD try to correct distance using TPC information on particle direct in CPV 
259   if (fESD != 0x0) {
260
261     Double_t rCPV = cpvGlobal.Pt() ;// Radius from IP to current point 
262
263     // Extrapolate the global track direction if any to CPV and find the closest track
264     Int_t nTracks = fESD->GetNumberOfTracks();
265     Int_t iClosestTrack = -1;
266     Double_t minDistance = 1.e6;
267     TVector3 inPHOS ; 
268
269     AliESDtrack *track;
270     Double_t xyz[3] ;
271     Double_t pxyz[3]; 
272     for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
273       track = fESD->GetTrack(iTrack);
274       if (!track->GetXYZAt(rCPV, fESD->GetMagneticField(), xyz))
275            continue; //track coord on the cylinder of PHOS radius
276       if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
277            continue;
278       //Check if this track hits PHOS
279       inPHOS.SetXYZ(xyz[0],xyz[1],xyz[2]);
280       distance2Track = inPHOS.Angle(cpvGlobal) ;
281       // Find the closest track to the CPV recpoint
282       if (distance2Track < minDistance) {
283          minDistance = distance2Track;
284          iClosestTrack = iTrack;
285       }
286     }
287
288     if (iClosestTrack != -1) {
289       track = fESD->GetTrack(iClosestTrack);
290       if (track->GetPxPyPzAt(rCPV, fESD->GetMagneticField(), pxyz)) { // track momentum ibid.
291         vecP.SetXYZ(pxyz[0],pxyz[1],pxyz[2]);
292         Int_t dummyMod ;
293         geom->ImpactOnEmc(vtxCPV,vecP.Theta(),vecP.Phi(),dummyMod,xCPV,zCPV) ;
294       }
295     }
296     
297     if(minDistance < fRtpc ){
298       trackindex = iClosestTrack ; 
299     }
300   }
301   if(trackindex!=-1){
302     dx=xCPV - vecEmc.X() ;
303     dz=zCPV - vecEmc.Z() ;
304   }
305   else{
306     dz=999. ;
307     dx=999. ;
308   }
309   return ;
310 }
311 //____________________________________________________________________________
312 void  AliPHOSTrackSegmentMakerv1::Init()
313 {
314   // Make all memory allocations that are not possible in default constructor
315   
316   AliPHOSGetter* gime = AliPHOSGetter::Instance();
317   if(!gime)
318     gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
319   
320   fLinkUpArray  = new TClonesArray("AliPHOSLink", 1000); 
321   if ( !gime->TrackSegmentMaker() ) {
322     gime->PostTrackSegmentMaker(this);
323   }
324 }
325
326 //____________________________________________________________________________
327 void  AliPHOSTrackSegmentMakerv1::InitParameters()
328 {
329   //Initializes parameters
330   fRcpv      = 10. ;
331   fRtpc      = 4. ;
332   fEmcFirst  = 0 ;    
333   fEmcLast   = 0 ;   
334   fCpvFirst  = 0 ;   
335   fCpvLast   = 0 ;   
336   fLinkUpArray = 0 ;
337   fWrite                   = kTRUE ;
338   fTrackSegmentsInRun       = 0 ; 
339   SetEventRange(0,-1) ;
340 }
341
342
343 //____________________________________________________________________________
344 void  AliPHOSTrackSegmentMakerv1::MakeLinks()const
345
346   // Finds distances (links) between all EMC and CPV clusters, 
347   // which are not further apart from each other than fRcpv 
348   // and sort them in accordance with this distance
349   
350   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
351   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
352   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
353
354   fLinkUpArray->Clear() ;    
355
356   AliPHOSCpvRecPoint * cpv ;
357   AliPHOSEmcRecPoint * emcclu ;
358
359   Int_t iLinkUp  = 0 ;
360   
361   Int_t iEmcRP;
362   for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
363     emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
364
365     //Bool_t toofar ;        
366     Int_t iCpv = 0 ;    
367     for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) { 
368       
369       cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
370       Int_t track = -1 ; 
371       Float_t dx,dz ;
372       GetDistanceInPHOSPlane(emcclu, cpv, track,dx,dz) ;     
373       if(TMath::Sqrt(dx*dx+dz*dz) < fRcpv ){ 
374         new ((*fLinkUpArray)[iLinkUp++])  AliPHOSLink(dx, dz, iEmcRP, iCpv, track) ;
375       }      
376     }
377   } 
378   
379   fLinkUpArray->Sort() ;  //first links with smallest distances
380 }
381
382 //____________________________________________________________________________
383 void  AliPHOSTrackSegmentMakerv1::MakePairs()
384
385   // Using the previously made list of "links", we found the smallest link - i.e. 
386   // link with the least distance between EMC and CPV and pointing to still 
387   // unassigned RecParticles. We assign these RecPoints to TrackSegment and 
388   // remove them from the list of "unassigned". 
389
390   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
391
392   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
393   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
394   TClonesArray * trackSegments = gime->TrackSegments();
395     
396   //Make arrays to mark clusters already chosen
397   Int_t * emcExist = 0;
398   if(fEmcLast > fEmcFirst)
399     emcExist = new Int_t[fEmcLast-fEmcFirst] ;
400   
401   Int_t index;
402   for(index = 0; index <fEmcLast-fEmcFirst; index ++)
403     emcExist[index] = 1 ;
404   
405   Bool_t * cpvExist = 0;
406   if(fCpvLast > fCpvFirst)
407     cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
408   for(index = 0; index <fCpvLast-fCpvFirst; index ++)
409     cpvExist[index] = kTRUE ;
410   
411   
412   // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance 
413   TIter nextUp(fLinkUpArray) ;
414   
415   AliPHOSLink * linkUp ;
416   
417   AliPHOSCpvRecPoint * nullpointer = 0 ;
418   
419   while ( (linkUp =  static_cast<AliPHOSLink *>(nextUp()) ) ){  
420
421     if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
422
423       if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
424          Float_t dx,dz ;
425          linkUp->GetXZ(dx,dz) ;
426          new ((* trackSegments)[fNTrackSegments]) 
427            AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) , 
428                                dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetCpv())) , 
429                                linkUp->GetTrack(),dx,dz) ;
430          
431        (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
432        fNTrackSegments++ ;
433        emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc  that Cpv was found 
434        //mark CPV recpoint as already used 
435        cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
436       } //if CpvUp still exist
437     } 
438   }        
439
440   //look through emc recPoints left without CPV
441   if(emcExist){ //if there is emc rec point
442     Int_t iEmcRP ;
443     for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst  ; iEmcRP++ ){
444       if(emcExist[iEmcRP] > 0 ){
445        new ((*trackSegments)[fNTrackSegments])  
446          AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)), 
447                            nullpointer) ;
448        (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
449        fNTrackSegments++;    
450       } 
451     }
452   }
453   delete [] emcExist ; 
454   delete [] cpvExist ; 
455 }
456
457 //____________________________________________________________________________
458 void  AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
459 {
460   // Steering method to perform track segment construction for events
461   // in the range from fFirstEvent to fLastEvent.
462   // This range is optionally set by SetEventRange().
463   // if fLastEvent=-1 (by default), then process events until the end.
464   
465   if(strstr(option,"tim"))
466     gBenchmark->Start("PHOSTSMaker");
467  
468   if(strstr(option,"print")) {
469     Print() ; 
470     return ; 
471   }
472   
473   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;  
474  
475   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
476
477   if (fLastEvent == -1) 
478     fLastEvent = gime->MaxEvent() - 1 ;
479   else 
480     fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
481   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
482
483   Int_t ievent ; 
484   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
485     gime->Event(ievent,"DR") ;
486    //Make some initializations 
487     fNTrackSegments = 0 ;
488     fEmcFirst = 0 ;    
489     fEmcLast  = 0 ;   
490     fCpvFirst = 0 ;   
491     fCpvLast  = 0 ;   
492     
493     gime->TrackSegments()->Clear();
494
495     //    if(!ReadRecPoints(ievent))   continue; //reads RecPoints for event ievent
496     
497     for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
498       FillOneModule() ; 
499       MakeLinks() ;
500       MakePairs() ;
501     }
502
503     WriteTrackSegments() ;
504
505     if(strstr(option,"deb"))
506       PrintTrackSegments(option);
507     
508     //increment the total number of track segments per run 
509     fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ; 
510   }
511   
512   if(strstr(option,"tim")){
513     gBenchmark->Stop("PHOSTSMaker");
514     Info("Exec", "took %f seconds for making TS %f seconds per event", 
515           gBenchmark->GetCpuTime("PHOSTSMaker"), 
516           gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
517    }
518   if(fWrite) //do not unload in "on flight" mode
519     Unload();
520 }
521 //____________________________________________________________________________
522 void AliPHOSTrackSegmentMakerv1::Unload() 
523 {
524   // Unloads the task from the folder
525   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;  
526   gime->PhosLoader()->UnloadRecPoints() ;
527   gime->PhosLoader()->UnloadTracks() ;
528 }
529
530 //____________________________________________________________________________
531 void AliPHOSTrackSegmentMakerv1::Print(const Option_t *)const
532 {
533   //  Print TrackSegmentMaker parameters
534
535   TString message("") ;
536   if( strcmp(GetName(), "") != 0 ) {
537     message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ; 
538     message += "Making Track segments\n" ;
539     message += "with parameters:\n" ; 
540     message += "     Maximal EMC - CPV distance (cm) %f\n" ;
541     message += "============================================\n" ;
542     Info("Print", message.Data(),fRcpv) ;
543   }
544   else
545     Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
546 }
547
548 //____________________________________________________________________________
549 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
550 {
551   // Writes found TrackSegments to TreeR. Creates branches 
552   // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
553   // In the former branch found TrackSegments are stored, while 
554   // in the latter all parameters, with which TS were made. 
555   // ROOT does not allow overwriting existing branches, therefore
556   // first we check, if branches with the same title already exist.
557   // If yes - exits without writing.
558
559   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
560
561   TClonesArray * trackSegments = gime->TrackSegments() ; 
562   trackSegments->Expand(trackSegments->GetEntriesFast()) ;
563
564   if(fWrite){ //We write TreeT
565     TTree * treeT = gime->TreeT();
566     
567     //First TS
568     Int_t bufferSize = 32000 ; 
569     TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
570     tsBranch->Fill() ;  
571     
572     gime->WriteTracks("OVERWRITE");
573     gime->WriteTrackSegmentMaker("OVERWRITE");
574   }
575 }
576
577
578 //____________________________________________________________________________
579 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
580 {
581   // option deb - prints # of found TrackSegments
582   // option deb all - prints as well indexed of found RecParticles assigned to the TS
583
584   TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ; 
585
586   Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ; 
587   printf("nevent: %d\n", gAlice->GetEvNumber()) ; 
588   printf("        Found %d TrackSegments\n", trackSegments->GetEntriesFast() ); 
589   
590   if(strstr(option,"all")) {  // printing found TS
591     printf("TrackSegment #  EMC RP#  CPV RP#\n") ; 
592     Int_t index;
593     for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
594       AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ; 
595       printf("   %d           %d        %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ; 
596     }   
597   }
598 }