]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSTrackSegmentMakerv1.cxx
e245695155eac01c997fa9ad534a13ca2084adb1
[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 /* History of cvs commits:
18  *
19  * $Log$
20  * Revision 1.74  2005/07/08 14:01:36  hristov
21  * Tracking in non-uniform nmagnetic field (Yu.Belikov)
22  *
23  * Revision 1.73  2005/05/28 14:19:05  schutz
24  * Compilation warnings fixed by T.P.
25  *
26  */
27
28 //_________________________________________________________________________
29 // Implementation version 1 of algorithm class to construct PHOS track segments
30 // Track segment for PHOS is list of 
31 //        EMC RecPoint + (possibly) CPV RecPoint
32 // To find TrackSegments we do the following: 
33 //  for each EMC RecPoints we look at
34 //   CPV RecPoints in the radious fRcpv. 
35 //  If there is such a CPV RecPoint, 
36 //   we make "Link" it is just indexes of EMC and CPV RecPoint and distance
37 //   between them in the PHOS plane. 
38 //  Then we sort "Links" and starting from the 
39 //   least "Link" pointing to the unassined EMC and CPV RecPoints assing them to 
40 //   new TrackSegment. 
41 // If there is no CPV RecPoint we make TrackSegment 
42 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
43 //// In principle this class should be called from AliPHOSReconstructor, but 
44 // one can use it as well in standalone mode.
45 // Use  case:
46 //  root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
47 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
48 //               // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
49 //               // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")                       
50 //  root [1] t->ExecuteTask()
51 //  root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
52 //  root [4] t->ExecuteTask("deb all time") 
53 //                 
54 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH) 
55 //
56
57 // --- ROOT system ---
58 #include "TTree.h"
59 #include "TBenchmark.h"
60
61 // --- Standard library ---
62 #include "Riostream.h"
63 // --- AliRoot header files ---
64 #include "AliPHOSGeometry.h"
65 #include "AliPHOSTrackSegmentMakerv1.h"
66 #include "AliPHOSTrackSegment.h"
67 #include "AliPHOSLink.h"
68 #include "AliPHOSGetter.h"
69 #include "AliESD.h"
70 #include "AliESDtrack.h"
71
72 ClassImp( AliPHOSTrackSegmentMakerv1) 
73
74
75 //____________________________________________________________________________
76   AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() : AliPHOSTrackSegmentMaker()
77 {
78   // default ctor (to be used mainly by Streamer)
79
80   InitParameters() ; 
81   fDefaultInit = kTRUE ; 
82 }
83
84 //____________________________________________________________________________
85  AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const TString alirunFileName, const TString eventFolderName)
86    :AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName)
87 {
88   // ctor
89
90   InitParameters() ; 
91   Init() ;
92   fDefaultInit = kFALSE ; 
93   fESD = 0;
94 }
95
96 //____________________________________________________________________________
97  AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
98
99   // dtor
100   // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
101   if (!fDefaultInit)  
102     delete fLinkUpArray ;
103 }
104
105
106 //____________________________________________________________________________
107 const TString AliPHOSTrackSegmentMakerv1::BranchName() const 
108 {  
109  
110   return GetName() ;
111 }
112
113 //____________________________________________________________________________
114 void  AliPHOSTrackSegmentMakerv1::FillOneModule()
115 {
116   // Finds first and last indexes between which 
117   // clusters from one PHOS module are
118
119   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
120   
121   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
122   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
123  
124   //First EMC clusters
125   Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
126   for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&  
127         ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule ); 
128       fEmcLast ++)  ;
129   
130   //Now CPV clusters
131   Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
132
133     for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) && 
134          ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule ); 
135        fCpvLast ++) ;
136       
137 }
138
139 //____________________________________________________________________________
140 Float_t  AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSCpvRecPoint * cpvClu, Int_t &trackindex) const
141 {
142   // Calculates the distance between the EMC RecPoint and the CPV RecPoint
143   // Clusters are sorted in "rows" and "columns" of width 1 cm
144
145   //Float_t delta = 1 ;  // Width of the rows in sorting of RecPoints (in cm)
146                        // if you change this value, change it as well in xxxRecPoint::Compare()
147   Float_t distance2Cpv   = fRcpv ;
148   Float_t distance2Track = fRtpc ; 
149
150   trackindex = -1 ; // closest track within fRCpv 
151
152   TVector3 vecEmc ;   // Local position of EMC recpoint
153   TVector3 vecCpv ;   // Local position of CPV recpoint propagated to EMC
154   TVector3 vecDist ;  // Distance between local positions of two points
155   
156   emcClu->GetLocalPosition(vecEmc) ;
157   cpvClu->GetLocalPosition(vecCpv) ;
158
159   //toofar = kTRUE ;
160   if(emcClu->GetPHOSMod() == cpvClu->GetPHOSMod()){ 
161
162     // Find EMC-CPV distance
163     distance2Cpv = (vecCpv - vecEmc).Mag() ;
164     
165     if (fESD != 0x0) {
166       // Extrapolate the global track direction if any to CPV and find the closest track
167       Int_t nTracks = fESD->GetNumberOfTracks();
168       Int_t iClosestTrack = -1;
169       Double_t minDistance = 1e6;
170       Double_t pxyz[3], xyz[3];
171       Double_t rPHOS = 
172         AliPHOSGetter::Instance()->PHOSGeometry()->GetIPtoCrystalSurface();
173       AliESDtrack *track;
174       for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
175         track = fESD->GetTrack(iTrack);
176         if (track->IsPHOS()) 
177           continue ; 
178         if (!track->GetXYZAt(rPHOS, fESD->GetMagneticField(), xyz))
179            continue; //track coord on the cylinder of PHOS radius
180         if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
181            continue;
182         if (!track->GetPxPyPzAt(rPHOS, fESD->GetMagneticField(), pxyz))
183            continue; // track momentum ibid.
184         PropagateToPlane(vecDist,xyz,pxyz,"CPV",cpvClu->GetPHOSMod());
185         //      Info("GetDistanceInPHOSPlane","Track %d propagation to CPV = (%f,%f,%f)",
186         //     iTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
187         vecDist -= vecCpv;
188         distance2Track = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
189         // Find the closest track to the EMC recpoint
190         if (distance2Track < minDistance) {
191           minDistance = distance2Track;
192           iClosestTrack = iTrack;
193         }
194       }
195
196       if (iClosestTrack != -1) {
197         track = fESD->GetTrack(iClosestTrack);
198         if (track->GetPxPyPzAt(rPHOS, fESD->GetMagneticField(), pxyz)) { // track momentum ibid.
199         TVector3 vecCpvGlobal; // Global position of the CPV recpoint
200         AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
201         const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
202         geom->GetGlobal((AliRecPoint*)cpvClu,vecCpvGlobal);
203         for (Int_t ixyz=0; ixyz<3; ixyz++)
204           xyz[ixyz] = vecCpvGlobal[ixyz];
205         PropagateToPlane(vecDist,xyz,pxyz,"EMC",cpvClu->GetPHOSMod());
206 //      Info("GetDistanceInPHOSPlane","Track %d propagation to EMC = (%f,%f,%f)",
207 //           iClosestTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
208         vecDist -= vecEmc;
209         distance2Track = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
210         }
211       }
212 //     } else {
213 //       // If no ESD exists, than simply find EMC-CPV distance
214 //       distance = (vecCpv - vecEmc).Mag() ;
215     
216       //if(distance2Track < fRcpv + 2*delta )
217       if(distance2Track < fRtpc )
218         trackindex = iClosestTrack ; 
219       //      toofar = kFALSE ;
220     }
221     //     Info("GetDistanceInPHOSPlane","cpv-emc distance is %f cm",
222     //   distance);
223   }
224   
225   return distance2Cpv ;
226 }
227
228 //____________________________________________________________________________
229 void AliPHOSTrackSegmentMakerv1::PropagateToPlane(TVector3& globalIntersection,
230                                                   Double_t *x,
231                                                   Double_t *p,
232                                                   const char *det,
233                                                   Int_t moduleNumber) const
234 {
235   // Propagate a straight-line track from the origin point x
236   // along the direction p to the CPV or EMC module moduleNumber
237   // Returns a local position of such a propagation
238
239   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
240   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
241   TVector3 moduleCenter;
242   geom->GetModuleCenter(moduleCenter,det,moduleNumber);
243   TVector3 vertex; vertex.SetXYZ(x[0],x[1],x[2]);
244   TVector3 direction; direction.SetXYZ(p[0],p[1],p[2]);
245
246 //   Info("PropagateToCPV","Center of the %s module %d is (%f,%f,%f)",
247 //        det,moduleNumber,moduleCenter[0],moduleCenter[1],moduleCenter[2]);
248
249   Double_t time = (moduleCenter.Mag2() - vertex.Dot(moduleCenter)) /
250     (direction.Dot(moduleCenter));
251   vertex += direction*time;
252   geom->Global2Local(globalIntersection,vertex,moduleNumber);
253 }
254
255 //____________________________________________________________________________
256 void  AliPHOSTrackSegmentMakerv1::Init()
257 {
258   // Make all memory allocations that are not possible in default constructor
259   
260   AliPHOSGetter* gime = AliPHOSGetter::Instance();
261   if(!gime)
262     gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
263   
264   fLinkUpArray  = new TClonesArray("AliPHOSLink", 1000); 
265   if ( !gime->TrackSegmentMaker() ) {
266     gime->PostTrackSegmentMaker(this);
267   }
268 }
269
270 //____________________________________________________________________________
271 void  AliPHOSTrackSegmentMakerv1::InitParameters()
272 {
273   //Initializes parameters
274   fRcpv      = 10. ;
275   fRtpc      = 4. ;
276   fEmcFirst  = 0 ;    
277   fEmcLast   = 0 ;   
278   fCpvFirst  = 0 ;   
279   fCpvLast   = 0 ;   
280   fLinkUpArray = 0 ;
281   fWrite                   = kTRUE ;
282   fTrackSegmentsInRun       = 0 ; 
283   SetEventRange(0,-1) ;
284 }
285
286
287 //____________________________________________________________________________
288 void  AliPHOSTrackSegmentMakerv1::MakeLinks()const
289
290   // Finds distances (links) between all EMC and CPV clusters, 
291   // which are not further apart from each other than fRcpv 
292   // and sort them in accordance with this distance
293   
294   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
295   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
296   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
297
298   fLinkUpArray->Clear() ;    
299
300   AliPHOSCpvRecPoint * cpv ;
301   AliPHOSEmcRecPoint * emcclu ;
302
303   Int_t iLinkUp  = 0 ;
304   
305   Int_t iEmcRP;
306   for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
307     emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
308
309     //Bool_t toofar ;        
310     Int_t iCpv = 0 ;    
311     for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) { 
312       
313       cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
314       Int_t track = -1 ; 
315       Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, track) ;     
316       //      if(toofar)
317       //        continue ;       
318       if(r < fRcpv) { 
319         new ((*fLinkUpArray)[iLinkUp++])  AliPHOSLink(r, iEmcRP, iCpv, track) ;
320       }      
321     }
322   } 
323   
324   fLinkUpArray->Sort() ;  //first links with smallest distances
325 }
326
327 //____________________________________________________________________________
328 void  AliPHOSTrackSegmentMakerv1::MakePairs()
329
330   // Using the previously made list of "links", we found the smallest link - i.e. 
331   // link with the least distance between EMC and CPV and pointing to still 
332   // unassigned RecParticles. We assign these RecPoints to TrackSegment and 
333   // remove them from the list of "unassigned". 
334
335   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
336
337   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
338   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
339   TClonesArray * trackSegments = gime->TrackSegments();
340     
341   //Make arrays to mark clusters already chosen
342   Int_t * emcExist = 0;
343   if(fEmcLast > fEmcFirst)
344     emcExist = new Int_t[fEmcLast-fEmcFirst] ;
345   
346   Int_t index;
347   for(index = 0; index <fEmcLast-fEmcFirst; index ++)
348     emcExist[index] = 1 ;
349   
350   Bool_t * cpvExist = 0;
351   if(fCpvLast > fCpvFirst)
352     cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
353   for(index = 0; index <fCpvLast-fCpvFirst; index ++)
354     cpvExist[index] = kTRUE ;
355   
356   
357   // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance 
358   TIter nextUp(fLinkUpArray) ;
359   
360   AliPHOSLink * linkUp ;
361   
362   AliPHOSCpvRecPoint * nullpointer = 0 ;
363   
364   while ( (linkUp =  static_cast<AliPHOSLink *>(nextUp()) ) ){  
365
366     if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
367
368       if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
369          new ((* trackSegments)[fNTrackSegments]) 
370            AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) , 
371                                dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetCpv())) , 
372                                linkUp->GetTrack()) ;
373          
374        (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
375        fNTrackSegments++ ;
376        emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc  that Cpv was found 
377        //mark CPV recpoint as already used 
378        cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
379       } //if CpvUp still exist
380     } 
381   }        
382
383   //look through emc recPoints left without CPV
384   if(emcExist){ //if there is emc rec point
385     Int_t iEmcRP ;
386     for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst  ; iEmcRP++ ){
387       if(emcExist[iEmcRP] > 0 ){
388        new ((*trackSegments)[fNTrackSegments])  
389          AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)), 
390                            nullpointer) ;
391        (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
392        fNTrackSegments++;    
393       } 
394     }
395   }
396   delete [] emcExist ; 
397   delete [] cpvExist ; 
398 }
399
400 //____________________________________________________________________________
401 void  AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
402 {
403   // Steering method to perform track segment construction for events
404   // in the range from fFirstEvent to fLastEvent.
405   // This range is optionally set by SetEventRange().
406   // if fLastEvent=-1 (by default), then process events until the end.
407   
408   if(strstr(option,"tim"))
409     gBenchmark->Start("PHOSTSMaker");
410  
411   if(strstr(option,"print")) {
412     Print() ; 
413     return ; 
414   }
415   
416   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;  
417  
418   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
419
420   if (fLastEvent == -1) 
421     fLastEvent = gime->MaxEvent() - 1 ;
422   else 
423     fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
424   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
425
426   Int_t ievent ; 
427   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
428     gime->Event(ievent,"R") ;
429    //Make some initializations 
430     fNTrackSegments = 0 ;
431     fEmcFirst = 0 ;    
432     fEmcLast  = 0 ;   
433     fCpvFirst = 0 ;   
434     fCpvLast  = 0 ;   
435     
436     gime->TrackSegments()->Clear();
437
438     //    if(!ReadRecPoints(ievent))   continue; //reads RecPoints for event ievent
439     
440     for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
441       FillOneModule() ; 
442       MakeLinks() ;
443       MakePairs() ;
444     }
445
446     WriteTrackSegments() ;
447
448     if(strstr(option,"deb"))
449       PrintTrackSegments(option);
450     
451     //increment the total number of track segments per run 
452     fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ; 
453   }
454   
455   if(strstr(option,"tim")){
456     gBenchmark->Stop("PHOSTSMaker");
457     Info("Exec", "took %f seconds for making TS %f seconds per event", 
458           gBenchmark->GetCpuTime("PHOSTSMaker"), 
459           gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
460    }
461   if(fWrite) //do not unload in "on flight" mode
462     Unload();
463 }
464
465 //____________________________________________________________________________
466 void AliPHOSTrackSegmentMakerv1::Unload() 
467 {
468   // Unloads the task from the folder
469   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;  
470   gime->PhosLoader()->UnloadRecPoints() ;
471   gime->PhosLoader()->UnloadTracks() ;
472 }
473
474 //____________________________________________________________________________
475 void AliPHOSTrackSegmentMakerv1::Print(const Option_t *)const
476 {
477   //  Print TrackSegmentMaker parameters
478
479   TString message("") ;
480   if( strcmp(GetName(), "") != 0 ) {
481     message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ; 
482     message += "Making Track segments\n" ;
483     message += "with parameters:\n" ; 
484     message += "     Maximal EMC - CPV distance (cm) %f\n" ;
485     message += "============================================\n" ;
486     Info("Print", message.Data(),fRcpv) ;
487   }
488   else
489     Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
490 }
491
492 //____________________________________________________________________________
493 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
494 {
495   // Writes found TrackSegments to TreeR. Creates branches 
496   // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
497   // In the former branch found TrackSegments are stored, while 
498   // in the latter all parameters, with which TS were made. 
499   // ROOT does not allow overwriting existing branches, therefore
500   // first we check, if branches with the same title already exist.
501   // If yes - exits without writing.
502
503   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
504
505   TClonesArray * trackSegments = gime->TrackSegments() ; 
506   trackSegments->Expand(trackSegments->GetEntriesFast()) ;
507
508   if(fWrite){ //We write TreeT
509     TTree * treeT = gime->TreeT();
510     
511     //First TS
512     Int_t bufferSize = 32000 ; 
513     TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
514     tsBranch->Fill() ;  
515     
516     gime->WriteTracks("OVERWRITE");
517     gime->WriteTrackSegmentMaker("OVERWRITE");
518   }
519 }
520
521
522 //____________________________________________________________________________
523 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
524 {
525   // option deb - prints # of found TrackSegments
526   // option deb all - prints as well indexed of found RecParticles assigned to the TS
527
528   TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ; 
529
530   Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ; 
531   printf("nevent: %d\n", gAlice->GetEvNumber()) ; 
532   printf("        Found %d TrackSegments\n", trackSegments->GetEntriesFast() ); 
533   
534   if(strstr(option,"all")) {  // printing found TS
535     printf("TrackSegment #  EMC RP#  CPV RP#\n") ; 
536     Int_t index;
537     for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
538       AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ; 
539       printf("   %d           %d        %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ; 
540     }   
541   }
542 }