]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSEmcRecPoint.cxx
overlapp between PUCP and PTCB corrected
[u/mrichter/AliRoot.git] / PHOS / AliPHOSEmcRecPoint.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
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 //  RecPoint implementation for PHOS-EMC 
20 //  An EmcRecPoint is a cluster of digits   
21 //           
22 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
23
24
25 // --- ROOT system ---
26 #include "TPad.h"
27 #include "TH2.h"
28 #include "TMath.h" 
29 #include "TCanvas.h" 
30
31 // --- Standard library ---
32
33 #include <iostream.h> 
34
35 // --- AliRoot header files ---
36
37  #include "AliGenerator.h"
38 #include "AliPHOSGeometry.h"
39 #include "AliPHOSEmcRecPoint.h"
40 #include "AliRun.h"
41
42 ClassImp(AliPHOSEmcRecPoint)
43
44 //____________________________________________________________________________
45 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
46 {
47   // ctor
48
49   fMulDigit   = 0 ;  
50   fAmp   = 0. ;   
51   fCoreEnergy = 0 ; 
52   fEnergyList = 0 ;
53   fGeom = AliPHOSGeometry::GetInstance() ;
54   fLocPos.SetX(1000000.)  ;      //Local position should be evaluated
55   
56 }
57
58 //____________________________________________________________________________
59 AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
60 {
61   // dtor
62   if ( fEnergyList )
63     delete[] fEnergyList ; 
64 }
65
66 //____________________________________________________________________________
67 void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
68 {
69   // Adds a digit to the RecPoint
70   //  and accumulates the total amplitude and the multiplicity 
71   
72   if(fEnergyList == 0)
73     fEnergyList =  new Float_t[fMaxDigit]; 
74
75   if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists 
76     fMaxDigit*=2 ; 
77     Int_t * tempo = new ( Int_t[fMaxDigit] ) ; 
78     Float_t * tempoE =  new ( Float_t[fMaxDigit] ) ;
79
80     Int_t index ;     
81     for ( index = 0 ; index < fMulDigit ; index++ ){
82       tempo[index]  = fDigitsList[index] ;
83       tempoE[index] = fEnergyList[index] ; 
84     }
85     
86     delete [] fDigitsList ; 
87     fDigitsList =  new ( Int_t[fMaxDigit] ) ;
88  
89     delete [] fEnergyList ;
90     fEnergyList =  new ( Float_t[fMaxDigit] ) ;
91
92     for ( index = 0 ; index < fMulDigit ; index++ ){
93       fDigitsList[index] = tempo[index] ;
94       fEnergyList[index] = tempoE[index] ; 
95     }
96  
97     delete [] tempo ;
98     delete [] tempoE ; 
99   } // if
100   
101   fDigitsList[fMulDigit]   = digit.GetIndexInList()  ; 
102   fEnergyList[fMulDigit]   = Energy ;
103   fMulDigit++ ; 
104   fAmp += Energy ; 
105
106   EvalPHOSMod(&digit) ;
107 }
108
109 //____________________________________________________________________________
110 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
111 {
112   // Tells if (true) or not (false) two digits are neighbors)
113   
114   Bool_t aren = kFALSE ;
115   
116   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
117   Int_t relid1[4] ; 
118   phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ; 
119
120   Int_t relid2[4] ; 
121   phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ; 
122   
123   Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
124   Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
125
126   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
127     aren = kTRUE ;
128   
129   return aren ;
130 }
131
132 //____________________________________________________________________________
133 Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
134 {
135   // Compares two RecPoints according to their position in the PHOS modules
136
137   Float_t delta = 1 ; //Width of "Sorting row". If you changibg this 
138                       //value (what is senseless) change as vell delta in
139                       //AliPHOSTrackSegmentMakerv* and other RecPoints...
140   Int_t rv ; 
141
142   AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ; 
143
144  
145   Int_t phosmod1 = GetPHOSMod() ;
146   Int_t phosmod2 = clu->GetPHOSMod() ;
147
148   TVector3 locpos1; 
149   GetLocalPosition(locpos1) ;
150   TVector3 locpos2;  
151   clu->GetLocalPosition(locpos2) ;  
152
153   if(phosmod1 == phosmod2 ) {
154     Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
155     if (rowdif> 0) 
156       rv = 1 ;
157      else if(rowdif < 0) 
158        rv = -1 ;
159     else if(locpos1.Z()>locpos2.Z()) 
160       rv = -1 ;
161     else 
162       rv = 1 ; 
163      }
164
165   else {
166     if(phosmod1 < phosmod2 ) 
167       rv = -1 ;
168     else 
169       rv = 1 ;
170   }
171
172   return rv ; 
173 }
174
175 //______________________________________________________________________________
176 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py)
177 {
178 //   Commented by Dmitri Peressounko: there is no possibility to ensure, 
179 //   that AliPHOSIndexToObject keeps the correct information.
180
181 //   // Execute action corresponding to one event
182 //   //  This member function is called when a AliPHOSRecPoint is clicked with the locator
183 //   //
184 //   //  If Left button is clicked on AliPHOSRecPoint, the digits are switched on    
185 //   //  and switched off when the mouse button is released.
186 //   //
187
188 //   //   static Int_t pxold, pyold;
189
190 //   AliPHOSIndexToObject * please =  AliPHOSIndexToObject::GetInstance() ; 
191   
192 //   static TGraph *  digitgraph = 0 ;
193   
194 //   if (!gPad->IsEditable()) return;
195   
196 //   TH2F * histo = 0 ;
197 //   TCanvas * histocanvas ; 
198   
199 //   switch (event) {
200     
201 //   case kButton1Down: {
202 //     AliPHOSDigit * digit ;
203 //     AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
204 //     Int_t iDigit;
205 //     Int_t relid[4] ;
206     
207 //     const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ; 
208 //     Float_t * xi = new Float_t[kMulDigit] ; 
209 //     Float_t * zi = new Float_t[kMulDigit] ; 
210     
211 //     // create the histogram for the single cluster 
212 //     // 1. gets histogram boundaries
213 //     Float_t ximax = -999. ; 
214 //     Float_t zimax = -999. ; 
215 //     Float_t ximin = 999. ; 
216 //     Float_t zimin = 999. ;
217     
218 //     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
219 //       digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
220 //       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
221 //       phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
222 //       if ( xi[iDigit] > ximax )
223 //      ximax = xi[iDigit] ; 
224 //       if ( xi[iDigit] < ximin )
225 //      ximin = xi[iDigit] ; 
226 //       if ( zi[iDigit] > zimax )
227 //      zimax = zi[iDigit] ; 
228 //       if ( zi[iDigit] < zimin )
229 //      zimin = zi[iDigit] ;     
230 //     }
231 //     ximax += phosgeom->GetCrystalSize(0) / 2. ;
232 //     zimax += phosgeom->GetCrystalSize(2) / 2. ;
233 //     ximin -= phosgeom->GetCrystalSize(0) / 2. ;
234 //     zimin -= phosgeom->GetCrystalSize(2) / 2. ;
235 //     Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5  ) ; 
236 //     Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
237     
238 //     // 2. gets the histogram title
239     
240 //     Text_t title[100] ; 
241 //     sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
242     
243 //     if (!histo) {
244 //       delete histo ; 
245 //       histo = 0 ; 
246 //     }
247 //     histo = new TH2F("cluster3D", title,  xdim, ximin, ximax, zdim, zimin, zimax)  ;
248     
249 //     Float_t x, z ; 
250 //     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
251 //       digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
252 //       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
253 //       phosgeom->RelPosInModule(relid, x, z);
254 //       histo->Fill(x, z, fEnergyList[iDigit] ) ;
255 //     }
256     
257 //     if (!digitgraph) {
258 //       digitgraph = new TGraph(kMulDigit,xi,zi);
259 //       digitgraph-> SetMarkerStyle(5) ; 
260 //       digitgraph-> SetMarkerSize(1.) ;
261 //       digitgraph-> SetMarkerColor(1) ;
262 //       digitgraph-> Paint("P") ;
263 //     }
264     
265 //     Print() ;
266 //     histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ; 
267 //     histocanvas->Draw() ; 
268 //     histo->Draw("lego1") ; 
269     
270 //     delete[] xi ; 
271 //     delete[] zi ; 
272     
273 //     break;
274 //   }
275   
276 //   case kButton1Up: 
277 //     if (digitgraph) {
278 //       delete digitgraph  ;
279 //       digitgraph = 0 ;
280 //     }
281 //     break;
282   
283 //    }
284 }
285
286 //____________________________________________________________________________
287 void  AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
288 {
289   // Calculates the dispersion of the shower at the origine of the RecPoint
290
291   Float_t d    = 0 ;
292   Float_t wtot = 0 ;
293
294   TVector3 locpos;
295   GetLocalPosition(locpos);
296   Float_t x = locpos.X() ;
297   Float_t z = locpos.Z() ;
298
299   AliPHOSDigit * digit ;
300   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
301   
302   Int_t iDigit;
303   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
304     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
305     Int_t relid[4] ;
306     Float_t xi ;
307     Float_t zi ;
308     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
309     phosgeom->RelPosInModule(relid, xi, zi);
310     Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
311     d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ; 
312     wtot+=w ;
313   }
314
315
316   d /= wtot ;
317
318   fDispersion = TMath::Sqrt(d) ;
319 }
320 //______________________________________________________________________________
321 void AliPHOSEmcRecPoint::EvalCoreEnergy(TClonesArray * digits)
322 {
323   //This function calculates energy in the core, 
324   //i.e. within radius rad = 3cm. Beyond this radius
325   //in accoradnce with shower profile energy deposition 
326   // should be less than 2%
327
328   Float_t coreRadius = 3 ;
329
330   TVector3 locpos;
331   GetLocalPosition(locpos);
332   Float_t x = locpos.X() ;
333   Float_t z = locpos.Z() ;
334
335   AliPHOSDigit * digit ;
336   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
337   
338   Int_t iDigit;
339   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
340     digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
341     Int_t relid[4] ;
342     Float_t xi ;
343     Float_t zi ;
344     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
345     phosgeom->RelPosInModule(relid, xi, zi);    
346     Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
347     if(distance < coreRadius)
348       fCoreEnergy += fEnergyList[iDigit] ;
349   }
350
351 }
352
353 //____________________________________________________________________________
354 void  AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
355 {
356   // Calculates the axis of the shower ellipsoid
357
358   Double_t wtot = 0. ;
359   Double_t x    = 0.;
360   Double_t z    = 0.;
361   Double_t dxx  = 0.;
362   Double_t dzz  = 0.;
363   Double_t dxz  = 0.;
364
365   AliPHOSDigit * digit ;
366   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
367   Int_t iDigit;
368
369   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
370     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
371     Int_t relid[4] ;
372     Float_t xi ;
373     Float_t zi ;
374     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
375     phosgeom->RelPosInModule(relid, xi, zi);
376     Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
377     dxx  += w * xi * xi ;
378     x    += w * xi ;
379     dzz  += w * zi * zi ;
380     z    += w * zi ; 
381     dxz  += w * xi * zi ; 
382     wtot += w ;
383   }
384   dxx /= wtot ;
385   x   /= wtot ;
386   dxx -= x * x ;
387   dzz /= wtot ;
388   z   /= wtot ;
389   dzz -= z * z ;
390   dxz /= wtot ;
391   dxz -= x * z ;
392
393 //   //Apply correction due to non-perpendicular incidence
394 //   Double_t CosX ;
395 //   Double_t CosZ ;
396 //   Double_t DistanceToIP= (Double_t ) ((AliPHOSGeometry *) fGeom)->GetIPtoCrystalSurface() ;
397   
398 //   CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
399 //   CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
400
401 //   dxx = dxx/(CosX*CosX) ;
402 //   dzz = dzz/(CosZ*CosZ) ;
403 //   dxz = dxz/(CosX*CosZ) ;
404
405
406   fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
407   if(fLambda[0] > 0)
408     fLambda[0] = TMath::Sqrt(fLambda[0]) ;
409
410   fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
411   if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
412     fLambda[1] = TMath::Sqrt(fLambda[1]) ;
413   else
414     fLambda[1]= 0. ;
415 }
416
417 //____________________________________________________________________________
418 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
419 {
420   AliPHOSRecPoint::EvalAll(logWeight,digits) ;
421   EvalLocalPosition(logWeight, digits) ;
422   EvalElipsAxis(logWeight, digits) ;
423   EvalDispersion(logWeight, digits) ;
424   EvalCoreEnergy(digits);
425 }
426 //____________________________________________________________________________
427 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
428 {
429   // Calculates the center of gravity in the local PHOS-module coordinates 
430
431   Float_t wtot = 0. ;
432  
433   Int_t relid[4] ;
434
435   Float_t x = 0. ;
436   Float_t z = 0. ;
437   
438   AliPHOSDigit * digit ;
439
440   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
441
442   Int_t iDigit;
443
444   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
445     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
446
447     Float_t xi ;
448     Float_t zi ;
449     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
450     phosgeom->RelPosInModule(relid, xi, zi);
451     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
452     x    += xi * w ;
453     z    += zi * w ;
454     wtot += w ;
455
456   }
457
458   x /= wtot ;
459   z /= wtot ;
460
461   // Correction for the depth of the shower starting point (TDR p 127)  
462   Float_t para = 0.925 ; 
463   Float_t parb = 6.52 ; 
464
465   Float_t xo,yo,zo ; //Coordinates of the origin
466   gAlice->Generator()->GetOrigin(xo,yo,zo) ;
467
468   Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
469
470   //Transform to the local ref.frame
471   Float_t xoL,yoL ;
472   xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
473   yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
474   
475   Float_t radius = TMath::Sqrt((xoL-x)*(xoL-x)+
476                                (phosgeom->GetIPtoCrystalSurface()-yoL)*(phosgeom->GetIPtoCrystalSurface()-yoL)+
477                                (zo-z)*(zo-z));
478  
479   Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ; 
480   Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
481  
482   Float_t depthx =  ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ; 
483   Float_t depthz =  ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ; 
484   
485
486   fLocPos.SetX(x - depthx)  ;
487   fLocPos.SetY(0.) ;
488   fLocPos.SetZ(z - depthz)  ;
489
490   fLocPosM = 0 ;
491 }
492
493 //____________________________________________________________________________
494 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
495 {
496   // Finds the maximum energy in the cluster
497   
498   Float_t menergy = 0. ;
499
500   Int_t iDigit;
501
502   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
503  
504     if(fEnergyList[iDigit] > menergy) 
505       menergy = fEnergyList[iDigit] ;
506   }
507   return menergy ;
508 }
509
510 //____________________________________________________________________________
511 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
512 {
513   // Calculates the multiplicity of digits with energy larger than H*energy 
514   
515   Int_t multipl   = 0 ;
516   Int_t iDigit ;
517   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
518
519     if(fEnergyList[iDigit] > H * fAmp) 
520       multipl++ ;
521   }
522   return multipl ;
523 }
524
525 //____________________________________________________________________________
526 Int_t  AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t *  maxAt, Float_t * maxAtEnergy,
527                                                Float_t locMaxCut,TClonesArray * digits) const
528
529   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
530   //  energy difference between two local maxima
531
532   AliPHOSDigit * digit ;
533   AliPHOSDigit * digitN ;
534   
535
536   Int_t iDigitN ;
537   Int_t iDigit ;
538
539   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
540     maxAt[iDigit] = (Int_t) digits->At(fDigitsList[iDigit])  ;
541
542   
543   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
544     if(maxAt[iDigit] != -1) {
545       digit = (AliPHOSDigit *) maxAt[iDigit] ;
546           
547       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
548         digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ; 
549         
550         if ( AreNeighbours(digit, digitN) ) {
551           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
552             maxAt[iDigitN] = -1 ;
553             // but may be digit too is not local max ?
554             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
555               maxAt[iDigit] = -1 ;
556           }
557           else {
558             maxAt[iDigit] = -1 ;
559             // but may be digitN too is not local max ?
560             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
561               maxAt[iDigitN] = -1 ; 
562           } 
563         } // if Areneighbours
564       } // while digitN
565     } // slot not empty
566   } // while digit
567   
568   iDigitN = 0 ;
569   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
570     if(maxAt[iDigit] != -1){
571       maxAt[iDigitN] = maxAt[iDigit] ;
572       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
573       iDigitN++ ; 
574     }
575   }
576   return iDigitN ;
577 }
578
579
580
581 //____________________________________________________________________________
582 void AliPHOSEmcRecPoint::Print(Option_t * option) 
583 {
584   // Print the list of digits belonging to the cluster
585   
586   cout << "AliPHOSEmcRecPoint: " << endl ;
587
588   Int_t iDigit;
589   cout << " digits # = " ;
590   for(iDigit=0; iDigit<fMulDigit; iDigit++)
591     cout << fDigitsList[iDigit] << "  " ;  
592   cout << endl ;
593   
594   cout << " Energies = " ;
595   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
596     cout  << fEnergyList[iDigit] << "  ";
597   cout << endl ;
598   
599   cout << " Primaries  " ;
600   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
601     cout << fTracksList[iDigit] << " " << endl ;
602         
603   cout << "       Multiplicity    = " << fMulDigit  << endl ;
604   cout << "       Cluster Energy  = " << fAmp << endl ;
605   cout << "       Number of primaries " << fMulTrack << endl ;
606   cout << "       Stored at position " << GetIndexInList() << endl ; 
607  
608 }
609  
610