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