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
56 //____________________________________________________________________________
57 AliEMCALTowerRecPoint::AliEMCALTowerRecPoint(const char * opt) : AliEMCALRecPoint(opt)
66 fLocPos.SetX(1000000.) ; //Local position should be evaluated
69 //____________________________________________________________________________
70 AliEMCALTowerRecPoint::~AliEMCALTowerRecPoint()
75 delete[] fEnergyList ;
78 //____________________________________________________________________________
79 void AliEMCALTowerRecPoint::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 ;
118 // EvalEMCALMod(&digit) ;
121 //____________________________________________________________________________
122 Bool_t AliEMCALTowerRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
124 // Tells if (true) or not (false) two digits are neighbors
126 Bool_t aren = kFALSE ;
128 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
129 AliEMCALGeometry * phosgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
132 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
135 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
137 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
138 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
140 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
146 //____________________________________________________________________________
147 Int_t AliEMCALTowerRecPoint::Compare(const TObject * obj) const
149 // Compares two RecPoints according to their position in the EMCAL modules
151 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
152 //value (what is senseless) change as vell delta in
153 //AliEMCALTrackSegmentMakerv* and other RecPoints...
156 AliEMCALTowerRecPoint * clu = (AliEMCALTowerRecPoint *)obj ;
159 Int_t phosmod1 = GetEMCALArm() ;
160 Int_t phosmod2 = clu->GetEMCALArm() ;
163 GetLocalPosition(locpos1) ;
165 clu->GetLocalPosition(locpos2) ;
167 if(phosmod1 == phosmod2 ) {
168 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
173 else if(locpos1.Z()>locpos2.Z())
180 if(phosmod1 < phosmod2 )
188 //______________________________________________________________________________
189 void AliEMCALTowerRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
192 // Execute action corresponding to one event
193 // This member function is called when a AliEMCALRecPoint is clicked with the locator
195 // If Left button is clicked on AliEMCALRecPoint, the digits are switched on
196 // and switched off when the mouse button is released.
199 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
200 // if(!gime) return ;
201 // AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
203 // static TGraph * digitgraph = 0 ;
205 // if (!gPad->IsEditable()) return;
207 // TH2F * histo = 0 ;
208 // TCanvas * histocanvas ;
210 // const TClonesArray * digits = gime->Digits() ;
214 // case kButton1Down: {
215 // AliEMCALDigit * digit ;
219 // const Int_t kMulDigit = AliEMCALTowerRecPoint::GetDigitsMultiplicity() ;
220 // Float_t * xi = new Float_t[kMulDigit] ;
221 // Float_t * zi = new Float_t[kMulDigit] ;
223 // // create the histogram for the single cluster
224 // // 1. gets histogram boundaries
225 // Float_t ximax = -999. ;
226 // Float_t zimax = -999. ;
227 // Float_t ximin = 999. ;
228 // Float_t zimin = 999. ;
230 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
231 // digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
232 // emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
233 // emcalgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
234 // if ( xi[iDigit] > ximax )
235 // ximax = xi[iDigit] ;
236 // if ( xi[iDigit] < ximin )
237 // ximin = xi[iDigit] ;
238 // if ( zi[iDigit] > zimax )
239 // zimax = zi[iDigit] ;
240 // if ( zi[iDigit] < zimin )
241 // zimin = zi[iDigit] ;
243 // ximax += emcalgeom->GetCrystalSize(0) / 2. ;
244 // zimax += emcalgeom->GetCrystalSize(2) / 2. ;
245 // ximin -= emcalgeom->GetCrystalSize(0) / 2. ;
246 // zimin -= emcalgeom->GetCrystalSize(2) / 2. ;
247 // Int_t xdim = (int)( (ximax - ximin ) / emcalgeom->GetCrystalSize(0) + 0.5 ) ;
248 // Int_t zdim = (int)( (zimax - zimin ) / emcalgeom->GetCrystalSize(2) + 0.5 ) ;
250 // // 2. gets the histogram title
252 // Text_t title[100] ;
253 // sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
259 // histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
262 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
263 // digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
264 // emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
265 // emcalgeom->RelPosInModule(relid, x, z);
266 // histo->Fill(x, z, fEnergyList[iDigit] ) ;
269 // if (!digitgraph) {
270 // digitgraph = new TGraph(kMulDigit,xi,zi);
271 // digitgraph-> SetMarkerStyle(5) ;
272 // digitgraph-> SetMarkerSize(1.) ;
273 // digitgraph-> SetMarkerColor(1) ;
274 // digitgraph-> Paint("P") ;
278 // histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
279 // histocanvas->Draw() ;
280 // histo->Draw("lego1") ;
290 // delete digitgraph ;
298 //____________________________________________________________________________
299 void AliEMCALTowerRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
301 // Calculates the dispersion of the shower at the origine of the RecPoint
306 AliEMCALDigit * digit ;
308 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
309 AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
312 // Calculates the center of gravity in the local EMCAL-module coordinates
317 if (!fTheta || !fPhi ) {
318 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
319 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
323 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
324 emcalgeom->PosInAlice(relid, thetai, phii);
325 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
326 fTheta = fTheta + thetai * w ;
341 const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ;
343 Float_t cyl_radius = emcalgeom->GetIP2Tower() ;
344 Float_t x = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
345 Float_t y = cyl_radius * TMath::Sin(fPhi * kDeg2Rad ) ;
346 Float_t z = cyl_radius / TMath::Tan(fTheta * kDeg2Rad ) ;
348 // Calculates the dispersion in coordinates
350 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
351 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
352 Float_t thetai = 0. ;
354 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
355 emcalgeom->PosInAlice(relid, thetai, phii);
356 Float_t xi = cyl_radius * TMath::Cos(phii * kDeg2Rad ) ;
357 Float_t yi = cyl_radius * TMath::Sin(phii * kDeg2Rad ) ;
358 Float_t zi = cyl_radius / TMath::Tan(thetai * kDeg2Rad ) ;
360 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
361 d += w*((xi-x)*(xi-x) + (yi-y)*(yi-y)+ (zi-z)*(zi-z) ) ;
370 fDispersion = TMath::Sqrt(d) ;
373 //______________________________________________________________________________
374 void AliEMCALTowerRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
376 // This function calculates energy in the core,
377 // i.e. within a radius rad = 3cm around the center. Beyond this radius
378 // in accordance with shower profile the energy deposition
379 // should be less than 2%
381 Float_t coreRadius = 10. ;
383 AliEMCALDigit * digit ;
387 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
388 AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
392 if (!fTheta || !fPhi ) {
393 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
394 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
398 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
399 emcalgeom->PosInAlice(relid, thetai, phii);
400 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
401 fTheta = fTheta + thetai * w ;
415 const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ;
417 Float_t cyl_radius = emcalgeom->GetIP2Tower();
418 Float_t x = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
419 Float_t y = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
420 Float_t z = cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ;
422 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
423 digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
425 Float_t thetai = 0. ;
427 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
428 emcalgeom->PosInAlice(relid, thetai, phii);
430 Float_t xi = cyl_radius * TMath::Cos(phii * kDeg2Rad ) ;
431 Float_t yi = cyl_radius * TMath::Sin(phii * kDeg2Rad ) ;
432 Float_t zi = cyl_radius * TMath::Tan(thetai * kDeg2Rad ) ;
434 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(yi-y)*(yi-y)+(zi-z)*(zi-z)) ;
435 if(distance < coreRadius)
436 fCoreEnergy += fEnergyList[iDigit] ;
441 //____________________________________________________________________________
442 void AliEMCALTowerRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
444 // Calculates the axis of the shower ellipsoid
453 AliEMCALDigit * digit ;
455 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
456 AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
459 const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ;
461 Float_t cyl_radius = emcalgeom->GetIP2Tower() ;
463 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
464 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
466 Float_t thetai = 0. ;
468 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
469 emcalgeom->PosInAlice(relid, thetai, phii);
470 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
471 Float_t xi = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
472 Float_t zi = cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ;
491 // //Apply correction due to non-perpendicular incidence
494 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
495 // AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
496 // Double_t DistanceToIP= (Double_t ) emcalgeom->GetIPDistance() ;
498 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
499 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
501 // dxx = dxx/(CosX*CosX) ;
502 // dzz = dzz/(CosZ*CosZ) ;
503 // dxz = dxz/(CosX*CosZ) ;
506 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
508 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
510 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
511 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
512 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
521 //____________________________________________________________________________
522 void AliEMCALTowerRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
524 // Evaluates all shower parameters
526 AliEMCALRecPoint::EvalAll(logWeight,digits) ;
527 EvalGlobalPosition(logWeight, digits) ;
528 EvalElipsAxis(logWeight, digits) ;
529 EvalDispersion(logWeight, digits) ;
530 EvalCoreEnergy(logWeight, digits);
534 //____________________________________________________________________________
535 void AliEMCALTowerRecPoint::EvalGlobalPosition(Float_t logWeight, TClonesArray * digits)
537 // Calculates the center of gravity in the local EMCAL-module coordinates
542 AliEMCALDigit * digit ;
544 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
545 AliEMCALGeometry * emcalgeom = static_cast<AliEMCALGeometry*>(gime->EMCALGeometry());
548 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
549 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
553 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
554 emcalgeom->PosInAlice(relid, thetai, phii);
555 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
556 fTheta = fTheta + thetai * w ;
576 //____________________________________________________________________________
577 Float_t AliEMCALTowerRecPoint::GetMaximalEnergy(void) const
579 // Finds the maximum energy in the cluster
581 Float_t menergy = 0. ;
585 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
587 if(fEnergyList[iDigit] > menergy)
588 menergy = fEnergyList[iDigit] ;
593 //____________________________________________________________________________
594 Int_t AliEMCALTowerRecPoint::GetMultiplicityAtLevel(const Float_t H) const
596 // Calculates the multiplicity of digits with energy larger than H*energy
600 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
602 if(fEnergyList[iDigit] > H * fAmp)
608 //____________________________________________________________________________
609 Int_t AliEMCALTowerRecPoint::GetNumberOfLocalMax(AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
610 Float_t locMaxCut,TClonesArray * digits) const
612 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
613 // energy difference between two local maxima
615 AliEMCALDigit * digit ;
616 AliEMCALDigit * digitN ;
622 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
623 maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit]) ;
626 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
628 digit = maxAt[iDigit] ;
630 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
631 digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ;
633 if ( AreNeighbours(digit, digitN) ) {
634 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
636 // but may be digit too is not local max ?
637 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
642 // but may be digitN too is not local max ?
643 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
646 } // if Areneighbours
652 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
654 maxAt[iDigitN] = maxAt[iDigit] ;
655 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
661 //____________________________________________________________________________
662 void AliEMCALTowerRecPoint::EvalTime(TClonesArray * digits){
666 for(Int_t idig=0; idig < fMulDigit; idig++){
667 if(fEnergyList[idig] > maxE){
668 maxE = fEnergyList[idig] ;
672 fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
675 //____________________________________________________________________________
676 void AliEMCALTowerRecPoint::Print(Option_t * option)
678 // Print the list of digits belonging to the cluster
680 TString message("\n") ;
683 message += "digits # = " ;
684 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
685 message += fDigitsList[iDigit] ;
689 message += "\nEnergies = " ;
690 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
691 message += fEnergyList[iDigit] ;
695 message += "\nPrimaries " ;
696 for(iDigit = 0;iDigit < fMulTrack; iDigit++) {
697 message += fTracksList[iDigit] ;
700 message += "\n Multiplicity = " ;
701 message += fMulDigit ;
702 message += "\n Cluster Energy = " ;
704 message += "\n Number of primaries " ;
705 message += fMulTrack ;
706 message += "\n Stored at position " ;
707 message += GetIndexInList() ;
709 Info("Print", message.Data() ) ;