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 **************************************************************************/
18 //_________________________________________________________________________
19 // RecPoint implementation for PHOS-EMC
20 // An EmcRecPoint is a cluster of digits
22 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
25 // --- ROOT system ---
31 // --- Standard library ---
33 // --- AliRoot header files ---
35 #include "AliPHOSLoader.h"
36 #include "AliGenerator.h"
37 #include "AliPHOSGeometry.h"
38 #include "AliPHOSDigit.h"
39 #include "AliPHOSEmcRecPoint.h"
41 ClassImp(AliPHOSEmcRecPoint)
43 //____________________________________________________________________________
44 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
52 fNExMax = 0 ; //Not unfolded yet
54 fLocPos.SetX(1000000.) ; //Local position should be evaluated
59 //____________________________________________________________________________
60 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const char * opt) : AliPHOSRecPoint(opt)
66 fNExMax = 0 ; //Not unfolded yet
70 fLocPos.SetX(1000000.) ; //Local position should be evaluated
75 //____________________________________________________________________________
76 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const AliPHOSEmcRecPoint & rp) : AliPHOSRecPoint(rp)
80 fMulDigit = rp.fMulDigit ;
82 fCoreEnergy = rp.fCoreEnergy ;
83 fEnergyList = new Float_t[rp.fMulDigit] ;
85 for(index = 0 ; index < fMulDigit ; index++)
86 fEnergyList[index] = rp.fEnergyList[index] ;
87 fNExMax = rp.fNExMax ;
91 //____________________________________________________________________________
92 AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
97 delete[] fEnergyList ;
100 //____________________________________________________________________________
101 void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
103 // Adds a digit to the RecPoint
104 // and accumulates the total amplitude and the multiplicity
107 fEnergyList = new Float_t[fMaxDigit];
109 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
111 Int_t * tempo = new Int_t[fMaxDigit];
112 Float_t * tempoE = new Float_t[fMaxDigit];
115 for ( index = 0 ; index < fMulDigit ; index++ ){
116 tempo[index] = fDigitsList[index] ;
117 tempoE[index] = fEnergyList[index] ;
120 delete [] fDigitsList ;
121 fDigitsList = new Int_t[fMaxDigit];
123 delete [] fEnergyList ;
124 fEnergyList = new Float_t[fMaxDigit];
126 for ( index = 0 ; index < fMulDigit ; index++ ){
127 fDigitsList[index] = tempo[index] ;
128 fEnergyList[index] = tempoE[index] ;
135 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
136 fEnergyList[fMulDigit] = Energy ;
140 EvalPHOSMod(&digit) ;
143 //____________________________________________________________________________
144 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
146 // Tells if (true) or not (false) two digits are neighbors
148 Bool_t aren = kFALSE ;
150 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
153 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
156 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
158 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
159 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
161 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
167 //____________________________________________________________________________
168 Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
170 // Compares two RecPoints according to their position in the PHOS modules
172 const Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
173 //value (what is senseless) change as vell delta in
174 //AliPHOSTrackSegmentMakerv* and other RecPoints...
177 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
180 Int_t phosmod1 = GetPHOSMod() ;
181 Int_t phosmod2 = clu->GetPHOSMod() ;
184 GetLocalPosition(locpos1) ;
186 clu->GetLocalPosition(locpos2) ;
188 if(phosmod1 == phosmod2 ) {
189 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
194 else if(locpos1.Z()>locpos2.Z())
201 if(phosmod1 < phosmod2 )
209 //______________________________________________________________________________
210 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t, Int_t) const
213 // Execute action corresponding to one event
214 // This member function is called when a AliPHOSRecPoint is clicked with the locator
216 // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
217 // and switched off when the mouse button is released.
220 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
222 static TGraph * digitgraph = 0 ;
224 if (!gPad->IsEditable()) return;
227 TCanvas * histocanvas ;
230 //try to get run loader from default event folder
231 AliRunLoader* rn = AliRunLoader::GetRunLoader(AliConfig::GetDefaultEventFolderName());
234 AliError(Form("Cannot find Run Loader in Default Event Folder"));
237 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(rn->GetLoader("PHOSLoader"));
240 AliError(Form("Cannot find PHOS Loader from Run Loader"));
245 const TClonesArray * digits = gime->Digits() ;
250 AliPHOSDigit * digit ;
254 const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
255 Float_t * xi = new Float_t[kMulDigit] ;
256 Float_t * zi = new Float_t[kMulDigit] ;
258 // create the histogram for the single cluster
259 // 1. gets histogram boundaries
260 Float_t ximax = -999. ;
261 Float_t zimax = -999. ;
262 Float_t ximin = 999. ;
263 Float_t zimin = 999. ;
265 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
266 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
267 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
268 phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
269 if ( xi[iDigit] > ximax )
271 if ( xi[iDigit] < ximin )
273 if ( zi[iDigit] > zimax )
275 if ( zi[iDigit] < zimin )
278 ximax += phosgeom->GetCrystalSize(0) / 2. ;
279 zimax += phosgeom->GetCrystalSize(2) / 2. ;
280 ximin -= phosgeom->GetCrystalSize(0) / 2. ;
281 zimin -= phosgeom->GetCrystalSize(2) / 2. ;
282 Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
283 Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
285 // 2. gets the histogram title
288 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
294 histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
297 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
298 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
299 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
300 phosgeom->RelPosInModule(relid, x, z);
301 histo->Fill(x, z, fEnergyList[iDigit] ) ;
305 digitgraph = new TGraph(kMulDigit,xi,zi);
306 digitgraph-> SetMarkerStyle(5) ;
307 digitgraph-> SetMarkerSize(1.) ;
308 digitgraph-> SetMarkerColor(1) ;
309 digitgraph-> Paint("P") ;
313 histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
314 histocanvas->Draw() ;
315 histo->Draw("lego1") ;
333 //____________________________________________________________________________
334 void AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
336 // Calculates the dispersion of the shower at the origine of the RecPoint
344 AliPHOSDigit * digit ;
346 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
348 // Calculates the center of gravity in the local PHOS-module coordinates
352 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
353 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
357 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
358 phosgeom->RelPosInModule(relid, xi, zi);
359 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
368 // Calculates the dispersion in coordinates
370 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
371 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
375 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
376 phosgeom->RelPosInModule(relid, xi, zi);
377 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
378 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
387 fDispersion = TMath::Sqrt(d) ;
390 //______________________________________________________________________________
391 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
393 // This function calculates energy in the core,
394 // i.e. within a radius rad = 3cm around the center. Beyond this radius
395 // in accordance with shower profile the energy deposition
396 // should be less than 2%
398 Float_t coreRadius = 3 ;
403 AliPHOSDigit * digit ;
405 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
409 // Calculates the center of gravity in the local PHOS-module coordinates
411 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
412 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
416 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
417 phosgeom->RelPosInModule(relid, xi, zi);
418 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
427 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
428 digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
432 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
433 phosgeom->RelPosInModule(relid, xi, zi);
434 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
435 if(distance < coreRadius)
436 fCoreEnergy += fEnergyList[iDigit] ;
441 //____________________________________________________________________________
442 void AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
444 // Calculates the axis of the shower ellipsoid
453 AliPHOSDigit * digit ;
455 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
460 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
461 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
465 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
466 phosgeom->RelPosInModule(relid, xi, zi);
467 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
484 // //Apply correction due to non-perpendicular incidence
487 // AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
488 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
489 // Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
491 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
492 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
494 // dxx = dxx/(CosX*CosX) ;
495 // dzz = dzz/(CosZ*CosZ) ;
496 // dxz = dxz/(CosX*CosZ) ;
499 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
501 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
503 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
504 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
505 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
510 //____________________________________________________________________________
511 void AliPHOSEmcRecPoint::EvalMoments(Float_t logWeight,TClonesArray * digits)
513 // Calculate the shower moments in the eigen reference system
514 // M2x, M2z, M3x, M4z
515 // Calculate the angle between the shower position vector and the eigen vector
523 Double_t lambda0=0, lambda1=0;
525 AliPHOSDigit * digit ;
527 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
531 // 1) Find covariance matrix elements:
539 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
540 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
541 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
542 phosgeom->RelPosInModule(relid, xi, zi);
543 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
560 // 2) Find covariance matrix eigen values lambda0 and lambda1
562 lambda0 = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
563 lambda1 = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
565 // 3) Find covariance matrix eigen vectors e0 and e1
569 e0.Set(1.,(lambda0-dxx)/dxz);
574 e1.Set(-e0.Y(),e0.X());
576 // 4) Rotate cluster tensor from (x,z) to (e0,e1) system
577 // and calculate moments M3x and M4z
579 Float_t cosPhi = e0.X();
580 Float_t sinPhi = e0.Y();
584 Double_t dx3, dz3, dz4;
594 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
595 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
596 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
597 phosgeom->RelPosInModule(relid, xiPHOS, ziPHOS);
598 xi = xiPHOS*cosPhi + ziPHOS*sinPhi;
599 zi = ziPHOS*cosPhi - xiPHOS*sinPhi;
600 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
606 dx3 += w * xi * xi * xi;
607 dz3 += w * zi * zi * zi ;
608 dz4 += w * zi * zi * zi * zi ;
619 dx3 += -3*dxx*x + 2*x*x*x;
620 dz4 += -4*dz3*z + 6*dzz*z*z -3*z*z*z*z;
625 // 5) Find an angle between cluster center vector and eigen vector e0
627 Float_t phi = TMath::ACos ((x*e0.X() + z*e0.Y()) / sqrt(x*x + z*z));
636 //____________________________________________________________________________
637 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
639 // Evaluates all shower parameters
640 EvalLocalPosition(logWeight, digits) ;
641 EvalElipsAxis(logWeight, digits) ;
642 EvalMoments(logWeight, digits) ;
643 EvalDispersion(logWeight, digits) ;
644 EvalCoreEnergy(logWeight, digits);
646 AliPHOSRecPoint::EvalAll(digits) ;
648 //____________________________________________________________________________
649 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
651 // Calculates the center of gravity in the local PHOS-module coordinates
659 AliPHOSDigit * digit ;
661 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
665 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
666 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
670 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
671 phosgeom->RelPosInModule(relid, xi, zi);
672 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
681 // Correction for the depth of the shower starting point (TDR p 127)
682 Float_t para = 0.925 ;
683 Float_t parb = 6.52 ;
685 Float_t xo,yo,zo ; //Coordinates of the origin
686 //We should check all 3 possibilities to avoid seg.v.
687 if(gAlice && gAlice->GetMCApp() && gAlice->Generator())
688 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
692 Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
694 //Transform to the local ref.frame
696 xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
697 yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
699 Float_t radius = phosgeom->GetIPtoCrystalSurface()-yoL;
701 Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ;
702 Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
704 Float_t depthx = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ;
705 Float_t depthz = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ;
707 fLocPos.SetX(x - depthx) ;
709 fLocPos.SetZ(z - depthz) ;
714 //____________________________________________________________________________
715 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
717 // Finds the maximum energy in the cluster
719 Float_t menergy = 0. ;
723 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
725 if(fEnergyList[iDigit] > menergy)
726 menergy = fEnergyList[iDigit] ;
731 //____________________________________________________________________________
732 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(Float_t H) const
734 // Calculates the multiplicity of digits with energy larger than H*energy
738 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
740 if(fEnergyList[iDigit] > H * fAmp)
746 //____________________________________________________________________________
747 Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax( AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
748 Float_t locMaxCut,TClonesArray * digits) const
750 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
751 // energy difference between two local maxima
753 AliPHOSDigit * digit ;
754 AliPHOSDigit * digitN ;
760 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
761 maxAt[iDigit] = (AliPHOSDigit*) digits->At(fDigitsList[iDigit]) ;
764 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
766 digit = maxAt[iDigit] ;
768 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
769 if(iDigit == iDigitN)
772 digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ;
774 if ( AreNeighbours(digit, digitN) ) {
775 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
777 // but may be digit too is not local max ?
778 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
783 // but may be digitN too is not local max ?
784 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
787 } // if Areneighbours
793 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
795 maxAt[iDigitN] = maxAt[iDigit] ;
796 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
802 //____________________________________________________________________________
803 void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits)
805 // Define a rec.point time as a time in the cell with the maximum energy
809 for(Int_t idig=0; idig < fMulDigit; idig++){
810 if(fEnergyList[idig] > maxE){
811 maxE = fEnergyList[idig] ;
815 fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
818 //____________________________________________________________________________
819 void AliPHOSEmcRecPoint::Purify(Float_t threshold){
820 //Removes digits below threshold
822 Int_t * tempo = new Int_t[fMaxDigit];
823 Float_t * tempoE = new Float_t[fMaxDigit];
826 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
827 if(fEnergyList[iDigit] > threshold){
828 tempo[mult] = fDigitsList[iDigit] ;
829 tempoE[mult] = fEnergyList[iDigit] ;
835 delete [] fDigitsList ;
836 delete [] fEnergyList ;
837 fDigitsList = new Int_t[fMulDigit];
838 fEnergyList = new Float_t[fMulDigit];
840 fAmp = 0.; //Recalculate total energy
841 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
842 fDigitsList[iDigit] = tempo[iDigit];
843 fEnergyList[iDigit] = tempoE[iDigit] ;
844 fAmp+=tempoE[iDigit];
851 //____________________________________________________________________________
852 void AliPHOSEmcRecPoint::Print(Option_t *) const
854 // Print the list of digits belonging to the cluster
857 message = "AliPHOSEmcRecPoint:\n" ;
858 message += " digits # = " ;
859 AliInfo(Form(message.Data())) ;
862 for(iDigit=0; iDigit<fMulDigit; iDigit++)
863 printf(" %d ", fDigitsList[iDigit] ) ;
865 printf(" Energies = ") ;
866 for(iDigit=0; iDigit<fMulDigit; iDigit++)
867 printf(" %f ", fEnergyList[iDigit] ) ;
869 printf(" Primaries ") ;
870 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
871 printf(" %d ", fTracksList[iDigit]) ;
873 message = " Multiplicity = %d" ;
874 message += " Cluster Energy = %f" ;
875 message += " Number of primaries %d" ;
876 message += " Stored at position %d" ;
878 printf(message.Data(), fMulDigit, fAmp, fMulTrack,GetIndexInList() ) ;