]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSEmcRecPoint.cxx
New modes only.
[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
baef0810 21//*--
b2a60966 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"
55fe9d13 30#include "TGraph.h"
d15a28e7 31
32// --- Standard library ---
33
d15a28e7 34// --- AliRoot header files ---
35
a6eedfad 36 #include "AliGenerator.h"
d15a28e7 37#include "AliPHOSGeometry.h"
38#include "AliPHOSEmcRecPoint.h"
39#include "AliRun.h"
7b7c1533 40#include "AliPHOSGetter.h"
d15a28e7 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 ;
94256c23 53 fNExMax = 0 ; //Not unfolded yet
9688c1dd 54 fTime = -1. ;
a6eedfad 55 fLocPos.SetX(1000000.) ; //Local position should be evaluated
7b7c1533 56
83974468 57}
58
73a68ccb 59//____________________________________________________________________________
60AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const char * opt) : AliPHOSRecPoint(opt)
61{
62 // ctor
63
64 fMulDigit = 0 ;
65 fAmp = 0. ;
94256c23 66 fNExMax = 0 ; //Not unfolded yet
73a68ccb 67 fCoreEnergy = 0 ;
68 fEnergyList = 0 ;
69 fTime = -1. ;
70 fLocPos.SetX(1000000.) ; //Local position should be evaluated
71
72}
73
55fe9d13 74//____________________________________________________________________________
75AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const AliPHOSEmcRecPoint & rp) : AliPHOSRecPoint(rp)
76{
77 // cpy ctor
78
79 fMulDigit = rp.fMulDigit ;
80 fAmp = rp.fAmp ;
81 fCoreEnergy = rp.fCoreEnergy ;
82 fEnergyList = new Float_t[rp.fMulDigit] ;
83 Int_t index ;
84 for(index = 0 ; index < fMulDigit ; index++)
85 fEnergyList[index] = rp.fEnergyList[index] ;
86 fNExMax = rp.fNExMax ;
87 fTime = rp.fTime ;
88}
89
83974468 90//____________________________________________________________________________
91AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
92{
88714635 93 // dtor
a4e98857 94
83974468 95 if ( fEnergyList )
96 delete[] fEnergyList ;
d15a28e7 97}
98
d15a28e7 99//____________________________________________________________________________
83974468 100void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
d15a28e7 101{
b2a60966 102 // Adds a digit to the RecPoint
a4e98857 103 // and accumulates the total amplitude and the multiplicity
d15a28e7 104
7932f811 105 if(fEnergyList == 0)
106 fEnergyList = new Float_t[fMaxDigit];
107
d15a28e7 108 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
9f616d61 109 fMaxDigit*=2 ;
83974468 110 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
9f616d61 111 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
112
113 Int_t index ;
d15a28e7 114 for ( index = 0 ; index < fMulDigit ; index++ ){
83974468 115 tempo[index] = fDigitsList[index] ;
d15a28e7 116 tempoE[index] = fEnergyList[index] ;
117 }
118
9f616d61 119 delete [] fDigitsList ;
83974468 120 fDigitsList = new ( Int_t[fMaxDigit] ) ;
9f616d61 121
122 delete [] fEnergyList ;
123 fEnergyList = new ( Float_t[fMaxDigit] ) ;
124
125 for ( index = 0 ; index < fMulDigit ; index++ ){
126 fDigitsList[index] = tempo[index] ;
127 fEnergyList[index] = tempoE[index] ;
128 }
129
130 delete [] tempo ;
131 delete [] tempoE ;
132 } // if
d15a28e7 133
83974468 134 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
135 fEnergyList[fMulDigit] = Energy ;
136 fMulDigit++ ;
d15a28e7 137 fAmp += Energy ;
7932f811 138
139 EvalPHOSMod(&digit) ;
d15a28e7 140}
141
142//____________________________________________________________________________
ad8cfaf4 143Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
d15a28e7 144{
a4e98857 145 // Tells if (true) or not (false) two digits are neighbors
d15a28e7 146
147 Bool_t aren = kFALSE ;
148
7b7c1533 149 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
150 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
151
d15a28e7 152 Int_t relid1[4] ;
92862013 153 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
d15a28e7 154
155 Int_t relid2[4] ;
92862013 156 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
d15a28e7 157
92862013 158 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
159 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
d15a28e7 160
92862013 161 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
d15a28e7 162 aren = kTRUE ;
163
164 return aren ;
165}
166
167//____________________________________________________________________________
2a941f4e 168Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
d15a28e7 169{
b2a60966 170 // Compares two RecPoints according to their position in the PHOS modules
171
7932f811 172 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
173 //value (what is senseless) change as vell delta in
174 //AliPHOSTrackSegmentMakerv* and other RecPoints...
d15a28e7 175 Int_t rv ;
176
177 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
178
179
ad8cfaf4 180 Int_t phosmod1 = GetPHOSMod() ;
92862013 181 Int_t phosmod2 = clu->GetPHOSMod() ;
d15a28e7 182
92862013 183 TVector3 locpos1;
7932f811 184 GetLocalPosition(locpos1) ;
92862013 185 TVector3 locpos2;
186 clu->GetLocalPosition(locpos2) ;
d15a28e7 187
92862013 188 if(phosmod1 == phosmod2 ) {
7932f811 189 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
d15a28e7 190 if (rowdif> 0)
7932f811 191 rv = 1 ;
d15a28e7 192 else if(rowdif < 0)
7932f811 193 rv = -1 ;
92862013 194 else if(locpos1.Z()>locpos2.Z())
d15a28e7 195 rv = -1 ;
196 else
197 rv = 1 ;
198 }
199
200 else {
92862013 201 if(phosmod1 < phosmod2 )
d15a28e7 202 rv = -1 ;
203 else
204 rv = 1 ;
205 }
206
207 return rv ;
208}
9f616d61 209//______________________________________________________________________________
baef0810 210void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
9f616d61 211{
9f616d61 212
9688c1dd 213 // Execute action corresponding to one event
214 // This member function is called when a AliPHOSRecPoint is clicked with the locator
215 //
216 // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
217 // and switched off when the mouse button is released.
218
219
220 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
221 if(!gime) return ;
222 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
83974468 223
9688c1dd 224 static TGraph * digitgraph = 0 ;
83974468 225
9688c1dd 226 if (!gPad->IsEditable()) return;
83974468 227
9688c1dd 228 TH2F * histo = 0 ;
229 TCanvas * histocanvas ;
230
afa51c4e 231 const TClonesArray * digits = gime->Digits() ;
9688c1dd 232
233 switch (event) {
83974468 234
9688c1dd 235 case kButton1Down: {
236 AliPHOSDigit * digit ;
237 Int_t iDigit;
238 Int_t relid[4] ;
83974468 239
9688c1dd 240 const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
241 Float_t * xi = new Float_t[kMulDigit] ;
242 Float_t * zi = new Float_t[kMulDigit] ;
83974468 243
9688c1dd 244 // create the histogram for the single cluster
245 // 1. gets histogram boundaries
246 Float_t ximax = -999. ;
247 Float_t zimax = -999. ;
248 Float_t ximin = 999. ;
249 Float_t zimin = 999. ;
83974468 250
9688c1dd 251 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
252 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
253 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
254 phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
255 if ( xi[iDigit] > ximax )
256 ximax = xi[iDigit] ;
257 if ( xi[iDigit] < ximin )
258 ximin = xi[iDigit] ;
259 if ( zi[iDigit] > zimax )
260 zimax = zi[iDigit] ;
261 if ( zi[iDigit] < zimin )
262 zimin = zi[iDigit] ;
263 }
264 ximax += phosgeom->GetCrystalSize(0) / 2. ;
265 zimax += phosgeom->GetCrystalSize(2) / 2. ;
266 ximin -= phosgeom->GetCrystalSize(0) / 2. ;
267 zimin -= phosgeom->GetCrystalSize(2) / 2. ;
268 Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
269 Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
83974468 270
9688c1dd 271 // 2. gets the histogram title
83974468 272
9688c1dd 273 Text_t title[100] ;
274 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
83974468 275
9688c1dd 276 if (!histo) {
277 delete histo ;
278 histo = 0 ;
279 }
280 histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
83974468 281
9688c1dd 282 Float_t x, z ;
283 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
284 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
285 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
286 phosgeom->RelPosInModule(relid, x, z);
287 histo->Fill(x, z, fEnergyList[iDigit] ) ;
288 }
83974468 289
9688c1dd 290 if (!digitgraph) {
291 digitgraph = new TGraph(kMulDigit,xi,zi);
292 digitgraph-> SetMarkerStyle(5) ;
293 digitgraph-> SetMarkerSize(1.) ;
294 digitgraph-> SetMarkerColor(1) ;
295 digitgraph-> Paint("P") ;
296 }
83974468 297
9688c1dd 298 // Print() ;
299 histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
300 histocanvas->Draw() ;
301 histo->Draw("lego1") ;
83974468 302
9688c1dd 303 delete[] xi ;
304 delete[] zi ;
83974468 305
9688c1dd 306 break;
307 }
83974468 308
9688c1dd 309 case kButton1Up:
310 if (digitgraph) {
311 delete digitgraph ;
312 digitgraph = 0 ;
313 }
314 break;
9f616d61 315
9688c1dd 316 }
9f616d61 317}
318
d15a28e7 319//____________________________________________________________________________
7932f811 320void AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
d15a28e7 321{
b2a60966 322 // Calculates the dispersion of the shower at the origine of the RecPoint
323
e5b16749 324 Float_t d = 0. ;
325 Float_t wtot = 0. ;
d15a28e7 326
d084d50d 327 Float_t x = 0.;
328 Float_t z = 0.;
d15a28e7 329
330 AliPHOSDigit * digit ;
7b7c1533 331
332 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
333 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
d15a28e7 334
e5b16749 335
336 // Calculates the center of gravity in the local PHOS-module coordinates
337
d15a28e7 338 Int_t iDigit;
e5b16749 339
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 w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
348 x += xi * w ;
349 z += zi * w ;
350 wtot += w ;
351 }
352 x /= wtot ;
353 z /= wtot ;
354
355
356// Calculates the dispersion in coordinates
357 wtot = 0.;
88714635 358 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
7932f811 359 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
d15a28e7 360 Int_t relid[4] ;
361 Float_t xi ;
362 Float_t zi ;
92862013 363 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
364 phosgeom->RelPosInModule(relid, xi, zi);
7932f811 365 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
92862013 366 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
d15a28e7 367 wtot+=w ;
d084d50d 368
369
d15a28e7 370 }
e5b16749 371
7932f811 372
92862013 373 d /= wtot ;
d15a28e7 374
7932f811 375 fDispersion = TMath::Sqrt(d) ;
c63c49e9 376
d15a28e7 377}
fad3e5b9 378//______________________________________________________________________________
e5b16749 379void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
fad3e5b9 380{
a4e98857 381 // This function calculates energy in the core,
382 // i.e. within a radius rad = 3cm around the center. Beyond this radius
383 // in accordance with shower profile the energy deposition
fad3e5b9 384 // should be less than 2%
385
fad3e5b9 386 Float_t coreRadius = 3 ;
387
e5b16749 388 Float_t x = 0 ;
389 Float_t z = 0 ;
fad3e5b9 390
391 AliPHOSDigit * digit ;
7b7c1533 392
393 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
394 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
395
fad3e5b9 396 Int_t iDigit;
e5b16749 397
398// Calculates the center of gravity in the local PHOS-module coordinates
399 Float_t wtot = 0;
400 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
401 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
402 Int_t relid[4] ;
403 Float_t xi ;
404 Float_t zi ;
405 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
406 phosgeom->RelPosInModule(relid, xi, zi);
407 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
408 x += xi * w ;
409 z += zi * w ;
410 wtot += w ;
411 }
412 x /= wtot ;
413 z /= wtot ;
414
415
fad3e5b9 416 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
7932f811 417 digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
fad3e5b9 418 Int_t relid[4] ;
419 Float_t xi ;
420 Float_t zi ;
421 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
422 phosgeom->RelPosInModule(relid, xi, zi);
423 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
424 if(distance < coreRadius)
7932f811 425 fCoreEnergy += fEnergyList[iDigit] ;
fad3e5b9 426 }
427
fad3e5b9 428}
d15a28e7 429
430//____________________________________________________________________________
7932f811 431void AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
d15a28e7 432{
b2a60966 433 // Calculates the axis of the shower ellipsoid
83974468 434
e8dbb96e 435 Double_t wtot = 0. ;
436 Double_t x = 0.;
437 Double_t z = 0.;
438 Double_t dxx = 0.;
439 Double_t dzz = 0.;
440 Double_t dxz = 0.;
d15a28e7 441
442 AliPHOSDigit * digit ;
7b7c1533 443
444 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
445 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
446
d15a28e7 447 Int_t iDigit;
448
e5b16749 449
d15a28e7 450 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
7932f811 451 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
d15a28e7 452 Int_t relid[4] ;
453 Float_t xi ;
454 Float_t zi ;
92862013 455 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
456 phosgeom->RelPosInModule(relid, xi, zi);
7932f811 457 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
92862013 458 dxx += w * xi * xi ;
d15a28e7 459 x += w * xi ;
92862013 460 dzz += w * zi * zi ;
d15a28e7 461 z += w * zi ;
92862013 462 dxz += w * xi * zi ;
d15a28e7 463 wtot += w ;
464 }
92862013 465 dxx /= wtot ;
d15a28e7 466 x /= wtot ;
92862013 467 dxx -= x * x ;
468 dzz /= wtot ;
d15a28e7 469 z /= wtot ;
92862013 470 dzz -= z * z ;
471 dxz /= wtot ;
472 dxz -= x * z ;
d15a28e7 473
69183710 474// //Apply correction due to non-perpendicular incidence
475// Double_t CosX ;
476// Double_t CosZ ;
7b7c1533 477// AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
478// AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
479 // Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
69183710 480
481// CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
482// CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
483
484// dxx = dxx/(CosX*CosX) ;
485// dzz = dzz/(CosZ*CosZ) ;
486// dxz = dxz/(CosX*CosZ) ;
487
488
7932f811 489 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
490 if(fLambda[0] > 0)
491 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
e8dbb96e 492
7932f811 493 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
494 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
495 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
e8dbb96e 496 else
7932f811 497 fLambda[1]= 0. ;
d15a28e7 498}
499
b2a60966 500//____________________________________________________________________________
7932f811 501void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
ad8cfaf4 502{
baef0810 503 // Evaluates all shower parameters
504
7932f811 505 EvalLocalPosition(logWeight, digits) ;
506 EvalElipsAxis(logWeight, digits) ;
507 EvalDispersion(logWeight, digits) ;
e5b16749 508 EvalCoreEnergy(logWeight, digits);
a6eedfad 509 EvalTime(digits) ;
bb53c80c 510 AliPHOSRecPoint::EvalAll(logWeight,digits) ;
ad8cfaf4 511}
512//____________________________________________________________________________
7932f811 513void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
b2a60966 514{
515 // Calculates the center of gravity in the local PHOS-module coordinates
b2a60966 516 Float_t wtot = 0. ;
517
518 Int_t relid[4] ;
519
520 Float_t x = 0. ;
521 Float_t z = 0. ;
522
523 AliPHOSDigit * digit ;
524
7b7c1533 525 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
526 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
b2a60966 527
528 Int_t iDigit;
529
b2a60966 530 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
7932f811 531 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
b2a60966 532
533 Float_t xi ;
534 Float_t zi ;
535 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
536 phosgeom->RelPosInModule(relid, xi, zi);
7932f811 537 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
b2a60966 538 x += xi * w ;
539 z += zi * w ;
540 wtot += w ;
b2a60966 541
9ce1a8d9 542 }
69183710 543
544 x /= wtot ;
545 z /= wtot ;
ad8cfaf4 546
547 // Correction for the depth of the shower starting point (TDR p 127)
548 Float_t para = 0.925 ;
549 Float_t parb = 6.52 ;
550
1b799736 551 Float_t xo,yo,zo ; //Coordinates of the origin
552 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
553
5830e1d9 554 Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
1b799736 555
556 //Transform to the local ref.frame
557 Float_t xoL,yoL ;
558 xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
559 yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
560
fd84cb73 561 Float_t radius = phosgeom->GetIPtoCrystalSurface()-yoL;
1b799736 562
563 Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ;
564 Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
ad8cfaf4 565
566 Float_t depthx = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ;
567 Float_t depthz = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ;
ad8cfaf4 568
569 fLocPos.SetX(x - depthx) ;
b2a60966 570 fLocPos.SetY(0.) ;
ad8cfaf4 571 fLocPos.SetZ(z - depthz) ;
b2a60966 572
a6eedfad 573 fLocPosM = 0 ;
b2a60966 574}
575
d15a28e7 576//____________________________________________________________________________
ad8cfaf4 577Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
d15a28e7 578{
b2a60966 579 // Finds the maximum energy in the cluster
580
d15a28e7 581 Float_t menergy = 0. ;
582
583 Int_t iDigit;
584
585 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
586
587 if(fEnergyList[iDigit] > menergy)
588 menergy = fEnergyList[iDigit] ;
589 }
590 return menergy ;
591}
592
593//____________________________________________________________________________
ad8cfaf4 594Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
d15a28e7 595{
b2a60966 596 // Calculates the multiplicity of digits with energy larger than H*energy
597
d15a28e7 598 Int_t multipl = 0 ;
599 Int_t iDigit ;
600 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
601
602 if(fEnergyList[iDigit] > H * fAmp)
603 multipl++ ;
604 }
605 return multipl ;
606}
607
608//____________________________________________________________________________
a0636361 609Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax( AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
7932f811 610 Float_t locMaxCut,TClonesArray * digits) const
d15a28e7 611{
b2a60966 612 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
a4e98857 613 // energy difference between two local maxima
b2a60966 614
d15a28e7 615 AliPHOSDigit * digit ;
616 AliPHOSDigit * digitN ;
617
618
619 Int_t iDigitN ;
620 Int_t iDigit ;
621
7932f811 622 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
a0636361 623 maxAt[iDigit] = (AliPHOSDigit*) digits->At(fDigitsList[iDigit]) ;
7932f811 624
d15a28e7 625
6ad0bfa0 626 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
a0636361 627 if(maxAt[iDigit]) {
628 digit = maxAt[iDigit] ;
83974468 629
6ad0bfa0 630 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
7932f811 631 digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ;
d15a28e7 632
9f616d61 633 if ( AreNeighbours(digit, digitN) ) {
d15a28e7 634 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
a0636361 635 maxAt[iDigitN] = 0 ;
6ad0bfa0 636 // but may be digit too is not local max ?
7932f811 637 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
a0636361 638 maxAt[iDigit] = 0 ;
d15a28e7 639 }
640 else {
a0636361 641 maxAt[iDigit] = 0 ;
6ad0bfa0 642 // but may be digitN too is not local max ?
7932f811 643 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
a0636361 644 maxAt[iDigitN] = 0 ;
d15a28e7 645 }
646 } // if Areneighbours
647 } // while digitN
648 } // slot not empty
649 } // while digit
650
651 iDigitN = 0 ;
6ad0bfa0 652 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
a0636361 653 if(maxAt[iDigit]){
d15a28e7 654 maxAt[iDigitN] = maxAt[iDigit] ;
9f616d61 655 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
656 iDigitN++ ;
d15a28e7 657 }
658 }
659 return iDigitN ;
660}
9688c1dd 661//____________________________________________________________________________
0bc3b8ed 662void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits)
663{
664 // Define a rec.point time as a time in the cell with the maximum energy
665
9688c1dd 666 Float_t maxE = 0;
667 Int_t maxAt = 0;
668 for(Int_t idig=0; idig < fMulDigit; idig++){
669 if(fEnergyList[idig] > maxE){
670 maxE = fEnergyList[idig] ;
671 maxAt = idig;
672 }
673 }
674 fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
675
676}
d15a28e7 677//____________________________________________________________________________
092b50ba 678void AliPHOSEmcRecPoint::Purify(Float_t threshold){
679 //Removes digits below threshold
680
681 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
682 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
683
684 Int_t mult = 0 ;
685 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
686 if(fEnergyList[iDigit] > threshold){
687 tempo[mult] = fDigitsList[iDigit] ;
688 tempoE[mult] = fEnergyList[iDigit] ;
689 mult++ ;
690 }
691 }
692
693 fMulDigit = mult ;
694 delete [] fDigitsList ;
695 delete [] fEnergyList ;
696 fDigitsList = new (Int_t[fMulDigit]) ;
697 fEnergyList = new ( Float_t[fMulDigit]) ;
698
699 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
700 fDigitsList[iDigit] = tempo[iDigit];
701 fEnergyList[iDigit] = tempoE[iDigit] ;
702 }
703
704 delete [] tempo ;
705 delete [] tempoE ;
706
707}
708//____________________________________________________________________________
709void AliPHOSEmcRecPoint::Print(Option_t * option) const
d15a28e7 710{
b2a60966 711 // Print the list of digits belonging to the cluster
712
21cd0c07 713 TString message ;
714 message = "AliPHOSEmcRecPoint:\n" ;
715 message += " digits # = " ;
716 Info("Print", message.Data()) ;
d15a28e7 717
d15a28e7 718 Int_t iDigit;
7932f811 719 for(iDigit=0; iDigit<fMulDigit; iDigit++)
53ad38d8 720 printf(" %d ", fDigitsList[iDigit] ) ;
7932f811 721
21cd0c07 722 Info("Print", " Energies = ") ;
7932f811 723 for(iDigit=0; iDigit<fMulDigit; iDigit++)
53ad38d8 724 printf(" %f ", fEnergyList[iDigit] ) ;
725 printf("\n") ;
21cd0c07 726 Info("Print", " Primaries ") ;
bc68d12c 727 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
53ad38d8 728 printf(" %d ", fTracksList[iDigit]) ;
729 printf("\n") ;
21cd0c07 730 message = " Multiplicity = %d" ;
731 message += " Cluster Energy = %f" ;
732 message += " Number of primaries %d" ;
733 message += " Stored at position %d" ;
83974468 734
21cd0c07 735 Info("Print", message.Data(), fMulDigit, fAmp, fMulTrack,GetIndexInList() ) ;
d15a28e7 736}
88714635 737
7932f811 738