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 ---
35 // --- AliRoot header files ---
37 #include "AliGenerator.h"
38 #include "AliEMCALGeometry.h"
39 #include "AliEMCALTowerRecPoint.h"
41 #include "AliEMCALGetter.h"
43 ClassImp(AliEMCALTowerRecPoint)
45 //____________________________________________________________________________
46 AliEMCALTowerRecPoint::AliEMCALTowerRecPoint() : AliEMCALRecPoint()
55 fLocPos.SetX(0.) ; //Local position should be evaluated
59 //____________________________________________________________________________
60 AliEMCALTowerRecPoint::AliEMCALTowerRecPoint(const char * opt) : AliEMCALRecPoint(opt)
69 fLocPos.SetX(1000000.) ; //Local position should be evaluated
73 //____________________________________________________________________________
74 AliEMCALTowerRecPoint::~AliEMCALTowerRecPoint()
79 delete[] fEnergyList ;
82 //____________________________________________________________________________
83 void AliEMCALTowerRecPoint::AddDigit(AliEMCALDigit & digit, Float_t Energy)
85 // Adds a digit to the RecPoint
86 // and accumulates the total amplitude and the multiplicity
89 fEnergyList = new Float_t[fMaxDigit];
91 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
93 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
94 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
97 for ( index = 0 ; index < fMulDigit ; index++ ){
98 tempo[index] = fDigitsList[index] ;
99 tempoE[index] = fEnergyList[index] ;
102 delete [] fDigitsList ;
103 fDigitsList = new ( Int_t[fMaxDigit] ) ;
105 delete [] fEnergyList ;
106 fEnergyList = new ( Float_t[fMaxDigit] ) ;
108 for ( index = 0 ; index < fMulDigit ; index++ ){
109 fDigitsList[index] = tempo[index] ;
110 fEnergyList[index] = tempoE[index] ;
117 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
118 fEnergyList[fMulDigit] = Energy ;
122 // EvalEMCALMod(&digit) ;
125 //____________________________________________________________________________
126 Bool_t AliEMCALTowerRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
128 // Tells if (true) or not (false) two digits are neighbors
130 Bool_t aren = kFALSE ;
132 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
133 AliEMCALGeometry * phosgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
136 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
139 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
141 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
142 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
144 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
150 //____________________________________________________________________________
151 Int_t AliEMCALTowerRecPoint::Compare(const TObject * obj) const
153 // Compares two RecPoints according to their position in the EMCAL modules
155 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
156 //value (what is senseless) change as vell delta in
157 //AliEMCALTrackSegmentMakerv* and other RecPoints...
160 AliEMCALTowerRecPoint * clu = (AliEMCALTowerRecPoint *)obj ;
163 Int_t phosmod1 = GetEMCALArm() ;
164 Int_t phosmod2 = clu->GetEMCALArm() ;
167 GetLocalPosition(locpos1) ;
169 clu->GetLocalPosition(locpos2) ;
171 if(phosmod1 == phosmod2 ) {
172 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
177 else if(locpos1.Z()>locpos2.Z())
184 if(phosmod1 < phosmod2 )
192 //______________________________________________________________________________
193 void AliEMCALTowerRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
196 // Execute action corresponding to one event
197 // This member function is called when a AliEMCALRecPoint is clicked with the locator
199 // If Left button is clicked on AliEMCALRecPoint, the digits are switched on
200 // and switched off when the mouse button is released.
203 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
204 // if(!gime) return ;
205 // AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
207 // static TGraph * digitgraph = 0 ;
209 // if (!gPad->IsEditable()) return;
211 // TH2F * histo = 0 ;
212 // TCanvas * histocanvas ;
214 // const TClonesArray * digits = gime->Digits() ;
218 // case kButton1Down: {
219 // AliEMCALDigit * digit ;
223 // const Int_t kMulDigit = AliEMCALTowerRecPoint::GetDigitsMultiplicity() ;
224 // Float_t * xi = new Float_t[kMulDigit] ;
225 // Float_t * zi = new Float_t[kMulDigit] ;
227 // // create the histogram for the single cluster
228 // // 1. gets histogram boundaries
229 // Float_t ximax = -999. ;
230 // Float_t zimax = -999. ;
231 // Float_t ximin = 999. ;
232 // Float_t zimin = 999. ;
234 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
235 // digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
236 // emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
237 // emcalgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
238 // if ( xi[iDigit] > ximax )
239 // ximax = xi[iDigit] ;
240 // if ( xi[iDigit] < ximin )
241 // ximin = xi[iDigit] ;
242 // if ( zi[iDigit] > zimax )
243 // zimax = zi[iDigit] ;
244 // if ( zi[iDigit] < zimin )
245 // zimin = zi[iDigit] ;
247 // ximax += emcalgeom->GetCrystalSize(0) / 2. ;
248 // zimax += emcalgeom->GetCrystalSize(2) / 2. ;
249 // ximin -= emcalgeom->GetCrystalSize(0) / 2. ;
250 // zimin -= emcalgeom->GetCrystalSize(2) / 2. ;
251 // Int_t xdim = (int)( (ximax - ximin ) / emcalgeom->GetCrystalSize(0) + 0.5 ) ;
252 // Int_t zdim = (int)( (zimax - zimin ) / emcalgeom->GetCrystalSize(2) + 0.5 ) ;
254 // // 2. gets the histogram title
256 // Text_t title[100] ;
257 // sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
263 // histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
266 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
267 // digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
268 // emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
269 // emcalgeom->RelPosInModule(relid, x, z);
270 // histo->Fill(x, z, fEnergyList[iDigit] ) ;
273 // if (!digitgraph) {
274 // digitgraph = new TGraph(kMulDigit,xi,zi);
275 // digitgraph-> SetMarkerStyle(5) ;
276 // digitgraph-> SetMarkerSize(1.) ;
277 // digitgraph-> SetMarkerColor(1) ;
278 // digitgraph-> Paint("P") ;
282 // histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
283 // histocanvas->Draw() ;
284 // histo->Draw("lego1") ;
294 // delete digitgraph ;
302 //____________________________________________________________________________
303 void AliEMCALTowerRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
305 // Calculates the dispersion of the shower at the origine of the RecPoint
310 AliEMCALDigit * digit ;
312 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
313 AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
316 // Calculates the center of gravity in the local EMCAL-module coordinates
321 if (!fTheta || !fPhi ) {
322 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
323 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
327 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
328 emcalgeom->PosInAlice(relid, thetai, phii);
329 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
330 fTheta = fTheta + thetai * w ;
345 const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ;
347 Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
348 Float_t x = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
349 Float_t y = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
350 Float_t z = cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ;
352 // Calculates the dispersion in coordinates
354 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
355 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
356 Float_t thetai = 0. ;
358 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
359 emcalgeom->PosInAlice(relid, thetai, phii);
360 Float_t xi = cyl_radius * TMath::Cos(phii * kDeg2Rad ) ;
361 Float_t yi = cyl_radius * TMath::Sin(phii * kDeg2Rad ) ;
362 Float_t zi = cyl_radius * TMath::Tan(thetai * kDeg2Rad ) ;
364 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
365 d += w*((xi-x)*(xi-x) + (yi-y)*(yi-y)+ (zi-z)*(zi-z) ) ;
374 fDispersion = TMath::Sqrt(d) ;
377 //______________________________________________________________________________
378 void AliEMCALTowerRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
380 // This function calculates energy in the core,
381 // i.e. within a radius rad = 3cm around the center. Beyond this radius
382 // in accordance with shower profile the energy deposition
383 // should be less than 2%
385 Float_t coreRadius = 10. ;
387 AliEMCALDigit * digit ;
391 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
392 AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
396 if (!fTheta || !fPhi ) {
397 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
398 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
402 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
403 emcalgeom->PosInAlice(relid, thetai, phii);
404 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
405 fTheta = fTheta + thetai * w ;
419 const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ;
421 Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
422 Float_t x = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
423 Float_t y = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
424 Float_t z = cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ;
426 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
427 digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
429 Float_t thetai = 0. ;
431 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
432 emcalgeom->PosInAlice(relid, thetai, phii);
434 Float_t xi = cyl_radius * TMath::Cos(phii * kDeg2Rad ) ;
435 Float_t yi = cyl_radius * TMath::Sin(phii * kDeg2Rad ) ;
436 Float_t zi = cyl_radius * TMath::Tan(thetai * kDeg2Rad ) ;
438 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(yi-y)*(yi-y)+(zi-z)*(zi-z)) ;
439 if(distance < coreRadius)
440 fCoreEnergy += fEnergyList[iDigit] ;
445 //____________________________________________________________________________
446 void AliEMCALTowerRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
448 // Calculates the axis of the shower ellipsoid
457 AliEMCALDigit * digit ;
459 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
460 AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
463 const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ;
465 Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
467 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
468 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
470 Float_t thetai = 0. ;
472 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
473 emcalgeom->PosInAlice(relid, thetai, phii);
474 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
475 Float_t xi = cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
476 Float_t zi = cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ;
495 // //Apply correction due to non-perpendicular incidence
498 // AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
499 // AliEMCALGeometry * emcalgeom = (AliEMCALGeometry*)gime->EMCALGeometry();
500 // Double_t DistanceToIP= (Double_t ) emcalgeom->GetIPDistance() ;
502 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
503 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
505 // dxx = dxx/(CosX*CosX) ;
506 // dzz = dzz/(CosZ*CosZ) ;
507 // dxz = dxz/(CosX*CosZ) ;
510 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
512 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
514 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
515 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
516 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
525 //____________________________________________________________________________
526 void AliEMCALTowerRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
528 // Evaluates all shower parameters
530 AliEMCALRecPoint::EvalAll(logWeight,digits) ;
531 EvalGlobalPosition(logWeight, digits) ;
532 EvalElipsAxis(logWeight, digits) ;
533 EvalDispersion(logWeight, digits) ;
534 EvalCoreEnergy(logWeight, digits);
538 //____________________________________________________________________________
539 void AliEMCALTowerRecPoint::EvalGlobalPosition(Float_t logWeight, TClonesArray * digits)
541 // Calculates the center of gravity in the local EMCAL-module coordinates
546 AliEMCALDigit * digit ;
548 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
549 AliEMCALGeometry * emcalgeom = static_cast<AliEMCALGeometry*>(gime->EMCALGeometry());
552 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
553 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
557 emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
558 emcalgeom->PosInAlice(relid, thetai, phii);
559 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
560 fTheta = fTheta + thetai * w ;
580 //____________________________________________________________________________
581 Float_t AliEMCALTowerRecPoint::GetMaximalEnergy(void) const
583 // Finds the maximum energy in the cluster
585 Float_t menergy = 0. ;
589 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
591 if(fEnergyList[iDigit] > menergy)
592 menergy = fEnergyList[iDigit] ;
597 //____________________________________________________________________________
598 Int_t AliEMCALTowerRecPoint::GetMultiplicityAtLevel(const Float_t H) const
600 // Calculates the multiplicity of digits with energy larger than H*energy
604 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
606 if(fEnergyList[iDigit] > H * fAmp)
612 //____________________________________________________________________________
613 Int_t AliEMCALTowerRecPoint::GetNumberOfLocalMax(Int_t * maxAt, Float_t * maxAtEnergy,
614 Float_t locMaxCut,TClonesArray * digits) const
616 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
617 // energy difference between two local maxima
619 AliEMCALDigit * digit ;
620 AliEMCALDigit * digitN ;
626 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
627 maxAt[iDigit] = (Int_t) digits->At(fDigitsList[iDigit]) ;
630 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
631 if(maxAt[iDigit] != -1) {
632 digit = (AliEMCALDigit *) maxAt[iDigit] ;
634 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
635 digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ;
637 if ( AreNeighbours(digit, digitN) ) {
638 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
639 maxAt[iDigitN] = -1 ;
640 // but may be digit too is not local max ?
641 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
646 // but may be digitN too is not local max ?
647 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
648 maxAt[iDigitN] = -1 ;
650 } // if Areneighbours
656 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
657 if(maxAt[iDigit] != -1){
658 maxAt[iDigitN] = maxAt[iDigit] ;
659 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
665 //____________________________________________________________________________
666 void AliEMCALTowerRecPoint::EvalTime(TClonesArray * digits){
670 for(Int_t idig=0; idig < fMulDigit; idig++){
671 if(fEnergyList[idig] > maxE){
672 maxE = fEnergyList[idig] ;
676 fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
679 //____________________________________________________________________________
680 void AliEMCALTowerRecPoint::Print(Option_t * option)
682 // Print the list of digits belonging to the cluster
684 cout << "AliEMCALTowerRecPoint: " << endl ;
687 cout << " digits # = " ;
688 for(iDigit=0; iDigit<fMulDigit; iDigit++)
689 cout << fDigitsList[iDigit] << " " ;
692 cout << " Energies = " ;
693 for(iDigit=0; iDigit<fMulDigit; iDigit++)
694 cout << fEnergyList[iDigit] << " ";
697 cout << " Primaries " ;
698 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
699 cout << fTracksList[iDigit] << " " << endl ;
701 cout << " Multiplicity = " << fMulDigit << endl ;
702 cout << " Cluster Energy = " << fAmp << endl ;
703 cout << " Number of primaries " << fMulTrack << endl ;
704 cout << " Stored at position " << GetIndexInList() << endl ;