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