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-CPV
20 // An CpvRecPoint is a cluster of digits
21 //*-- Author: Yuri Kharlov
22 // (after Dmitri Peressounko (RRC KI & SUBATECH))
25 // --- ROOT system ---
31 // --- Standard library ---
35 // --- AliRoot header files ---
37 #include "AliPHOSGeometry.h"
38 #include "AliPHOSCpvRecPoint.h"
39 #include "AliPHOSPpsdRecPoint.h"
41 #include "AliPHOSIndexToObject.h"
43 ClassImp(AliPHOSCpvRecPoint)
45 //____________________________________________________________________________
46 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint(Float_t W0, Float_t LocMaxCut)
53 fEnergyList = new Float_t[fMaxDigit];
54 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
55 fDelta = phosgeom->GetCrystalSize(0) ;
57 fLocMaxCut = LocMaxCut ;
58 fLocPos.SetX(1000000.) ; //Local position should be evaluated
63 //____________________________________________________________________________
64 AliPHOSCpvRecPoint::~AliPHOSCpvRecPoint()
67 if ( fEnergyList ) delete[] fEnergyList ;
70 //____________________________________________________________________________
71 void AliPHOSCpvRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
73 // Adds a digit to the RecPoint
74 // and accumulates the total amplitude and the multiplicity
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 ;
108 //____________________________________________________________________________
109 Bool_t AliPHOSCpvRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 )
111 // Tells if (true) or not (false) two digits are neighbors)
113 Bool_t aren = kFALSE ;
115 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
117 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
120 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
122 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
123 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
125 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
131 //____________________________________________________________________________
132 Int_t AliPHOSCpvRecPoint::Compare(TObject * obj)
134 // Compares two RecPoints according to their position in the PHOS modules
138 if( (strcmp(obj->ClassName() , "AliPHOSPpsdRecPoint" )) == 0) // PPSD Rec Point
140 AliPHOSPpsdRecPoint * clu = (AliPHOSPpsdRecPoint *)obj ;
141 if(this->GetPHOSMod() < clu->GetPHOSMod() )
149 AliPHOSCpvRecPoint * clu = (AliPHOSCpvRecPoint *) obj ;
151 Int_t phosmod1 = this->GetPHOSMod() ;
152 Int_t phosmod2 = clu->GetPHOSMod() ;
155 this->GetLocalPosition(locpos1) ;
157 clu->GetLocalPosition(locpos2) ;
159 if(phosmod1 == phosmod2 ) {
160 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/fDelta)-(Int_t)TMath::Ceil(locpos2.X()/fDelta) ;
165 else if(locpos1.Z()>locpos2.Z())
172 if(phosmod1 < phosmod2 )
182 //______________________________________________________________________________
183 void AliPHOSCpvRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py)
185 // Execute action corresponding to one event
186 // This member function is called when a AliPHOSRecPoint is clicked with the locator
188 // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
189 // and switched off when the mouse button is released.
192 // static Int_t pxold, pyold;
194 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
196 static TGraph * digitgraph = 0 ;
198 if (!gPad->IsEditable()) return;
201 TCanvas * histocanvas ;
206 AliPHOSDigit * digit ;
207 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
211 const Int_t kMulDigit = AliPHOSCpvRecPoint::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 *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
224 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
225 phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
226 if ( xi[iDigit] > ximax )
228 if ( xi[iDigit] < ximin )
230 if ( zi[iDigit] > zimax )
232 if ( zi[iDigit] < zimin )
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
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 *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
256 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
257 phosgeom->RelPosInModule(relid, x, z);
258 histo->Fill(x, z, fEnergyList[iDigit] ) ;
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("cluser", "a single cluster", 600, 500) ;
271 histocanvas->Draw() ;
272 histo->Draw("lego1") ;
290 //____________________________________________________________________________
291 Float_t AliPHOSCpvRecPoint::GetDispersion()
293 // Calculates the dispersion of the shower at the origine of the RecPoint
295 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
301 GetLocalPosition(locpos);
302 Float_t x = locpos.X() ;
303 Float_t z = locpos.Z() ;
305 AliPHOSDigit * digit ;
306 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
309 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
310 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
314 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
315 phosgeom->RelPosInModule(relid, xi, zi);
316 Float_t w = TMath::Max(0.,fW0+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
317 d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
323 return TMath::Sqrt(d) ;
326 //____________________________________________________________________________
327 void AliPHOSCpvRecPoint::GetElipsAxis(Float_t * lambda)
329 // Calculates the axis of the shower ellipsoid
331 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
340 AliPHOSDigit * digit ;
341 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
344 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
345 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
349 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
350 phosgeom->RelPosInModule(relid, xi, zi);
351 Float_t w = TMath::Max(0.,fW0+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
369 Float_t tmp0 = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz );
370 Float_t tmp1 = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz );
371 if (tmp0>=0) lambda[0] = TMath::Sqrt(tmp0);
372 else lambda[0] = -TMath::Sqrt(-tmp0);
373 if (tmp1>=0) lambda[1] = TMath::Sqrt(tmp1);
374 else lambda[1] = -TMath::Sqrt(-tmp1);
375 // lambda[0] = TMath::Sqrt( 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ) ;
376 // lambda[1] = TMath::Sqrt( 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ) ;
379 //____________________________________________________________________________
380 void AliPHOSCpvRecPoint::GetLocalPosition(TVector3 &LPos)
382 // Calculates the center of gravity in the local PHOS-module coordinates
384 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
386 if( fLocPos.X() < 1000000.) { // already evaluated
398 AliPHOSDigit * digit ;
400 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
404 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
405 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) );
409 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
410 phosgeom->RelPosInModule(relid, xi, zi);
411 Float_t w = TMath::Max( 0., fW0 + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
423 if (fMulDigit != 0) cout << "AliPHOSCpvRecPoint: too low log weight factor "
424 << "to evaluate cluster's center\n";
433 //____________________________________________________________________________
434 void AliPHOSCpvRecPoint::GetClusterLengths(Int_t &lengX, Int_t &lengZ)
436 // Calculates the cluster lengths along X and Z axes
437 // These characteristics are needed for CPV to tune
438 // digitization+reconstruction to experimental data
439 // Yuri Kharlov. 24 October 2000
441 if( fLengX != -1 && fLengZ != -1 ) { // already evaluated
447 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
451 AliPHOSDigit * digit ;
453 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
455 const Int_t kMaxLeng=20;
456 Int_t idX[kMaxLeng], idZ[kMaxLeng];
461 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
462 digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) );
463 Int_t absId = digit->GetId();
464 phosgeom->AbsToRelNumbering(absId, relid) ;
468 for (i=0; i<lengX; i++) if (relid[2]==idX[i]) { dejavu=kTRUE; break; }
472 lengX = TMath::Min(lengX,kMaxLeng);
476 for (i=0; i<lengZ; i++) if (relid[3]==idZ[i]) { dejavu=kTRUE; break; }
480 lengZ = TMath::Min(lengZ,kMaxLeng);
487 //____________________________________________________________________________
488 Float_t AliPHOSCpvRecPoint::GetMaximalEnergy(void)
490 // Finds the maximum energy in the cluster
492 Float_t menergy = 0. ;
496 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
498 if(fEnergyList[iDigit] > menergy)
499 menergy = fEnergyList[iDigit] ;
504 //____________________________________________________________________________
505 Int_t AliPHOSCpvRecPoint::GetMultiplicityAtLevel(const Float_t H)
507 // Calculates the multiplicity of digits with energy larger than H*energy
511 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
513 if(fEnergyList[iDigit] > H * fAmp)
519 //____________________________________________________________________________
520 Int_t AliPHOSCpvRecPoint::GetNumberOfLocalMax(Int_t * maxAt, Float_t * maxAtEnergy)
522 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
523 // energy difference between two local maxima
525 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
527 AliPHOSDigit * digit ;
528 AliPHOSDigit * digitN ;
534 for(iDigit = 0; iDigit < fMulDigit; iDigit++){
535 maxAt[iDigit] = (Int_t) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
538 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
539 if(maxAt[iDigit] != -1) {
540 digit = (AliPHOSDigit *) maxAt[iDigit] ;
542 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
543 digitN = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigitN]) ) ;
545 if ( AreNeighbours(digit, digitN) ) {
546 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
547 maxAt[iDigitN] = -1 ;
548 // but may be digit too is not local max ?
549 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + fLocMaxCut)
554 // but may be digitN too is not local max ?
555 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - fLocMaxCut)
556 maxAt[iDigitN] = -1 ;
558 } // if Areneighbours
564 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
565 if(maxAt[iDigit] != -1){
566 maxAt[iDigitN] = maxAt[iDigit] ;
567 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
575 //____________________________________________________________________________
576 void AliPHOSCpvRecPoint::Print(Option_t * option)
578 // Print the list of digits belonging to the cluster
580 cout << "AliPHOSCpvRecPoint: " << endl ;
582 AliPHOSDigit * digit ;
584 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
590 AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
592 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
593 digit = please->GimeDigit( fDigitsList[iDigit] ) ;
594 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
595 phosgeom->RelPosInModule(relid, xi, zi);
596 cout << " Id = " << digit->GetId() ;
597 cout << " module = " << relid[0] ;
598 cout << " x = " << xi ;
599 cout << " z = " << zi ;
600 cout << " Energy = " << fEnergyList[iDigit] << endl ;
602 cout << " Multiplicity = " << fMulDigit << endl ;
603 cout << " Cluster Energy = " << fAmp << endl ;
604 cout << " Stored at position " << GetIndexInList() << endl ;