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 "AliPHOSGeometry.h"
38 #include "AliPHOSEmcRecPoint.h"
40 #include "AliPHOSIndexToObject.h"
42 ClassImp(AliPHOSEmcRecPoint)
44 //____________________________________________________________________________
45 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(Float_t W0, Float_t LocMaxCut)
52 fEnergyList = new Float_t[fMaxDigit];
53 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
54 fDelta = phosgeom->GetCrystalSize(0) ;
56 fLocMaxCut = LocMaxCut ;
57 fLocPos.SetX(1000000.) ; //Local position should be evaluated
61 //____________________________________________________________________________
62 AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
66 delete[] fEnergyList ;
69 //____________________________________________________________________________
70 void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
72 // Adds a digit to the RecPoint
73 // and accumulates the total amplitude and the multiplicity
75 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
77 Int_t * tempo = new ( Int_t[fMaxDigit] ) ;
78 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
81 for ( index = 0 ; index < fMulDigit ; index++ ){
82 tempo[index] = fDigitsList[index] ;
83 tempoE[index] = fEnergyList[index] ;
86 delete [] fDigitsList ;
87 fDigitsList = new ( Int_t[fMaxDigit] ) ;
89 delete [] fEnergyList ;
90 fEnergyList = new ( Float_t[fMaxDigit] ) ;
92 for ( index = 0 ; index < fMulDigit ; index++ ){
93 fDigitsList[index] = tempo[index] ;
94 fEnergyList[index] = tempoE[index] ;
101 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
102 fEnergyList[fMulDigit] = Energy ;
107 //____________________________________________________________________________
108 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 )
110 // Tells if (true) or not (false) two digits are neighbors)
112 Bool_t aren = kFALSE ;
114 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
116 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
119 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
121 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
122 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
124 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
130 //____________________________________________________________________________
131 Int_t AliPHOSEmcRecPoint::Compare(TObject * obj)
133 // Compares two RecPoints according to their position in the PHOS modules
137 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
140 Int_t phosmod1 = this->GetPHOSMod() ;
141 Int_t phosmod2 = clu->GetPHOSMod() ;
144 this->GetLocalPosition(locpos1) ;
146 clu->GetLocalPosition(locpos2) ;
148 if(phosmod1 == phosmod2 ) {
149 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/fDelta)-(Int_t)TMath::Ceil(locpos2.X()/fDelta) ;
154 else if(locpos1.Z()>locpos2.Z())
161 if(phosmod1 < phosmod2 )
170 //______________________________________________________________________________
171 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py)
173 // Execute action corresponding to one event
174 // This member function is called when a AliPHOSRecPoint is clicked with the locator
176 // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
177 // and switched off when the mouse button is released.
180 // static Int_t pxold, pyold;
182 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
184 static TGraph * digitgraph = 0 ;
186 if (!gPad->IsEditable()) return;
189 TCanvas * histocanvas ;
194 AliPHOSDigit * digit ;
195 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
199 const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
200 Float_t * xi = new Float_t[kMulDigit] ;
201 Float_t * zi = new Float_t[kMulDigit] ;
203 // create the histogram for the single cluster
204 // 1. gets histogram boundaries
205 Float_t ximax = -999. ;
206 Float_t zimax = -999. ;
207 Float_t ximin = 999. ;
208 Float_t zimin = 999. ;
210 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
211 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
212 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
213 phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
214 if ( xi[iDigit] > ximax )
216 if ( xi[iDigit] < ximin )
218 if ( zi[iDigit] > zimax )
220 if ( zi[iDigit] < zimin )
223 ximax += phosgeom->GetCrystalSize(0) / 2. ;
224 zimax += phosgeom->GetCrystalSize(2) / 2. ;
225 ximin -= phosgeom->GetCrystalSize(0) / 2. ;
226 zimin -= phosgeom->GetCrystalSize(2) / 2. ;
227 Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
228 Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
230 // 2. gets the histogram title
233 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
239 histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
242 for(iDigit=0; iDigit<kMulDigit; iDigit++) {
243 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
244 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
245 phosgeom->RelPosInModule(relid, x, z);
246 histo->Fill(x, z, fEnergyList[iDigit] ) ;
250 digitgraph = new TGraph(kMulDigit,xi,zi);
251 digitgraph-> SetMarkerStyle(5) ;
252 digitgraph-> SetMarkerSize(1.) ;
253 digitgraph-> SetMarkerColor(1) ;
254 digitgraph-> Paint("P") ;
258 histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ;
259 histocanvas->Draw() ;
260 histo->Draw("lego1") ;
278 //____________________________________________________________________________
279 Float_t AliPHOSEmcRecPoint::GetDispersion()
281 // Calculates the dispersion of the shower at the origine of the RecPoint
283 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
289 GetLocalPosition(locpos);
290 Float_t x = locpos.X() ;
291 Float_t z = locpos.Z() ;
293 AliPHOSDigit * digit ;
294 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
297 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
298 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
302 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
303 phosgeom->RelPosInModule(relid, xi, zi);
304 Float_t w = TMath::Max(0.,fW0+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
305 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
311 return TMath::Sqrt(d) ;
314 //____________________________________________________________________________
315 void AliPHOSEmcRecPoint::GetElipsAxis(Float_t * lambda)
317 // Calculates the axis of the shower ellipsoid
319 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
328 AliPHOSDigit * digit ;
329 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
332 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
333 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
337 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
338 phosgeom->RelPosInModule(relid, xi, zi);
339 Double_t w = TMath::Max(0.,fW0+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
356 // //Apply correction due to non-perpendicular incidence
359 // Double_t DistanceToIP= (Double_t ) ((AliPHOSGeometry *) fGeom)->GetIPtoCrystalSurface() ;
361 // CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
362 // CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
364 // dxx = dxx/(CosX*CosX) ;
365 // dzz = dzz/(CosZ*CosZ) ;
366 // dxz = dxz/(CosX*CosZ) ;
369 lambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
371 lambda[0] = TMath::Sqrt(lambda[0]) ;
373 lambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
374 if(lambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
375 lambda[1] = TMath::Sqrt(lambda[1]) ;
380 //____________________________________________________________________________
381 void AliPHOSEmcRecPoint::GetLocalPosition(TVector3 &LPos)
383 // Calculates the center of gravity in the local PHOS-module coordinates
385 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
387 if( fLocPos.X() < 1000000.) { // already evaluated
399 AliPHOSDigit * digit ;
401 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
406 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
407 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) );
411 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
412 phosgeom->RelPosInModule(relid, xi, zi);
413 Float_t w = TMath::Max( 0., fW0 + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
429 //____________________________________________________________________________
430 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void)
432 // Finds the maximum energy in the cluster
434 Float_t menergy = 0. ;
438 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
440 if(fEnergyList[iDigit] > menergy)
441 menergy = fEnergyList[iDigit] ;
446 //____________________________________________________________________________
447 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H)
449 // Calculates the multiplicity of digits with energy larger than H*energy
453 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
455 if(fEnergyList[iDigit] > H * fAmp)
461 //____________________________________________________________________________
462 Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t * maxAt, Float_t * maxAtEnergy)
464 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
465 // energy difference between two local maxima
467 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
469 AliPHOSDigit * digit ;
470 AliPHOSDigit * digitN ;
476 for(iDigit = 0; iDigit < fMulDigit; iDigit++){
477 maxAt[iDigit] = (Int_t) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
480 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
481 if(maxAt[iDigit] != -1) {
482 digit = (AliPHOSDigit *) maxAt[iDigit] ;
484 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
485 digitN = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigitN]) ) ;
487 if ( AreNeighbours(digit, digitN) ) {
488 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
489 maxAt[iDigitN] = -1 ;
490 // but may be digit too is not local max ?
491 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + fLocMaxCut)
496 // but may be digitN too is not local max ?
497 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - fLocMaxCut)
498 maxAt[iDigitN] = -1 ;
500 } // if Areneighbours
506 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
507 if(maxAt[iDigit] != -1){
508 maxAt[iDigitN] = maxAt[iDigit] ;
509 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
517 // //____________________________________________________________________________
518 // AliPHOSEmcRecPoint& AliPHOSEmcRecPoint::operator = (AliPHOSEmcRecPoint Clu)
520 // int * dl = Clu.GetDigitsList() ;
523 // delete fDigitsList ;
525 // AliPHOSDigit * digit ;
529 // for(iDigit=0; iDigit<fMulDigit; iDigit++) {
530 // digit = (AliPHOSDigit *) dl[iDigit];
531 // AddDigit(*digit) ;
534 // fAmp = Clu.GetTotalEnergy() ;
535 // fGeom = Clu.GetGeom() ;
537 // Clu.GetLocalPosition(locpos) ;
539 // fMulDigit = Clu.GetMultiplicity() ;
540 // fMaxDigit = Clu.GetMaximumMultiplicity() ;
541 // fPHOSMod = Clu.GetPHOSMod() ;
542 // fW0 = Clu.GetLogWeightCut() ;
543 // fDelta = Clu.GetDelta() ;
544 // fLocMaxCut = Clu.GetLocMaxCut() ;
551 //____________________________________________________________________________
552 void AliPHOSEmcRecPoint::Print(Option_t * option)
554 // Print the list of digits belonging to the cluster
556 cout << "AliPHOSEmcRecPoint: " << endl ;
558 AliPHOSDigit * digit ;
560 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
566 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
568 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
569 digit = please->GimeDigit( fDigitsList[iDigit] ) ;
570 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
571 phosgeom->RelPosInModule(relid, xi, zi);
572 cout << " Id = " << digit->GetId() ;
573 cout << " module = " << relid[0] ;
574 cout << " x = " << xi ;
575 cout << " z = " << zi ;
576 cout << " Energy = " << fEnergyList[iDigit] << endl ;
578 cout << " Multiplicity = " << fMulDigit << endl ;
579 cout << " Cluster Energy = " << fAmp << endl ;
580 cout << " Stored at position " << GetIndexInList() << endl ;
583 //______________________________________________________________________________
584 // void AliPHOSEmcRecPoint::Streamer(TBuffer &R__b)
586 // // Stream an object of class AliPHOSEmcRecPoint.
588 // if (R__b.IsReading()) {
589 // Version_t R__v = R__b.ReadVersion(); if (R__v) { }
590 // AliPHOSRecPoint::Streamer(R__b);
592 // fEnergyList = new Float_t[fMulDigit] ;
593 // R__b.ReadFastArray(fEnergyList, fMulDigit);
594 // R__b >> fLocMaxCut;
597 // R__b.WriteVersion(AliPHOSEmcRecPoint::IsA());
598 // AliPHOSRecPoint::Streamer(R__b);
600 // R__b.WriteFastArray(fEnergyList, fMulDigit);
601 // R__b << fLocMaxCut;