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