]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSTrackSegmentMakerv1.cxx
merge conflicts solved
[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 // Implementation version 1 of algorithm class to construct PHOS track segments
18 // Track segment for PHOS is list of 
19 //        EMC RecPoint + (possibly) CPV RecPoint + (possibly) PPSD RecPoint
20 // To find TrackSegments we do the following: 
21 //  for each EMC RecPoints we look at
22 //   CPV/PPSD RecPoints in the radious fR0. 
23 //  If there is such a CPV RecPoint, 
24 //   we make "Link" it is just indexes of EMC and CPV/PPSD RecPoint and distance
25 //   between them in the PHOS plane. 
26 //  Then we sort "Links" and starting from the 
27 //   least "Link" pointing to the unassined EMC and CPV RecPoints assing them to 
28 //   new TrackSegment. 
29 // If there is no CPV/PPSD RecPoint we make TrackSegment 
30 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
31 //
32 // In principle this class should be called from AliPHOSReconstructioner, but 
33 // one can use it as well in standalone mode.
34 // Use  case:
35 //  root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root")
36 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
37 //  root [1] t->ExecuteTask()
38 //  root [2] t->SetMaxEmcPpsdDistance(5)
39 //  root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
40 //  root [4] t->ExecuteTask("deb all time") 
41 //                 
42 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH)
43 //
44
45 // --- ROOT system ---
46 #include "TROOT.h"
47 #include "TFile.h"
48 #include "TTree.h"
49 #include "TSystem.h"
50 #include "TBenchmark.h"
51 // --- Standard library ---
52
53 #include <iostream.h>
54 #include <iomanip.h>
55
56 // --- AliRoot header files ---
57
58 #include "AliPHOSTrackSegmentMakerv1.h"
59 #include "AliPHOSClusterizerv1.h"
60 #include "AliPHOSTrackSegment.h"
61 #include "AliPHOSCpvRecPoint.h"
62 #include "AliPHOSPpsdRecPoint.h"
63 #include "AliPHOSLink.h"
64 #include "AliPHOS.h"
65 #include "AliRun.h"
66
67 ClassImp( AliPHOSTrackSegmentMakerv1) 
68
69
70 //____________________________________________________________________________
71  AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() : AliPHOSTrackSegmentMaker()
72 {
73   // ctor
74   SetTitle("version 1") ;
75   SetName("AliPHOSTrackSegmentMaker") ;
76   fR0 = 10. ;   
77   fEmcFirst = 0 ;    
78   fEmcLast  = 0 ;   
79   fCpvFirst = 0 ;   
80   fCpvLast  = 0 ;   
81   fPpsdFirst= 0 ;   
82   fPpsdLast = 0 ;   
83   fLinkLowArray = 0 ;
84   fLinkUpArray  = 0 ;
85   fIsInitialized = kFALSE ;
86 }
87 //____________________________________________________________________________
88  AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const char* headerFile, const char* branchTitle): 
89 AliPHOSTrackSegmentMaker()
90 {
91   // ctor
92   SetTitle("version 1") ;
93   SetName("AliPHOSTrackSegmentMaker") ;
94   fR0 = 10. ;   
95   fEmcFirst = 0 ;    
96   fEmcLast  = 0 ;   
97   fCpvFirst = 0 ;   
98   fCpvLast  = 0 ;   
99   fPpsdFirst= 0 ;   
100   fPpsdLast = 0 ;   
101
102   fHeaderFileName = headerFile ;
103   fRecPointsBranchTitle = branchTitle ;
104     
105   TFile * file = (TFile*) gROOT->GetFile(fHeaderFileName.Data() ) ;
106   
107   if(file == 0){
108     file = new TFile(fHeaderFileName.Data(),"update") ;
109     gAlice = (AliRun *) file->Get("gAlice") ;
110   }
111   
112   AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;    
113   fGeom  = AliPHOSGeometry::GetInstance(phos->GetGeometry()->GetName(),phos->GetGeometry()->GetTitle() );
114   
115   fEmcRecPoints = new TObjArray(200) ;
116   fCpvRecPoints = new TObjArray(200) ;
117   fClusterizer  = new AliPHOSClusterizerv1() ;
118   
119   fTrackSegments = new TClonesArray("AliPHOSTrackSegment",200) ;
120   
121   fLinkLowArray = new TClonesArray("AliPHOSLink", 1000);
122   fLinkUpArray  = new TClonesArray("AliPHOSLink", 1000); 
123   
124   fIsInitialized = kTRUE ;
125
126 }
127 //____________________________________________________________________________
128 void  AliPHOSTrackSegmentMakerv1::Init()
129 {
130   // Make all memory allocations that are not possible in default constructor
131
132   if(!fIsInitialized){
133     if(fHeaderFileName.IsNull())
134       fHeaderFileName = "galice.root" ;
135     
136     
137     TFile * file = (TFile*) gROOT->GetFile(fHeaderFileName.Data() ) ;
138     
139     if(file == 0){
140       file = new TFile(fHeaderFileName.Data(),"update") ;
141       gAlice = (AliRun *) file->Get("gAlice") ;
142     }
143     
144     AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;    
145     fGeom  = AliPHOSGeometry::GetInstance(phos->GetGeometry()->GetName(),phos->GetGeometry()->GetTitle() );
146
147
148     fEmcRecPoints = new TObjArray(200) ;
149     fCpvRecPoints = new TObjArray(200) ;
150     fClusterizer  = new AliPHOSClusterizerv1() ;
151
152     
153     fTrackSegments = new TClonesArray("AliPHOSTrackSegment",200) ;
154
155     fLinkLowArray = new TClonesArray("AliPHOSLink", 1000);
156     fLinkUpArray  = new TClonesArray("AliPHOSLink", 1000); 
157     
158     fIsInitialized = kTRUE ;
159    }
160 }
161
162 //____________________________________________________________________________
163  AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
164
165   // dtor
166   if(fLinkLowArray) delete fLinkLowArray ;
167   if(fLinkUpArray)  delete fLinkUpArray  ;
168 }
169
170 //____________________________________________________________________________
171 void  AliPHOSTrackSegmentMakerv1::FillOneModule()
172 {
173   // Finds first and last indexes between which 
174   // clusters from one PHOS module are
175  
176
177   //First EMC clusters
178   Int_t totalEmc = fEmcRecPoints->GetEntriesFast() ;
179   for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&  
180         (((AliPHOSRecPoint *) fEmcRecPoints->At(fEmcLast))->GetPHOSMod() == fModule ); 
181       fEmcLast ++)  ;
182   
183   
184   //Now CPV clusters
185   Int_t totalCpv = fCpvRecPoints->GetEntriesFast() ;
186
187   if(fModule <= fGeom->GetNCPVModules()){ // in CPV geometry
188     
189     for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) && 
190           (((AliPHOSRecPoint *) fCpvRecPoints->At(fCpvLast))->GetPHOSMod() == fModule ); 
191         fCpvLast ++) ;
192     
193     fPpsdFirst = fCpvLast ; //To avoid scanning RecPoints between fPpsdFirst and fPpsdLast
194     fPpsdLast  = fCpvLast ; //and to be ready to switch to mixed geometry 
195   }
196   else{  //in PPSD geometry    
197     fCpvLast = fPpsdLast ;
198     //Upper layer first
199     for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&  
200           (((AliPHOSPpsdRecPoint *) fCpvRecPoints->At(fCpvLast))->GetPHOSMod() == fModule ) &&
201           (((AliPHOSPpsdRecPoint *) fCpvRecPoints->At(fCpvLast))->GetUp()) ; 
202         fCpvLast ++)  ;
203     
204     fPpsdLast= fCpvLast ;
205     for(fPpsdFirst = fPpsdLast; (fPpsdLast < totalCpv)  &&
206           (((AliPHOSPpsdRecPoint *) fCpvRecPoints->At(fPpsdLast))->GetPHOSMod() == fModule ) &&
207           (!((AliPHOSPpsdRecPoint *) fCpvRecPoints->At(fPpsdLast))->GetUp()) ; 
208         fPpsdLast ++) ;
209   }
210     
211 }
212 //____________________________________________________________________________
213 Float_t  AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSRecPoint * cpvClu, Bool_t &toofar)const
214 {
215   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
216   // Clusters are sorted in "rows" and "columns" of width 1 cm
217
218   Float_t delta = 1 ;  // Width of the rows in sorting of RecPoints (in cm)
219                        // if you change this value, change it as well in xxxRecPoint::Compare()
220   Float_t r = fR0 ;
221  
222   TVector3 vecEmc ;
223   TVector3 vecCpv ;
224   
225   emcClu->GetLocalPosition(vecEmc) ;
226   cpvClu->GetLocalPosition(vecCpv)  ; 
227
228   if(emcClu->GetPHOSMod() == cpvClu->GetPHOSMod()){ 
229     if(vecCpv.X() <= vecEmc.X() + fR0 + 2*delta ){ 
230
231       vecCpv = vecCpv  - vecEmc ; 
232       r = vecCpv.Mag() ;
233       toofar = kFALSE ;
234
235     } // if  xPpsd >= xEmc + ...
236     else 
237       toofar = kTRUE ;
238   } 
239   else 
240     toofar = kTRUE ;
241
242   //toofar = kFALSE ;
243  
244   
245   return r ;
246 }
247
248 //____________________________________________________________________________
249 void  AliPHOSTrackSegmentMakerv1::MakeLinks()const
250
251   // Finds distances (links) between all EMC and PPSD clusters, 
252   // which are not further apart from each other than fR0 
253   // and sort them in accordance with this distance
254   
255   fLinkUpArray->Clear() ;    
256   fLinkLowArray->Clear() ;
257
258   AliPHOSRecPoint * ppsd ; 
259   AliPHOSRecPoint * cpv ;
260   AliPHOSEmcRecPoint * emcclu ;
261
262   Int_t iLinkLow = 0 ;
263   Int_t iLinkUp  = 0 ;
264   
265   Int_t iEmcRP;
266   for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
267     emcclu = (AliPHOSEmcRecPoint *) fEmcRecPoints->At(iEmcRP) ;
268
269     Bool_t toofar ;    
270     Int_t iPpsd ;
271     for(iPpsd = fPpsdFirst; iPpsd < fPpsdLast;iPpsd++ ) {
272       
273       ppsd = (AliPHOSRecPoint *) fCpvRecPoints->At(iPpsd) ;
274       Float_t r = GetDistanceInPHOSPlane(emcclu, ppsd, toofar) ;
275
276       if(toofar) 
277         break ;  
278       if(r < fR0)
279         new ((*fLinkLowArray)[iLinkLow++])  AliPHOSLink(r, iEmcRP, iPpsd) ;
280     }
281     
282     Int_t iCpv = 0 ;    
283     for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) { 
284       
285       cpv = (AliPHOSRecPoint *) fCpvRecPoints->At(iCpv) ;
286       Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, toofar) ;
287       
288       if(toofar)
289         break ;  
290       if(r < fR0) { 
291         new ((*fLinkUpArray)[iLinkUp++])  AliPHOSLink(r, iEmcRP, iCpv) ;
292       }      
293     }
294   } 
295   
296   fLinkLowArray->Sort() ; //first links with smallest distances
297   fLinkUpArray->Sort() ;
298 }
299
300 //____________________________________________________________________________
301 void  AliPHOSTrackSegmentMakerv1::MakePairs()
302
303   // Using the previously made list of "links", we found the smallest link - i.e. 
304   // link with the least distance between EMC and CPV and pointing to still 
305   // unassigned RecParticles. We assign these RecPoints to TrackSegment and 
306   // remove them from the list of "unassigned". 
307   
308   //Make arrays to mark clusters already chousen
309   Int_t * emcExist = 0;
310   if(fEmcLast > fEmcFirst)
311     emcExist = new Int_t[fEmcLast-fEmcFirst] ;
312   
313   Int_t index;
314   for(index = 0; index <fEmcLast-fEmcFirst; index ++)
315     emcExist[index] = 1 ;
316   
317   Bool_t * cpvExist = 0;
318   if(fCpvLast > fCpvFirst)
319     cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
320   for(index = 0; index <fCpvLast-fCpvFirst; index ++)
321     cpvExist[index] = kTRUE ;
322   
323   Bool_t * ppsdExist = 0;
324   if(fPpsdLast > fPpsdFirst)
325     ppsdExist = new Bool_t[fPpsdLast-fPpsdFirst] ;
326   for(index = 0; index <fPpsdLast-fPpsdFirst; index ++)
327     ppsdExist[index] = kTRUE ;
328   
329   // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance 
330   TIter nextLow(fLinkLowArray) ;
331   TIter nextUp(fLinkUpArray) ;
332   
333   AliPHOSLink * linkLow ;
334   AliPHOSLink * linkUp ;
335
336
337   AliPHOSRecPoint * nullpointer = 0 ;
338
339   while ( (linkLow =  (AliPHOSLink *)nextLow() ) ){
340   
341     if( (emcExist[linkLow->GetEmc()-fEmcFirst]> 0) && 
342         ppsdExist[linkLow->GetPpsd()-fPpsdFirst]  ){ // RecPoints not removed yet 
343       new ((*fTrackSegments)[fNTrackSegments]) AliPHOSTrackSegment((AliPHOSEmcRecPoint *) fEmcRecPoints->At(linkLow->GetEmc()), 
344                                                  nullpointer, 
345                                                 (AliPHOSRecPoint *)fCpvRecPoints->At(linkLow->GetPpsd()) ) ;
346          
347       ((AliPHOSTrackSegment* )fTrackSegments->At(fNTrackSegments))->SetIndexInList(fNTrackSegments);    
348       //replace index of emc to negative and shifted index of TS      
349       emcExist[linkLow->GetEmc()-fEmcFirst] = -2 - fNTrackSegments ;  
350       //mark ppsd as used
351       ppsdExist[linkLow->GetPpsd()-fPpsdFirst] = kFALSE ; 
352       fNTrackSegments++ ;
353     } 
354   } 
355          
356
357   while ( (linkUp =  (AliPHOSLink *)nextUp() ) ){  
358     if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){ //without ppsd Up yet 
359
360       if(cpvExist[linkUp->GetPpsd()-fCpvFirst]){ //CPV still exist
361         
362         if(emcExist[linkUp->GetEmc()-fEmcFirst] > 0){ //without ppsd Low => create new TS
363
364           new ((* fTrackSegments)[fNTrackSegments]) AliPHOSTrackSegment((AliPHOSEmcRecPoint *) fEmcRecPoints->At(linkUp->GetEmc()) , 
365                                                                       (AliPHOSRecPoint *)fCpvRecPoints->At(linkUp->GetPpsd()), 
366                                                                       nullpointer) ;
367           ((AliPHOSTrackSegment *) fTrackSegments->At(fNTrackSegments))->SetIndexInList(fNTrackSegments);
368           fNTrackSegments++ ;
369         }
370         else{ // append ppsd Up to existing TS
371           ((AliPHOSTrackSegment *)fTrackSegments->At(-2-emcExist[linkUp->GetEmc()-fEmcFirst]))->SetCpvRecPoint((AliPHOSCpvRecPoint *)fCpvRecPoints->At(linkUp->GetPpsd()));
372         }
373
374         emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc  that Cpv was found 
375         //mark CPV recpoint as already used 
376         cpvExist[linkUp->GetPpsd()-fCpvFirst] = kFALSE ;
377       } //if ppsdUp still exist
378     } 
379   }      
380
381   //look through emc recPoints left without CPV/PPSD
382   if(emcExist){ //if there is emc rec point
383     Int_t iEmcRP ;
384     for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst  ; iEmcRP++ ){
385       if(emcExist[iEmcRP] > 0 ){
386         new ((*fTrackSegments)[fNTrackSegments])  AliPHOSTrackSegment((AliPHOSEmcRecPoint *)fEmcRecPoints->At(iEmcRP+fEmcFirst), 
387                                                                     nullpointer, 
388                                                                     nullpointer ) ;
389         ((AliPHOSTrackSegment *) fTrackSegments->At(fNTrackSegments))->SetIndexInList(fNTrackSegments);
390         fNTrackSegments++;    
391       } 
392     }
393   }
394   
395 }
396
397 //____________________________________________________________________________
398 void  AliPHOSTrackSegmentMakerv1::Exec(Option_t * option)
399 {
400   // STEERing method
401
402   if(! fIsInitialized) Init() ;
403
404   if(strstr(option,"tim"))
405     gBenchmark->Start("PHOSTSMaker");  
406
407   Int_t nEvents = (Int_t) gAlice->TreeE()->GetEntries() ;
408   
409   for(fEvent = 0;fEvent< nEvents; fEvent++){
410     if(!ReadRecPoints())  //reads RecPoints for event fEvent
411       return;
412     
413     for(fModule = 1; fModule <= fGeom->GetNModules() ; fModule++ ){
414       
415       FillOneModule() ; 
416       
417       MakeLinks() ;
418       
419       MakePairs() ;
420       
421     }
422
423     WriteTrackSegments() ;
424     if(strstr(option,"deb"))
425       PrintTrackSegments(option) ;
426   }
427
428   if(strstr(option,"tim")){
429     gBenchmark->Stop("PHOSTSMaker");
430     cout << "AliPHOSTSMaker:" << endl ;
431     cout << "  took " << gBenchmark->GetCpuTime("PHOSTSMaker") << " seconds for making TS " 
432          <<  gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents << " seconds per event " << endl ;
433     cout << endl ;
434   }
435
436
437 }
438 //____________________________________________________________________________
439 void AliPHOSTrackSegmentMakerv1::Print(Option_t * option)const
440 {
441   //  Print TrackSegmentMaker parameters
442
443   if(fIsInitialized){
444     cout <<  "======== AliPHOSTrackSegmentMakerv1 ========" << endl ;
445     cout <<  "Making Track segments "<< endl ;
446     cout <<  "    Headers file:                   " << fHeaderFileName.Data() << endl ;
447     cout <<  "    RecPoints branch file name:     " << fRecPointsBranchTitle.Data() << endl ;
448     cout <<  "    TrackSegments Branch file name: " << fTSBranchTitle.Data() << endl ;
449     cout <<  "with parameters: " << endl ;
450     cout <<  "    Maximal EMC - CPV (PPSD) distance (cm)" << fR0 << endl ;
451     cout <<  "============================================" << endl ;
452   }
453   else
454     cout << "AliPHOSTrackSegmentMakerv1 not initialized " << endl ;
455 }
456 //____________________________________________________________________________
457 Bool_t AliPHOSTrackSegmentMakerv1::ReadRecPoints()
458 {
459   // Reads Emc and CPV recPoints with given title (fRecPointsBranchTitle) 
460   // made previously with Clusterizer.
461
462
463   //Make some initializations 
464   fEmcRecPoints->Clear() ;
465   fCpvRecPoints->Clear() ;
466   fTrackSegments->Clear() ;
467   fNTrackSegments = 0 ;
468   fEmcFirst = 0 ;    
469   fEmcLast  = 0 ;   
470   fCpvFirst = 0 ;   
471   fCpvLast  = 0 ;   
472   fPpsdFirst= 0 ;   
473   fPpsdLast = 0 ;   
474
475
476   gAlice->GetEvent(fEvent) ;
477
478   // Get TreeR header from file
479   if(gAlice->TreeR()==0){
480     char treeName[20]; 
481     sprintf(treeName,"TreeR%d",fEvent);
482     cout << "Error in AliPHOSTrackSegmentMakerv1 : no "<<treeName << endl  ;
483     cout << "   Do nothing " << endl ;
484     return kFALSE ;
485   }
486
487
488   //Find RecPoints with title fRecPointsBranchTitle
489   TBranch * emcBranch = 0;
490   TBranch * cpvBranch = 0;
491   TBranch * clusterizerBranch = 0;
492
493   TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
494   Int_t ibranch;
495   Bool_t emcNotFound = kTRUE ;
496   Bool_t cpvNotFound = kTRUE ;  
497   Bool_t clusterizerNotFound = kTRUE ;
498   
499   for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
500
501     if(emcNotFound){
502       emcBranch=(TBranch *) branches->At(ibranch) ;
503       if( fRecPointsBranchTitle.CompareTo(emcBranch->GetTitle())==0 )
504         if( strcmp(emcBranch->GetName(),"PHOSEmcRP") == 0) {
505           emcNotFound = kFALSE ;
506         }
507     }
508     
509     if(cpvNotFound){
510       cpvBranch=(TBranch *) branches->At(ibranch) ;
511       if( fRecPointsBranchTitle.CompareTo(cpvBranch->GetTitle())==0 )
512         if( strcmp(cpvBranch->GetName(),"PHOSCpvRP") == 0) 
513           cpvNotFound = kFALSE ;
514     }
515     
516     if(clusterizerNotFound){
517       clusterizerBranch = (TBranch *) branches->At(ibranch) ;
518       if( fRecPointsBranchTitle.CompareTo(clusterizerBranch->GetTitle()) == 0)
519         if( strcmp(clusterizerBranch->GetName(),"AliPHOSClusterizer") == 0) 
520           clusterizerNotFound = kFALSE ;
521     }
522     
523   }
524
525   if(clusterizerNotFound || emcNotFound || cpvNotFound){
526     cout << "AliPHOSTrackSegmentMakerv1: " << endl ;
527     cout << "    Can't find Branch with RecPoints or Clusterizer " ;
528     cout << "    Do nothing" <<endl  ;
529     return kFALSE ;
530   }
531   
532   emcBranch->SetAddress(&fEmcRecPoints) ;
533   cpvBranch->SetAddress(&fCpvRecPoints) ;
534   clusterizerBranch->SetAddress(&fClusterizer) ;
535
536   emcBranch->GetEntry(0) ;
537   cpvBranch->GetEntry(0) ;
538   clusterizerBranch->GetEntry(0) ;
539   
540   return kTRUE ;
541   
542 }
543 //____________________________________________________________________________
544 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
545 {
546   // Writes found TrackSegments to TreeR. Creates branches 
547   // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
548   // In the former branch found TrackSegments are stored, while 
549   // in the latter all parameters, with which TS were made. 
550   // ROOT does not allow overwriting existing branches, therefore
551   // first we check, if branches with the same title already exist.
552   // If yes - exits without writing.
553   
554   //First, check, if branches already exist
555   TBranch * tsMakerBranch = 0;
556   TBranch * tsBranch = 0;
557   
558   TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
559   Int_t ibranch;
560   Bool_t tsMakerNotFound = kTRUE ;
561   Bool_t tsNotFound = kTRUE ;
562   
563   for(ibranch = 0;(ibranch <branches->GetEntries())&&(tsMakerNotFound||tsNotFound);ibranch++){
564     if(tsMakerNotFound){
565       tsMakerBranch=(TBranch *) branches->At(ibranch) ;
566       if( (strcmp(tsMakerBranch->GetName(),"AliPHOSTrackSegmentMaker") == 0) &&
567           (fTSBranchTitle.CompareTo( tsMakerBranch->GetTitle())==0 ))
568         tsMakerNotFound = kFALSE ;
569     }
570     if(tsNotFound){
571       tsBranch=(TBranch *) branches->At(ibranch) ;
572       if( (strcmp(tsBranch->GetName(),"PHOSTS") == 0)  &&
573           (fTSBranchTitle.CompareTo( tsBranch->GetTitle())==0 ))
574         tsNotFound = kFALSE ;
575     }
576   }
577
578   if(!(tsMakerNotFound && tsNotFound )){ 
579     cout << "AliPHOSTrackSegmentMakerv1 error:"<< endl ;
580     cout << "       Branches PHOSTS and AliPHOSTrackSegementMaker " << endl ;
581     cout << "       with title '"<<fTSBranchTitle.Data() << "' already exist " << endl ;
582     cout << "       can not overwrite " << endl ;
583     return ;
584   }
585
586   //Make branch in TreeR for TrackSegments 
587   char * filename = 0;
588   if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){   //generating file name
589     filename = new char[strlen(gAlice->GetBaseFile())+20] ;
590     sprintf(filename,"%s/PHOS.Reco.root",gAlice->GetBaseFile()) ; 
591   }
592
593   TDirectory *cwd = gDirectory;
594   
595   //First TS
596   Int_t bufferSize = 32000 ;    
597   tsBranch = gAlice->TreeR()->Branch("PHOSTS",&fTrackSegments,bufferSize);
598   tsBranch->SetTitle(fTSBranchTitle.Data());
599   if (filename) {
600     tsBranch->SetFile(filename);
601     TIter next( tsBranch->GetListOfBranches());
602     TBranch * sb ;
603     while ((sb=(TBranch*)next())) {
604       sb->SetFile(filename);
605     }   
606     cwd->cd();
607   } 
608   
609   //Second -TSMaker
610   Int_t splitlevel = 0 ;
611   AliPHOSTrackSegmentMakerv1 * ts = this ;
612   tsMakerBranch = gAlice->TreeR()->Branch("AliPHOSTrackSegmentMaker","AliPHOSTrackSegmentMakerv1",
613                                           &ts,bufferSize,splitlevel);
614   tsMakerBranch->SetTitle(fTSBranchTitle.Data());
615   if (filename) {
616     tsMakerBranch->SetFile(filename);
617     TIter next( tsMakerBranch->GetListOfBranches());
618     TBranch * sb;
619     while ((sb=(TBranch*)next())) {
620       sb->SetFile(filename);
621     }   
622     cwd->cd();
623   } 
624   
625   tsBranch->Fill() ;  
626   tsMakerBranch->Fill() ;
627
628   gAlice->TreeR()->Write(0,kOverwrite) ;  
629   
630 }
631
632
633 //____________________________________________________________________________
634 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
635 {
636   // option deb - prints # of found TrackSegments
637   // option deb all - prints as well indexed of found RecParticles assigned to the TS
638
639   
640   cout << "AliPHOSTrackSegmentMakerv1: " << endl ;
641   cout << "       Found " << fTrackSegments->GetEntriesFast() << "  trackSegments " << endl ;
642   
643   if(strstr(option,"all")) {  // printing found TS
644     cout << "TrackSegment # " << "    EMC RP#    " << "    CPV RP#    " << "     PPSD RP#" << endl ; 
645     
646     Int_t index;
647     for (index = 0 ; index <fTrackSegments->GetEntriesFast() ; index++) {
648       AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )fTrackSegments->At(index) ; 
649       cout<<"   "<< setw(4) << ts->GetIndexInList() << "            " 
650           <<setw(4) << ts->GetEmcIndex()<< "            " 
651           <<setw(4) << ts->GetCpvIndex()<< "            " 
652           <<setw(4) << ts->GetPpsdIndex()<< endl ;
653     }   
654     
655     cout << "-------------------------------------------------------"<< endl ;
656   }
657 }