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