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