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