]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - PHOS/AliPHOSEmcRecPoint.cxx
can now only write mixed digits into one of the file that contains summable digits...
[u/mrichter/AliRoot.git] / PHOS / AliPHOSEmcRecPoint.cxx
... / ...
CommitLineData
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#include "AliPHOSIndexToObject.h"
42
43ClassImp(AliPHOSEmcRecPoint)
44
45//____________________________________________________________________________
46AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(Float_t W0, Float_t LocMaxCut)
47 : AliPHOSRecPoint()
48{
49 // ctor
50
51 fMulDigit = 0 ;
52 fAmp = 0. ;
53 fEnergyList = new Float_t[fMaxDigit];
54 fGeom = AliPHOSGeometry::GetInstance() ;
55 fDelta = ((AliPHOSGeometry *) fGeom)->GetCrystalSize(0) ;
56 fW0 = W0 ;
57 fLocMaxCut = LocMaxCut ;
58 fLocPos.SetX(1000000.) ; //Local position should be evaluated
59
60}
61
62//____________________________________________________________________________
63AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
64{
65 // dtor
66 if ( fEnergyList )
67 delete[] fEnergyList ;
68}
69
70//____________________________________________________________________________
71void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
72{
73 // Adds a digit to the RecPoint
74 // and accumulates the total amplitude and the multiplicity
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
108//____________________________________________________________________________
109Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
110{
111 // Tells if (true) or not (false) two digits are neighbors)
112
113 Bool_t aren = kFALSE ;
114
115 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
116 Int_t relid1[4] ;
117 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
118
119 Int_t relid2[4] ;
120 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
121
122 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
123 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
124
125 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
126 aren = kTRUE ;
127
128 return aren ;
129}
130
131//____________________________________________________________________________
132Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
133{
134 // Compares two RecPoints according to their position in the PHOS modules
135
136 Int_t rv ;
137
138 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
139
140
141 Int_t phosmod1 = GetPHOSMod() ;
142 Int_t phosmod2 = clu->GetPHOSMod() ;
143
144 TVector3 locpos1;
145 this->GetLocalPosition(locpos1) ;
146 TVector3 locpos2;
147 clu->GetLocalPosition(locpos2) ;
148
149 if(phosmod1 == phosmod2 ) {
150 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/fDelta)-(Int_t)TMath::Ceil(locpos2.X()/fDelta) ;
151 if (rowdif> 0)
152 rv = -1 ;
153 else if(rowdif < 0)
154 rv = 1 ;
155 else if(locpos1.Z()>locpos2.Z())
156 rv = -1 ;
157 else
158 rv = 1 ;
159 }
160
161 else {
162 if(phosmod1 < phosmod2 )
163 rv = -1 ;
164 else
165 rv = 1 ;
166 }
167
168 return rv ;
169}
170
171//______________________________________________________________________________
172void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py)
173{
174 // Execute action corresponding to one event
175 // This member function is called when a AliPHOSRecPoint is clicked with the locator
176 //
177 // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
178 // and switched off when the mouse button is released.
179 //
180
181 // static Int_t pxold, pyold;
182
183 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
184
185 static TGraph * digitgraph = 0 ;
186
187 if (!gPad->IsEditable()) return;
188
189 TH2F * histo = 0 ;
190 TCanvas * histocanvas ;
191
192 switch (event) {
193
194 case kButton1Down: {
195 AliPHOSDigit * digit ;
196 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
197 Int_t iDigit;
198 Int_t relid[4] ;
199
200 const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
201 Float_t * xi = new Float_t[kMulDigit] ;
202 Float_t * zi = new Float_t[kMulDigit] ;
203
204 // create the histogram for the single cluster
205 // 1. gets histogram boundaries
206 Float_t ximax = -999. ;
207 Float_t zimax = -999. ;
208 Float_t ximin = 999. ;
209 Float_t zimin = 999. ;
210
211 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
212 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
213 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
214 phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
215 if ( xi[iDigit] > ximax )
216 ximax = xi[iDigit] ;
217 if ( xi[iDigit] < ximin )
218 ximin = xi[iDigit] ;
219 if ( zi[iDigit] > zimax )
220 zimax = zi[iDigit] ;
221 if ( zi[iDigit] < zimin )
222 zimin = zi[iDigit] ;
223 }
224 ximax += phosgeom->GetCrystalSize(0) / 2. ;
225 zimax += phosgeom->GetCrystalSize(2) / 2. ;
226 ximin -= phosgeom->GetCrystalSize(0) / 2. ;
227 zimin -= phosgeom->GetCrystalSize(2) / 2. ;
228 Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
229 Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
230
231 // 2. gets the histogram title
232
233 Text_t title[100] ;
234 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
235
236 if (!histo) {
237 delete histo ;
238 histo = 0 ;
239 }
240 histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
241
242 Float_t x, z ;
243 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
244 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
245 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
246 phosgeom->RelPosInModule(relid, x, z);
247 histo->Fill(x, z, fEnergyList[iDigit] ) ;
248 }
249
250 if (!digitgraph) {
251 digitgraph = new TGraph(kMulDigit,xi,zi);
252 digitgraph-> SetMarkerStyle(5) ;
253 digitgraph-> SetMarkerSize(1.) ;
254 digitgraph-> SetMarkerColor(1) ;
255 digitgraph-> Paint("P") ;
256 }
257
258 Print() ;
259 histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ;
260 histocanvas->Draw() ;
261 histo->Draw("lego1") ;
262
263 delete[] xi ;
264 delete[] zi ;
265
266 break;
267 }
268
269 case kButton1Up:
270 if (digitgraph) {
271 delete digitgraph ;
272 digitgraph = 0 ;
273 }
274 break;
275
276 }
277}
278
279//____________________________________________________________________________
280Float_t AliPHOSEmcRecPoint::GetDispersion() const
281{
282 // Calculates the dispersion of the shower at the origine of the RecPoint
283
284 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
285
286 Float_t d = 0 ;
287 Float_t wtot = 0 ;
288
289 TVector3 locpos;
290 GetLocalPosition(locpos);
291 Float_t x = locpos.X() ;
292 Float_t z = locpos.Z() ;
293
294 AliPHOSDigit * digit ;
295 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
296
297 Int_t iDigit;
298 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
299 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
300 Int_t relid[4] ;
301 Float_t xi ;
302 Float_t zi ;
303 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
304 phosgeom->RelPosInModule(relid, xi, zi);
305 Float_t w = TMath::Max(0.,fW0+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
306 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
307 wtot+=w ;
308 }
309
310 d /= wtot ;
311
312 return TMath::Sqrt(d) ;
313}
314//______________________________________________________________________________
315Float_t AliPHOSEmcRecPoint::CoreEnergy()
316{
317 //This function calculates energy in the core,
318 //i.e. within radius rad = 3cm. Beyond this radius
319 //in accoradnce with shower profile energy deposition
320 // should be less than 2%
321
322 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
323
324 Float_t eCore = 0 ;
325 Float_t coreRadius = 3 ;
326
327 TVector3 locpos;
328 GetLocalPosition(locpos);
329 Float_t x = locpos.X() ;
330 Float_t z = locpos.Z() ;
331
332 AliPHOSDigit * digit ;
333 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
334
335 Int_t iDigit;
336 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
337 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
338 Int_t relid[4] ;
339 Float_t xi ;
340 Float_t zi ;
341 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
342 phosgeom->RelPosInModule(relid, xi, zi);
343 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
344 if(distance < coreRadius)
345 eCore += fEnergyList[iDigit] ;
346 }
347
348return eCore ;
349}
350
351//____________________________________________________________________________
352void AliPHOSEmcRecPoint::GetElipsAxis(Float_t * lambda)
353{
354 // Calculates the axis of the shower ellipsoid
355
356 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
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 *) ( please->GimeDigit(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.,fW0+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 lambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
407 if(lambda[0] > 0)
408 lambda[0] = TMath::Sqrt(lambda[0]) ;
409
410 lambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
411 if(lambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
412 lambda[1] = TMath::Sqrt(lambda[1]) ;
413 else
414 lambda[1]= 0. ;
415}
416
417//____________________________________________________________________________
418void AliPHOSEmcRecPoint::EvalAll(void )
419{
420 AliPHOSRecPoint::EvalAll() ;
421 EvalLocalPosition() ;
422}
423//____________________________________________________________________________
424void AliPHOSEmcRecPoint::EvalLocalPosition(void )
425{
426 // Calculates the center of gravity in the local PHOS-module coordinates
427
428 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
429
430 Float_t wtot = 0. ;
431
432 Int_t relid[4] ;
433
434 Float_t x = 0. ;
435 Float_t z = 0. ;
436
437 AliPHOSDigit * digit ;
438
439 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
440
441 Int_t iDigit;
442
443 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
444 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) );
445
446 Float_t xi ;
447 Float_t zi ;
448 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
449 phosgeom->RelPosInModule(relid, xi, zi);
450 Float_t w = TMath::Max( 0., fW0 + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
451 x += xi * w ;
452 z += zi * w ;
453 wtot += w ;
454
455 }
456
457 x /= wtot ;
458 z /= wtot ;
459
460 // Correction for the depth of the shower starting point (TDR p 127)
461 Float_t para = 0.925 ;
462 Float_t parb = 6.52 ;
463
464 Float_t xo,yo,zo ; //Coordinates of the origin
465 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
466
467 Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
468
469 //Transform to the local ref.frame
470 Float_t xoL,yoL ;
471 xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
472 yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
473
474 Float_t radius = TMath::Sqrt((xoL-x)*(xoL-x)+
475 (phosgeom->GetIPtoCrystalSurface()-yoL)*(phosgeom->GetIPtoCrystalSurface()-yoL)+
476 (zo-z)*(zo-z));
477
478 Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ;
479 Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
480
481 Float_t depthx = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ;
482 Float_t depthz = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ;
483
484
485 fLocPos.SetX(x - depthx) ;
486 fLocPos.SetY(0.) ;
487 fLocPos.SetZ(z - depthz) ;
488
489}
490
491//____________________________________________________________________________
492Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
493{
494 // Finds the maximum energy in the cluster
495
496 Float_t menergy = 0. ;
497
498 Int_t iDigit;
499
500 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
501
502 if(fEnergyList[iDigit] > menergy)
503 menergy = fEnergyList[iDigit] ;
504 }
505 return menergy ;
506}
507
508//____________________________________________________________________________
509Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
510{
511 // Calculates the multiplicity of digits with energy larger than H*energy
512
513 Int_t multipl = 0 ;
514 Int_t iDigit ;
515 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
516
517 if(fEnergyList[iDigit] > H * fAmp)
518 multipl++ ;
519 }
520 return multipl ;
521}
522
523//____________________________________________________________________________
524Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t * maxAt, Float_t * maxAtEnergy) const
525{
526 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
527 // energy difference between two local maxima
528
529 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
530
531 AliPHOSDigit * digit ;
532 AliPHOSDigit * digitN ;
533
534
535 Int_t iDigitN ;
536 Int_t iDigit ;
537
538 for(iDigit = 0; iDigit < fMulDigit; iDigit++){
539 maxAt[iDigit] = (Int_t) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
540 }
541
542 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
543 if(maxAt[iDigit] != -1) {
544 digit = (AliPHOSDigit *) maxAt[iDigit] ;
545
546 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
547 digitN = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigitN]) ) ;
548
549 if ( AreNeighbours(digit, digitN) ) {
550 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
551 maxAt[iDigitN] = -1 ;
552 // but may be digit too is not local max ?
553 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + fLocMaxCut)
554 maxAt[iDigit] = -1 ;
555 }
556 else {
557 maxAt[iDigit] = -1 ;
558 // but may be digitN too is not local max ?
559 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - fLocMaxCut)
560 maxAt[iDigitN] = -1 ;
561 }
562 } // if Areneighbours
563 } // while digitN
564 } // slot not empty
565 } // while digit
566
567 iDigitN = 0 ;
568 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
569 if(maxAt[iDigit] != -1){
570 maxAt[iDigitN] = maxAt[iDigit] ;
571 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
572 iDigitN++ ;
573 }
574 }
575 return iDigitN ;
576}
577
578
579// //____________________________________________________________________________
580// AliPHOSEmcRecPoint& AliPHOSEmcRecPoint::operator = (AliPHOSEmcRecPoint Clu)
581// {
582// int * dl = Clu.GetDigitsList() ;
583
584// if(fDigitsList)
585// delete fDigitsList ;
586
587// AliPHOSDigit * digit ;
588
589// Int_t iDigit;
590
591// for(iDigit=0; iDigit<fMulDigit; iDigit++) {
592// digit = (AliPHOSDigit *) dl[iDigit];
593// AddDigit(*digit) ;
594// }
595
596// fAmp = Clu.GetEnergy() ;
597// fGeom = Clu.GetGeom() ;
598// TVector3 locpos;
599// Clu.GetLocalPosition(locpos) ;
600// fLocPos = locpos;
601// fMulDigit = Clu.GetMultiplicity() ;
602// fMaxDigit = Clu.GetMaximumMultiplicity() ;
603// fPHOSMod = Clu.GetPHOSMod() ;
604// fW0 = Clu.GetLogWeightCut() ;
605// fDelta = Clu.GetDelta() ;
606// fLocMaxCut = Clu.GetLocMaxCut() ;
607
608// delete dl ;
609
610// return *this ;
611// }
612
613//____________________________________________________________________________
614void AliPHOSEmcRecPoint::Print(Option_t * option)
615{
616 // Print the list of digits belonging to the cluster
617
618 cout << "AliPHOSEmcRecPoint: " << endl ;
619
620 AliPHOSDigit * digit ;
621 Int_t iDigit;
622 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
623
624 Float_t xi ;
625 Float_t zi ;
626 Int_t relid[4] ;
627
628 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
629
630 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
631 digit = please->GimeDigit( fDigitsList[iDigit] ) ;
632 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
633 phosgeom->RelPosInModule(relid, xi, zi);
634 cout << " Id = " << digit->GetId() ;
635 cout << " module = " << relid[0] ;
636 cout << " x = " << xi ;
637 cout << " z = " << zi ;
638 cout << " Energy = " << fEnergyList[iDigit] << endl ;
639 }
640 cout << " Multiplicity = " << fMulDigit << endl ;
641 cout << " Cluster Energy = " << fAmp << endl ;
642 cout << " Stored at position " << GetIndexInList() << endl ;
643
644}
645//______________________________________________________________________________
646// void AliPHOSEmcRecPoint::Streamer(TBuffer &R__b)
647// {
648// // Stream an object of class AliPHOSEmcRecPoint.
649
650// if (R__b.IsReading()) {
651// Version_t R__v = R__b.ReadVersion(); if (R__v) { }
652// AliPHOSRecPoint::Streamer(R__b);
653// R__b >> fDelta;
654// fEnergyList = new Float_t[fMulDigit] ;
655// R__b.ReadFastArray(fEnergyList, fMulDigit);
656// R__b >> fLocMaxCut;
657// R__b >> fW0;
658// } else {
659// R__b.WriteVersion(AliPHOSEmcRecPoint::IsA());
660// AliPHOSRecPoint::Streamer(R__b);
661// R__b << fDelta;
662// R__b.WriteFastArray(fEnergyList, fMulDigit);
663// R__b << fLocMaxCut;
664// R__b << fW0;
665// }
666// }
667