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 **************************************************************************/
18 //_________________________________________________________________________
19 // RecPoint implementation for EMCAL-EMC
20 // An TowerRecPoint is a cluster of digits
22 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
25 // --- ROOT system ---
31 // --- Standard library ---
33 // --- AliRoot header files ---
35 #include "AliGenerator.h"
36 #include "AliEMCALGeometry.h"
37 #include "AliEMCALTowerRecPoint.h"
39 #include "AliEMCALGetter.h"
41 ClassImp(AliEMCALTowerRecPoint)
43 //____________________________________________________________________________
44 AliEMCALTowerRecPoint::AliEMCALTowerRecPoint() : AliEMCALRecPoint()
53 fLocPos.SetX(0.) ; //Local position should be evaluated
57 //____________________________________________________________________________
58 AliEMCALTowerRecPoint::AliEMCALTowerRecPoint(const char * opt) : AliEMCALRecPoint(opt)
67 fLocPos.SetX(1000000.) ; //Local position should be evaluated
71 //____________________________________________________________________________
72 AliEMCALTowerRecPoint::~AliEMCALTowerRecPoint()
77 delete[] fEnergyList ;
80 //____________________________________________________________________________
81 void AliEMCALTowerRecPoint::AddDigit(AliEMCALDigit & digit, Float_t Energy)
83 // Adds a digit to the RecPoint
84 // and accumulates the total amplitude and the multiplicity
87 fEnergyList = new Float_t[fMaxDigit];
89 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
91 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
92 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
95 for ( index = 0 ; index < fMulDigit ; index++ ){
96 tempo[index] = fDigitsList[index] ;
97 tempoE[index] = fEnergyList[index] ;
100 delete [] fDigitsList ;
101 fDigitsList = new ( Int_t[fMaxDigit] ) ;
103 delete [] fEnergyList ;
104 fEnergyList = new ( Float_t[fMaxDigit] ) ;
106 for ( index = 0 ; index < fMulDigit ; index++ ){
107 fDigitsList[index] = tempo[index] ;
108 fEnergyList[index] = tempoE[index] ;
115 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
116 fEnergyList[fMulDigit] = Energy ;
120 // EvalEMCALMod(&digit) ;
123 //____________________________________________________________________________
124 Bool_t AliEMCALTowerRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
126 // Tells if (true) or not (false) two digits are neighbors
128 Bool_t aren = kFALSE ;
130 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
131 AliEMCALGeometry * phosgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
134 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
137 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
139 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
140 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
142 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
148 //____________________________________________________________________________
149 Int_t AliEMCALTowerRecPoint::Compare(const TObject * obj) const
151 // Compares two RecPoints according to their position in the EMCAL modules
153 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
154 //value (what is senseless) change as vell delta in
155 //AliEMCALTrackSegmentMakerv* and other RecPoints...
158 AliEMCALTowerRecPoint * clu = (AliEMCALTowerRecPoint *)obj ;
161 Int_t phosmod1 = GetEMCALArm() ;
162 Int_t phosmod2 = clu->GetEMCALArm() ;
165 GetLocalPosition(locpos1) ;
167 clu->GetLocalPosition(locpos2) ;
169 if(phosmod1 == phosmod2 ) {
170 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
175 else if(locpos1.Z()>locpos2.Z())
182 if(phosmod1 < phosmod2 )
190 //______________________________________________________________________________
191 void AliEMCALTowerRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
194 // Execute action corresponding to one event
195 // This member function is called when a AliEMCALRecPoint is clicked with the locator
197 // If Left button is clicked on AliEMCALRecPoint, the digits are switched on
198 // and switched off when the mouse button is released.
201 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
202 // if(!gime) return ;
203 // AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
205 // static TGraph * digitgraph = 0 ;
207 // if (!gPad->IsEditable()) return;
209 // TH2F * histo = 0 ;
210 // TCanvas * histocanvas ;
212 // const TClonesArray * digits = gime->Digits() ;
216 // case kButton1Down: {
217 // AliEMCALDigit * digit ;
221 // const Int_t kMulDigit = AliEMCALTowerRecPoint::GetDigitsMultiplicity() ;
222 // Float_t * xi = new Float_t[kMulDigit] ;
223 // Float_t * zi = new Float_t[kMulDigit] ;
225 // // create the histogram for the single cluster
226 // // 1. gets histogram boundaries
227 // Float_t ximax = -999. ;
228 // Float_t zimax = -999. ;
229 // Float_t ximin = 999. ;
230 // Float_t zimin = 999. ;
232 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
233 // digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
234 // emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
235 // emcalgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
236 // if ( xi[iDigit] > ximax )
237 // ximax = xi[iDigit] ;
238 // if ( xi[iDigit] < ximin )
239 // ximin = xi[iDigit] ;
240 // if ( zi[iDigit] > zimax )
241 // zimax = zi[iDigit] ;
242 // if ( zi[iDigit] < zimin )
243 // zimin = zi[iDigit] ;
245 // ximax += emcalgeom->GetCrystalSize(0) / 2. ;
246 // zimax += emcalgeom->GetCrystalSize(2) / 2. ;
247 // ximin -= emcalgeom->GetCrystalSize(0) / 2. ;
248 // zimin -= emcalgeom->GetCrystalSize(2) / 2. ;
249 // Int_t xdim = (int)( (ximax - ximin ) / emcalgeom->GetCrystalSize(0) + 0.5 ) ;
250 // Int_t zdim = (int)( (zimax - zimin ) / emcalgeom->GetCrystalSize(2) + 0.5 ) ;
252 // // 2. gets the histogram title
254 // Text_t title[100] ;
255 // sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
261 // histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
264 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
265 // digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
266 // emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
267 // emcalgeom->RelPosInModule(relid, x, z);
268 // histo->Fill(x, z, fEnergyList[iDigit] ) ;
271 // if (!digitgraph) {
272 // digitgraph = new TGraph(kMulDigit,xi,zi);
273 // digitgraph-> SetMarkerStyle(5) ;
274 // digitgraph-> SetMarkerSize(1.) ;
275 // digitgraph-> SetMarkerColor(1) ;
276 // digitgraph-> Paint("P") ;
280 // histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
281 // histocanvas->Draw() ;
282 // histo->Draw("lego1") ;
292 // delete digitgraph ;
300 //____________________________________________________________________________
301 void AliEMCALTowerRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
303 // Calculates the dispersion of the shower at the origine of the RecPoint
308 AliEMCALDigit * digit ;
310 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
311 AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
314 // Calculates the center of gravity in the local EMCAL-module coordinates
319 if (!fTheta || !fPhi ) {
320 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
321 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
325 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
326 emcalgeom->PosInAlice(relid, thetai, phii);
327 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
328 fTheta = fTheta + thetai * w ;
343 const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ;
345 Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
346 Float_t x = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
347 Float_t y = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
348 Float_t z = cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ;
350 // Calculates the dispersion in coordinates
352 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
353 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
354 Float_t thetai = 0. ;
356 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
357 emcalgeom->PosInAlice(relid, thetai, phii);
358 Float_t xi = cyl_radius * TMath::Cos(phii * kDeg2Rad ) ;
359 Float_t yi = cyl_radius * TMath::Sin(phii * kDeg2Rad ) ;
360 Float_t zi = cyl_radius * TMath::Tan(thetai * kDeg2Rad ) ;
362 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
363 d += w*((xi-x)*(xi-x) + (yi-y)*(yi-y)+ (zi-z)*(zi-z) ) ;
372 fDispersion = TMath::Sqrt(d) ;
375 //______________________________________________________________________________
376 void AliEMCALTowerRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
378 // This function calculates energy in the core,
379 // i.e. within a radius rad = 3cm around the center. Beyond this radius
380 // in accordance with shower profile the energy deposition
381 // should be less than 2%
383 Float_t coreRadius = 10. ;
385 AliEMCALDigit * digit ;
389 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
390 AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
394 if (!fTheta || !fPhi ) {
395 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
396 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
400 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
401 emcalgeom->PosInAlice(relid, thetai, phii);
402 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
403 fTheta = fTheta + thetai * w ;
417 const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ;
419 Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
420 Float_t x = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
421 Float_t y = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
422 Float_t z = cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ;
424 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
425 digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
427 Float_t thetai = 0. ;
429 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
430 emcalgeom->PosInAlice(relid, thetai, phii);
432 Float_t xi = cyl_radius * TMath::Cos(phii * kDeg2Rad ) ;
433 Float_t yi = cyl_radius * TMath::Sin(phii * kDeg2Rad ) ;
434 Float_t zi = cyl_radius * TMath::Tan(thetai * kDeg2Rad ) ;
436 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(yi-y)*(yi-y)+(zi-z)*(zi-z)) ;
437 if(distance < coreRadius)
438 fCoreEnergy += fEnergyList[iDigit] ;
443 //____________________________________________________________________________
444 void AliEMCALTowerRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
446 // Calculates the axis of the shower ellipsoid
455 AliEMCALDigit * digit ;
457 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
458 AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
461 const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ;
463 Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
465 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
466 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
468 Float_t thetai = 0. ;
470 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
471 emcalgeom->PosInAlice(relid, thetai, phii);
472 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
473 Float_t xi = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
474 Float_t zi = cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ;
493 // //Apply correction due to non-perpendicular incidence
496 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
497 // AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
498 // Double_t DistanceToIP= (Double_t ) emcalgeom->GetIPDistance() ;
500 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
501 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
503 // dxx = dxx/(CosX*CosX) ;
504 // dzz = dzz/(CosZ*CosZ) ;
505 // dxz = dxz/(CosX*CosZ) ;
508 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
510 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
512 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
513 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
514 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
523 //____________________________________________________________________________
524 void AliEMCALTowerRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
526 // Evaluates all shower parameters
528 AliEMCALRecPoint::EvalAll(logWeight,digits) ;
529 EvalGlobalPosition(logWeight, digits) ;
530 EvalElipsAxis(logWeight, digits) ;
531 EvalDispersion(logWeight, digits) ;
532 EvalCoreEnergy(logWeight, digits);
536 //____________________________________________________________________________
537 void AliEMCALTowerRecPoint::EvalGlobalPosition(Float_t logWeight, TClonesArray * digits)
539 // Calculates the center of gravity in the local EMCAL-module coordinates
544 AliEMCALDigit * digit ;
546 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
547 AliEMCALGeometry * emcalgeom = static_cast<AliEMCALGeometry*>(gime->EMCALGeometry());
550 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
551 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
555 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
556 emcalgeom->PosInAlice(relid, thetai, phii);
557 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
558 fTheta = fTheta + thetai * w ;
578 //____________________________________________________________________________
579 Float_t AliEMCALTowerRecPoint::GetMaximalEnergy(void) const
581 // Finds the maximum energy in the cluster
583 Float_t menergy = 0. ;
587 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
589 if(fEnergyList[iDigit] > menergy)
590 menergy = fEnergyList[iDigit] ;
595 //____________________________________________________________________________
596 Int_t AliEMCALTowerRecPoint::GetMultiplicityAtLevel(const Float_t H) const
598 // Calculates the multiplicity of digits with energy larger than H*energy
602 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
604 if(fEnergyList[iDigit] > H * fAmp)
610 //____________________________________________________________________________
611 Int_t AliEMCALTowerRecPoint::GetNumberOfLocalMax(AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
612 Float_t locMaxCut,TClonesArray * digits) const
614 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
615 // energy difference between two local maxima
617 AliEMCALDigit * digit ;
618 AliEMCALDigit * digitN ;
624 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
625 maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit]) ;
628 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
630 digit = maxAt[iDigit] ;
632 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
633 digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ;
635 if ( AreNeighbours(digit, digitN) ) {
636 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
638 // but may be digit too is not local max ?
639 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
644 // but may be digitN too is not local max ?
645 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
648 } // if Areneighbours
654 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
656 maxAt[iDigitN] = maxAt[iDigit] ;
657 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
663 //____________________________________________________________________________
664 void AliEMCALTowerRecPoint::EvalTime(TClonesArray * digits){
668 for(Int_t idig=0; idig < fMulDigit; idig++){
669 if(fEnergyList[idig] > maxE){
670 maxE = fEnergyList[idig] ;
674 fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
677 //____________________________________________________________________________
678 void AliEMCALTowerRecPoint::Print(Option_t * option)
680 // Print the list of digits belonging to the cluster
682 TString message("\n") ;
685 message += "digits # = " ;
686 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
687 message += fDigitsList[iDigit] ;
691 message += "\nEnergies = " ;
692 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
693 message += fEnergyList[iDigit] ;
697 message += "\nPrimaries " ;
698 for(iDigit = 0;iDigit < fMulTrack; iDigit++) {
699 message += fTracksList[iDigit] ;
702 message += "\n Multiplicity = " ;
703 message += fMulDigit ;
704 message += "\n Cluster Energy = " ;
706 message += "\n Number of primaries " ;
707 message += fMulTrack ;
708 message += "\n Stored at position " ;
709 message += GetIndexInList() ;
711 Info("Print", message.Data() ) ;