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 PHOS-EMC
20 // An EmcRecPoint 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 "AliPHOSGeometry.h"
39 #include "AliPHOSEmcRecPoint.h"
41 #include "AliPHOSGetter.h"
43 ClassImp(AliPHOSEmcRecPoint)
45 //____________________________________________________________________________
46 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
55 fLocPos.SetX(1000000.) ; //Local position should be evaluated
59 //____________________________________________________________________________
60 AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
65 delete[] fEnergyList ;
68 //____________________________________________________________________________
69 void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
71 // Adds a digit to the RecPoint
72 // and accumulates the total amplitude and the multiplicity
75 fEnergyList = new Float_t[fMaxDigit];
77 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
79 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
80 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
83 for ( index = 0 ; index < fMulDigit ; index++ ){
84 tempo[index] = fDigitsList[index] ;
85 tempoE[index] = fEnergyList[index] ;
88 delete [] fDigitsList ;
89 fDigitsList = new ( Int_t[fMaxDigit] ) ;
91 delete [] fEnergyList ;
92 fEnergyList = new ( Float_t[fMaxDigit] ) ;
94 for ( index = 0 ; index < fMulDigit ; index++ ){
95 fDigitsList[index] = tempo[index] ;
96 fEnergyList[index] = tempoE[index] ;
103 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
104 fEnergyList[fMulDigit] = Energy ;
108 EvalPHOSMod(&digit) ;
111 //____________________________________________________________________________
112 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
114 // Tells if (true) or not (false) two digits are neighbors
116 Bool_t aren = kFALSE ;
118 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
119 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
122 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
125 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
127 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
128 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
130 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
136 //____________________________________________________________________________
137 Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
139 // Compares two RecPoints according to their position in the PHOS modules
141 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
142 //value (what is senseless) change as vell delta in
143 //AliPHOSTrackSegmentMakerv* and other RecPoints...
146 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
149 Int_t phosmod1 = GetPHOSMod() ;
150 Int_t phosmod2 = clu->GetPHOSMod() ;
153 GetLocalPosition(locpos1) ;
155 clu->GetLocalPosition(locpos2) ;
157 if(phosmod1 == phosmod2 ) {
158 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
163 else if(locpos1.Z()>locpos2.Z())
170 if(phosmod1 < phosmod2 )
178 //______________________________________________________________________________
179 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
182 // Execute action corresponding to one event
183 // This member function is called when a AliPHOSRecPoint is clicked with the locator
185 // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
186 // and switched off when the mouse button is released.
189 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
191 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
193 static TGraph * digitgraph = 0 ;
195 if (!gPad->IsEditable()) return;
198 TCanvas * histocanvas ;
200 TClonesArray * digits = gime->Digits() ;
205 AliPHOSDigit * digit ;
209 const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
210 Float_t * xi = new Float_t[kMulDigit] ;
211 Float_t * zi = new Float_t[kMulDigit] ;
213 // create the histogram for the single cluster
214 // 1. gets histogram boundaries
215 Float_t ximax = -999. ;
216 Float_t zimax = -999. ;
217 Float_t ximin = 999. ;
218 Float_t zimin = 999. ;
220 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
221 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
222 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
223 phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
224 if ( xi[iDigit] > ximax )
226 if ( xi[iDigit] < ximin )
228 if ( zi[iDigit] > zimax )
230 if ( zi[iDigit] < zimin )
233 ximax += phosgeom->GetCrystalSize(0) / 2. ;
234 zimax += phosgeom->GetCrystalSize(2) / 2. ;
235 ximin -= phosgeom->GetCrystalSize(0) / 2. ;
236 zimin -= phosgeom->GetCrystalSize(2) / 2. ;
237 Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
238 Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
240 // 2. gets the histogram title
243 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
249 histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
252 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
253 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
254 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
255 phosgeom->RelPosInModule(relid, x, z);
256 histo->Fill(x, z, fEnergyList[iDigit] ) ;
260 digitgraph = new TGraph(kMulDigit,xi,zi);
261 digitgraph-> SetMarkerStyle(5) ;
262 digitgraph-> SetMarkerSize(1.) ;
263 digitgraph-> SetMarkerColor(1) ;
264 digitgraph-> Paint("P") ;
268 histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
269 histocanvas->Draw() ;
270 histo->Draw("lego1") ;
288 //____________________________________________________________________________
289 void AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
291 // Calculates the dispersion of the shower at the origine of the RecPoint
299 AliPHOSDigit * digit ;
301 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
302 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
305 // Calculates the center of gravity in the local PHOS-module coordinates
309 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
310 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
314 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
315 phosgeom->RelPosInModule(relid, xi, zi);
316 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
325 // Calculates the dispersion in coordinates
327 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
328 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
332 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
333 phosgeom->RelPosInModule(relid, xi, zi);
334 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
335 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
344 fDispersion = TMath::Sqrt(d) ;
347 //______________________________________________________________________________
348 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
350 // This function calculates energy in the core,
351 // i.e. within a radius rad = 3cm around the center. Beyond this radius
352 // in accordance with shower profile the energy deposition
353 // should be less than 2%
355 Float_t coreRadius = 3 ;
360 AliPHOSDigit * digit ;
362 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
363 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
367 // Calculates the center of gravity in the local PHOS-module coordinates
369 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
370 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
374 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
375 phosgeom->RelPosInModule(relid, xi, zi);
376 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
385 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
386 digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
390 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
391 phosgeom->RelPosInModule(relid, xi, zi);
392 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
393 if(distance < coreRadius)
394 fCoreEnergy += fEnergyList[iDigit] ;
399 //____________________________________________________________________________
400 void AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
402 // Calculates the axis of the shower ellipsoid
411 AliPHOSDigit * digit ;
413 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
414 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
419 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
420 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
424 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
425 phosgeom->RelPosInModule(relid, xi, zi);
426 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
443 // //Apply correction due to non-perpendicular incidence
446 // AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
447 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
448 // Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
450 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
451 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
453 // dxx = dxx/(CosX*CosX) ;
454 // dzz = dzz/(CosZ*CosZ) ;
455 // dxz = dxz/(CosX*CosZ) ;
458 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
460 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
462 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
463 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
464 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
469 //____________________________________________________________________________
470 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
472 // Evaluates all shower parameters
474 AliPHOSRecPoint::EvalAll(logWeight,digits) ;
475 EvalLocalPosition(logWeight, digits) ;
476 EvalElipsAxis(logWeight, digits) ;
477 EvalDispersion(logWeight, digits) ;
478 EvalCoreEnergy(logWeight, digits);
480 //____________________________________________________________________________
481 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
483 // Calculates the center of gravity in the local PHOS-module coordinates
491 AliPHOSDigit * digit ;
493 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
494 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
498 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
499 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
503 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
504 phosgeom->RelPosInModule(relid, xi, zi);
505 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
515 // Correction for the depth of the shower starting point (TDR p 127)
516 Float_t para = 0.925 ;
517 Float_t parb = 6.52 ;
519 Float_t xo,yo,zo ; //Coordinates of the origin
520 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
522 Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
524 //Transform to the local ref.frame
526 xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
527 yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
529 Float_t radius = TMath::Sqrt((xoL-x)*(xoL-x)+
530 (phosgeom->GetIPtoCrystalSurface()-yoL)*(phosgeom->GetIPtoCrystalSurface()-yoL)+
533 Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ;
534 Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
536 Float_t depthx = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ;
537 Float_t depthz = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ;
540 fLocPos.SetX(x - depthx) ;
542 fLocPos.SetZ(z - depthz) ;
547 //____________________________________________________________________________
548 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
550 // Finds the maximum energy in the cluster
552 Float_t menergy = 0. ;
556 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
558 if(fEnergyList[iDigit] > menergy)
559 menergy = fEnergyList[iDigit] ;
564 //____________________________________________________________________________
565 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
567 // Calculates the multiplicity of digits with energy larger than H*energy
571 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
573 if(fEnergyList[iDigit] > H * fAmp)
579 //____________________________________________________________________________
580 Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t * maxAt, Float_t * maxAtEnergy,
581 Float_t locMaxCut,TClonesArray * digits) const
583 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
584 // energy difference between two local maxima
586 AliPHOSDigit * digit ;
587 AliPHOSDigit * digitN ;
593 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
594 maxAt[iDigit] = (Int_t) digits->At(fDigitsList[iDigit]) ;
597 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
598 if(maxAt[iDigit] != -1) {
599 digit = (AliPHOSDigit *) maxAt[iDigit] ;
601 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
602 digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ;
604 if ( AreNeighbours(digit, digitN) ) {
605 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
606 maxAt[iDigitN] = -1 ;
607 // but may be digit too is not local max ?
608 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
613 // but may be digitN too is not local max ?
614 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
615 maxAt[iDigitN] = -1 ;
617 } // if Areneighbours
623 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
624 if(maxAt[iDigit] != -1){
625 maxAt[iDigitN] = maxAt[iDigit] ;
626 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
632 //____________________________________________________________________________
633 void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits){
637 for(Int_t idig=0; idig < fMulDigit; idig++){
638 if(fEnergyList[idig] > maxE){
639 maxE = fEnergyList[idig] ;
643 fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
646 //____________________________________________________________________________
647 void AliPHOSEmcRecPoint::Print(Option_t * option)
649 // Print the list of digits belonging to the cluster
651 cout << "AliPHOSEmcRecPoint: " << endl ;
654 cout << " digits # = " ;
655 for(iDigit=0; iDigit<fMulDigit; iDigit++)
656 cout << fDigitsList[iDigit] << " " ;
659 cout << " Energies = " ;
660 for(iDigit=0; iDigit<fMulDigit; iDigit++)
661 cout << fEnergyList[iDigit] << " ";
664 cout << " Primaries " ;
665 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
666 cout << fTracksList[iDigit] << " " << endl ;
668 cout << " Multiplicity = " << fMulDigit << endl ;
669 cout << " Cluster Energy = " << fAmp << endl ;
670 cout << " Number of primaries " << fMulTrack << endl ;
671 cout << " Stored at position " << GetIndexInList() << endl ;