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