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 ---
33 // --- AliRoot header files ---
35 #include "AliGenerator.h"
36 #include "AliPHOSGeometry.h"
37 #include "AliPHOSEmcRecPoint.h"
39 #include "AliPHOSGetter.h"
41 ClassImp(AliPHOSEmcRecPoint)
43 //____________________________________________________________________________
44 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
53 fLocPos.SetX(1000000.) ; //Local position should be evaluated
57 //____________________________________________________________________________
58 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const char * opt) : AliPHOSRecPoint(opt)
67 fLocPos.SetX(1000000.) ; //Local position should be evaluated
71 //____________________________________________________________________________
72 AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
77 delete[] fEnergyList ;
80 //____________________________________________________________________________
81 void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
83 // Adds a digit to the RecPoint
84 // and accumulates the total amplitude and the multiplicity
87 fEnergyList = new Float_t[fMaxDigit];
89 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
91 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
92 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
95 for ( index = 0 ; index < fMulDigit ; index++ ){
96 tempo[index] = fDigitsList[index] ;
97 tempoE[index] = fEnergyList[index] ;
100 delete [] fDigitsList ;
101 fDigitsList = new ( Int_t[fMaxDigit] ) ;
103 delete [] fEnergyList ;
104 fEnergyList = new ( Float_t[fMaxDigit] ) ;
106 for ( index = 0 ; index < fMulDigit ; index++ ){
107 fDigitsList[index] = tempo[index] ;
108 fEnergyList[index] = tempoE[index] ;
115 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
116 fEnergyList[fMulDigit] = Energy ;
120 EvalPHOSMod(&digit) ;
123 //____________________________________________________________________________
124 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
126 // Tells if (true) or not (false) two digits are neighbors
128 Bool_t aren = kFALSE ;
130 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
131 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
134 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
137 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
139 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
140 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
142 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
148 //____________________________________________________________________________
149 Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
151 // Compares two RecPoints according to their position in the PHOS modules
153 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
154 //value (what is senseless) change as vell delta in
155 //AliPHOSTrackSegmentMakerv* and other RecPoints...
158 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
161 Int_t phosmod1 = GetPHOSMod() ;
162 Int_t phosmod2 = clu->GetPHOSMod() ;
165 GetLocalPosition(locpos1) ;
167 clu->GetLocalPosition(locpos2) ;
169 if(phosmod1 == phosmod2 ) {
170 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
175 else if(locpos1.Z()>locpos2.Z())
182 if(phosmod1 < phosmod2 )
190 //______________________________________________________________________________
191 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
194 // Execute action corresponding to one event
195 // This member function is called when a AliPHOSRecPoint is clicked with the locator
197 // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
198 // and switched off when the mouse button is released.
201 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
203 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
205 static TGraph * digitgraph = 0 ;
207 if (!gPad->IsEditable()) return;
210 TCanvas * histocanvas ;
212 const TClonesArray * digits = gime->Digits() ;
217 AliPHOSDigit * digit ;
221 const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
222 Float_t * xi = new Float_t[kMulDigit] ;
223 Float_t * zi = new Float_t[kMulDigit] ;
225 // create the histogram for the single cluster
226 // 1. gets histogram boundaries
227 Float_t ximax = -999. ;
228 Float_t zimax = -999. ;
229 Float_t ximin = 999. ;
230 Float_t zimin = 999. ;
232 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
233 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
234 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
235 phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
236 if ( xi[iDigit] > ximax )
238 if ( xi[iDigit] < ximin )
240 if ( zi[iDigit] > zimax )
242 if ( zi[iDigit] < zimin )
245 ximax += phosgeom->GetCrystalSize(0) / 2. ;
246 zimax += phosgeom->GetCrystalSize(2) / 2. ;
247 ximin -= phosgeom->GetCrystalSize(0) / 2. ;
248 zimin -= phosgeom->GetCrystalSize(2) / 2. ;
249 Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
250 Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
252 // 2. gets the histogram title
255 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
261 histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
264 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
265 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
266 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
267 phosgeom->RelPosInModule(relid, x, z);
268 histo->Fill(x, z, fEnergyList[iDigit] ) ;
272 digitgraph = new TGraph(kMulDigit,xi,zi);
273 digitgraph-> SetMarkerStyle(5) ;
274 digitgraph-> SetMarkerSize(1.) ;
275 digitgraph-> SetMarkerColor(1) ;
276 digitgraph-> Paint("P") ;
280 histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
281 histocanvas->Draw() ;
282 histo->Draw("lego1") ;
300 //____________________________________________________________________________
301 void AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
303 // Calculates the dispersion of the shower at the origine of the RecPoint
311 AliPHOSDigit * digit ;
313 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
314 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
317 // Calculates the center of gravity in the local PHOS-module coordinates
321 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
322 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
326 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
327 phosgeom->RelPosInModule(relid, xi, zi);
328 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
337 // Calculates the dispersion in coordinates
339 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
340 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
344 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
345 phosgeom->RelPosInModule(relid, xi, zi);
346 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
347 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
356 fDispersion = TMath::Sqrt(d) ;
359 //______________________________________________________________________________
360 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
362 // This function calculates energy in the core,
363 // i.e. within a radius rad = 3cm around the center. Beyond this radius
364 // in accordance with shower profile the energy deposition
365 // should be less than 2%
367 Float_t coreRadius = 3 ;
372 AliPHOSDigit * digit ;
374 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
375 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
379 // Calculates the center of gravity in the local PHOS-module coordinates
381 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
382 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
386 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
387 phosgeom->RelPosInModule(relid, xi, zi);
388 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
397 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
398 digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
402 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
403 phosgeom->RelPosInModule(relid, xi, zi);
404 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
405 if(distance < coreRadius)
406 fCoreEnergy += fEnergyList[iDigit] ;
411 //____________________________________________________________________________
412 void AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
414 // Calculates the axis of the shower ellipsoid
423 AliPHOSDigit * digit ;
425 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
426 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
431 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
432 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
436 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
437 phosgeom->RelPosInModule(relid, xi, zi);
438 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
455 // //Apply correction due to non-perpendicular incidence
458 // AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
459 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
460 // Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
462 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
463 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
465 // dxx = dxx/(CosX*CosX) ;
466 // dzz = dzz/(CosZ*CosZ) ;
467 // dxz = dxz/(CosX*CosZ) ;
470 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
472 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
474 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
475 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
476 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
481 //____________________________________________________________________________
482 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
484 // Evaluates all shower parameters
486 AliPHOSRecPoint::EvalAll(logWeight,digits) ;
487 EvalLocalPosition(logWeight, digits) ;
488 EvalElipsAxis(logWeight, digits) ;
489 EvalDispersion(logWeight, digits) ;
490 EvalCoreEnergy(logWeight, digits);
493 //____________________________________________________________________________
494 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
496 // Calculates the center of gravity in the local PHOS-module coordinates
504 AliPHOSDigit * digit ;
506 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
507 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
511 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
512 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
516 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
517 phosgeom->RelPosInModule(relid, xi, zi);
518 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
528 // Correction for the depth of the shower starting point (TDR p 127)
529 Float_t para = 0.925 ;
530 Float_t parb = 6.52 ;
532 Float_t xo,yo,zo ; //Coordinates of the origin
533 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
535 Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
537 //Transform to the local ref.frame
539 xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
540 yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
542 Float_t radius = phosgeom->GetIPtoCrystalSurface()-yoL;
544 Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ;
545 Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
547 Float_t depthx = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ;
548 Float_t depthz = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ;
550 fLocPos.SetX(x - depthx) ;
552 fLocPos.SetZ(z - depthz) ;
557 //____________________________________________________________________________
558 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
560 // Finds the maximum energy in the cluster
562 Float_t menergy = 0. ;
566 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
568 if(fEnergyList[iDigit] > menergy)
569 menergy = fEnergyList[iDigit] ;
574 //____________________________________________________________________________
575 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
577 // Calculates the multiplicity of digits with energy larger than H*energy
581 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
583 if(fEnergyList[iDigit] > H * fAmp)
589 //____________________________________________________________________________
590 Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax( AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
591 Float_t locMaxCut,TClonesArray * digits) const
593 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
594 // energy difference between two local maxima
596 AliPHOSDigit * digit ;
597 AliPHOSDigit * digitN ;
603 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
604 maxAt[iDigit] = (AliPHOSDigit*) digits->At(fDigitsList[iDigit]) ;
607 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
609 digit = maxAt[iDigit] ;
611 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
612 digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ;
614 if ( AreNeighbours(digit, digitN) ) {
615 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
617 // but may be digit too is not local max ?
618 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
623 // but may be digitN too is not local max ?
624 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
627 } // if Areneighbours
633 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
635 maxAt[iDigitN] = maxAt[iDigit] ;
636 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
642 //____________________________________________________________________________
643 void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits){
647 for(Int_t idig=0; idig < fMulDigit; idig++){
648 if(fEnergyList[idig] > maxE){
649 maxE = fEnergyList[idig] ;
653 fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
656 //____________________________________________________________________________
657 void AliPHOSEmcRecPoint::Purify(Float_t threshold){
658 //Removes digits below threshold
660 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
661 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
664 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
665 if(fEnergyList[iDigit] > threshold){
666 tempo[mult] = fDigitsList[iDigit] ;
667 tempoE[mult] = fEnergyList[iDigit] ;
673 delete [] fDigitsList ;
674 delete [] fEnergyList ;
675 fDigitsList = new (Int_t[fMulDigit]) ;
676 fEnergyList = new ( Float_t[fMulDigit]) ;
678 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
679 fDigitsList[iDigit] = tempo[iDigit];
680 fEnergyList[iDigit] = tempoE[iDigit] ;
687 //____________________________________________________________________________
688 void AliPHOSEmcRecPoint::Print(Option_t * option) const
690 // Print the list of digits belonging to the cluster
693 message = "AliPHOSEmcRecPoint:\n" ;
694 message += " digits # = " ;
695 Info("Print", message.Data()) ;
698 for(iDigit=0; iDigit<fMulDigit; iDigit++)
699 Info("Print", " %d ", fDigitsList[iDigit] ) ;
701 Info("Print", " Energies = ") ;
702 for(iDigit=0; iDigit<fMulDigit; iDigit++)
703 Info("Print", " %f ", fEnergyList[iDigit] ) ;
705 Info("Print", " Primaries ") ;
706 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
707 Info("Print", " %d ", fTracksList[iDigit]) ;
709 message = " Multiplicity = %d" ;
710 message += " Cluster Energy = %f" ;
711 message += " Number of primaries %d" ;
712 message += " Stored at position %d" ;
714 Info("Print", message.Data(), fMulDigit, fAmp, fMulTrack,GetIndexInList() ) ;