]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSEmcRecPoint.cxx
Simplification for Fast Simulator
[u/mrichter/AliRoot.git] / PHOS / AliPHOSEmcRecPoint.cxx
CommitLineData
d15a28e7 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
b2a60966 16/* $Id$ */
17
d15a28e7 18//_________________________________________________________________________
b2a60966 19// RecPoint implementation for PHOS-EMC
20// An EmcRecPoint is a cluster of digits
21//
22//*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
23
d15a28e7 24
25// --- ROOT system ---
9f616d61 26#include "TPad.h"
27#include "TH2.h"
d15a28e7 28#include "TMath.h"
9f616d61 29#include "TCanvas.h"
d15a28e7 30
31// --- Standard library ---
32
de9ec31b 33#include <iostream.h>
d15a28e7 34
35// --- AliRoot header files ---
36
7932f811 37 #include "AliGenerator.h"
d15a28e7 38#include "AliPHOSGeometry.h"
39#include "AliPHOSEmcRecPoint.h"
40#include "AliRun.h"
41
42ClassImp(AliPHOSEmcRecPoint)
43
44//____________________________________________________________________________
7932f811 45AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
d15a28e7 46{
47 // ctor
48
49 fMulDigit = 0 ;
50 fAmp = 0. ;
7932f811 51 fCoreEnergy = 0 ;
52 fEnergyList = 0 ;
ad8cfaf4 53 fGeom = AliPHOSGeometry::GetInstance() ;
d15a28e7 54 fLocPos.SetX(1000000.) ; //Local position should be evaluated
69183710 55
83974468 56}
57
58//____________________________________________________________________________
59AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
60{
88714635 61 // dtor
83974468 62 if ( fEnergyList )
63 delete[] fEnergyList ;
d15a28e7 64}
65
d15a28e7 66//____________________________________________________________________________
83974468 67void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
d15a28e7 68{
b2a60966 69 // Adds a digit to the RecPoint
70 // and accumulates the total amplitude and the multiplicity
d15a28e7 71
7932f811 72 if(fEnergyList == 0)
73 fEnergyList = new Float_t[fMaxDigit];
74
d15a28e7 75 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
9f616d61 76 fMaxDigit*=2 ;
83974468 77 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
9f616d61 78 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
79
80 Int_t index ;
d15a28e7 81 for ( index = 0 ; index < fMulDigit ; index++ ){
83974468 82 tempo[index] = fDigitsList[index] ;
d15a28e7 83 tempoE[index] = fEnergyList[index] ;
84 }
85
9f616d61 86 delete [] fDigitsList ;
83974468 87 fDigitsList = new ( Int_t[fMaxDigit] ) ;
9f616d61 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
d15a28e7 100
83974468 101 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
102 fEnergyList[fMulDigit] = Energy ;
103 fMulDigit++ ;
d15a28e7 104 fAmp += Energy ;
7932f811 105
106 EvalPHOSMod(&digit) ;
d15a28e7 107}
108
109//____________________________________________________________________________
ad8cfaf4 110Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
d15a28e7 111{
b2a60966 112 // Tells if (true) or not (false) two digits are neighbors)
d15a28e7 113
114 Bool_t aren = kFALSE ;
115
92862013 116 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
d15a28e7 117 Int_t relid1[4] ;
92862013 118 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
d15a28e7 119
120 Int_t relid2[4] ;
92862013 121 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
d15a28e7 122
92862013 123 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
124 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
d15a28e7 125
92862013 126 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
d15a28e7 127 aren = kTRUE ;
128
129 return aren ;
130}
131
132//____________________________________________________________________________
2a941f4e 133Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
d15a28e7 134{
b2a60966 135 // Compares two RecPoints according to their position in the PHOS modules
136
7932f811 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...
d15a28e7 140 Int_t rv ;
141
142 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
143
144
ad8cfaf4 145 Int_t phosmod1 = GetPHOSMod() ;
92862013 146 Int_t phosmod2 = clu->GetPHOSMod() ;
d15a28e7 147
92862013 148 TVector3 locpos1;
7932f811 149 GetLocalPosition(locpos1) ;
92862013 150 TVector3 locpos2;
151 clu->GetLocalPosition(locpos2) ;
d15a28e7 152
92862013 153 if(phosmod1 == phosmod2 ) {
7932f811 154 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
d15a28e7 155 if (rowdif> 0)
7932f811 156 rv = 1 ;
d15a28e7 157 else if(rowdif < 0)
7932f811 158 rv = -1 ;
92862013 159 else if(locpos1.Z()>locpos2.Z())
d15a28e7 160 rv = -1 ;
161 else
162 rv = 1 ;
163 }
164
165 else {
92862013 166 if(phosmod1 < phosmod2 )
d15a28e7 167 rv = -1 ;
168 else
169 rv = 1 ;
170 }
171
172 return rv ;
173}
174
9f616d61 175//______________________________________________________________________________
176void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py)
177{
7932f811 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// //
9f616d61 187
7932f811 188// // static Int_t pxold, pyold;
9f616d61 189
7932f811 190// AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
9f616d61 191
7932f811 192// static TGraph * digitgraph = 0 ;
83974468 193
7932f811 194// if (!gPad->IsEditable()) return;
83974468 195
7932f811 196// TH2F * histo = 0 ;
197// TCanvas * histocanvas ;
83974468 198
7932f811 199// switch (event) {
83974468 200
7932f811 201// case kButton1Down: {
202// AliPHOSDigit * digit ;
203// AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
204// Int_t iDigit;
205// Int_t relid[4] ;
83974468 206
7932f811 207// const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
208// Float_t * xi = new Float_t[kMulDigit] ;
209// Float_t * zi = new Float_t[kMulDigit] ;
83974468 210
7932f811 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. ;
83974468 217
7932f811 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 ) ;
83974468 237
7932f811 238// // 2. gets the histogram title
83974468 239
7932f811 240// Text_t title[100] ;
241// sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
83974468 242
7932f811 243// if (!histo) {
244// delete histo ;
245// histo = 0 ;
246// }
247// histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
83974468 248
7932f811 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// }
83974468 256
7932f811 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// }
83974468 264
7932f811 265// Print() ;
266// histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
267// histocanvas->Draw() ;
268// histo->Draw("lego1") ;
83974468 269
7932f811 270// delete[] xi ;
271// delete[] zi ;
83974468 272
7932f811 273// break;
274// }
83974468 275
7932f811 276// case kButton1Up:
277// if (digitgraph) {
278// delete digitgraph ;
279// digitgraph = 0 ;
280// }
281// break;
9f616d61 282
7932f811 283// }
9f616d61 284}
285
d15a28e7 286//____________________________________________________________________________
7932f811 287void AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
d15a28e7 288{
b2a60966 289 // Calculates the dispersion of the shower at the origine of the RecPoint
290
92862013 291 Float_t d = 0 ;
d15a28e7 292 Float_t wtot = 0 ;
293
92862013 294 TVector3 locpos;
295 GetLocalPosition(locpos);
296 Float_t x = locpos.X() ;
297 Float_t z = locpos.Z() ;
d15a28e7 298
299 AliPHOSDigit * digit ;
92862013 300 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
d15a28e7 301
302 Int_t iDigit;
88714635 303 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
7932f811 304 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
d15a28e7 305 Int_t relid[4] ;
306 Float_t xi ;
307 Float_t zi ;
92862013 308 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
309 phosgeom->RelPosInModule(relid, xi, zi);
7932f811 310 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
92862013 311 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
d15a28e7 312 wtot+=w ;
313 }
314
7932f811 315
92862013 316 d /= wtot ;
d15a28e7 317
7932f811 318 fDispersion = TMath::Sqrt(d) ;
d15a28e7 319}
fad3e5b9 320//______________________________________________________________________________
7932f811 321void AliPHOSEmcRecPoint::EvalCoreEnergy(TClonesArray * digits)
fad3e5b9 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
fad3e5b9 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++) {
7932f811 340 digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
fad3e5b9 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)
7932f811 348 fCoreEnergy += fEnergyList[iDigit] ;
fad3e5b9 349 }
350
fad3e5b9 351}
d15a28e7 352
353//____________________________________________________________________________
7932f811 354void AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
d15a28e7 355{
b2a60966 356 // Calculates the axis of the shower ellipsoid
83974468 357
e8dbb96e 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.;
d15a28e7 364
365 AliPHOSDigit * digit ;
92862013 366 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
d15a28e7 367 Int_t iDigit;
368
369 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
7932f811 370 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
d15a28e7 371 Int_t relid[4] ;
372 Float_t xi ;
373 Float_t zi ;
92862013 374 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
375 phosgeom->RelPosInModule(relid, xi, zi);
7932f811 376 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
92862013 377 dxx += w * xi * xi ;
d15a28e7 378 x += w * xi ;
92862013 379 dzz += w * zi * zi ;
d15a28e7 380 z += w * zi ;
92862013 381 dxz += w * xi * zi ;
d15a28e7 382 wtot += w ;
383 }
92862013 384 dxx /= wtot ;
d15a28e7 385 x /= wtot ;
92862013 386 dxx -= x * x ;
387 dzz /= wtot ;
d15a28e7 388 z /= wtot ;
92862013 389 dzz -= z * z ;
390 dxz /= wtot ;
391 dxz -= x * z ;
d15a28e7 392
69183710 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
7932f811 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]) ;
e8dbb96e 409
7932f811 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]) ;
e8dbb96e 413 else
7932f811 414 fLambda[1]= 0. ;
d15a28e7 415}
416
b2a60966 417//____________________________________________________________________________
7932f811 418void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
ad8cfaf4 419{
7932f811 420 AliPHOSRecPoint::EvalAll(logWeight,digits) ;
421 EvalLocalPosition(logWeight, digits) ;
422 EvalElipsAxis(logWeight, digits) ;
423 EvalDispersion(logWeight, digits) ;
424 EvalCoreEnergy(digits);
ad8cfaf4 425}
426//____________________________________________________________________________
7932f811 427void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
b2a60966 428{
429 // Calculates the center of gravity in the local PHOS-module coordinates
83974468 430
b2a60966 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
b2a60966 444 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
7932f811 445 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
b2a60966 446
447 Float_t xi ;
448 Float_t zi ;
449 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
450 phosgeom->RelPosInModule(relid, xi, zi);
7932f811 451 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
b2a60966 452 x += xi * w ;
453 z += zi * w ;
454 wtot += w ;
b2a60966 455
9ce1a8d9 456 }
69183710 457
458 x /= wtot ;
459 z /= wtot ;
ad8cfaf4 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
1b799736 465 Float_t xo,yo,zo ; //Coordinates of the origin
466 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
467
5830e1d9 468 Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
1b799736 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)+
5830e1d9 476 (phosgeom->GetIPtoCrystalSurface()-yoL)*(phosgeom->GetIPtoCrystalSurface()-yoL)+
477 (zo-z)*(zo-z));
1b799736 478
479 Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ;
480 Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
ad8cfaf4 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) ;
b2a60966 487 fLocPos.SetY(0.) ;
ad8cfaf4 488 fLocPos.SetZ(z - depthz) ;
b2a60966 489
b2a60966 490}
491
d15a28e7 492//____________________________________________________________________________
ad8cfaf4 493Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
d15a28e7 494{
b2a60966 495 // Finds the maximum energy in the cluster
496
d15a28e7 497 Float_t menergy = 0. ;
498
499 Int_t iDigit;
500
501 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
502
503 if(fEnergyList[iDigit] > menergy)
504 menergy = fEnergyList[iDigit] ;
505 }
506 return menergy ;
507}
508
509//____________________________________________________________________________
ad8cfaf4 510Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
d15a28e7 511{
b2a60966 512 // Calculates the multiplicity of digits with energy larger than H*energy
513
d15a28e7 514 Int_t multipl = 0 ;
515 Int_t iDigit ;
516 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
517
518 if(fEnergyList[iDigit] > H * fAmp)
519 multipl++ ;
520 }
521 return multipl ;
522}
523
524//____________________________________________________________________________
7932f811 525Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t * maxAt, Float_t * maxAtEnergy,
526 Float_t locMaxCut,TClonesArray * digits) const
d15a28e7 527{
b2a60966 528 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
529 // energy difference between two local maxima
530
d15a28e7 531 AliPHOSDigit * digit ;
532 AliPHOSDigit * digitN ;
533
534
535 Int_t iDigitN ;
536 Int_t iDigit ;
537
7932f811 538 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
539 maxAt[iDigit] = (Int_t) digits->At(fDigitsList[iDigit]) ;
540
d15a28e7 541
6ad0bfa0 542 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
9f616d61 543 if(maxAt[iDigit] != -1) {
544 digit = (AliPHOSDigit *) maxAt[iDigit] ;
83974468 545
6ad0bfa0 546 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
7932f811 547 digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ;
d15a28e7 548
9f616d61 549 if ( AreNeighbours(digit, digitN) ) {
d15a28e7 550 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
551 maxAt[iDigitN] = -1 ;
6ad0bfa0 552 // but may be digit too is not local max ?
7932f811 553 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
d15a28e7 554 maxAt[iDigit] = -1 ;
555 }
556 else {
557 maxAt[iDigit] = -1 ;
6ad0bfa0 558 // but may be digitN too is not local max ?
7932f811 559 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
d15a28e7 560 maxAt[iDigitN] = -1 ;
561 }
562 } // if Areneighbours
563 } // while digitN
564 } // slot not empty
565 } // while digit
566
567 iDigitN = 0 ;
6ad0bfa0 568 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
d15a28e7 569 if(maxAt[iDigit] != -1){
570 maxAt[iDigitN] = maxAt[iDigit] ;
9f616d61 571 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
572 iDigitN++ ;
d15a28e7 573 }
574 }
575 return iDigitN ;
576}
577
d15a28e7 578
d15a28e7 579
580//____________________________________________________________________________
581void AliPHOSEmcRecPoint::Print(Option_t * option)
582{
b2a60966 583 // Print the list of digits belonging to the cluster
584
d15a28e7 585 cout << "AliPHOSEmcRecPoint: " << endl ;
586
d15a28e7 587 Int_t iDigit;
7932f811 588 cout << " digits # = " ;
589 for(iDigit=0; iDigit<fMulDigit; iDigit++)
590 cout << fDigitsList[iDigit] << " " ;
591 cout << endl ;
592
593 cout << " Energies = " ;
594 for(iDigit=0; iDigit<fMulDigit; iDigit++)
595 cout << fEnergyList[iDigit] << " ";
596 cout << endl ;
83974468 597
d15a28e7 598 cout << " Multiplicity = " << fMulDigit << endl ;
599 cout << " Cluster Energy = " << fAmp << endl ;
83974468 600 cout << " Stored at position " << GetIndexInList() << endl ;
601
d15a28e7 602}
88714635 603
7932f811 604