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 = "" ;
91 //____________________________________________________________________________
92 AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const char * headerFile, const char * name) : AliPHOSTrackSegmentMaker(headerFile, name)
104 fHeaderFileName = GetTitle() ;
105 fRecPointsBranchTitle = GetName() ;
106 fTrackSegmentsBranchTitle = GetName() ;
108 TString tempo(GetName()) ;
109 tempo.Append(Version()) ;
110 SetName(tempo.Data()) ;
116 //____________________________________________________________________________
117 AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
120 if(fLinkLowArray) delete fLinkLowArray ;
121 if(fLinkUpArray) delete fLinkUpArray ;
124 //____________________________________________________________________________
125 void AliPHOSTrackSegmentMakerv1::FillOneModule()
127 // Finds first and last indexes between which
128 // clusters from one PHOS module are
130 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
131 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
132 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
133 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
136 Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
137 for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&
138 (((AliPHOSRecPoint *) emcRecPoints->At(fEmcLast))->GetPHOSMod() == fModule );
142 Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
144 if(fModule <= geom->GetNCPVModules()){ // in CPV geometry
146 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
147 (((AliPHOSRecPoint *) cpvRecPoints->At(fCpvLast))->GetPHOSMod() == fModule );
150 fPpsdFirst = fCpvLast ; //To avoid scanning RecPoints between fPpsdFirst and fPpsdLast
151 fPpsdLast = fCpvLast ; //and to be ready to switch to mixed geometry
153 else{ //in PPSD geometry
154 fCpvLast = fPpsdLast ;
156 for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) &&
157 (((AliPHOSPpsdRecPoint *) cpvRecPoints->At(fCpvLast))->GetPHOSMod() == fModule ) &&
158 (((AliPHOSPpsdRecPoint *) cpvRecPoints->At(fCpvLast))->GetUp()) ;
161 fPpsdLast= fCpvLast ;
162 for(fPpsdFirst = fPpsdLast; (fPpsdLast < totalCpv) &&
163 (((AliPHOSPpsdRecPoint *) cpvRecPoints->At(fPpsdLast))->GetPHOSMod() == fModule ) &&
164 (!((AliPHOSPpsdRecPoint *) cpvRecPoints->At(fPpsdLast))->GetUp()) ;
170 //____________________________________________________________________________
171 Float_t AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSRecPoint * cpvClu, Bool_t &toofar)const
173 // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
174 // Clusters are sorted in "rows" and "columns" of width 1 cm
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()
183 emcClu->GetLocalPosition(vecEmc) ;
184 cpvClu->GetLocalPosition(vecCpv) ;
186 if(emcClu->GetPHOSMod() == cpvClu->GetPHOSMod()){
187 if(vecCpv.X() <= vecEmc.X() + fR0 + 2*delta ){
189 vecCpv = vecCpv - vecEmc ;
193 } // if xPpsd >= xEmc + ...
206 //____________________________________________________________________________
207 void AliPHOSTrackSegmentMakerv1::Init()
209 // Make all memory allocations that are not possible in default constructor
211 if ( strcmp(GetTitle(), "") == 0 )
212 SetTitle("galice.root") ;
214 TString taskName(GetName()) ;
215 taskName.ReplaceAll(Version(), "") ;
217 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), taskName.Data()) ;
219 cerr << "ERROR: AliPHOSTrackSegmentMakerv1::Init -> Could not obtain the Getter object !" << endl ;
223 fLinkLowArray = new TClonesArray("AliPHOSLink", 1000);
224 fLinkUpArray = new TClonesArray("AliPHOSLink", 1000);
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") ;
230 // create a folder on the white board //YSAlice/WhiteBoard/RecPoints/PHOS/trackSegmentsName
232 gime->Post(GetTitle(), "T", taskName.Data() ) ;
236 //____________________________________________________________________________
237 void AliPHOSTrackSegmentMakerv1::MakeLinks()const
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
243 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
244 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
245 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
247 fLinkUpArray->Clear() ;
248 fLinkLowArray->Clear() ;
250 AliPHOSRecPoint * ppsd ;
251 AliPHOSRecPoint * cpv ;
252 AliPHOSEmcRecPoint * emcclu ;
258 for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
259 emcclu = (AliPHOSEmcRecPoint *) emcRecPoints->At(iEmcRP) ;
263 for(iPpsd = fPpsdFirst; iPpsd < fPpsdLast;iPpsd++ ) {
265 ppsd = (AliPHOSRecPoint *) cpvRecPoints->At(iPpsd) ;
266 Float_t r = GetDistanceInPHOSPlane(emcclu, ppsd, toofar) ;
271 new ((*fLinkLowArray)[iLinkLow++]) AliPHOSLink(r, iEmcRP, iPpsd) ;
275 for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) {
277 cpv = (AliPHOSRecPoint *) cpvRecPoints->At(iCpv) ;
278 Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, toofar) ;
283 new ((*fLinkUpArray)[iLinkUp++]) AliPHOSLink(r, iEmcRP, iCpv) ;
288 fLinkLowArray->Sort() ; //first links with smallest distances
289 fLinkUpArray->Sort() ;
292 //____________________________________________________________________________
293 void AliPHOSTrackSegmentMakerv1::MakePairs()
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".
300 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
301 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
302 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
303 TClonesArray * trackSegments = gime->TrackSegments() ;
305 //Make arrays to mark clusters already chosen
306 Int_t * emcExist = 0;
307 if(fEmcLast > fEmcFirst)
308 emcExist = new Int_t[fEmcLast-fEmcFirst] ;
311 for(index = 0; index <fEmcLast-fEmcFirst; index ++)
312 emcExist[index] = 1 ;
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 ;
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 ;
326 // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance
327 TIter nextLow(fLinkLowArray) ;
328 TIter nextUp(fLinkUpArray) ;
330 AliPHOSLink * linkLow ;
331 AliPHOSLink * linkUp ;
334 AliPHOSRecPoint * nullpointer = 0 ;
336 while ( (linkLow = (AliPHOSLink *)nextLow() ) ){
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()),
342 (AliPHOSRecPoint *)cpvRecPoints->At(linkLow->GetPpsd()) ) ;
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 ;
348 ppsdExist[linkLow->GetPpsd()-fPpsdFirst] = kFALSE ;
354 while ( (linkUp = (AliPHOSLink *)nextUp() ) ){
355 if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){ //without ppsd Up yet
357 if(cpvExist[linkUp->GetPpsd()-fCpvFirst]){ //CPV still exist
359 if(emcExist[linkUp->GetEmc()-fEmcFirst] > 0){ //without ppsd Low => create new TS
361 new ((* trackSegments)[fNTrackSegments]) AliPHOSTrackSegment((AliPHOSEmcRecPoint *) emcRecPoints->At(linkUp->GetEmc()) ,
362 (AliPHOSRecPoint *)cpvRecPoints->At(linkUp->GetPpsd()),
364 ((AliPHOSTrackSegment *) trackSegments->At(fNTrackSegments))->SetIndexInList(fNTrackSegments);
367 else{ // append ppsd Up to existing TS
368 ((AliPHOSTrackSegment *)trackSegments->At(-2-emcExist[linkUp->GetEmc()-fEmcFirst]))->SetCpvRecPoint((AliPHOSCpvRecPoint *)cpvRecPoints->At(linkUp->GetPpsd()));
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
378 //look through emc recPoints left without CPV/PPSD
379 if(emcExist){ //if there is emc rec point
381 for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst ; iEmcRP++ ){
382 if(emcExist[iEmcRP] > 0 ){
383 new ((*trackSegments)[fNTrackSegments]) AliPHOSTrackSegment((AliPHOSEmcRecPoint *)emcRecPoints->At(iEmcRP+fEmcFirst),
386 ((AliPHOSTrackSegment *) trackSegments->At(fNTrackSegments))->SetIndexInList(fNTrackSegments);
394 //____________________________________________________________________________
395 void AliPHOSTrackSegmentMakerv1::Exec(Option_t * option)
399 if( strcmp(GetName(), "")== 0 )
402 if(strstr(option,"tim"))
403 gBenchmark->Start("PHOSTrackSegmentMakerv1");
405 if(strstr(option,"print")) {
410 //check, if the branch with name of this" already exits?
411 TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
413 TBranch * branch = 0 ;
414 Bool_t phostsfound = kFALSE, tracksegmentmakerfound = kFALSE ;
416 TString taskName(GetName()) ;
417 taskName.ReplaceAll(Version(), "") ;
419 while ( (branch = (TBranch*)next()) && (!phostsfound || !tracksegmentmakerfound) ) {
420 if ( (strcmp(branch->GetName(), "PHOSTS")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
421 phostsfound = kTRUE ;
423 else if ( (strcmp(branch->GetName(), "AliPHOSTrackSegmentMaker")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
424 tracksegmentmakerfound = kTRUE ;
427 if ( phostsfound || tracksegmentmakerfound ) {
428 cerr << "WARNING: AliPHOSTrackSegmentMakerv1::Exec -> TrackSegments and/or TrackSegmentMaker branch with name "
429 << taskName.Data() << " already exits" << endl ;
433 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
434 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
437 for(ievent = 0; ievent < nevents; ievent++){
438 gAlice->GetEvent(ievent) ;
439 gAlice->SetEvent(ievent) ;
441 if(!ReadRecPoints(ievent)) //reads RecPoints for event ievent
444 for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ){
454 WriteTrackSegments(ievent) ;
456 if(strstr(option,"deb"))
457 PrintTrackSegments(option) ;
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 ;
471 //____________________________________________________________________________
472 void AliPHOSTrackSegmentMakerv1::Print(Option_t * option)const
474 // Print TrackSegmentMaker parameters
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 ;
487 cout << "AliPHOSTrackSegmentMakerv1 not initialized " << endl ;
489 //____________________________________________________________________________
490 Bool_t AliPHOSTrackSegmentMakerv1::ReadRecPoints(Int_t event)
492 // Reads Emc and CPV recPoints
493 // made previously with Clusterizer.
496 //Make some initializations
498 fNTrackSegments = 0 ;
506 // Get TreeR header from file
507 if(gAlice->TreeR()==0){
508 cerr << "ERROR: AliPHOSTrackSegmentMakerv1::ReadRecPoints -> There is no Reconstruction Tree" << endl;
514 TBranch * emcbranch = 0;
515 TBranch * cpvbranch = 0;
516 TBranch * clusterizerbranch = 0;
517 TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
519 TBranch * branch = 0 ;
520 Bool_t phosemcfound = kFALSE, phoscpvfound = kFALSE, clusterizerfound = kFALSE ;
522 TString taskName(GetName()) ;
523 taskName.ReplaceAll(Version(), "") ;
525 while ( (branch = (TBranch*)next()) && (!phosemcfound || !phoscpvfound || !clusterizerfound) ) {
526 if ( (strcmp(branch->GetName(), "PHOSEmcRP")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
527 phosemcfound = kTRUE ;
531 else if ( (strcmp(branch->GetName(), "PHOSCpvRP")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
532 phoscpvfound = kTRUE ;
535 } else if ( (strcmp(branch->GetName(), "AliPHOSClusterizer")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
536 clusterizerfound = kTRUE ;
537 clusterizerbranch = branch ;
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 ;
546 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
548 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
549 emcRecPoints->Clear() ;
550 emcbranch->SetAddress(&emcRecPoints) ;
552 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
553 cpvRecPoints->Clear() ;
554 cpvbranch->SetAddress(&cpvRecPoints) ;
556 TClonesArray * trackSegments = gime->TrackSegments() ;
557 trackSegments->Clear() ;
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) ;
566 emcbranch->GetEntry(0) ;
567 cpvbranch->GetEntry(0) ;
569 clusterizerbranch->GetEntry(0) ;
571 printf("ReadRecPoint: EMC=%d, CPV=%d\n",emcRecPoints->GetEntries(),cpvRecPoints->GetEntries());
577 //____________________________________________________________________________
578 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments(Int_t event)
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.
588 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
589 TClonesArray * trackSegments = gime->TrackSegments() ;
591 //Make branch in TreeR for TrackSegments
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()) ;
598 TDirectory *cwd = gDirectory;
601 Int_t bufferSize = 32000 ;
602 TBranch * tsBranch = gAlice->TreeR()->Branch("PHOSTS",&trackSegments,bufferSize);
603 tsBranch->SetTitle(fTrackSegmentsBranchTitle.Data());
605 tsBranch->SetFile(filename);
606 TIter next( tsBranch->GetListOfBranches());
608 while ((sb=(TBranch*)next())) {
609 sb->SetFile(filename);
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());
621 tsMakerBranch->SetFile(filename);
622 TIter next( tsMakerBranch->GetListOfBranches());
624 while ((sb=(TBranch*)next())) {
625 sb->SetFile(filename);
631 tsMakerBranch->Fill() ;
633 gAlice->TreeR()->Write(0,kOverwrite) ;
638 //____________________________________________________________________________
639 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
641 // option deb - prints # of found TrackSegments
642 // option deb all - prints as well indexed of found RecParticles assigned to the TS
644 TClonesArray * trackSegments = AliPHOSGetter::GetInstance()->TrackSegments() ;
646 cout << "AliPHOSTrackSegmentMakerv1: event "<<gAlice->GetEvNumber() << endl ;
647 cout << " Found " << trackSegments->GetEntriesFast() << " trackSegments " << endl ;
649 if(strstr(option,"all")) { // printing found TS
650 cout << "TrackSegment # " << " EMC RP# " << " CPV RP# " << " PPSD RP#" << endl ;
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 ;
661 cout << "-------------------------------------------------------"<< endl ;