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