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()
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 )
177 //______________________________________________________________________________
178 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
181 // Execute action corresponding to one event
182 // This member function is called when a AliPHOSRecPoint is clicked with the locator
184 // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
185 // and switched off when the mouse button is released.
188 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
190 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
192 static TGraph * digitgraph = 0 ;
194 if (!gPad->IsEditable()) return;
197 TCanvas * histocanvas ;
199 TClonesArray * digits = gime->Digits() ;
204 AliPHOSDigit * digit ;
208 const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
209 Float_t * xi = new Float_t[kMulDigit] ;
210 Float_t * zi = new Float_t[kMulDigit] ;
212 // create the histogram for the single cluster
213 // 1. gets histogram boundaries
214 Float_t ximax = -999. ;
215 Float_t zimax = -999. ;
216 Float_t ximin = 999. ;
217 Float_t zimin = 999. ;
219 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
220 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
221 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
222 phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
223 if ( xi[iDigit] > ximax )
225 if ( xi[iDigit] < ximin )
227 if ( zi[iDigit] > zimax )
229 if ( zi[iDigit] < zimin )
232 ximax += phosgeom->GetCrystalSize(0) / 2. ;
233 zimax += phosgeom->GetCrystalSize(2) / 2. ;
234 ximin -= phosgeom->GetCrystalSize(0) / 2. ;
235 zimin -= phosgeom->GetCrystalSize(2) / 2. ;
236 Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
237 Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
239 // 2. gets the histogram title
242 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
248 histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
251 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
252 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
253 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
254 phosgeom->RelPosInModule(relid, x, z);
255 histo->Fill(x, z, fEnergyList[iDigit] ) ;
259 digitgraph = new TGraph(kMulDigit,xi,zi);
260 digitgraph-> SetMarkerStyle(5) ;
261 digitgraph-> SetMarkerSize(1.) ;
262 digitgraph-> SetMarkerColor(1) ;
263 digitgraph-> Paint("P") ;
267 histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
268 histocanvas->Draw() ;
269 histo->Draw("lego1") ;
287 //____________________________________________________________________________
288 void AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
290 // Calculates the dispersion of the shower at the origine of the RecPoint
298 AliPHOSDigit * digit ;
300 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
301 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
304 // Calculates the center of gravity in the local PHOS-module coordinates
308 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
309 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
313 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
314 phosgeom->RelPosInModule(relid, xi, zi);
315 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
324 // Calculates the dispersion in coordinates
326 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
327 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
331 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
332 phosgeom->RelPosInModule(relid, xi, zi);
333 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
334 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
343 fDispersion = TMath::Sqrt(d) ;
346 //______________________________________________________________________________
347 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
349 // This function calculates energy in the core,
350 // i.e. within a radius rad = 3cm around the center. Beyond this radius
351 // in accordance with shower profile the energy deposition
352 // should be less than 2%
354 Float_t coreRadius = 3 ;
359 AliPHOSDigit * digit ;
361 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
362 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
366 // Calculates the center of gravity in the local PHOS-module coordinates
368 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
369 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
373 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
374 phosgeom->RelPosInModule(relid, xi, zi);
375 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
384 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
385 digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
389 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
390 phosgeom->RelPosInModule(relid, xi, zi);
391 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
392 if(distance < coreRadius)
393 fCoreEnergy += fEnergyList[iDigit] ;
398 //____________________________________________________________________________
399 void AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
401 // Calculates the axis of the shower ellipsoid
410 AliPHOSDigit * digit ;
412 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
413 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
418 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
419 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
423 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
424 phosgeom->RelPosInModule(relid, xi, zi);
425 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
442 // //Apply correction due to non-perpendicular incidence
445 // AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
446 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
447 // Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
449 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
450 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
452 // dxx = dxx/(CosX*CosX) ;
453 // dzz = dzz/(CosZ*CosZ) ;
454 // dxz = dxz/(CosX*CosZ) ;
457 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
459 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
461 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
462 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
463 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
468 //____________________________________________________________________________
469 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
471 // Evaluates all shower parameters
473 AliPHOSRecPoint::EvalAll(logWeight,digits) ;
474 EvalLocalPosition(logWeight, digits) ;
475 EvalElipsAxis(logWeight, digits) ;
476 EvalDispersion(logWeight, digits) ;
477 EvalCoreEnergy(logWeight, digits);
479 //____________________________________________________________________________
480 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
482 // Calculates the center of gravity in the local PHOS-module coordinates
490 AliPHOSDigit * digit ;
492 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
493 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
497 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
498 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
502 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
503 phosgeom->RelPosInModule(relid, xi, zi);
504 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
514 // Correction for the depth of the shower starting point (TDR p 127)
515 Float_t para = 0.925 ;
516 Float_t parb = 6.52 ;
518 Float_t xo,yo,zo ; //Coordinates of the origin
519 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
521 Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
523 //Transform to the local ref.frame
525 xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
526 yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
528 Float_t radius = TMath::Sqrt((xoL-x)*(xoL-x)+
529 (phosgeom->GetIPtoCrystalSurface()-yoL)*(phosgeom->GetIPtoCrystalSurface()-yoL)+
532 Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ;
533 Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
535 Float_t depthx = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ;
536 Float_t depthz = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ;
539 fLocPos.SetX(x - depthx) ;
541 fLocPos.SetZ(z - depthz) ;
545 //____________________________________________________________________________
546 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
548 // Finds the maximum energy in the cluster
550 Float_t menergy = 0. ;
554 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
556 if(fEnergyList[iDigit] > menergy)
557 menergy = fEnergyList[iDigit] ;
562 //____________________________________________________________________________
563 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
565 // Calculates the multiplicity of digits with energy larger than H*energy
569 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
571 if(fEnergyList[iDigit] > H * fAmp)
577 //____________________________________________________________________________
578 Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t * maxAt, Float_t * maxAtEnergy,
579 Float_t locMaxCut,TClonesArray * digits) const
581 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
582 // energy difference between two local maxima
584 AliPHOSDigit * digit ;
585 AliPHOSDigit * digitN ;
591 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
592 maxAt[iDigit] = (Int_t) digits->At(fDigitsList[iDigit]) ;
595 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
596 if(maxAt[iDigit] != -1) {
597 digit = (AliPHOSDigit *) maxAt[iDigit] ;
599 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
600 digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ;
602 if ( AreNeighbours(digit, digitN) ) {
603 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
604 maxAt[iDigitN] = -1 ;
605 // but may be digit too is not local max ?
606 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
611 // but may be digitN too is not local max ?
612 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
613 maxAt[iDigitN] = -1 ;
615 } // if Areneighbours
621 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
622 if(maxAt[iDigit] != -1){
623 maxAt[iDigitN] = maxAt[iDigit] ;
624 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
630 //____________________________________________________________________________
631 void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits){
635 for(Int_t idig=0; idig < fMulDigit; idig++){
636 if(fEnergyList[idig] > maxE){
637 maxE = fEnergyList[idig] ;
641 fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
644 //____________________________________________________________________________
645 void AliPHOSEmcRecPoint::Print(Option_t * option)
647 // Print the list of digits belonging to the cluster
649 cout << "AliPHOSEmcRecPoint: " << endl ;
652 cout << " digits # = " ;
653 for(iDigit=0; iDigit<fMulDigit; iDigit++)
654 cout << fDigitsList[iDigit] << " " ;
657 cout << " Energies = " ;
658 for(iDigit=0; iDigit<fMulDigit; iDigit++)
659 cout << fEnergyList[iDigit] << " ";
662 cout << " Primaries " ;
663 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
664 cout << fTracksList[iDigit] << " " << endl ;
666 cout << " Multiplicity = " << fMulDigit << endl ;
667 cout << " Cluster Energy = " << fAmp << endl ;
668 cout << " Number of primaries " << fMulTrack << endl ;
669 cout << " Stored at position " << GetIndexInList() << endl ;