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