1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
29 // If there is no CPV/PPSD RecPoint we make TrackSegment
30 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
32 // In principle this class should be called from AliPHOSReconstructioner, but
33 // one can use it as well in standalone mode.
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")
42 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH)
45 // --- ROOT system ---
51 #include "TBenchmark.h"
52 // --- Standard library ---
57 // --- AliRoot header files ---
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"
69 ClassImp( AliPHOSTrackSegmentMakerv1)
72 //____________________________________________________________________________
73 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() : AliPHOSTrackSegmentMaker()
75 // default ctor (to be used mainly by Streamer)
86 fHeaderFileName = "" ;
87 fRecPointsBranchTitle = "" ;
88 fTrackSegmentsBranchTitle = "" ;
89 fTrackSegmentsInRun = 0 ;
92 //____________________________________________________________________________
93 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const char * headerFile, const char * name) : AliPHOSTrackSegmentMaker(headerFile, name)
105 fHeaderFileName = GetTitle() ;
106 fRecPointsBranchTitle = GetName() ;
107 fTrackSegmentsBranchTitle = GetName() ;
108 fTrackSegmentsInRun = 0 ;
110 TString tempo(GetName()) ;
112 tempo.Append(Version()) ;
113 SetName(tempo.Data()) ;
119 //____________________________________________________________________________
120 AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
123 if(fLinkLowArray) delete fLinkLowArray ;
124 if(fLinkUpArray) delete fLinkUpArray ;
127 //____________________________________________________________________________
128 void AliPHOSTrackSegmentMakerv1::FillOneModule()
130 // Finds first and last indexes between which
131 // clusters from one PHOS module are
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) ;
141 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
142 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
143 (((AliPHOSRecPoint *) emcRecPoints->At(fEmcLast))->GetPHOSMod() == fModule );
147 Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
149 if(fModule <= geom->GetNCPVModules()){ // in CPV geometry
151 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
152 (((AliPHOSRecPoint *) cpvRecPoints->At(fCpvLast))->GetPHOSMod() == fModule );
155 fPpsdFirst = fCpvLast ; //To avoid scanning RecPoints between fPpsdFirst and fPpsdLast
156 fPpsdLast = fCpvLast ; //and to be ready to switch to mixed geometry
158 else{ //in PPSD geometry
159 fCpvLast = fPpsdLast ;
161 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
162 (((AliPHOSPpsdRecPoint *) cpvRecPoints->At(fCpvLast))->GetPHOSMod() == fModule ) &&
163 (((AliPHOSPpsdRecPoint *) cpvRecPoints->At(fCpvLast))->GetUp()) ;
166 fPpsdLast= fCpvLast ;
167 for(fPpsdFirst = fPpsdLast; (fPpsdLast < totalCpv) &&
168 (((AliPHOSPpsdRecPoint *) cpvRecPoints->At(fPpsdLast))->GetPHOSMod() == fModule ) &&
169 (!((AliPHOSPpsdRecPoint *) cpvRecPoints->At(fPpsdLast))->GetUp()) ;
175 //____________________________________________________________________________
176 Float_t AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSRecPoint * cpvClu, Bool_t &toofar)const
178 // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
179 // Clusters are sorted in "rows" and "columns" of width 1 cm
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()
188 emcClu->GetLocalPosition(vecEmc) ;
189 cpvClu->GetLocalPosition(vecCpv) ;
191 if(emcClu->GetPHOSMod() == cpvClu->GetPHOSMod()){
192 if(vecCpv.X() <= vecEmc.X() + fR0 + 2*delta ){
194 vecCpv = vecCpv - vecEmc ;
198 } // if xPpsd >= xEmc + ...
211 //____________________________________________________________________________
212 void AliPHOSTrackSegmentMakerv1::Init()
214 // Make all memory allocations that are not possible in default constructor
216 if ( strcmp(GetTitle(), "") == 0 )
217 SetTitle("galice.root") ;
219 TString taskName(GetName()) ;
220 taskName.Remove(taskName.Index(Version())-1) ;
222 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), taskName) ;
224 cerr << "ERROR: AliPHOSTrackSegmentMakerv1::Init -> Could not obtain the Getter object !" << endl ;
228 fLinkLowArray = new TClonesArray("AliPHOSLink", 1000);
229 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
231 //add Task to //YSAlice/tasks/Reconstructioner/PHOS
232 gime->PostTrackSegmentMaker(this) ;
234 // create a folder on the white board //YSAlice/WhiteBoard/RecPoints/PHOS/trackSegmentsName
235 gime->PostTrackSegments(taskName) ;
239 //____________________________________________________________________________
240 void AliPHOSTrackSegmentMakerv1::MakeLinks()const
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
246 TString taskName(GetName()) ;
247 taskName.Remove(taskName.Index(Version())-1) ;
249 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
250 TObjArray * emcRecPoints = gime->EmcRecPoints(taskName) ;
251 TObjArray * cpvRecPoints = gime->CpvRecPoints(taskName) ;
253 fLinkUpArray->Clear() ;
254 fLinkLowArray->Clear() ;
256 AliPHOSRecPoint * ppsd ;
257 AliPHOSRecPoint * cpv ;
258 AliPHOSEmcRecPoint * emcclu ;
264 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
265 emcclu = (AliPHOSEmcRecPoint *) emcRecPoints->At(iEmcRP) ;
269 for(iPpsd = fPpsdFirst; iPpsd < fPpsdLast;iPpsd++ ) {
271 ppsd = (AliPHOSRecPoint *) cpvRecPoints->At(iPpsd) ;
272 Float_t r = GetDistanceInPHOSPlane(emcclu, ppsd, toofar) ;
277 new ((*fLinkLowArray)[iLinkLow++]) AliPHOSLink(r, iEmcRP, iPpsd) ;
281 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
283 cpv = (AliPHOSRecPoint *) cpvRecPoints->At(iCpv) ;
284 Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, toofar) ;
289 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(r, iEmcRP, iCpv) ;
294 fLinkLowArray->Sort() ; //first links with smallest distances
295 fLinkUpArray->Sort() ;
298 //____________________________________________________________________________
299 void AliPHOSTrackSegmentMakerv1::MakePairs()
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".
306 TString taskName(GetName()) ;
307 taskName.Remove(taskName.Index(Version())-1) ;
309 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
310 TObjArray * emcRecPoints = gime->EmcRecPoints(taskName) ;
311 TObjArray * cpvRecPoints = gime->CpvRecPoints(taskName) ;
312 TClonesArray * trackSegments = gime->TrackSegments(taskName) ;
314 //Make arrays to mark clusters already chosen
315 Int_t * emcExist = 0;
316 if(fEmcLast > fEmcFirst)
317 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
320 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
321 emcExist[index] = 1 ;
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 ;
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 ;
335 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
336 TIter nextLow(fLinkLowArray) ;
337 TIter nextUp(fLinkUpArray) ;
339 AliPHOSLink * linkLow ;
340 AliPHOSLink * linkUp ;
343 AliPHOSRecPoint * nullpointer = 0 ;
345 while ( (linkLow = (AliPHOSLink *)nextLow() ) ){
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()),
351 (AliPHOSRecPoint *)cpvRecPoints->At(linkLow->GetPpsd()) ) ;
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 ;
357 ppsdExist[linkLow->GetPpsd()-fPpsdFirst] = kFALSE ;
363 while ( (linkUp = (AliPHOSLink *)nextUp() ) ){
364 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){ //without ppsd Up yet
366 if(cpvExist[linkUp->GetPpsd()-fCpvFirst]){ //CPV still exist
368 if(emcExist[linkUp->GetEmc()-fEmcFirst] > 0){ //without ppsd Low => create new TS
370 new ((* trackSegments)[fNTrackSegments]) AliPHOSTrackSegment((AliPHOSEmcRecPoint *) emcRecPoints->At(linkUp->GetEmc()) ,
371 (AliPHOSRecPoint *)cpvRecPoints->At(linkUp->GetPpsd()),
373 ((AliPHOSTrackSegment *) trackSegments->At(fNTrackSegments))->SetIndexInList(fNTrackSegments);
376 else{ // append ppsd Up to existing TS
377 ((AliPHOSTrackSegment *)trackSegments->At(-2-emcExist[linkUp->GetEmc()-fEmcFirst]))->SetCpvRecPoint((AliPHOSCpvRecPoint *)cpvRecPoints->At(linkUp->GetPpsd()));
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
387 //look through emc recPoints left without CPV/PPSD
388 if(emcExist){ //if there is emc rec point
390 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
391 if(emcExist[iEmcRP] > 0 ){
392 new ((*trackSegments)[fNTrackSegments]) AliPHOSTrackSegment((AliPHOSEmcRecPoint *)emcRecPoints->At(iEmcRP+fEmcFirst),
395 ((AliPHOSTrackSegment *) trackSegments->At(fNTrackSegments))->SetIndexInList(fNTrackSegments);
403 //____________________________________________________________________________
404 void AliPHOSTrackSegmentMakerv1::Exec(Option_t * option)
408 if( strcmp(GetName(), "")== 0 )
411 if(strstr(option,"tim"))
412 gBenchmark->Start("PHOSTSMaker");
414 if(strstr(option,"print")) {
419 //check, if the branch with name of this" already exits?
420 gAlice->GetEvent(0) ;
421 TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
423 TBranch * branch = 0 ;
424 Bool_t phostsfound = kFALSE, tracksegmentmakerfound = kFALSE ;
426 TString taskName(GetName()) ;
427 taskName.Remove(taskName.Index(Version())-1) ;
429 while ( (branch = (TBranch*)next()) && (!phostsfound || !tracksegmentmakerfound) ) {
430 if ( (strcmp(branch->GetName(), "PHOSTS")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
431 phostsfound = kTRUE ;
433 else if ( (strcmp(branch->GetName(), "AliPHOSTrackSegmentMaker")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
434 tracksegmentmakerfound = kTRUE ;
437 if ( phostsfound || tracksegmentmakerfound ) {
438 cerr << "WARNING: AliPHOSTrackSegmentMakerv1::Exec -> TrackSegments and/or TrackSegmentMaker branch with name "
439 << taskName.Data() << " already exits" << endl ;
443 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
444 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
445 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
448 for(ievent = 0; ievent < nevents; ievent++){
450 // if(!ReadRecPoints(ievent)) //reads RecPoints for event ievent
452 gime->Event(ievent,"R") ;
454 //Make some initializations
455 fNTrackSegments = 0 ;
462 gime->TrackSegments(taskName)->Clear() ;
464 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ){
474 WriteTrackSegments(ievent) ;
476 if(strstr(option,"deb"))
477 PrintTrackSegments(option) ;
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 ;
491 //____________________________________________________________________________
492 void AliPHOSTrackSegmentMakerv1::Print(Option_t * option)const
494 // Print TrackSegmentMaker parameters
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 ;
507 cout << "AliPHOSTrackSegmentMakerv1 not initialized " << endl ;
509 //____________________________________________________________________________
510 Bool_t AliPHOSTrackSegmentMakerv1::ReadRecPoints(Int_t event)
512 // Reads Emc and CPV recPoints
513 // made previously with Clusterizer.
516 //Make some initializations
518 fNTrackSegments = 0 ;
526 // Get TreeR header from file
527 if(gAlice->TreeR()==0){
528 cerr << "ERROR: AliPHOSTrackSegmentMakerv1::ReadRecPoints -> There is no Reconstruction Tree" << endl;
534 TBranch * emcbranch = 0;
535 TBranch * cpvbranch = 0;
536 TBranch * clusterizerbranch = 0;
537 TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
539 TBranch * branch = 0 ;
540 Bool_t phosemcfound = kFALSE, phoscpvfound = kFALSE, clusterizerfound = kFALSE ;
542 TString taskName(GetName()) ;
543 taskName.ReplaceAll(Version(), "") ;
545 while ( (branch = (TBranch*)next()) && (!phosemcfound || !phoscpvfound || !clusterizerfound) ) {
546 if ( (strcmp(branch->GetName(), "PHOSEmcRP")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
547 phosemcfound = kTRUE ;
551 else if ( (strcmp(branch->GetName(), "PHOSCpvRP")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
552 phoscpvfound = kTRUE ;
555 } else if ( (strcmp(branch->GetName(), "AliPHOSClusterizer")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
556 clusterizerfound = kTRUE ;
557 clusterizerbranch = branch ;
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 ;
566 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
568 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
569 emcRecPoints->Clear() ;
570 emcbranch->SetAddress(&emcRecPoints) ;
572 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
573 cpvRecPoints->Clear() ;
574 cpvbranch->SetAddress(&cpvRecPoints) ;
576 TClonesArray * trackSegments = gime->TrackSegments() ;
577 trackSegments->Clear() ;
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) ;
586 emcbranch->GetEntry(0) ;
587 cpvbranch->GetEntry(0) ;
589 clusterizerbranch->GetEntry(0) ;
591 printf("ReadRecPoint: EMC=%d, CPV=%d\n",emcRecPoints->GetEntries(),cpvRecPoints->GetEntries());
597 //____________________________________________________________________________
598 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments(Int_t event)
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.
608 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
609 TString taskName(GetName()) ;
610 taskName.Remove(taskName.Index(Version())-1) ;
612 TClonesArray * trackSegments = gime->TrackSegments(taskName) ;
613 trackSegments->Expand(trackSegments->GetEntriesFast()) ;
615 //Make branch in TreeR for TrackSegments
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()) ;
622 TDirectory *cwd = gDirectory;
625 Int_t bufferSize = 32000 ;
626 TBranch * tsBranch = gAlice->TreeR()->Branch("PHOSTS",&trackSegments,bufferSize);
627 tsBranch->SetTitle(fTrackSegmentsBranchTitle);
629 tsBranch->SetFile(filename);
630 TIter next( tsBranch->GetListOfBranches());
632 while ((sb=(TBranch*)next())) {
633 sb->SetFile(filename);
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);
645 tsMakerBranch->SetFile(filename);
646 TIter next( tsMakerBranch->GetListOfBranches());
648 while ((sb=(TBranch*)next())) {
649 sb->SetFile(filename);
655 tsMakerBranch->Fill() ;
657 gAlice->TreeR()->Write(0,kOverwrite) ;
662 //____________________________________________________________________________
663 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
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) ;
670 TClonesArray * trackSegments = AliPHOSGetter::GetInstance()->TrackSegments(taskName) ;
672 cout << "AliPHOSTrackSegmentMakerv1: event "<<gAlice->GetEvNumber() << endl ;
673 cout << " Found " << trackSegments->GetEntriesFast() << " trackSegments " << endl ;
675 fTrackSegmentsInRun += trackSegments->GetEntriesFast() ;
677 if(strstr(option,"all")) { // printing found TS
678 cout << "TrackSegment # " << " EMC RP# " << " CPV RP# " << " PPSD RP#" << endl ;
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 ;
689 cout << "-------------------------------------------------------"<< endl ;