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