]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSTrackSegmentMakerv1.cxx
restyled Print
[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 fRcpv. 
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 //// In principle this class should be called from AliPHOSReconstructioner, but 
32 // one can use it as well in standalone mode.
33 // Use  case:
34 //  root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
35 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
36 //               // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
37 //               // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")                       
38 //  root [1] t->ExecuteTask()
39 //  root [2] t->SetMaxEmcPpsdDistance(5)
40 //  root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
41 //  root [4] t->ExecuteTask("deb all time") 
42 //                 
43 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH) 
44 //
45
46 // --- ROOT system ---
47 #include "TROOT.h"
48 #include "TFile.h"
49 #include "TFolder.h"
50 #include "TTree.h"
51 #include "TSystem.h"
52 #include "TBenchmark.h"
53
54 // --- Standard library ---
55
56 // --- AliRoot header files ---
57
58 #include "AliPHOSTrackSegmentMakerv1.h"
59 #include "AliPHOSClusterizerv1.h"
60 #include "AliPHOSTrackSegment.h"
61 #include "AliPHOSCpvRecPoint.h"
62 #include "AliPHOSLink.h"
63 #include "AliPHOSGetter.h"
64 #include "AliPHOS.h"
65 #include "AliRun.h"
66
67 ClassImp( AliPHOSTrackSegmentMakerv1) 
68
69
70 //____________________________________________________________________________
71   AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() : AliPHOSTrackSegmentMaker()
72 {
73   // default ctor (to be used mainly by Streamer)
74
75   InitParameters() ; 
76 //   fHeaderFileName           = "" ;
77 //   fRecPointsBranchTitle     = "" ;
78 //   fTrackSegmentsBranchTitle = "" ; 
79 //   fFrom                     = "" ; 
80
81   fTrackSegmentsInRun       = 0 ; 
82
83   fDefaultInit = kTRUE ; 
84 }
85
86 //____________________________________________________________________________
87  AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const char * headerFile, const char * name, const Bool_t toSplit) : AliPHOSTrackSegmentMaker(headerFile, name, toSplit)
88 {
89   // ctor
90
91   InitParameters() ; 
92 //   fHeaderFileName           = GetTitle() ;
93 //   fRecPointsBranchTitle     = GetName() ;
94 //   fTrackSegmentsBranchTitle = GetName() ; 
95   fTrackSegmentsInRun       = 0 ; 
96
97 //   if ( from == 0 ) 
98 //     fFrom = name ; 
99 //   else
100 //     fFrom = from ; 
101   Init() ;
102
103   fDefaultInit = kFALSE ; 
104
105 }
106
107 //____________________________________________________________________________
108  AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
109
110   // dtor
111   // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
112   
113   if (!fDefaultInit) {
114     delete fLinkUpArray  ;
115     
116 //     AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
117     
118 //     // remove the task from the folder list
119 //     gime->RemoveTask("T",GetName()) ;
120 //     TString name(GetName()) ; 
121 //     name.ReplaceAll("tsm", "clu") ; 
122 //     gime->RemoveTask("C",name) ;
123     
124 //     // remove the data from the folder list
125 //     name = GetName() ; 
126 //     name.Remove(name.Index(":")) ; 
127 //     gime->RemoveObjects("RE", name) ; // EMCARecPoints
128 //     gime->RemoveObjects("RC", name) ; // CPVRecPoints
129 //     gime->RemoveObjects("T", name) ;  // TrackSegments
130     
131 //     // Delete gAlice
132 //     gime->CloseFile() ; 
133     
134     fSplitFile = 0 ; 
135   }
136 }
137
138
139 //____________________________________________________________________________
140 const TString AliPHOSTrackSegmentMakerv1::BranchName() const 
141 {  
142   TString branchName(GetName() ) ;
143   branchName.Remove(branchName.Index(Version())-1) ;
144   return branchName ;
145 }
146
147 //____________________________________________________________________________
148 void  AliPHOSTrackSegmentMakerv1::FillOneModule()
149 {
150   // Finds first and last indexes between which 
151   // clusters from one PHOS module are
152   
153   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
154   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
155   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
156  
157   //First EMC clusters
158   Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
159   for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&  
160         ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule ); 
161       fEmcLast ++)  ;
162   
163   //Now CPV clusters
164   Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
165
166     for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) && 
167           ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule ); 
168         fCpvLast ++) ;
169       
170 }
171
172 //____________________________________________________________________________
173 Float_t  AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSRecPoint * cpvClu, Bool_t &toofar)const
174 {
175   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
176   // Clusters are sorted in "rows" and "columns" of width 1 cm
177
178   Float_t delta = 1 ;  // Width of the rows in sorting of RecPoints (in cm)
179                        // if you change this value, change it as well in xxxRecPoint::Compare()
180   Float_t r = fRcpv ;
181  
182   TVector3 vecEmc ;
183   TVector3 vecCpv ;
184   
185   emcClu->GetLocalPosition(vecEmc) ;
186   cpvClu->GetLocalPosition(vecCpv)  ; 
187
188   if(emcClu->GetPHOSMod() == cpvClu->GetPHOSMod()){ 
189     if(vecCpv.X() <= vecEmc.X() + fRcpv + 2*delta ){ 
190
191       vecCpv = vecCpv  - vecEmc ; 
192       r = vecCpv.Mag() ;
193       toofar = kFALSE ;
194
195     } // if  xPpsd >= xEmc + ...
196     else 
197       toofar = kTRUE ;
198   } 
199   else 
200     toofar = kTRUE ;
201
202   //toofar = kFALSE ;
203  
204   
205   return r ;
206 }
207
208 //____________________________________________________________________________
209 void  AliPHOSTrackSegmentMakerv1::Init()
210 {
211   // Make all memory allocations that are not possible in default constructor
212   
213   if ( strcmp(GetTitle(), "") == 0 )
214     SetTitle("galice.root") ;
215     
216   TString branchname = GetName() ;
217   branchname.Remove(branchname.Index(Version())-1) ;
218   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(),branchname.Data(), fToSplit ) ; 
219   if ( gime == 0 ) {
220     Error("Init", "Could not obtain the Getter object !") ; 
221     return ;
222   } 
223   
224   fSplitFile = 0 ;
225   if(fToSplit){
226     //First - extract full path if necessary
227     TString fileName(GetTitle()) ;
228     Ssiz_t islash = fileName.Last('/') ;
229     if(islash<fileName.Length())
230       fileName.Remove(islash+1,fileName.Length()) ;
231     else
232       fileName="" ;
233     fileName+="PHOS.RecData." ;
234     if((strcmp(branchname.Data(),"Default")!=0)&&(strcmp(branchname.Data(),"")!=0)){
235       fileName+=branchname ;
236       fileName+="." ;
237     }
238     fileName+="root" ;
239     fSplitFile = static_cast<TFile*>(gROOT->GetFile(fileName.Data()));   
240     if(!fSplitFile)
241       fSplitFile =  TFile::Open(fileName.Data(),"update") ;
242   }
243   
244   fLinkUpArray  = new TClonesArray("AliPHOSLink", 1000); 
245   
246
247   gime->PostTrackSegmentMaker(this) ;
248   gime->PostTrackSegments(BranchName()) ; 
249
250 }
251
252 //____________________________________________________________________________
253 void  AliPHOSTrackSegmentMakerv1::InitParameters()
254 {
255   fRcpv      = 10. ;   
256   fEmcFirst  = 0 ;    
257   fEmcLast   = 0 ;   
258   fCpvFirst  = 0 ;   
259   fCpvLast   = 0 ;   
260   fLinkUpArray = 0 ;
261   TString tsmName( GetName()) ; 
262   if (tsmName.IsNull() ) 
263     tsmName = "Default" ; 
264   tsmName.Append(":") ; 
265   tsmName.Append(Version()) ; 
266   SetName(tsmName) ;
267 }
268
269
270 //____________________________________________________________________________
271 void  AliPHOSTrackSegmentMakerv1::MakeLinks()const
272
273   // Finds distances (links) between all EMC and PPSD clusters, 
274   // which are not further apart from each other than fRcpv 
275   // and sort them in accordance with this distance
276   
277   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
278   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
279   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
280
281   fLinkUpArray->Clear() ;    
282
283   AliPHOSRecPoint * cpv ;
284   AliPHOSEmcRecPoint * emcclu ;
285
286   Int_t iLinkUp  = 0 ;
287   
288   Int_t iEmcRP;
289   for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
290     emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
291
292     Bool_t toofar ;        
293     Int_t iCpv = 0 ;    
294     for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) { 
295       
296       cpv = dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(iCpv)) ;
297       Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, toofar) ;
298       
299       if(toofar)
300         break ;  
301       if(r < fRcpv) { 
302         new ((*fLinkUpArray)[iLinkUp++])  AliPHOSLink(r, iEmcRP, iCpv) ;
303       }      
304     }
305   } 
306   
307   fLinkUpArray->Sort() ;  //first links with smallest distances
308 }
309
310 //____________________________________________________________________________
311 void  AliPHOSTrackSegmentMakerv1::MakePairs()
312
313   // Using the previously made list of "links", we found the smallest link - i.e. 
314   // link with the least distance between EMC and CPV and pointing to still 
315   // unassigned RecParticles. We assign these RecPoints to TrackSegment and 
316   // remove them from the list of "unassigned". 
317   
318   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
319   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
320   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
321   TClonesArray * trackSegments = gime->TrackSegments(BranchName()) ;   
322     
323   //Make arrays to mark clusters already chosen
324   Int_t * emcExist = 0;
325   if(fEmcLast > fEmcFirst)
326     emcExist = new Int_t[fEmcLast-fEmcFirst] ;
327   
328   Int_t index;
329   for(index = 0; index <fEmcLast-fEmcFirst; index ++)
330     emcExist[index] = 1 ;
331   
332   Bool_t * cpvExist = 0;
333   if(fCpvLast > fCpvFirst)
334     cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
335   for(index = 0; index <fCpvLast-fCpvFirst; index ++)
336     cpvExist[index] = kTRUE ;
337   
338   
339   // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance 
340   TIter nextUp(fLinkUpArray) ;
341   
342   AliPHOSLink * linkUp ;
343   
344   AliPHOSRecPoint * nullpointer = 0 ;
345   
346   while ( (linkUp =  static_cast<AliPHOSLink *>(nextUp()) ) ){  
347
348     if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){ //without ppsd Up yet 
349
350       if(cpvExist[linkUp->GetPpsd()-fCpvFirst]){ //CPV still exist
351         
352         new ((* trackSegments)[fNTrackSegments]) 
353           AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) , 
354                               dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(linkUp->GetPpsd()))) ;
355         (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
356         fNTrackSegments++ ;
357         
358         emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc  that Cpv was found 
359         //mark CPV recpoint as already used 
360         cpvExist[linkUp->GetPpsd()-fCpvFirst] = kFALSE ;
361       } //if ppsdUp still exist
362     } 
363   }      
364
365   //look through emc recPoints left without CPV/PPSD
366   if(emcExist){ //if there is emc rec point
367     Int_t iEmcRP ;
368     for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst  ; iEmcRP++ ){
369       if(emcExist[iEmcRP] > 0 ){
370         new ((*trackSegments)[fNTrackSegments])  
371           AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)), 
372                               nullpointer) ;
373         (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
374         fNTrackSegments++;    
375       } 
376     }
377   }
378   delete [] emcExist ; 
379   delete [] cpvExist ; 
380 }
381
382 //____________________________________________________________________________
383 void  AliPHOSTrackSegmentMakerv1::Exec(Option_t * option)
384 {
385   // STEERing method
386
387   if( strcmp(GetName(), "")== 0 ) 
388     Init() ;
389
390   if(strstr(option,"tim"))
391     gBenchmark->Start("PHOSTSMaker");
392  
393   if(strstr(option,"print")) {
394     Print("") ; 
395     return ; 
396   }
397
398 //   gAlice->GetEvent(0) ;
399 //   //check, if the branch with name of this" already exits?
400 //   if (gAlice->TreeR()) { 
401 //     TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeR()->GetListOfBranches()) ;
402 //     TIter next(lob) ; 
403 //     TBranch * branch = 0 ;  
404 //     Bool_t phostsfound = kFALSE, tracksegmentmakerfound = kFALSE ; 
405     
406 //     TString branchname = GetName() ;
407 //     branchname.Remove(branchname.Index(Version())-1) ;
408     
409 //     while ( (branch = static_cast<TBranch*>(next())) && (!phostsfound || !tracksegmentmakerfound) ) {
410 //       if ( (strcmp(branch->GetName(), "PHOSTS")==0) && (strcmp(branch->GetTitle(), branchname.Data())==0) ) 
411 //      phostsfound = kTRUE ;
412       
413 //       else if ( (strcmp(branch->GetName(), "AliPHOSTrackSegmentMaker")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
414 //      tracksegmentmakerfound = kTRUE ; 
415 //     }
416     
417 //     if ( phostsfound || tracksegmentmakerfound ) {
418 //       Warning("Exec", "TrackSegments and/or TrackSegmentMaker branch with name %s already exists", branchname.Data() ) ;
419 //       return ; 
420 //     }       
421 //   }
422
423   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
424   if(gime->BranchExists("TrackSegments") )
425     return ;
426   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
427   Int_t nevents = gime->MaxEvent() ;       //(Int_t) gAlice->TreeE()->GetEntries() ;
428   Int_t ievent ;
429
430   for(ievent = 0; ievent < nevents; ievent++){
431
432     gime->Event(ievent,"R") ;
433     //Make some initializations 
434     fNTrackSegments = 0 ;
435     fEmcFirst = 0 ;    
436     fEmcLast  = 0 ;   
437     fCpvFirst = 0 ;   
438     fCpvLast  = 0 ;   
439     gime->TrackSegments(BranchName())->Clear() ; 
440
441     //    if(!ReadRecPoints(ievent))   continue; //reads RecPoints for event ievent
442     
443     for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ){
444       
445       FillOneModule() ; 
446       
447       MakeLinks() ;
448       
449       MakePairs() ;
450       
451     }
452
453     WriteTrackSegments(ievent) ;
454
455     if(strstr(option,"deb"))
456       PrintTrackSegments(option) ;
457     
458     //increment the total number of track segments per run 
459     fTrackSegmentsInRun += gime->TrackSegments(BranchName())->GetEntriesFast() ; 
460
461   }
462
463   if(strstr(option,"tim")){
464     gBenchmark->Stop("PHOSTSMaker");
465     Info("Exec", "took %f seconds for making TS %f seconds per event", 
466          gBenchmark->GetCpuTime("PHOSTSMaker"), gBenchmark->GetCpuTime("PHOSTSMaker")/nevents) ;
467   }
468     
469 }
470
471 //____________________________________________________________________________
472 void AliPHOSTrackSegmentMakerv1::Print(Option_t * option)const
473 {
474   //  Print TrackSegmentMaker parameters
475
476   TString message("") ;
477   if( strcmp(GetName(), "") != 0 ) {
478     message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ; 
479     message += "Making Track segments\n" ;
480     message += "with parameters:\n" ; 
481     message += "     Maximal EMC - CPV (PPSD) distance (cm) %f\n" ;
482     message += "============================================\n" ;
483     Info("Print", message.Data(),fRcpv) ;
484   }
485   else
486     Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
487 }
488
489 //____________________________________________________________________________
490 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments(Int_t event)
491 {
492   // Writes found TrackSegments to TreeR. Creates branches 
493   // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
494   // In the former branch found TrackSegments are stored, while 
495   // in the latter all parameters, with which TS were made. 
496   // ROOT does not allow overwriting existing branches, therefore
497   // first we check, if branches with the same title already exist.
498   // If yes - exits without writing.
499   
500   AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ; 
501
502   TClonesArray * trackSegments = gime->TrackSegments() ; 
503   trackSegments->Expand(trackSegments->GetEntriesFast()) ;
504   TTree * treeR ;
505   
506   if(fToSplit){
507     if(!fSplitFile)
508       return ;
509     fSplitFile->cd() ;
510     char name[10] ;
511     sprintf(name,"%s%d", "TreeR",event) ;
512     treeR = dynamic_cast<TTree*>(fSplitFile->Get(name)); 
513   }
514   else{
515     treeR = gAlice->TreeR();
516   }
517   
518   if(!treeR){
519     gAlice->MakeTree("R", fSplitFile);
520     treeR = gAlice->TreeR() ;
521   }
522   
523   //First TS
524   Int_t bufferSize = 32000 ;    
525   TBranch * tsBranch = treeR->Branch("PHOSTS",&trackSegments,bufferSize);
526   tsBranch->SetTitle(BranchName());
527
528   //Second -TSMaker
529   Int_t splitlevel = 0 ;
530   AliPHOSTrackSegmentMakerv1 * ts = this ;
531   TBranch * tsMakerBranch = treeR->Branch("AliPHOSTrackSegmentMaker","AliPHOSTrackSegmentMakerv1",
532                                           &ts,bufferSize,splitlevel);
533   tsMakerBranch->SetTitle(BranchName());
534
535   tsBranch->Fill() ;  
536   tsMakerBranch->Fill() ;
537
538   treeR->AutoSave() ; //Write(0,kOverwrite) ;  
539   if(gAlice->TreeR()!=treeR)
540     treeR->Delete();
541 }
542
543
544 //____________________________________________________________________________
545 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
546 {
547   // option deb - prints # of found TrackSegments
548   // option deb all - prints as well indexed of found RecParticles assigned to the TS
549   TString taskName(GetName()) ; 
550   taskName.Remove(taskName.Index(Version())-1) ;
551   
552   TClonesArray * trackSegments = AliPHOSGetter::GetInstance()->TrackSegments(taskName) ; 
553
554   TString message ; 
555   message  = "\nevent " ;
556   message += gAlice->GetEvNumber() ; 
557   message += "\n      Found " ;
558   message += trackSegments->GetEntriesFast() ; 
559   message += " TrackSegments\n" ; 
560
561   if(strstr(option,"all")) {  // printing found TS
562     message += "TrackSegment #  EMC RP#  CPV RP#\n" ; 
563     Int_t index;
564     for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
565       AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ; 
566       message += "\n" ; 
567       message += ts->GetIndexInList() ; 
568       message += " " ; 
569       message += ts->GetEmcIndex() ; 
570       message += " " ; 
571       message += ts->GetCpvIndex() ; 
572     }   
573   }
574   Info("Print", message.Data() ) ; 
575 }