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 ---
34 #include "AliPHOSLoader.h"
35 #include "AliGenerator.h"
36 #include "AliPHOSGeometry.h"
37 #include "AliPHOSDigit.h"
38 #include "AliPHOSEmcRecPoint.h"
40 ClassImp(AliPHOSEmcRecPoint)
42 //____________________________________________________________________________
43 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
51 fNExMax = 0 ; //Not unfolded yet
53 fLocPos.SetX(1000000.) ; //Local position should be evaluated
58 //____________________________________________________________________________
59 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const char * opt) : AliPHOSRecPoint(opt)
65 fNExMax = 0 ; //Not unfolded yet
69 fLocPos.SetX(1000000.) ; //Local position should be evaluated
74 //____________________________________________________________________________
75 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const AliPHOSEmcRecPoint & rp) : AliPHOSRecPoint(rp)
79 fMulDigit = rp.fMulDigit ;
81 fCoreEnergy = rp.fCoreEnergy ;
82 fEnergyList = new Float_t[rp.fMulDigit] ;
84 for(index = 0 ; index < fMulDigit ; index++)
85 fEnergyList[index] = rp.fEnergyList[index] ;
86 fNExMax = rp.fNExMax ;
90 //____________________________________________________________________________
91 AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
96 delete[] fEnergyList ;
99 //____________________________________________________________________________
100 void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
102 // Adds a digit to the RecPoint
103 // and accumulates the total amplitude and the multiplicity
106 fEnergyList = new Float_t[fMaxDigit];
108 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
110 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
111 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
114 for ( index = 0 ; index < fMulDigit ; index++ ){
115 tempo[index] = fDigitsList[index] ;
116 tempoE[index] = fEnergyList[index] ;
119 delete [] fDigitsList ;
120 fDigitsList = new ( Int_t[fMaxDigit] ) ;
122 delete [] fEnergyList ;
123 fEnergyList = new ( Float_t[fMaxDigit] ) ;
125 for ( index = 0 ; index < fMulDigit ; index++ ){
126 fDigitsList[index] = tempo[index] ;
127 fEnergyList[index] = tempoE[index] ;
134 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
135 fEnergyList[fMulDigit] = Energy ;
139 EvalPHOSMod(&digit) ;
142 //____________________________________________________________________________
143 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
145 // Tells if (true) or not (false) two digits are neighbors
147 Bool_t aren = kFALSE ;
149 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
152 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
155 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
157 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
158 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
160 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
166 //____________________________________________________________________________
167 Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
169 // Compares two RecPoints according to their position in the PHOS modules
171 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
172 //value (what is senseless) change as vell delta in
173 //AliPHOSTrackSegmentMakerv* and other RecPoints...
176 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
179 Int_t phosmod1 = GetPHOSMod() ;
180 Int_t phosmod2 = clu->GetPHOSMod() ;
183 GetLocalPosition(locpos1) ;
185 clu->GetLocalPosition(locpos2) ;
187 if(phosmod1 == phosmod2 ) {
188 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
193 else if(locpos1.Z()>locpos2.Z())
200 if(phosmod1 < phosmod2 )
208 //______________________________________________________________________________
209 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t, Int_t) const
212 // Execute action corresponding to one event
213 // This member function is called when a AliPHOSRecPoint is clicked with the locator
215 // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
216 // and switched off when the mouse button is released.
219 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
221 static TGraph * digitgraph = 0 ;
223 if (!gPad->IsEditable()) return;
226 TCanvas * histocanvas ;
229 //try to get run loader from default event folder
230 AliRunLoader* rn = AliRunLoader::GetRunLoader(AliConfig::fgkDefaultEventFolderName);
233 Error("ExecuteEvent","Can not find Run Loader in Default Event Folder");
236 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(rn->GetLoader("PHOSLoader"));
239 Error("ExecuteEvent","Can not find PHOS Loader from Run Loader");
244 const TClonesArray * digits = gime->Digits() ;
249 AliPHOSDigit * digit ;
253 const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
254 Float_t * xi = new Float_t[kMulDigit] ;
255 Float_t * zi = new Float_t[kMulDigit] ;
257 // create the histogram for the single cluster
258 // 1. gets histogram boundaries
259 Float_t ximax = -999. ;
260 Float_t zimax = -999. ;
261 Float_t ximin = 999. ;
262 Float_t zimin = 999. ;
264 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
265 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
266 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
267 phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
268 if ( xi[iDigit] > ximax )
270 if ( xi[iDigit] < ximin )
272 if ( zi[iDigit] > zimax )
274 if ( zi[iDigit] < zimin )
277 ximax += phosgeom->GetCrystalSize(0) / 2. ;
278 zimax += phosgeom->GetCrystalSize(2) / 2. ;
279 ximin -= phosgeom->GetCrystalSize(0) / 2. ;
280 zimin -= phosgeom->GetCrystalSize(2) / 2. ;
281 Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
282 Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
284 // 2. gets the histogram title
287 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
293 histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
296 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
297 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
298 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
299 phosgeom->RelPosInModule(relid, x, z);
300 histo->Fill(x, z, fEnergyList[iDigit] ) ;
304 digitgraph = new TGraph(kMulDigit,xi,zi);
305 digitgraph-> SetMarkerStyle(5) ;
306 digitgraph-> SetMarkerSize(1.) ;
307 digitgraph-> SetMarkerColor(1) ;
308 digitgraph-> Paint("P") ;
312 histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
313 histocanvas->Draw() ;
314 histo->Draw("lego1") ;
332 //____________________________________________________________________________
333 void AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
335 // Calculates the dispersion of the shower at the origine of the RecPoint
343 AliPHOSDigit * digit ;
345 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
347 // Calculates the center of gravity in the local PHOS-module coordinates
351 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
352 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
356 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
357 phosgeom->RelPosInModule(relid, xi, zi);
358 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
367 // Calculates the dispersion in coordinates
369 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
370 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
374 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
375 phosgeom->RelPosInModule(relid, xi, zi);
376 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
377 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
386 fDispersion = TMath::Sqrt(d) ;
389 //______________________________________________________________________________
390 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
392 // This function calculates energy in the core,
393 // i.e. within a radius rad = 3cm around the center. Beyond this radius
394 // in accordance with shower profile the energy deposition
395 // should be less than 2%
397 Float_t coreRadius = 3 ;
402 AliPHOSDigit * digit ;
404 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
408 // Calculates the center of gravity in the local PHOS-module coordinates
410 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
411 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
415 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
416 phosgeom->RelPosInModule(relid, xi, zi);
417 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
426 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
427 digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
431 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
432 phosgeom->RelPosInModule(relid, xi, zi);
433 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
434 if(distance < coreRadius)
435 fCoreEnergy += fEnergyList[iDigit] ;
440 //____________________________________________________________________________
441 void AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
443 // Calculates the axis of the shower ellipsoid
452 AliPHOSDigit * digit ;
454 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
459 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
460 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
464 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
465 phosgeom->RelPosInModule(relid, xi, zi);
466 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
483 // //Apply correction due to non-perpendicular incidence
486 // AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
487 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
488 // Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
490 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
491 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
493 // dxx = dxx/(CosX*CosX) ;
494 // dzz = dzz/(CosZ*CosZ) ;
495 // dxz = dxz/(CosX*CosZ) ;
498 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
500 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
502 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
503 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
504 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
509 //____________________________________________________________________________
510 void AliPHOSEmcRecPoint::EvalMoments(Float_t logWeight,TClonesArray * digits)
512 // Calculate the shower moments in the eigen reference system
513 // M2x, M2z, M3x, M4z
514 // Calculate the angle between the shower position vector and the eigen vector
522 Double_t lambda0=0, lambda1=0;
524 AliPHOSDigit * digit ;
526 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
530 // 1) Find covariance matrix elements:
538 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
539 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
540 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
541 phosgeom->RelPosInModule(relid, xi, zi);
542 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
559 // 2) Find covariance matrix eigen values lambda0 and lambda1
561 lambda0 = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
562 lambda1 = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
564 // 3) Find covariance matrix eigen vectors e0 and e1
568 e0.Set(1.,(lambda0-dxx)/dxz);
573 e1.Set(-e0.Y(),e0.X());
575 // 4) Rotate cluster tensor from (x,z) to (e0,e1) system
576 // and calculate moments M3x and M4z
578 Float_t cosPhi = e0.X();
579 Float_t sinPhi = e0.Y();
583 Double_t dx3, dz3, dz4;
593 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
594 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
595 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
596 phosgeom->RelPosInModule(relid, xiPHOS, ziPHOS);
597 xi = xiPHOS*cosPhi + ziPHOS*sinPhi;
598 zi = ziPHOS*cosPhi - xiPHOS*sinPhi;
599 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
605 dx3 += w * xi * xi * xi;
606 dz3 += w * zi * zi * zi ;
607 dz4 += w * zi * zi * zi * zi ;
618 dx3 += -3*dxx*x + 2*x*x*x;
619 dz4 += -4*dz3*z + 6*dzz*z*z -3*z*z*z*z;
624 // 5) Find an angle between cluster center vector and eigen vector e0
626 Float_t phi = TMath::ACos ((x*e0.X() + z*e0.Y()) / sqrt(x*x + z*z));
635 //____________________________________________________________________________
636 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
638 // 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 = AliPHOSLoader::GetPHOSGeometry();
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 ) ) ;
682 // Correction for the depth of the shower starting point (TDR p 127)
683 Float_t para = 0.925 ;
684 Float_t parb = 6.52 ;
686 Float_t xo,yo,zo ; //Coordinates of the origin
687 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
689 Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
691 //Transform to the local ref.frame
693 xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
694 yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
696 Float_t radius = phosgeom->GetIPtoCrystalSurface()-yoL;
698 Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ;
699 Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
701 Float_t depthx = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ;
702 Float_t depthz = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ;
704 fLocPos.SetX(x - depthx) ;
706 fLocPos.SetZ(z - depthz) ;
711 //____________________________________________________________________________
712 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
714 // Finds the maximum energy in the cluster
716 Float_t menergy = 0. ;
720 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
722 if(fEnergyList[iDigit] > menergy)
723 menergy = fEnergyList[iDigit] ;
728 //____________________________________________________________________________
729 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
731 // Calculates the multiplicity of digits with energy larger than H*energy
735 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
737 if(fEnergyList[iDigit] > H * fAmp)
743 //____________________________________________________________________________
744 Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax( AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
745 Float_t locMaxCut,TClonesArray * digits) const
747 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
748 // energy difference between two local maxima
750 AliPHOSDigit * digit ;
751 AliPHOSDigit * digitN ;
757 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
758 maxAt[iDigit] = (AliPHOSDigit*) digits->At(fDigitsList[iDigit]) ;
761 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
763 digit = maxAt[iDigit] ;
765 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
766 digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ;
768 if ( AreNeighbours(digit, digitN) ) {
769 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
771 // but may be digit too is not local max ?
772 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
777 // but may be digitN too is not local max ?
778 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
781 } // if Areneighbours
787 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
789 maxAt[iDigitN] = maxAt[iDigit] ;
790 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
796 //____________________________________________________________________________
797 void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits)
799 // Define a rec.point time as a time in the cell with the maximum energy
803 for(Int_t idig=0; idig < fMulDigit; idig++){
804 if(fEnergyList[idig] > maxE){
805 maxE = fEnergyList[idig] ;
809 fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
812 //____________________________________________________________________________
813 void AliPHOSEmcRecPoint::Purify(Float_t threshold){
814 //Removes digits below threshold
816 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
817 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
820 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
821 if(fEnergyList[iDigit] > threshold){
822 tempo[mult] = fDigitsList[iDigit] ;
823 tempoE[mult] = fEnergyList[iDigit] ;
829 delete [] fDigitsList ;
830 delete [] fEnergyList ;
831 fDigitsList = new (Int_t[fMulDigit]) ;
832 fEnergyList = new ( Float_t[fMulDigit]) ;
834 fAmp = 0.; //Recalculate total energy
835 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
836 fDigitsList[iDigit] = tempo[iDigit];
837 fEnergyList[iDigit] = tempoE[iDigit] ;
838 fAmp+=tempoE[iDigit];
845 //____________________________________________________________________________
846 void AliPHOSEmcRecPoint::Print(Option_t *) const
848 // Print the list of digits belonging to the cluster
851 message = "AliPHOSEmcRecPoint:\n" ;
852 message += " digits # = " ;
853 Info("Print", message.Data()) ;
856 for(iDigit=0; iDigit<fMulDigit; iDigit++)
857 printf(" %d ", fDigitsList[iDigit] ) ;
859 Info("Print", " Energies = ") ;
860 for(iDigit=0; iDigit<fMulDigit; iDigit++)
861 printf(" %f ", fEnergyList[iDigit] ) ;
863 Info("Print", " Primaries ") ;
864 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
865 printf(" %d ", fTracksList[iDigit]) ;
867 message = " Multiplicity = %d" ;
868 message += " Cluster Energy = %f" ;
869 message += " Number of primaries %d" ;
870 message += " Stored at position %d" ;
872 Info("Print", message.Data(), fMulDigit, fAmp, fMulTrack,GetIndexInList() ) ;