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 ---
32 // --- Standard library ---
34 // --- AliRoot header files ---
36 #include "AliGenerator.h"
37 #include "AliPHOSGeometry.h"
38 #include "AliPHOSEmcRecPoint.h"
40 #include "AliPHOSGetter.h"
42 ClassImp(AliPHOSEmcRecPoint)
44 //____________________________________________________________________________
45 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
53 fNExMax = 0 ; //Not unfolded yet
55 fLocPos.SetX(1000000.) ; //Local position should be evaluated
59 //____________________________________________________________________________
60 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const char * opt) : AliPHOSRecPoint(opt)
66 fNExMax = 0 ; //Not unfolded yet
70 fLocPos.SetX(1000000.) ; //Local position should be evaluated
74 //____________________________________________________________________________
75 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const AliPHOSEmcRecPoint & rp) : AliPHOSRecPoint(rp)
79 fMulDigit = rp.fMulDigit ;
81 fCoreEnergy = rp.fCoreEnergy ;
82 fEnergyList = new Float_t[rp.fMulDigit] ;
84 for(index = 0 ; index < fMulDigit ; index++)
85 fEnergyList[index] = rp.fEnergyList[index] ;
86 fNExMax = rp.fNExMax ;
90 //____________________________________________________________________________
91 AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
96 delete[] fEnergyList ;
99 //____________________________________________________________________________
100 void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
102 // Adds a digit to the RecPoint
103 // and accumulates the total amplitude and the multiplicity
106 fEnergyList = new Float_t[fMaxDigit];
108 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
110 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
111 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
114 for ( index = 0 ; index < fMulDigit ; index++ ){
115 tempo[index] = fDigitsList[index] ;
116 tempoE[index] = fEnergyList[index] ;
119 delete [] fDigitsList ;
120 fDigitsList = new ( Int_t[fMaxDigit] ) ;
122 delete [] fEnergyList ;
123 fEnergyList = new ( Float_t[fMaxDigit] ) ;
125 for ( index = 0 ; index < fMulDigit ; index++ ){
126 fDigitsList[index] = tempo[index] ;
127 fEnergyList[index] = tempoE[index] ;
134 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
135 fEnergyList[fMulDigit] = Energy ;
139 EvalPHOSMod(&digit) ;
142 //____________________________________________________________________________
143 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
145 // Tells if (true) or not (false) two digits are neighbors
147 Bool_t aren = kFALSE ;
149 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
150 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
153 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
156 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
158 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
159 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
161 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
167 //____________________________________________________________________________
168 Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
170 // Compares two RecPoints according to their position in the PHOS modules
172 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
173 //value (what is senseless) change as vell delta in
174 //AliPHOSTrackSegmentMakerv* and other RecPoints...
177 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
180 Int_t phosmod1 = GetPHOSMod() ;
181 Int_t phosmod2 = clu->GetPHOSMod() ;
184 GetLocalPosition(locpos1) ;
186 clu->GetLocalPosition(locpos2) ;
188 if(phosmod1 == phosmod2 ) {
189 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
194 else if(locpos1.Z()>locpos2.Z())
201 if(phosmod1 < phosmod2 )
209 //______________________________________________________________________________
210 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
213 // Execute action corresponding to one event
214 // This member function is called when a AliPHOSRecPoint is clicked with the locator
216 // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
217 // and switched off when the mouse button is released.
220 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
222 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
224 static TGraph * digitgraph = 0 ;
226 if (!gPad->IsEditable()) return;
229 TCanvas * histocanvas ;
231 const TClonesArray * digits = gime->Digits() ;
236 AliPHOSDigit * digit ;
240 const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
241 Float_t * xi = new Float_t[kMulDigit] ;
242 Float_t * zi = new Float_t[kMulDigit] ;
244 // create the histogram for the single cluster
245 // 1. gets histogram boundaries
246 Float_t ximax = -999. ;
247 Float_t zimax = -999. ;
248 Float_t ximin = 999. ;
249 Float_t zimin = 999. ;
251 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
252 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
253 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
254 phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
255 if ( xi[iDigit] > ximax )
257 if ( xi[iDigit] < ximin )
259 if ( zi[iDigit] > zimax )
261 if ( zi[iDigit] < zimin )
264 ximax += phosgeom->GetCrystalSize(0) / 2. ;
265 zimax += phosgeom->GetCrystalSize(2) / 2. ;
266 ximin -= phosgeom->GetCrystalSize(0) / 2. ;
267 zimin -= phosgeom->GetCrystalSize(2) / 2. ;
268 Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
269 Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
271 // 2. gets the histogram title
274 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
280 histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
283 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
284 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
285 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
286 phosgeom->RelPosInModule(relid, x, z);
287 histo->Fill(x, z, fEnergyList[iDigit] ) ;
291 digitgraph = new TGraph(kMulDigit,xi,zi);
292 digitgraph-> SetMarkerStyle(5) ;
293 digitgraph-> SetMarkerSize(1.) ;
294 digitgraph-> SetMarkerColor(1) ;
295 digitgraph-> Paint("P") ;
299 histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
300 histocanvas->Draw() ;
301 histo->Draw("lego1") ;
319 //____________________________________________________________________________
320 void AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
322 // Calculates the dispersion of the shower at the origine of the RecPoint
330 AliPHOSDigit * digit ;
332 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
333 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
336 // Calculates the center of gravity in the local PHOS-module coordinates
340 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
341 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
345 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
346 phosgeom->RelPosInModule(relid, xi, zi);
347 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
356 // Calculates the dispersion in coordinates
358 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
359 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
363 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
364 phosgeom->RelPosInModule(relid, xi, zi);
365 Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
366 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
375 fDispersion = TMath::Sqrt(d) ;
378 //______________________________________________________________________________
379 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
381 // This function calculates energy in the core,
382 // i.e. within a radius rad = 3cm around the center. Beyond this radius
383 // in accordance with shower profile the energy deposition
384 // should be less than 2%
386 Float_t coreRadius = 3 ;
391 AliPHOSDigit * digit ;
393 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
394 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
398 // Calculates the center of gravity in the local PHOS-module coordinates
400 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
401 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
405 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
406 phosgeom->RelPosInModule(relid, xi, zi);
407 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
416 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
417 digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
421 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
422 phosgeom->RelPosInModule(relid, xi, zi);
423 Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
424 if(distance < coreRadius)
425 fCoreEnergy += fEnergyList[iDigit] ;
430 //____________________________________________________________________________
431 void AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
433 // Calculates the axis of the shower ellipsoid
442 AliPHOSDigit * digit ;
444 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
445 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
450 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
451 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
455 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
456 phosgeom->RelPosInModule(relid, xi, zi);
457 Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
474 // //Apply correction due to non-perpendicular incidence
477 // AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
478 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
479 // Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
481 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
482 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
484 // dxx = dxx/(CosX*CosX) ;
485 // dzz = dzz/(CosZ*CosZ) ;
486 // dxz = dxz/(CosX*CosZ) ;
489 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
491 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
493 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
494 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
495 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
500 //____________________________________________________________________________
501 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
503 // Evaluates all shower parameters
505 EvalLocalPosition(logWeight, digits) ;
506 EvalElipsAxis(logWeight, digits) ;
507 EvalDispersion(logWeight, digits) ;
508 EvalCoreEnergy(logWeight, digits);
510 AliPHOSRecPoint::EvalAll(logWeight,digits) ;
512 //____________________________________________________________________________
513 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
515 // Calculates the center of gravity in the local PHOS-module coordinates
523 AliPHOSDigit * digit ;
525 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
526 AliPHOSGeometry * phosgeom = (AliPHOSGeometry*)gime->PHOSGeometry();
530 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
531 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
535 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
536 phosgeom->RelPosInModule(relid, xi, zi);
537 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
547 // Correction for the depth of the shower starting point (TDR p 127)
548 Float_t para = 0.925 ;
549 Float_t parb = 6.52 ;
551 Float_t xo,yo,zo ; //Coordinates of the origin
552 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
554 Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
556 //Transform to the local ref.frame
558 xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
559 yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
561 Float_t radius = phosgeom->GetIPtoCrystalSurface()-yoL;
563 Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ;
564 Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
566 Float_t depthx = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ;
567 Float_t depthz = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ;
569 fLocPos.SetX(x - depthx) ;
571 fLocPos.SetZ(z - depthz) ;
576 //____________________________________________________________________________
577 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
579 // Finds the maximum energy in the cluster
581 Float_t menergy = 0. ;
585 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
587 if(fEnergyList[iDigit] > menergy)
588 menergy = fEnergyList[iDigit] ;
593 //____________________________________________________________________________
594 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
596 // Calculates the multiplicity of digits with energy larger than H*energy
600 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
602 if(fEnergyList[iDigit] > H * fAmp)
608 //____________________________________________________________________________
609 Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax( AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
610 Float_t locMaxCut,TClonesArray * digits) const
612 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
613 // energy difference between two local maxima
615 AliPHOSDigit * digit ;
616 AliPHOSDigit * digitN ;
622 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
623 maxAt[iDigit] = (AliPHOSDigit*) digits->At(fDigitsList[iDigit]) ;
626 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
628 digit = maxAt[iDigit] ;
630 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
631 digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ;
633 if ( AreNeighbours(digit, digitN) ) {
634 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
636 // but may be digit too is not local max ?
637 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
642 // but may be digitN too is not local max ?
643 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
646 } // if Areneighbours
652 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
654 maxAt[iDigitN] = maxAt[iDigit] ;
655 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
661 //____________________________________________________________________________
662 void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits)
664 // Define a rec.point time as a time in the cell with the maximum energy
668 for(Int_t idig=0; idig < fMulDigit; idig++){
669 if(fEnergyList[idig] > maxE){
670 maxE = fEnergyList[idig] ;
674 fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
677 //____________________________________________________________________________
678 void AliPHOSEmcRecPoint::Purify(Float_t threshold){
679 //Removes digits below threshold
681 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
682 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
685 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
686 if(fEnergyList[iDigit] > threshold){
687 tempo[mult] = fDigitsList[iDigit] ;
688 tempoE[mult] = fEnergyList[iDigit] ;
694 delete [] fDigitsList ;
695 delete [] fEnergyList ;
696 fDigitsList = new (Int_t[fMulDigit]) ;
697 fEnergyList = new ( Float_t[fMulDigit]) ;
699 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
700 fDigitsList[iDigit] = tempo[iDigit];
701 fEnergyList[iDigit] = tempoE[iDigit] ;
708 //____________________________________________________________________________
709 void AliPHOSEmcRecPoint::Print(Option_t * option) const
711 // Print the list of digits belonging to the cluster
714 message = "AliPHOSEmcRecPoint:\n" ;
715 message += " digits # = " ;
716 Info("Print", message.Data()) ;
719 for(iDigit=0; iDigit<fMulDigit; iDigit++)
720 printf(" %d ", fDigitsList[iDigit] ) ;
722 Info("Print", " Energies = ") ;
723 for(iDigit=0; iDigit<fMulDigit; iDigit++)
724 printf(" %f ", fEnergyList[iDigit] ) ;
726 Info("Print", " Primaries ") ;
727 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
728 printf(" %d ", fTracksList[iDigit]) ;
730 message = " Multiplicity = %d" ;
731 message += " Cluster Energy = %f" ;
732 message += " Number of primaries %d" ;
733 message += " Stored at position %d" ;
735 Info("Print", message.Data(), fMulDigit, fAmp, fMulTrack,GetIndexInList() ) ;