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()
54 fLocPos.SetX(1000000.) ; //Local position should be evaluated
58 //____________________________________________________________________________
59 AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
64 delete[] fEnergyList ;
67 //____________________________________________________________________________
68 void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
70 // Adds a digit to the RecPoint
71 // and accumulates the total amplitude and the multiplicity
74 fEnergyList = new Float_t[fMaxDigit];
76 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
78 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
79 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
82 for ( index = 0 ; index < fMulDigit ; index++ ){
83 tempo[index] = fDigitsList[index] ;
84 tempoE[index] = fEnergyList[index] ;
87 delete [] fDigitsList ;
88 fDigitsList = new ( Int_t[fMaxDigit] ) ;
90 delete [] fEnergyList ;
91 fEnergyList = new ( Float_t[fMaxDigit] ) ;
93 for ( index = 0 ; index < fMulDigit ; index++ ){
94 fDigitsList[index] = tempo[index] ;
95 fEnergyList[index] = tempoE[index] ;
102 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
103 fEnergyList[fMulDigit] = Energy ;
107 EvalPHOSMod(&digit) ;
110 //____________________________________________________________________________
111 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
113 // Tells if (true) or not (false) two digits are neighbors
115 Bool_t aren = kFALSE ;
117 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
118 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
121 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
124 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
126 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
127 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
129 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
135 //____________________________________________________________________________
136 Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
138 // Compares two RecPoints according to their position in the PHOS modules
140 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
141 //value (what is senseless) change as vell delta in
142 //AliPHOSTrackSegmentMakerv* and other RecPoints...
145 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
148 Int_t phosmod1 = GetPHOSMod() ;
149 Int_t phosmod2 = clu->GetPHOSMod() ;
152 GetLocalPosition(locpos1) ;
154 clu->GetLocalPosition(locpos2) ;
156 if(phosmod1 == phosmod2 ) {
157 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
162 else if(locpos1.Z()>locpos2.Z())
169 if(phosmod1 < phosmod2 )
178 //______________________________________________________________________________
179 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
181 // Commented by Dmitri Peressounko: there is no possibility to ensure,
182 // that AliPHOSGetter keeps the correct information.
184 // // Execute action corresponding to one event
185 // // This member function is called when a AliPHOSRecPoint is clicked with the locator
187 // // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
188 // // and switched off when the mouse button is released.
191 // // static Int_t pxold, pyold;
193 // AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
195 // static TGraph * digitgraph = 0 ;
197 // if (!gPad->IsEditable()) return;
199 // TH2F * histo = 0 ;
200 // TCanvas * histocanvas ;
204 // case kButton1Down: {
205 // AliPHOSDigit * digit ;
206 // AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
207 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
211 // const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
212 // Float_t * xi = new Float_t[kMulDigit] ;
213 // Float_t * zi = new Float_t[kMulDigit] ;
215 // // create the histogram for the single cluster
216 // // 1. gets histogram boundaries
217 // Float_t ximax = -999. ;
218 // Float_t zimax = -999. ;
219 // Float_t ximin = 999. ;
220 // Float_t zimin = 999. ;
222 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
223 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
224 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
225 // phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
226 // if ( xi[iDigit] > ximax )
227 // ximax = xi[iDigit] ;
228 // if ( xi[iDigit] < ximin )
229 // ximin = xi[iDigit] ;
230 // if ( zi[iDigit] > zimax )
231 // zimax = zi[iDigit] ;
232 // if ( zi[iDigit] < zimin )
233 // zimin = zi[iDigit] ;
235 // ximax += phosgeom->GetCrystalSize(0) / 2. ;
236 // zimax += phosgeom->GetCrystalSize(2) / 2. ;
237 // ximin -= phosgeom->GetCrystalSize(0) / 2. ;
238 // zimin -= phosgeom->GetCrystalSize(2) / 2. ;
239 // Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
240 // Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
242 // // 2. gets the histogram title
244 // Text_t title[100] ;
245 // sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
251 // histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
254 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
255 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
256 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
257 // phosgeom->RelPosInModule(relid, x, z);
258 // histo->Fill(x, z, fEnergyList[iDigit] ) ;
261 // if (!digitgraph) {
262 // digitgraph = new TGraph(kMulDigit,xi,zi);
263 // digitgraph-> SetMarkerStyle(5) ;
264 // digitgraph-> SetMarkerSize(1.) ;
265 // digitgraph-> SetMarkerColor(1) ;
266 // digitgraph-> Paint("P") ;
270 // histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
271 // histocanvas->Draw() ;
272 // histo->Draw("lego1") ;
282 // delete digitgraph ;
290 //____________________________________________________________________________
291 void AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
293 // Calculates the dispersion of the shower at the origine of the RecPoint
301 AliPHOSDigit * digit ;
303 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
304 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
307 // Calculates the center of gravity in the local PHOS-module coordinates
311 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
312 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
316 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
317 phosgeom->RelPosInModule(relid, xi, zi);
318 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
327 // Calculates the dispersion in coordinates
329 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
330 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
334 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
335 phosgeom->RelPosInModule(relid, xi, zi);
336 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
337 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
346 fDispersion = TMath::Sqrt(d) ;
349 //______________________________________________________________________________
350 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
352 // This function calculates energy in the core,
353 // i.e. within a radius rad = 3cm around the center. Beyond this radius
354 // in accordance with shower profile the energy deposition
355 // should be less than 2%
357 Float_t coreRadius = 3 ;
362 AliPHOSDigit * digit ;
364 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
365 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
369 // Calculates the center of gravity in the local PHOS-module coordinates
371 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
372 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
376 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
377 phosgeom->RelPosInModule(relid, xi, zi);
378 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
387 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
388 digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
392 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
393 phosgeom->RelPosInModule(relid, xi, zi);
394 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
395 if(distance < coreRadius)
396 fCoreEnergy += fEnergyList[iDigit] ;
401 //____________________________________________________________________________
402 void AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
404 // Calculates the axis of the shower ellipsoid
413 AliPHOSDigit * digit ;
415 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
416 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
421 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
422 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
426 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
427 phosgeom->RelPosInModule(relid, xi, zi);
428 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
445 // //Apply correction due to non-perpendicular incidence
448 // AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
449 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
450 // Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
452 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
453 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
455 // dxx = dxx/(CosX*CosX) ;
456 // dzz = dzz/(CosZ*CosZ) ;
457 // dxz = dxz/(CosX*CosZ) ;
460 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
462 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
464 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
465 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
466 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
471 //____________________________________________________________________________
472 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
474 // Evaluates all shower parameters
476 AliPHOSRecPoint::EvalAll(logWeight,digits) ;
477 EvalLocalPosition(logWeight, digits) ;
478 EvalElipsAxis(logWeight, digits) ;
479 EvalDispersion(logWeight, digits) ;
480 EvalCoreEnergy(logWeight, digits);
482 //____________________________________________________________________________
483 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
485 // Calculates the center of gravity in the local PHOS-module coordinates
493 AliPHOSDigit * digit ;
495 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
496 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
500 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
501 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
505 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
506 phosgeom->RelPosInModule(relid, xi, zi);
507 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
517 // Correction for the depth of the shower starting point (TDR p 127)
518 Float_t para = 0.925 ;
519 Float_t parb = 6.52 ;
521 Float_t xo,yo,zo ; //Coordinates of the origin
522 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
524 Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
526 //Transform to the local ref.frame
528 xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
529 yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
531 Float_t radius = TMath::Sqrt((xoL-x)*(xoL-x)+
532 (phosgeom->GetIPtoCrystalSurface()-yoL)*(phosgeom->GetIPtoCrystalSurface()-yoL)+
535 Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ;
536 Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
538 Float_t depthx = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ;
539 Float_t depthz = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ;
542 fLocPos.SetX(x - depthx) ;
544 fLocPos.SetZ(z - depthz) ;
549 //____________________________________________________________________________
550 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
552 // Finds the maximum energy in the cluster
554 Float_t menergy = 0. ;
558 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
560 if(fEnergyList[iDigit] > menergy)
561 menergy = fEnergyList[iDigit] ;
566 //____________________________________________________________________________
567 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
569 // Calculates the multiplicity of digits with energy larger than H*energy
573 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
575 if(fEnergyList[iDigit] > H * fAmp)
581 //____________________________________________________________________________
582 Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t * maxAt, Float_t * maxAtEnergy,
583 Float_t locMaxCut,TClonesArray * digits) const
585 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
586 // energy difference between two local maxima
588 AliPHOSDigit * digit ;
589 AliPHOSDigit * digitN ;
595 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
596 maxAt[iDigit] = (Int_t) digits->At(fDigitsList[iDigit]) ;
599 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
600 if(maxAt[iDigit] != -1) {
601 digit = (AliPHOSDigit *) maxAt[iDigit] ;
603 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
604 digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ;
606 if ( AreNeighbours(digit, digitN) ) {
607 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
608 maxAt[iDigitN] = -1 ;
609 // but may be digit too is not local max ?
610 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
615 // but may be digitN too is not local max ?
616 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
617 maxAt[iDigitN] = -1 ;
619 } // if Areneighbours
625 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
626 if(maxAt[iDigit] != -1){
627 maxAt[iDigitN] = maxAt[iDigit] ;
628 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
637 //____________________________________________________________________________
638 void AliPHOSEmcRecPoint::Print(Option_t * option)
640 // Print the list of digits belonging to the cluster
642 cout << "AliPHOSEmcRecPoint: " << endl ;
645 cout << " digits # = " ;
646 for(iDigit=0; iDigit<fMulDigit; iDigit++)
647 cout << fDigitsList[iDigit] << " " ;
650 cout << " Energies = " ;
651 for(iDigit=0; iDigit<fMulDigit; iDigit++)
652 cout << fEnergyList[iDigit] << " ";
655 cout << " Primaries " ;
656 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
657 cout << fTracksList[iDigit] << " " << endl ;
659 cout << " Multiplicity = " << fMulDigit << endl ;
660 cout << " Cluster Energy = " << fAmp << endl ;
661 cout << " Number of primaries " << fMulTrack << endl ;
662 cout << " Stored at position " << GetIndexInList() << endl ;