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