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