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