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