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 /* History of cvs commits:
23 //_________________________________________________________________________
24 // RecPoint implementation for PHOS-EMC
25 // An EmcRecPoint is a cluster of digits
27 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
30 // --- ROOT system ---
36 // --- Standard library ---
38 // --- AliRoot header files ---
40 #include "AliPHOSLoader.h"
41 #include "AliGenerator.h"
42 #include "AliPHOSGeometry.h"
43 #include "AliPHOSDigit.h"
44 #include "AliPHOSEmcRecPoint.h"
46 ClassImp(AliPHOSEmcRecPoint)
48 //____________________________________________________________________________
49 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
57 fNExMax = 0 ; //Not unfolded yet
59 fLocPos.SetX(1000000.) ; //Local position should be evaluated
64 //____________________________________________________________________________
65 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const char * opt) : AliPHOSRecPoint(opt)
71 fNExMax = 0 ; //Not unfolded yet
75 fLocPos.SetX(1000000.) ; //Local position should be evaluated
80 //____________________________________________________________________________
81 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const AliPHOSEmcRecPoint & rp) : AliPHOSRecPoint(rp)
85 fMulDigit = rp.fMulDigit ;
87 fCoreEnergy = rp.fCoreEnergy ;
88 fEnergyList = new Float_t[rp.fMulDigit] ;
90 for(index = 0 ; index < fMulDigit ; index++)
91 fEnergyList[index] = rp.fEnergyList[index] ;
92 fNExMax = rp.fNExMax ;
96 //____________________________________________________________________________
97 AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
102 delete[] fEnergyList ;
105 //____________________________________________________________________________
106 void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
108 // Adds a digit to the RecPoint
109 // and accumulates the total amplitude and the multiplicity
112 fEnergyList = new Float_t[fMaxDigit];
114 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
116 Int_t * tempo = new Int_t[fMaxDigit];
117 Float_t * tempoE = new Float_t[fMaxDigit];
120 for ( index = 0 ; index < fMulDigit ; index++ ){
121 tempo[index] = fDigitsList[index] ;
122 tempoE[index] = fEnergyList[index] ;
125 delete [] fDigitsList ;
126 fDigitsList = new Int_t[fMaxDigit];
128 delete [] fEnergyList ;
129 fEnergyList = new Float_t[fMaxDigit];
131 for ( index = 0 ; index < fMulDigit ; index++ ){
132 fDigitsList[index] = tempo[index] ;
133 fEnergyList[index] = tempoE[index] ;
140 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
141 fEnergyList[fMulDigit] = Energy ;
145 EvalPHOSMod(&digit) ;
148 //____________________________________________________________________________
149 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
151 // Tells if (true) or not (false) two digits are neighbors
153 Bool_t aren = kFALSE ;
155 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
158 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
161 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
163 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
164 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
166 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
172 //____________________________________________________________________________
173 Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
175 // Compares two RecPoints according to their position in the PHOS modules
177 const Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
178 //value (what is senseless) change as vell delta in
179 //AliPHOSTrackSegmentMakerv* and other RecPoints...
182 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
185 Int_t phosmod1 = GetPHOSMod() ;
186 Int_t phosmod2 = clu->GetPHOSMod() ;
189 GetLocalPosition(locpos1) ;
191 clu->GetLocalPosition(locpos2) ;
193 if(phosmod1 == phosmod2 ) {
194 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
199 else if(locpos1.Z()>locpos2.Z())
206 if(phosmod1 < phosmod2 )
214 //______________________________________________________________________________
215 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t, Int_t) /*const*/
218 // Execute action corresponding to one event
219 // This member function is called when a AliPHOSRecPoint is clicked with the locator
221 // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
222 // and switched off when the mouse button is released.
225 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
227 static TGraph * digitgraph = 0 ;
229 if (!gPad->IsEditable()) return;
232 TCanvas * histocanvas ;
235 //try to get run loader from default event folder
236 AliRunLoader* rn = AliRunLoader::GetRunLoader(AliConfig::GetDefaultEventFolderName());
239 AliError(Form("Cannot find Run Loader in Default Event Folder"));
242 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(rn->GetLoader("PHOSLoader"));
245 AliError(Form("Cannot find PHOS Loader from Run Loader"));
250 const TClonesArray * digits = gime->Digits() ;
255 AliPHOSDigit * digit ;
259 const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
260 Float_t * xi = new Float_t[kMulDigit] ;
261 Float_t * zi = new Float_t[kMulDigit] ;
263 // create the histogram for the single cluster
264 // 1. gets histogram boundaries
265 Float_t ximax = -999. ;
266 Float_t zimax = -999. ;
267 Float_t ximin = 999. ;
268 Float_t zimin = 999. ;
270 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
271 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
272 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
273 phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
274 if ( xi[iDigit] > ximax )
276 if ( xi[iDigit] < ximin )
278 if ( zi[iDigit] > zimax )
280 if ( zi[iDigit] < zimin )
283 ximax += phosgeom->GetCrystalSize(0) / 2. ;
284 zimax += phosgeom->GetCrystalSize(2) / 2. ;
285 ximin -= phosgeom->GetCrystalSize(0) / 2. ;
286 zimin -= phosgeom->GetCrystalSize(2) / 2. ;
287 Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
288 Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
290 // 2. gets the histogram title
293 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
299 histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
302 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
303 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
304 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
305 phosgeom->RelPosInModule(relid, x, z);
306 histo->Fill(x, z, fEnergyList[iDigit] ) ;
310 digitgraph = new TGraph(kMulDigit,xi,zi);
311 digitgraph-> SetMarkerStyle(5) ;
312 digitgraph-> SetMarkerSize(1.) ;
313 digitgraph-> SetMarkerColor(1) ;
314 digitgraph-> Paint("P") ;
318 histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
319 histocanvas->Draw() ;
320 histo->Draw("lego1") ;
338 //____________________________________________________________________________
339 void AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
341 // Calculates the dispersion of the shower at the origine of the RecPoint
349 AliPHOSDigit * digit ;
351 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
353 // Calculates the center of gravity in the local PHOS-module coordinates
357 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
358 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
362 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
363 phosgeom->RelPosInModule(relid, xi, zi);
364 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
373 // Calculates the dispersion in coordinates
375 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
376 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
380 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
381 phosgeom->RelPosInModule(relid, xi, zi);
382 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
383 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
392 fDispersion = TMath::Sqrt(d) ;
395 //______________________________________________________________________________
396 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
398 // This function calculates energy in the core,
399 // i.e. within a radius rad = 3cm around the center. Beyond this radius
400 // in accordance with shower profile the energy deposition
401 // should be less than 2%
403 Float_t coreRadius = 3 ;
408 AliPHOSDigit * digit ;
410 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
414 // Calculates the center of gravity in the local PHOS-module coordinates
416 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
417 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
421 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
422 phosgeom->RelPosInModule(relid, xi, zi);
423 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
432 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
433 digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
437 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
438 phosgeom->RelPosInModule(relid, xi, zi);
439 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
440 if(distance < coreRadius)
441 fCoreEnergy += fEnergyList[iDigit] ;
446 //____________________________________________________________________________
447 void AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
449 // Calculates the axis of the shower ellipsoid
458 AliPHOSDigit * digit ;
460 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
465 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
466 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
470 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
471 phosgeom->RelPosInModule(relid, xi, zi);
472 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
489 // //Apply correction due to non-perpendicular incidence
492 // AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
493 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
494 // Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
496 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
497 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
499 // dxx = dxx/(CosX*CosX) ;
500 // dzz = dzz/(CosZ*CosZ) ;
501 // dxz = dxz/(CosX*CosZ) ;
504 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
506 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
508 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
509 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
510 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
515 //____________________________________________________________________________
516 void AliPHOSEmcRecPoint::EvalMoments(Float_t logWeight,TClonesArray * digits)
518 // Calculate the shower moments in the eigen reference system
519 // M2x, M2z, M3x, M4z
520 // Calculate the angle between the shower position vector and the eigen vector
528 Double_t lambda0=0, lambda1=0;
530 AliPHOSDigit * digit ;
532 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
536 // 1) Find covariance matrix elements:
544 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
545 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
546 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
547 phosgeom->RelPosInModule(relid, xi, zi);
548 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
565 // 2) Find covariance matrix eigen values lambda0 and lambda1
567 lambda0 = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
568 lambda1 = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
570 // 3) Find covariance matrix eigen vectors e0 and e1
574 e0.Set(1.,(lambda0-dxx)/dxz);
579 e1.Set(-e0.Y(),e0.X());
581 // 4) Rotate cluster tensor from (x,z) to (e0,e1) system
582 // and calculate moments M3x and M4z
584 Float_t cosPhi = e0.X();
585 Float_t sinPhi = e0.Y();
589 Double_t dx3, dz3, dz4;
599 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
600 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
601 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
602 phosgeom->RelPosInModule(relid, xiPHOS, ziPHOS);
603 xi = xiPHOS*cosPhi + ziPHOS*sinPhi;
604 zi = ziPHOS*cosPhi - xiPHOS*sinPhi;
605 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
611 dx3 += w * xi * xi * xi;
612 dz3 += w * zi * zi * zi ;
613 dz4 += w * zi * zi * zi * zi ;
624 dx3 += -3*dxx*x + 2*x*x*x;
625 dz4 += -4*dz3*z + 6*dzz*z*z -3*z*z*z*z;
630 // 5) Find an angle between cluster center vector and eigen vector e0
632 Float_t phi = TMath::ACos ((x*e0.X() + z*e0.Y()) / sqrt(x*x + z*z));
641 //____________________________________________________________________________
642 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
644 // Evaluates all shower parameters
645 EvalLocalPosition(logWeight, digits) ;
646 EvalElipsAxis(logWeight, digits) ;
647 EvalMoments(logWeight, digits) ;
648 EvalDispersion(logWeight, digits) ;
649 EvalCoreEnergy(logWeight, digits);
651 AliPHOSRecPoint::EvalAll(digits) ;
653 //____________________________________________________________________________
654 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
656 // Calculates the center of gravity in the local PHOS-module coordinates
664 AliPHOSDigit * digit ;
666 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
670 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
671 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
675 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
676 phosgeom->RelPosInModule(relid, xi, zi);
677 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
686 // Correction for the depth of the shower starting point (TDR p 127)
687 Float_t para = 0.925 ;
688 Float_t parb = 6.52 ;
690 Float_t xo,yo,zo ; //Coordinates of the origin
691 //We should check all 3 possibilities to avoid seg.v.
692 if(gAlice && gAlice->GetMCApp() && gAlice->Generator())
693 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
697 Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
699 //Transform to the local ref.frame
701 xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
702 yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
704 Float_t radius = phosgeom->GetIPtoCrystalSurface()-yoL;
706 Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ;
707 Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
709 Float_t depthx = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ;
710 Float_t depthz = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ;
712 fLocPos.SetX(x - depthx) ;
714 fLocPos.SetZ(z - depthz) ;
719 //____________________________________________________________________________
720 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
722 // Finds the maximum energy in the cluster
724 Float_t menergy = 0. ;
728 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
730 if(fEnergyList[iDigit] > menergy)
731 menergy = fEnergyList[iDigit] ;
736 //____________________________________________________________________________
737 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(Float_t H) const
739 // Calculates the multiplicity of digits with energy larger than H*energy
743 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
745 if(fEnergyList[iDigit] > H * fAmp)
751 //____________________________________________________________________________
752 Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax( AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
753 Float_t locMaxCut,TClonesArray * digits) const
755 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
756 // energy difference between two local maxima
758 AliPHOSDigit * digit ;
759 AliPHOSDigit * digitN ;
765 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
766 maxAt[iDigit] = (AliPHOSDigit*) digits->At(fDigitsList[iDigit]) ;
769 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
771 digit = maxAt[iDigit] ;
773 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
774 if(iDigit == iDigitN)
777 digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ;
779 if ( AreNeighbours(digit, digitN) ) {
780 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
782 // but may be digit too is not local max ?
783 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
788 // but may be digitN too is not local max ?
789 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
792 } // if Areneighbours
798 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
800 maxAt[iDigitN] = maxAt[iDigit] ;
801 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
807 //____________________________________________________________________________
808 void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits)
810 // Define a rec.point time as a time in the cell with the maximum energy
814 for(Int_t idig=0; idig < fMulDigit; idig++){
815 if(fEnergyList[idig] > maxE){
816 maxE = fEnergyList[idig] ;
820 fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
823 //____________________________________________________________________________
824 void AliPHOSEmcRecPoint::Purify(Float_t threshold){
825 //Removes digits below threshold
827 Int_t * tempo = new Int_t[fMaxDigit];
828 Float_t * tempoE = new Float_t[fMaxDigit];
831 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
832 if(fEnergyList[iDigit] > threshold){
833 tempo[mult] = fDigitsList[iDigit] ;
834 tempoE[mult] = fEnergyList[iDigit] ;
840 delete [] fDigitsList ;
841 delete [] fEnergyList ;
842 fDigitsList = new Int_t[fMulDigit];
843 fEnergyList = new Float_t[fMulDigit];
845 fAmp = 0.; //Recalculate total energy
846 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
847 fDigitsList[iDigit] = tempo[iDigit];
848 fEnergyList[iDigit] = tempoE[iDigit] ;
849 fAmp+=tempoE[iDigit];
856 //____________________________________________________________________________
857 void AliPHOSEmcRecPoint::Print(Option_t *) const
859 // Print the list of digits belonging to the cluster
862 message = "AliPHOSEmcRecPoint:\n" ;
863 message += " digits # = " ;
864 AliInfo(Form(message.Data())) ;
867 for(iDigit=0; iDigit<fMulDigit; iDigit++)
868 printf(" %d ", fDigitsList[iDigit] ) ;
870 printf(" Energies = ") ;
871 for(iDigit=0; iDigit<fMulDigit; iDigit++)
872 printf(" %f ", fEnergyList[iDigit] ) ;
874 printf(" Primaries ") ;
875 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
876 printf(" %d ", fTracksList[iDigit]) ;
878 message = " Multiplicity = %d" ;
879 message += " Cluster Energy = %f" ;
880 message += " Number of primaries %d" ;
881 message += " Stored at position %d" ;
883 printf(message.Data(), fMulDigit, fAmp, fMulTrack,GetIndexInList() ) ;