1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //_________________________________________________________________________
17 // Reconstructed Points for the EMCAL
18 // A RecPoint is a cluster of digits
19 //*-- Author: Yves Schutz (SUBATECH)
20 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
21 //*-- Author: Heather Gray (LBL) merged AliEMCALRecPoint and AliEMCALTowerRecPoint 02/04
23 // --- ROOT system ---
26 #include "TPaveText.h"
27 #include "TClonesArray.h"
30 // --- Standard library ---
32 // --- AliRoot header files ---
33 #include "AliGenerator.h"
34 #include "AliEMCALGeometry.h"
35 #include "AliEMCALDigit.h"
36 #include "AliEMCALRecPoint.h"
37 #include "AliEMCALGetter.h"
39 ClassImp(AliEMCALRecPoint)
42 //____________________________________________________________________________
43 AliEMCALRecPoint::AliEMCALRecPoint()
53 fLocPos.SetX(0.) ; //Local position should be evaluated
54 fCoreRadius = 10; //HG Check this
57 //____________________________________________________________________________
58 AliEMCALRecPoint::AliEMCALRecPoint(const char * opt) : AliRecPoint(opt)
67 fLocPos.SetX(1000000.) ; //Local position should be evaluated
68 fCoreRadius = 10; //HG Check this
70 //____________________________________________________________________________
71 AliEMCALRecPoint::~AliEMCALRecPoint()
75 delete[] fEnergyList ;
78 //____________________________________________________________________________
79 void AliEMCALRecPoint::AddDigit(AliEMCALDigit & digit, Float_t Energy)
81 // Adds a digit to the RecPoint
82 // and accumulates the total amplitude and the multiplicity
85 fEnergyList = new Float_t[fMaxDigit];
87 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
89 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
90 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
93 for ( index = 0 ; index < fMulDigit ; index++ ){
94 tempo[index] = fDigitsList[index] ;
95 tempoE[index] = fEnergyList[index] ;
98 delete [] fDigitsList ;
99 fDigitsList = new ( Int_t[fMaxDigit] ) ;
101 delete [] fEnergyList ;
102 fEnergyList = new ( Float_t[fMaxDigit] ) ;
104 for ( index = 0 ; index < fMulDigit ; index++ ){
105 fDigitsList[index] = tempo[index] ;
106 fEnergyList[index] = tempoE[index] ;
113 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
114 fEnergyList[fMulDigit] = Energy ;
119 //____________________________________________________________________________
120 Bool_t AliEMCALRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
122 // Tells if (true) or not (false) two digits are neighbours
123 // A neighbour is defined as being two digits which share a corner
125 Bool_t areNeighbours = kFALSE ;
127 AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
130 geom->AbsToRelNumbering(digit1->GetId(), relid1) ;
133 geom->AbsToRelNumbering(digit2->GetId(), relid2) ;
135 Int_t rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;
136 Int_t coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;
138 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
139 areNeighbours = kTRUE ;
141 return areNeighbours;
144 //____________________________________________________________________________
145 Int_t AliEMCALRecPoint::Compare(const TObject * obj) const
147 // Compares two RecPoints according to their position in the EMCAL modules
149 Float_t delta = 1 ; //Width of "Sorting row". If you change this
150 //value (what is senseless) change as well delta in
151 //AliEMCALTrackSegmentMakerv* and other RecPoints...
154 AliEMCALRecPoint * clu = (AliEMCALRecPoint *)obj ;
157 GetLocalPosition(locpos1);
159 clu->GetLocalPosition(locpos2);
161 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
166 else if(locpos1.Y()>locpos2.Y())
174 //____________________________________________________________________________
175 Int_t AliEMCALRecPoint::DistancetoPrimitive(Int_t px, Int_t py)
177 // Compute distance from point px,py to a AliEMCALRecPoint considered as a Tmarker
178 // Compute the closest distance of approach from point px,py to this marker.
179 // The distance is computed in pixels units.
180 // HG Still need to update -> Not sure what this should achieve
182 TVector3 pos(0.,0.,0.) ;
183 GetLocalPosition(pos) ;
184 Float_t x = pos.X() ;
185 Float_t y = pos.Y() ;
186 const Int_t kMaxDiff = 10;
187 Int_t pxm = gPad->XtoAbsPixel(x);
188 Int_t pym = gPad->YtoAbsPixel(y);
189 Int_t dist = (px-pxm)*(px-pxm) + (py-pym)*(py-pym);
191 if (dist > kMaxDiff) return 9999;
195 //___________________________________________________________________________
196 void AliEMCALRecPoint::Draw(Option_t *option)
198 // Draw this AliEMCALRecPoint with its current attributes
203 //______________________________________________________________________________
204 void AliEMCALRecPoint::ExecuteEvent(Int_t /*event*/, Int_t, Int_t)
206 // Execute action corresponding to one event
207 // This member function is called when a AliEMCALRecPoint is clicked with the locator
209 // If Left button is clicked on AliEMCALRecPoint, the digits are switched on
210 // and switched off when the mouse button is released.
212 // static Int_t pxold, pyold;
214 /* static TGraph * digitgraph = 0 ;
215 static TPaveText* clustertext = 0 ;
217 if (!gPad->IsEditable()) return;
223 AliEMCALDigit * digit ;
224 AliEMCALGeometry * emcalgeom = (AliEMCALGetter::Instance())->EMCALGeometry() ;
229 const Int_t kMulDigit=AliEMCALRecPoint::GetDigitsMultiplicity() ;
230 Float_t * xi = new Float_t [kMulDigit] ;
231 Float_t * zi = new Float_t [kMulDigit] ;
233 for(iDigit = 0; iDigit < kMulDigit; iDigit++) {
234 Fatal("AliEMCALRecPoint::ExecuteEvent", " -> Something wrong with the code");
235 digit = 0 ; //dynamic_cast<AliEMCALDigit *>((fDigitsList)[iDigit]);
236 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
237 emcalgeom->PosInAlice(relid, xi[iDigit], zi[iDigit]) ;
241 digitgraph = new TGraph(fMulDigit,xi,zi);
242 digitgraph-> SetMarkerStyle(5) ;
243 digitgraph-> SetMarkerSize(1.) ;
244 digitgraph-> SetMarkerColor(1) ;
245 digitgraph-> Draw("P") ;
249 TVector3 pos(0.,0.,0.) ;
250 GetLocalPosition(pos) ;
251 clustertext = new TPaveText(pos.X()-10,pos.Z()+10,pos.X()+50,pos.Z()+35,"") ;
254 sprintf(line1,"Energy=%1.2f GeV",GetEnergy()) ;
255 sprintf(line2,"%d Digits",GetDigitsMultiplicity()) ;
256 clustertext ->AddText(line1) ;
257 clustertext ->AddText(line2) ;
258 clustertext ->Draw("");
282 //____________________________________________________________________________
283 void AliEMCALRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits)
285 // Evaluates all shower parameters
287 EvalLocalPosition(logWeight, digits) ;
288 EvalElipsAxis(logWeight, digits) ;
289 EvalDispersion(logWeight, digits) ;
290 EvalCoreEnergy(logWeight, digits);
293 //EvalPrimaries(digits) ;
296 //____________________________________________________________________________
297 void AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits)
299 // Calculates the dispersion of the shower at the origin of the RecPoint
304 AliEMCALDigit * digit ;
306 AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
308 // Calculates the centre of gravity in the local EMCAL-module coordinates
311 if (!fLocPos.X() || !fLocPos.Y())
312 EvalLocalPosition(logWeight, digits) ;
314 const Float_t kDeg2Rad = TMath::DegToRad() ;
316 Float_t cluEta = fLocPos.X() ;
317 Float_t cluPhi = fLocPos.Y() ;
318 Float_t cluR = fLocPos.Z() ;
321 printf("EvalDispersion: eta,phi,r = %f,%f,%f", cluEta, cluPhi, cluR) ;
323 // Calculates the dispersion in coordinates
325 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
326 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
329 geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
330 phii = phii * kDeg2Rad;
332 printf("EvalDispersion: id = %d, etai,phii = %f,%f", digit->GetId(), etai, phii) ;
334 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
335 d += w * ( (etai-cluEta)*(etai-cluEta) + (phii-cluPhi)*(phii-cluPhi) ) ;
344 fDispersion = TMath::Sqrt(d) ;
348 //____________________________________________________________________________
349 void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
351 // Calculates the center of gravity in the local EMCAL-module coordinates
356 AliEMCALDigit * digit ;
357 AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
361 const Float_t kDeg2Rad = TMath::DegToRad();
363 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
364 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
368 geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
369 phii = phii * kDeg2Rad;
370 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
371 cluEta += (etai * w) ;
372 cluPhi += (phii * w );
384 fLocPos.SetX(cluEta);
385 fLocPos.SetY(cluPhi);
386 fLocPos.SetZ(geom->GetIP2ECASection());
389 printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ;
393 //______________________________________________________________________________
394 void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
396 // This function calculates energy in the core,
397 // i.e. within a radius rad = 3cm around the center. Beyond this radius
398 // in accordance with shower profile the energy deposition
399 // should be less than 2%
401 AliEMCALDigit * digit ;
402 const Float_t kDeg2Rad = TMath::DegToRad() ;
403 AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
406 if (!fLocPos.X() || !fLocPos.Y() ) {
407 EvalLocalPosition(logWeight, digits);
410 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
411 digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
414 geom->PosInAlice(digit->GetId(), etai, phii);
415 phii = phii * kDeg2Rad;
417 Float_t distance = TMath::Sqrt((etai-fLocPos.X())*(etai-fLocPos.X())+(phii-fLocPos.Y())*(phii-fLocPos.Y())) ;
418 if(distance < fCoreRadius)
419 fCoreEnergy += fEnergyList[iDigit] ;
423 //____________________________________________________________________________
424 void AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
426 // Calculates the axis of the shower ellipsoid in eta and phi
435 AliEMCALDigit * digit ;
437 AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
441 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
442 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
445 geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
446 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
447 dxx += w * etai * etai ;
449 dzz += w * phii * phii ;
451 dxz += w * etai * etai ;
464 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
466 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
470 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
471 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
472 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
481 //______________________________________________________________________________
482 void AliEMCALRecPoint::EvalPrimaries(TClonesArray * digits)
484 // Constructs the list of primary particles (tracks) which have contributed to this RecPoint
486 AliEMCALDigit * digit ;
487 Int_t * tempo = new Int_t[fMaxTrack] ;
490 for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
491 digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ;
492 Int_t nprimaries = digit->GetNprimary() ;
493 Int_t * newprimaryarray = new Int_t[nprimaries] ;
495 for ( ii = 0 ; ii < nprimaries ; ii++)
496 newprimaryarray[ii] = digit->GetPrimary(ii+1) ;
499 for ( jndex = 0 ; jndex < nprimaries ; jndex++ ) { // all primaries in digit
500 if ( fMulTrack > fMaxTrack ) {
502 Error("GetNprimaries", "increase fMaxTrack ") ;
505 Int_t newprimary = newprimaryarray[jndex] ;
507 Bool_t already = kFALSE ;
508 for ( kndex = 0 ; kndex < fMulTrack ; kndex++ ) { //check if not already stored
509 if ( newprimary == tempo[kndex] ){
514 if ( !already) { // store it
515 tempo[fMulTrack] = newprimary ;
518 } // all primaries in digit
519 delete newprimaryarray ;
523 fTracksList = new Int_t[fMulTrack] ;
524 for(index = 0; index < fMulTrack; index++)
525 fTracksList[index] = tempo[index] ;
531 //____________________________________________________________________________
532 void AliEMCALRecPoint::GetLocalPosition(TVector3 & lpos) const
534 // returns the position of the cluster in the local reference system of ALICE
535 // X = eta, Y = phi, Z = r (a constant for the EMCAL)
537 lpos.SetX(fLocPos.X()) ;
538 lpos.SetY(fLocPos.Y()) ;
539 lpos.SetZ(fLocPos.Z()) ;
542 //____________________________________________________________________________
543 void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const
545 // returns the position of the cluster in the global reference system of ALICE
546 // These are now the Cartesian X, Y and Z
548 AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
549 Int_t absid = geom->TowerIndexFromEtaPhi(fLocPos.X(), TMath::RadToDeg()*fLocPos.Y());
550 geom->XYZFromIndex(absid, gpos);
553 //____________________________________________________________________________
554 Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const
556 // Finds the maximum energy in the cluster
558 Float_t menergy = 0. ;
562 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
564 if(fEnergyList[iDigit] > menergy)
565 menergy = fEnergyList[iDigit] ;
570 //____________________________________________________________________________
571 Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const
573 // Calculates the multiplicity of digits with energy larger than H*energy
577 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
579 if(fEnergyList[iDigit] > H * fAmp)
585 //____________________________________________________________________________
586 Int_t AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
587 Float_t locMaxCut,TClonesArray * digits) const
589 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
590 // energy difference between two local maxima
592 AliEMCALDigit * digit ;
593 AliEMCALDigit * digitN ;
598 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
599 maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit]) ;
601 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
603 digit = maxAt[iDigit] ;
605 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
606 digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ;
608 if ( AreNeighbours(digit, digitN) ) {
609 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
611 // but may be digit too is not local max ?
612 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
617 // but may be digitN too is not local max ?
618 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
621 } // if Areneighbours
627 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
629 maxAt[iDigitN] = maxAt[iDigit] ;
630 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
636 //____________________________________________________________________________
637 void AliEMCALRecPoint::EvalTime(TClonesArray * digits){
638 // time is set to the time of the digit with the maximum energy
642 for(Int_t idig=0; idig < fMulDigit; idig++){
643 if(fEnergyList[idig] > maxE){
644 maxE = fEnergyList[idig] ;
648 fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
652 //______________________________________________________________________________
653 void AliEMCALRecPoint::Paint(Option_t *)
655 // Paint this ALiRecPoint as a TMarker with its current attributes
657 TVector3 pos(0.,0.,0.) ;
658 GetLocalPosition(pos) ;
659 Coord_t x = pos.X() ;
660 Coord_t y = pos.Z() ;
661 Color_t markercolor = 1 ;
662 Size_t markersize = 1. ;
663 Style_t markerstyle = 5 ;
665 if (!gPad->IsBatch()) {
666 gVirtualX->SetMarkerColor(markercolor) ;
667 gVirtualX->SetMarkerSize (markersize) ;
668 gVirtualX->SetMarkerStyle(markerstyle) ;
670 gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ;
671 gPad->PaintPolyMarker(1,&x,&y,"") ;
674 //______________________________________________________________________________
675 Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const
677 //Converts Theta (Radians) to Eta(Radians)
678 return (2.*TMath::ATan(TMath::Exp(-arg)));
681 //______________________________________________________________________________
682 Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const
684 //Converts Eta (Radians) to Theta(Radians)
685 return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));