]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EMCAL/AliEMCALRecPoint.cxx
Using the new (TGeo) geometry of TPC
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALRecPoint.cxx
CommitLineData
ab48128d 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/* $Id$ */
16//_________________________________________________________________________
70a93198 17// Reconstructed Points for the EMCAL
18// A RecPoint is a cluster of digits
d64c959b 19//*-- Author: Yves Schutz (SUBATECH)
70a93198 20//*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
21//*-- Author: Heather Gray (LBL) merged AliEMCALRecPoint and AliEMCALTowerRecPoint 02/04
ab48128d 22
23// --- ROOT system ---
24#include "TPad.h"
d64c959b 25#include "TGraph.h"
26#include "TPaveText.h"
ab48128d 27#include "TClonesArray.h"
ff1e7e2f 28#include "TMath.h"
ab48128d 29
30// --- Standard library ---
ab48128d 31
32// --- AliRoot header files ---
70a93198 33#include "AliGenerator.h"
ab48128d 34#include "AliEMCALGeometry.h"
35#include "AliEMCALDigit.h"
36#include "AliEMCALRecPoint.h"
ab48128d 37
38ClassImp(AliEMCALRecPoint)
39
ab48128d 40//____________________________________________________________________________
41AliEMCALRecPoint::AliEMCALRecPoint()
42 : AliRecPoint()
43{
44 // ctor
692088ae 45 fMaxTrack = 0 ;
70a93198 46 fMulDigit = 0 ;
87cdc3be 47 fMaxParent = 0;
48 fMulParent = 0;
70a93198 49 fAmp = 0. ;
50 fCoreEnergy = 0 ;
51 fEnergyList = 0 ;
87cdc3be 52 fParentsList = 0;
70a93198 53 fTime = 0. ;
54 fLocPos.SetX(0.) ; //Local position should be evaluated
55 fCoreRadius = 10; //HG Check this
ab48128d 56}
57
58//____________________________________________________________________________
59AliEMCALRecPoint::AliEMCALRecPoint(const char * opt) : AliRecPoint(opt)
60{
61 // ctor
ff1e7e2f 62 fMaxTrack = 1000 ;
63 fMaxParent = 1000;
87cdc3be 64 fMulDigit = 0 ;
65 fMulParent = 0;
70a93198 66 fAmp = 0. ;
67 fCoreEnergy = 0 ;
68 fEnergyList = 0 ;
87cdc3be 69 fParentsList = new Int_t[fMaxParent];
70a93198 70 fTime = -1. ;
71 fLocPos.SetX(1000000.) ; //Local position should be evaluated
72 fCoreRadius = 10; //HG Check this
73}
74//____________________________________________________________________________
75AliEMCALRecPoint::~AliEMCALRecPoint()
76{
77 // dtor
78 if ( fEnergyList )
79 delete[] fEnergyList ;
87cdc3be 80 if ( fParentsList)
81 delete[] fParentsList;
70a93198 82}
83
84//____________________________________________________________________________
85void AliEMCALRecPoint::AddDigit(AliEMCALDigit & digit, Float_t Energy)
86{
87 // Adds a digit to the RecPoint
88 // and accumulates the total amplitude and the multiplicity
89
90 if(fEnergyList == 0)
91 fEnergyList = new Float_t[fMaxDigit];
92
93 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
94 fMaxDigit*=2 ;
a64a06d6 95 Int_t * tempo = new Int_t[fMaxDigit];
96 Float_t * tempoE = new Float_t[fMaxDigit];
70a93198 97
98 Int_t index ;
99 for ( index = 0 ; index < fMulDigit ; index++ ){
100 tempo[index] = fDigitsList[index] ;
101 tempoE[index] = fEnergyList[index] ;
102 }
103
104 delete [] fDigitsList ;
a64a06d6 105 fDigitsList = new Int_t[fMaxDigit];
70a93198 106
107 delete [] fEnergyList ;
a64a06d6 108 fEnergyList = new Float_t[fMaxDigit];
70a93198 109
110 for ( index = 0 ; index < fMulDigit ; index++ ){
111 fDigitsList[index] = tempo[index] ;
112 fEnergyList[index] = tempoE[index] ;
113 }
114
115 delete [] tempo ;
116 delete [] tempoE ;
117 } // if
118
119 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
120 fEnergyList[fMulDigit] = Energy ;
121 fMulDigit++ ;
122 fAmp += Energy ;
123
124}
125//____________________________________________________________________________
126Bool_t AliEMCALRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
127{
128 // Tells if (true) or not (false) two digits are neighbours
129 // A neighbour is defined as being two digits which share a corner
130
131 Bool_t areNeighbours = kFALSE ;
132
5dee926e 133 //AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
134 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
70a93198 135
136 Int_t relid1[2] ;
1963b290 137 //copied for shish-kebab geometry, ieta,iphi is cast as float as eta,phi conversion
138 // for this geometry does not exist
139 int nSupMod=0, nTower=0, nIphi=0, nIeta=0;
140 int iphi=0, ieta=0;
141 geom->GetCellIndex(digit1->GetId(), nSupMod,nTower,nIphi,nIeta);
d87bd045 142 geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
1963b290 143 relid1[0]=ieta;
144 relid1[1]=iphi;
145// geom->AbsToRelNumbering(digit1->GetId(), relid1) ;
70a93198 146
147 Int_t relid2[2] ;
1963b290 148 //copied for shish-kebab geometry, ieta,iphi is cast as float as eta,phi conversion
149 // for this geometry does not exist
150 int nSupMod1=0, nTower1=0, nIphi1=0, nIeta1=0;
151 int iphi1=0, ieta1=0;
152 geom->GetCellIndex(digit2->GetId(), nSupMod1,nTower1,nIphi1,nIeta1);
d87bd045 153 geom->GetCellPhiEtaIndexInSModule(nSupMod1,nTower1,nIphi1,nIeta1, iphi1,ieta1);
1963b290 154 relid2[0]=ieta1;
155 relid2[1]=iphi1;
156// geom->AbsToRelNumbering(digit2->GetId(), relid2) ;
70a93198 157
158 Int_t rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;
159 Int_t coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;
160
161 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
162 areNeighbours = kTRUE ;
ab48128d 163
70a93198 164 return areNeighbours;
165}
166
167//____________________________________________________________________________
168Int_t AliEMCALRecPoint::Compare(const TObject * obj) const
169{
170 // Compares two RecPoints according to their position in the EMCAL modules
171
172 Float_t delta = 1 ; //Width of "Sorting row". If you change this
173 //value (what is senseless) change as well delta in
174 //AliEMCALTrackSegmentMakerv* and other RecPoints...
175 Int_t rv ;
176
177 AliEMCALRecPoint * clu = (AliEMCALRecPoint *)obj ;
178
179 TVector3 locpos1;
180 GetLocalPosition(locpos1);
181 TVector3 locpos2;
182 clu->GetLocalPosition(locpos2);
183
184 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
185 if (rowdif> 0)
186 rv = 1 ;
187 else if(rowdif < 0)
188 rv = -1 ;
189 else if(locpos1.Y()>locpos2.Y())
190 rv = -1 ;
191 else
192 rv = 1 ;
193
194 return rv ;
ab48128d 195}
196
197//____________________________________________________________________________
198Int_t AliEMCALRecPoint::DistancetoPrimitive(Int_t px, Int_t py)
199{
200 // Compute distance from point px,py to a AliEMCALRecPoint considered as a Tmarker
201 // Compute the closest distance of approach from point px,py to this marker.
202 // The distance is computed in pixels units.
70a93198 203 // HG Still need to update -> Not sure what this should achieve
ab48128d 204
205 TVector3 pos(0.,0.,0.) ;
70a93198 206 GetLocalPosition(pos) ;
ab48128d 207 Float_t x = pos.X() ;
70a93198 208 Float_t y = pos.Y() ;
ab48128d 209 const Int_t kMaxDiff = 10;
210 Int_t pxm = gPad->XtoAbsPixel(x);
211 Int_t pym = gPad->YtoAbsPixel(y);
212 Int_t dist = (px-pxm)*(px-pxm) + (py-pym)*(py-pym);
213
214 if (dist > kMaxDiff) return 9999;
215 return dist;
216}
217
218//___________________________________________________________________________
219 void AliEMCALRecPoint::Draw(Option_t *option)
220 {
221 // Draw this AliEMCALRecPoint with its current attributes
222
223 AppendPad(option);
224 }
225
226//______________________________________________________________________________
70a93198 227void AliEMCALRecPoint::ExecuteEvent(Int_t /*event*/, Int_t, Int_t)
ab48128d 228{
229 // Execute action corresponding to one event
230 // This member function is called when a AliEMCALRecPoint is clicked with the locator
231 //
232 // If Left button is clicked on AliEMCALRecPoint, the digits are switched on
233 // and switched off when the mouse button is released.
234
235 // static Int_t pxold, pyold;
236
70a93198 237 /* static TGraph * digitgraph = 0 ;
ab48128d 238 static TPaveText* clustertext = 0 ;
239
240 if (!gPad->IsEditable()) return;
241
242 switch (event) {
243
244
245 case kButton1Down:{
246 AliEMCALDigit * digit ;
88cb7938 247 AliEMCALGeometry * emcalgeom = (AliEMCALGetter::Instance())->EMCALGeometry() ;
ab48128d 248
249 Int_t iDigit;
70a93198 250 Int_t relid[2] ;
ab48128d 251
252 const Int_t kMulDigit=AliEMCALRecPoint::GetDigitsMultiplicity() ;
253 Float_t * xi = new Float_t [kMulDigit] ;
254 Float_t * zi = new Float_t [kMulDigit] ;
255
256 for(iDigit = 0; iDigit < kMulDigit; iDigit++) {
257 Fatal("AliEMCALRecPoint::ExecuteEvent", " -> Something wrong with the code");
258 digit = 0 ; //dynamic_cast<AliEMCALDigit *>((fDigitsList)[iDigit]);
259 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
260 emcalgeom->PosInAlice(relid, xi[iDigit], zi[iDigit]) ;
261 }
262
263 if (!digitgraph) {
264 digitgraph = new TGraph(fMulDigit,xi,zi);
265 digitgraph-> SetMarkerStyle(5) ;
266 digitgraph-> SetMarkerSize(1.) ;
267 digitgraph-> SetMarkerColor(1) ;
268 digitgraph-> Draw("P") ;
269 }
270 if (!clustertext) {
271
272 TVector3 pos(0.,0.,0.) ;
273 GetLocalPosition(pos) ;
274 clustertext = new TPaveText(pos.X()-10,pos.Z()+10,pos.X()+50,pos.Z()+35,"") ;
275 Text_t line1[40] ;
276 Text_t line2[40] ;
277 sprintf(line1,"Energy=%1.2f GeV",GetEnergy()) ;
278 sprintf(line2,"%d Digits",GetDigitsMultiplicity()) ;
279 clustertext ->AddText(line1) ;
280 clustertext ->AddText(line2) ;
281 clustertext ->Draw("");
282 }
283 gPad->Update() ;
9e5d2067 284 Print("") ;
ab48128d 285 delete[] xi ;
286 delete[] zi ;
287 }
288
289break;
290
291 case kButton1Up:
292 if (digitgraph) {
293 delete digitgraph ;
294 digitgraph = 0 ;
295 }
296 if (clustertext) {
297 delete clustertext ;
298 clustertext = 0 ;
299 }
300
301 break;
302
70a93198 303 }*/
304}
305//____________________________________________________________________________
306void AliEMCALRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits)
307{
308 // Evaluates all shower parameters
309
310 EvalLocalPosition(logWeight, digits) ;
500aeccc 311 // printf("eval position done\n");
70a93198 312 EvalElipsAxis(logWeight, digits) ;
500aeccc 313 // printf("eval axis done\n");
70a93198 314 EvalDispersion(logWeight, digits) ;
500aeccc 315 // printf("eval dispersion done\n");
1963b290 316 // EvalCoreEnergy(logWeight, digits);
317 // printf("eval energy done\n");
70a93198 318 EvalTime(digits) ;
500aeccc 319 // printf("eval time done\n");
70a93198 320
87cdc3be 321 EvalPrimaries(digits) ;
500aeccc 322 // printf("eval pri done\n");
87cdc3be 323 EvalParents(digits);
500aeccc 324 // printf("eval parent done\n");
70a93198 325}
326
327//____________________________________________________________________________
328void AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits)
329{
330 // Calculates the dispersion of the shower at the origin of the RecPoint
331
332 Float_t d = 0. ;
333 Float_t wtot = 0. ;
334
335 AliEMCALDigit * digit ;
336
5dee926e 337 //AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
338 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
339
70a93198 340 // Calculates the centre of gravity in the local EMCAL-module coordinates
341 Int_t iDigit;
342
343 if (!fLocPos.X() || !fLocPos.Y())
344 EvalLocalPosition(logWeight, digits) ;
345
1963b290 346 //Sub const Float_t kDeg2Rad = TMath::DegToRad() ;
70a93198 347
348 Float_t cluEta = fLocPos.X() ;
349 Float_t cluPhi = fLocPos.Y() ;
350 Float_t cluR = fLocPos.Z() ;
351
352 if (gDebug == 2)
353 printf("EvalDispersion: eta,phi,r = %f,%f,%f", cluEta, cluPhi, cluR) ;
354
355 // Calculates the dispersion in coordinates
356 wtot = 0.;
357 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
358 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
359 Float_t etai = 0.;
360 Float_t phii = 0.;
1963b290 361 //copied for shish-kebab geometry, ieta,iphi is cast as float as eta,phi conversion
362 // for this geometry does not exist
363 int nSupMod=0, nTower=0, nIphi=0, nIeta=0;
364 int iphi=0, ieta=0;
365 geom->GetCellIndex(digit->GetId(), nSupMod,nTower,nIphi,nIeta);
d87bd045 366 geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
1963b290 367 etai=(Float_t)ieta;
368 phii=(Float_t)iphi;
500aeccc 369 // printf("%f,%d,%d \n", fAmp, ieta, iphi) ;
1963b290 370
371/* Sub
372 geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
70a93198 373 phii = phii * kDeg2Rad;
1963b290 374*/
375///////////////////////////
500aeccc 376// if(fAmp>0)printf("%f %d %f", fAmp,iDigit,fEnergyList[iDigit]) ;
1963b290 377/////////////////////////
378
70a93198 379 if (gDebug == 2)
380 printf("EvalDispersion: id = %d, etai,phii = %f,%f", digit->GetId(), etai, phii) ;
381
382 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
383 d += w * ( (etai-cluEta)*(etai-cluEta) + (phii-cluPhi)*(phii-cluPhi) ) ;
384 wtot+=w ;
ab48128d 385 }
70a93198 386
387 if ( wtot > 0 )
388 d /= wtot ;
389 else
390 d = 0. ;
391
392 fDispersion = TMath::Sqrt(d) ;
500aeccc 393 // printf("Dispersion: = %f", fDispersion) ;
ab48128d 394}
70a93198 395
ab48128d 396//____________________________________________________________________________
70a93198 397void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
88cb7938 398{
70a93198 399 // Calculates the center of gravity in the local EMCAL-module coordinates
400 Float_t wtot = 0. ;
401
402 // Int_t relid[3] ;
403
404 AliEMCALDigit * digit ;
5dee926e 405 //AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
406 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
70a93198 407 Int_t iDigit;
408 Float_t cluEta = 0;
409 Float_t cluPhi = 0;
1963b290 410 //Sub const Float_t kDeg2Rad = TMath::DegToRad();
70a93198 411
412 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
413 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
414
415 Float_t etai ;
416 Float_t phii ;
1963b290 417 //copied for shish-kebab geometry, ieta,iphi is cast as float as eta,phi conversion
418 // for this geometry does not exist
419 int nSupMod=0, nTower=0, nIphi=0, nIeta=0;
420 int iphi=0, ieta=0;
421 geom->GetCellIndex(digit->GetId(), nSupMod,nTower,nIphi,nIeta);
d87bd045 422 geom->GetCellPhiEtaIndexInSModule(nSupMod, nTower,nIphi,nIeta, iphi,ieta); //19-oct-05
1963b290 423 etai=(Float_t)ieta;
424 phii=(Float_t)iphi;
425//Sub geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
426//Sub phii = phii * kDeg2Rad;
70a93198 427 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
428 cluEta += (etai * w) ;
429 cluPhi += (phii * w );
430 wtot += w ;
431 }
ab48128d 432
70a93198 433 if ( wtot > 0 ) {
434 cluEta /= wtot ;
435 cluPhi /= wtot ;
436 } else {
437 cluEta = -1 ;
438 cluPhi = -1.;
439 }
440
441 fLocPos.SetX(cluEta);
442 fLocPos.SetY(cluPhi);
443 fLocPos.SetZ(geom->GetIP2ECASection());
444
1963b290 445// if (gDebug==2)
500aeccc 446// printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ;
70a93198 447 fLocPosM = 0 ;
ab48128d 448}
449
70a93198 450//______________________________________________________________________________
451void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
452{
453 // This function calculates energy in the core,
454 // i.e. within a radius rad = 3cm around the center. Beyond this radius
455 // in accordance with shower profile the energy deposition
456 // should be less than 2%
457
458 AliEMCALDigit * digit ;
459 const Float_t kDeg2Rad = TMath::DegToRad() ;
5dee926e 460 //AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
461 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
462
70a93198 463 Int_t iDigit;
464
465 if (!fLocPos.X() || !fLocPos.Y() ) {
466 EvalLocalPosition(logWeight, digits);
467 }
468
469 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
470 digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
471 Float_t etai = 0. ;
472 Float_t phii = 0. ;
473 geom->PosInAlice(digit->GetId(), etai, phii);
474 phii = phii * kDeg2Rad;
475
476 Float_t distance = TMath::Sqrt((etai-fLocPos.X())*(etai-fLocPos.X())+(phii-fLocPos.Y())*(phii-fLocPos.Y())) ;
477 if(distance < fCoreRadius)
478 fCoreEnergy += fEnergyList[iDigit] ;
479 }
480
481}
ab48128d 482//____________________________________________________________________________
70a93198 483void AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
ab48128d 484{
70a93198 485 // Calculates the axis of the shower ellipsoid in eta and phi
ab48128d 486
70a93198 487 Double_t wtot = 0. ;
488 Double_t x = 0.;
489 Double_t z = 0.;
490 Double_t dxx = 0.;
491 Double_t dzz = 0.;
492 Double_t dxz = 0.;
ab48128d 493
ff1e7e2f 494 const Float_t kDeg2Rad = TMath::DegToRad();
70a93198 495 AliEMCALDigit * digit ;
496
5dee926e 497 //AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
498 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
1963b290 499 TString gn(geom->GetName());
70a93198 500
501 Int_t iDigit;
ab48128d 502
70a93198 503 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
504 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
505 Float_t etai = 0. ;
506 Float_t phii = 0. ;
1963b290 507 if(gn.Contains("SHISH")) {
508 //copied for shish-kebab geometry, ieta,iphi is cast as float as eta,phi conversion
509 // for this geometry does not exist
510 int nSupMod=0, nTower=0, nIphi=0, nIeta=0;
511 int iphi=0, ieta=0;
512 geom->GetCellIndex(digit->GetId(), nSupMod,nTower,nIphi,nIeta);
d87bd045 513 geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
1963b290 514 etai=(Float_t)ieta;
515 phii=(Float_t)iphi;
516 } else {
517 geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
518 phii = phii * kDeg2Rad;
519 }
520
70a93198 521 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
ff1e7e2f 522
70a93198 523 dxx += w * etai * etai ;
524 x += w * etai ;
525 dzz += w * phii * phii ;
526 z += w * phii ;
1963b290 527
ff1e7e2f 528 dxz += w * etai * phii ;
1963b290 529
70a93198 530 wtot += w ;
531 }
ff1e7e2f 532
70a93198 533 if ( wtot > 0 ) {
534 dxx /= wtot ;
535 x /= wtot ;
536 dxx -= x * x ;
537 dzz /= wtot ;
538 z /= wtot ;
539 dzz -= z * z ;
540 dxz /= wtot ;
541 dxz -= x * z ;
ab48128d 542
70a93198 543 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
544 if(fLambda[0] > 0)
545 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
546 else
547 fLambda[0] = 0;
548
549 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
ff1e7e2f 550
70a93198 551 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
552 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
553 else
554 fLambda[1]= 0. ;
555 } else {
556 fLambda[0]= 0. ;
557 fLambda[1]= 0. ;
ab48128d 558 }
ff1e7e2f 559
1963b290 560 // printf("Evalaxis: lambdas = %f,%f", fLambda[0],fLambda[1]) ;
ff1e7e2f 561
ab48128d 562}
563
564//______________________________________________________________________________
565void AliEMCALRecPoint::EvalPrimaries(TClonesArray * digits)
566{
567 // Constructs the list of primary particles (tracks) which have contributed to this RecPoint
568
569 AliEMCALDigit * digit ;
570 Int_t * tempo = new Int_t[fMaxTrack] ;
571
572 Int_t index ;
573 for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
574 digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ;
575 Int_t nprimaries = digit->GetNprimary() ;
5c0368b8 576 if ( nprimaries == 0 ) continue ;
ab48128d 577 Int_t * newprimaryarray = new Int_t[nprimaries] ;
578 Int_t ii ;
579 for ( ii = 0 ; ii < nprimaries ; ii++)
580 newprimaryarray[ii] = digit->GetPrimary(ii+1) ;
581
582 Int_t jndex ;
583 for ( jndex = 0 ; jndex < nprimaries ; jndex++ ) { // all primaries in digit
584 if ( fMulTrack > fMaxTrack ) {
f792c312 585 fMulTrack = fMaxTrack ;
9859bfc0 586 Error("GetNprimaries", "increase fMaxTrack ") ;
ab48128d 587 break ;
588 }
589 Int_t newprimary = newprimaryarray[jndex] ;
590 Int_t kndex ;
591 Bool_t already = kFALSE ;
592 for ( kndex = 0 ; kndex < fMulTrack ; kndex++ ) { //check if not already stored
593 if ( newprimary == tempo[kndex] ){
594 already = kTRUE ;
595 break ;
596 }
597 } // end of check
5c0368b8 598 if ( !already && (fMulTrack < fMaxTrack)) { // store it
ab48128d 599 tempo[fMulTrack] = newprimary ;
600 fMulTrack++ ;
601 } // store it
602 } // all primaries in digit
f792c312 603 delete [] newprimaryarray ;
ab48128d 604 } // all digits
605
606
607 fTracksList = new Int_t[fMulTrack] ;
608 for(index = 0; index < fMulTrack; index++)
609 fTracksList[index] = tempo[index] ;
610
f792c312 611 delete [] tempo ;
ab48128d 612
613}
7ee5c5be 614
87cdc3be 615//______________________________________________________________________________
616void AliEMCALRecPoint::EvalParents(TClonesArray * digits)
617{
618 // Constructs the list of parent particles (tracks) which have contributed to this RecPoint
619
620 AliEMCALDigit * digit ;
621 Int_t * tempo = new Int_t[fMaxParent] ;
622
623 Int_t index ;
624 for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
625 digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ;
626 Int_t nparents = digit->GetNiparent() ;
5c0368b8 627 if ( nparents == 0 ) continue ;
87cdc3be 628 Int_t * newparentarray = new Int_t[nparents] ;
629 Int_t ii ;
630 for ( ii = 0 ; ii < nparents ; ii++)
631 newparentarray[ii] = digit->GetIparent(ii+1) ;
632
633 Int_t jndex ;
634 for ( jndex = 0 ; jndex < nparents ; jndex++ ) { // all primaries in digit
635 if ( fMulParent > fMaxParent ) {
636 fMulTrack = - 1 ;
637 Error("GetNiparent", "increase fMaxParent") ;
638 break ;
639 }
640 Int_t newparent = newparentarray[jndex] ;
641 Int_t kndex ;
642 Bool_t already = kFALSE ;
f1d429fd 643 for ( kndex = 0 ; kndex < fMulParent ; kndex++ ) { //check if not already stored
87cdc3be 644 if ( newparent == tempo[kndex] ){
645 already = kTRUE ;
646 break ;
647 }
648 } // end of check
5c0368b8 649 if ( !already && (fMulTrack < fMaxTrack)) { // store it
87cdc3be 650 tempo[fMulParent] = newparent ;
651 fMulParent++ ;
652 } // store it
653 } // all parents in digit
27e2a47c 654 delete [] newparentarray ;
87cdc3be 655 } // all digits
656
27e2a47c 657 if (fMulParent>0) {
658 fParentsList = new Int_t[fMulParent] ;
659 for(index = 0; index < fMulParent; index++)
660 fParentsList[index] = tempo[index] ;
661 }
87cdc3be 662
27e2a47c 663 delete [] tempo ;
87cdc3be 664
665}
666
70a93198 667//____________________________________________________________________________
668void AliEMCALRecPoint::GetLocalPosition(TVector3 & lpos) const
669{
670 // returns the position of the cluster in the local reference system of ALICE
671 // X = eta, Y = phi, Z = r (a constant for the EMCAL)
672
673 lpos.SetX(fLocPos.X()) ;
674 lpos.SetY(fLocPos.Y()) ;
675 lpos.SetZ(fLocPos.Z()) ;
676}
677
ab48128d 678//____________________________________________________________________________
679void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const
680{
681 // returns the position of the cluster in the global reference system of ALICE
70a93198 682 // These are now the Cartesian X, Y and Z
683
5dee926e 684 //AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
685 AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
70a93198 686 Int_t absid = geom->TowerIndexFromEtaPhi(fLocPos.X(), TMath::RadToDeg()*fLocPos.Y());
687 geom->XYZFromIndex(absid, gpos);
688}
689
690//____________________________________________________________________________
691Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const
692{
693 // Finds the maximum energy in the cluster
ab48128d 694
70a93198 695 Float_t menergy = 0. ;
696
697 Int_t iDigit;
698
699 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
700
701 if(fEnergyList[iDigit] > menergy)
702 menergy = fEnergyList[iDigit] ;
703 }
704 return menergy ;
ab48128d 705}
706
aad8e277 707//____________________________________________________________________________
70a93198 708Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const
aad8e277 709{
70a93198 710 // Calculates the multiplicity of digits with energy larger than H*energy
711
712 Int_t multipl = 0 ;
713 Int_t iDigit ;
714 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
715
716 if(fEnergyList[iDigit] > H * fAmp)
717 multipl++ ;
718 }
719 return multipl ;
720}
721
722//____________________________________________________________________________
723Int_t AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
724 Float_t locMaxCut,TClonesArray * digits) const
725{
726 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
727 // energy difference between two local maxima
728
729 AliEMCALDigit * digit ;
730 AliEMCALDigit * digitN ;
731
732 Int_t iDigitN ;
733 Int_t iDigit ;
734
735 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
736 maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit]) ;
737
738 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
739 if(maxAt[iDigit]) {
740 digit = maxAt[iDigit] ;
741
742 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
743 digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ;
744
745 if ( AreNeighbours(digit, digitN) ) {
746 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
747 maxAt[iDigitN] = 0 ;
748 // but may be digit too is not local max ?
749 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
750 maxAt[iDigit] = 0 ;
751 }
752 else {
753 maxAt[iDigit] = 0 ;
754 // but may be digitN too is not local max ?
755 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
756 maxAt[iDigitN] = 0 ;
757 }
758 } // if Areneighbours
759 } // while digitN
760 } // slot not empty
761 } // while digit
762
763 iDigitN = 0 ;
764 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
765 if(maxAt[iDigit] ){
766 maxAt[iDigitN] = maxAt[iDigit] ;
767 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
768 iDigitN++ ;
769 }
770 }
771 return iDigitN ;
772}
773//____________________________________________________________________________
774void AliEMCALRecPoint::EvalTime(TClonesArray * digits){
775 // time is set to the time of the digit with the maximum energy
776
777 Float_t maxE = 0;
778 Int_t maxAt = 0;
779 for(Int_t idig=0; idig < fMulDigit; idig++){
780 if(fEnergyList[idig] > maxE){
781 maxE = fEnergyList[idig] ;
782 maxAt = idig;
783 }
784 }
785 fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
aad8e277 786
aad8e277 787}
ab48128d 788
789//______________________________________________________________________________
790void AliEMCALRecPoint::Paint(Option_t *)
791{
792 // Paint this ALiRecPoint as a TMarker with its current attributes
793
794 TVector3 pos(0.,0.,0.) ;
795 GetLocalPosition(pos) ;
796 Coord_t x = pos.X() ;
797 Coord_t y = pos.Z() ;
798 Color_t markercolor = 1 ;
799 Size_t markersize = 1. ;
800 Style_t markerstyle = 5 ;
801
802 if (!gPad->IsBatch()) {
803 gVirtualX->SetMarkerColor(markercolor) ;
804 gVirtualX->SetMarkerSize (markersize) ;
805 gVirtualX->SetMarkerStyle(markerstyle) ;
806 }
807 gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ;
808 gPad->PaintPolyMarker(1,&x,&y,"") ;
809}
70a93198 810
811//______________________________________________________________________________
812Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const
813{
814 //Converts Theta (Radians) to Eta(Radians)
815 return (2.*TMath::ATan(TMath::Exp(-arg)));
816}
817
818//______________________________________________________________________________
819Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const
820{
821 //Converts Eta (Radians) to Theta(Radians)
822 return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));
823}
261b1065 824
825//____________________________________________________________________________
826void AliEMCALRecPoint::Print(Option_t *) const
827{
828 // Print the list of digits belonging to the cluster
829
830 TString message ;
831 message = "AliPHOSEmcRecPoint:\n" ;
832 message += " digits # = " ;
833 Info("Print", message.Data()) ;
834
835 Int_t iDigit;
836 for(iDigit=0; iDigit<fMulDigit; iDigit++)
837 printf(" %d ", fDigitsList[iDigit] ) ;
838
839 Info("Print", " Energies = ") ;
840 for(iDigit=0; iDigit<fMulDigit; iDigit++)
841 printf(" %f ", fEnergyList[iDigit] ) ;
842 printf("\n") ;
843 Info("Print", " Primaries ") ;
844 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
845 printf(" %d ", fTracksList[iDigit]) ;
846 printf("\n") ;
847 message = " Multiplicity = %d" ;
848 message += " Cluster Energy = %f" ;
849 message += " Core energy = %f" ;
850 message += " Core radius = %f" ;
851 message += " Number of primaries %d" ;
852 message += " Stored at position %d" ;
853
854 Info("Print", message.Data(), fMulDigit, fAmp, fCoreEnergy, fCoreRadius, fMulTrack, GetIndexInList() ) ;
855}