]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSEmcRecPoint.cxx
Improve documentation
[u/mrichter/AliRoot.git] / PHOS / AliPHOSEmcRecPoint.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 //  RecPoint implementation for PHOS-EMC 
20 //  An EmcRecPoint is a cluster of digits   
21 //           
22 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
23
24
25 // --- ROOT system ---
26 #include "TPad.h"
27 #include "TH2.h"
28 #include "TMath.h" 
29 #include "TCanvas.h" 
30
31 // --- Standard library ---
32
33 #include <iostream>
34
35 // --- AliRoot header files ---
36
37 #include "AliPHOSGeometry.h"
38 #include "AliPHOSEmcRecPoint.h"
39 #include "AliRun.h"
40
41 ClassImp(AliPHOSEmcRecPoint)
42
43 //____________________________________________________________________________
44 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(Float_t W0, Float_t LocMaxCut)
45   : AliPHOSRecPoint()
46 {
47   // ctor
48
49   fMulDigit   = 0 ;  
50   fAmp   = 0. ;   
51   fEnergyList = new Float_t[fMaxDigit]; 
52   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
53   fDelta     =  phosgeom->GetCrystalSize(0) ; 
54   fW0        = W0 ;          
55   fLocMaxCut = LocMaxCut ; 
56   fLocPos.SetX(1000000.)  ;      //Local position should be evaluated
57 }
58
59 //____________________________________________________________________________
60 void AliPHOSEmcRecPoint::AddDigit(AliDigitNew & digit, Float_t Energy)
61 {
62   // Adds a digit to the RecPoint
63   //  and accumulates the total amplitude and the multiplicity 
64   
65   if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists 
66     fMaxDigit*=2 ; 
67     int * tempo = new ( int[fMaxDigit] ) ; 
68     Float_t * tempoE =  new ( Float_t[fMaxDigit] ) ;
69
70     Int_t index ;     
71     for ( index = 0 ; index < fMulDigit ; index++ ){
72       tempo[index] = fDigitsList[index] ;
73       tempoE[index] = fEnergyList[index] ; 
74     }
75     
76     delete [] fDigitsList ; 
77     fDigitsList =  new ( int[fMaxDigit] ) ;
78  
79     delete [] fEnergyList ;
80     fEnergyList =  new ( Float_t[fMaxDigit] ) ;
81
82     for ( index = 0 ; index < fMulDigit ; index++ ){
83       fDigitsList[index] = tempo[index] ;
84       fEnergyList[index] = tempoE[index] ; 
85     }
86  
87     delete [] tempo ;
88     delete [] tempoE ; 
89   } // if
90   
91   fDigitsList[fMulDigit]   =  (int) &digit  ; 
92   fEnergyList[fMulDigit++] = Energy ;
93   fAmp += Energy ; 
94 }
95
96 //____________________________________________________________________________
97 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) 
98 {
99   // Tells if (true) or not (false) two digits are neighbors)
100   
101   Bool_t aren = kFALSE ;
102   
103   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
104   Int_t relid1[4] ; 
105   phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ; 
106
107   Int_t relid2[4] ; 
108   phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ; 
109   
110   Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
111   Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
112
113   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
114     aren = kTRUE ;
115   
116   return aren ;
117 }
118
119 //____________________________________________________________________________
120 Int_t AliPHOSEmcRecPoint::Compare(TObject * obj)
121 {
122   // Compares two RecPoints according to their position in the PHOS modules
123
124   Int_t rv ; 
125
126   AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ; 
127
128  
129   Int_t phosmod1 = this->GetPHOSMod() ;
130   Int_t phosmod2 = clu->GetPHOSMod() ;
131
132   TVector3 locpos1; 
133   this->GetLocalPosition(locpos1) ;
134   TVector3 locpos2;  
135   clu->GetLocalPosition(locpos2) ;  
136
137   if(phosmod1 == phosmod2 ) {
138     Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/fDelta)-(Int_t)TMath::Ceil(locpos2.X()/fDelta) ;
139     if (rowdif> 0) 
140       rv = -1 ;
141      else if(rowdif < 0) 
142        rv = 1 ;
143     else if(locpos1.Z()>locpos2.Z()) 
144       rv = -1 ;
145     else 
146       rv = 1 ; 
147      }
148
149   else {
150     if(phosmod1 < phosmod2 ) 
151       rv = -1 ;
152     else 
153       rv = 1 ;
154   }
155
156   return rv ; 
157 }
158
159 //______________________________________________________________________________
160 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py)
161 {
162   // Execute action corresponding to one event
163   //  This member function is called when a AliPHOSRecPoint is clicked with the locator
164   //
165   //  If Left button is clicked on AliPHOSRecPoint, the digits are switched on    
166   //  and switched off when the mouse button is released.
167   //
168
169   //   static Int_t pxold, pyold;
170
171    static TGraph *  digitgraph = 0 ;
172
173    if (!gPad->IsEditable()) return;
174
175    TH2F * histo = 0 ;
176    TCanvas * histocanvas ; 
177    
178    switch (event) {
179    
180    case kButton1Down: {
181      AliPHOSDigit * digit ;
182      AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
183      Int_t iDigit;
184      Int_t relid[4] ;
185      Float_t xi[fMulDigit] ;
186      Float_t zi[fMulDigit] ;
187
188      // create the histogram for the single cluster 
189      // 1. gets histogram boundaries
190      Float_t ximax = -999. ; 
191      Float_t zimax = -999. ; 
192      Float_t ximin = 999. ; 
193      Float_t zimin = 999. ;
194  
195      for(iDigit=0; iDigit<fMulDigit; iDigit++) {
196        digit = (AliPHOSDigit *) fDigitsList[iDigit];
197        phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
198        phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
199        if ( xi[iDigit] > ximax )
200          ximax = xi[iDigit] ; 
201        if ( xi[iDigit] < ximin )
202          ximin = xi[iDigit] ; 
203        if ( zi[iDigit] > zimax )
204          zimax = zi[iDigit] ; 
205        if ( zi[iDigit] < zimin )
206          zimin = zi[iDigit] ;     
207      }
208      ximax += phosgeom->GetCrystalSize(0) / 2. ;
209      zimax += phosgeom->GetCrystalSize(2) / 2. ;
210      ximin -= phosgeom->GetCrystalSize(0) / 2. ;
211      zimin -= phosgeom->GetCrystalSize(2) / 2. ;
212      Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5  ) ; 
213      Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
214  
215      // 2. gets the histogram title
216
217      Text_t title[100] ; 
218      sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
219   
220      if (!histo) {
221        delete histo ; 
222        histo = 0 ; 
223      }
224      histo = new TH2F("cluster3D", title,  xdim, ximin, ximax, zdim, zimin, zimax)  ;
225
226      Float_t x, z ; 
227      for(iDigit=0; iDigit<fMulDigit; iDigit++) {
228        digit = (AliPHOSDigit *) fDigitsList[iDigit];
229        phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
230        phosgeom->RelPosInModule(relid, x, z);
231        histo->Fill(x, z, fEnergyList[iDigit] ) ;
232      }
233
234      if (!digitgraph) {
235        digitgraph = new TGraph(fMulDigit,xi,zi);
236        digitgraph-> SetMarkerStyle(5) ; 
237        digitgraph-> SetMarkerSize(1.) ;
238        digitgraph-> SetMarkerColor(1) ;
239        digitgraph-> Paint("P") ;
240      }
241
242      Print() ;
243      histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ; 
244      histocanvas->Draw() ; 
245      histo->Draw("lego1") ; 
246
247      break;
248    }
249
250    case kButton1Up: 
251      if (digitgraph) {
252        delete digitgraph  ;
253        digitgraph = 0 ;
254      }
255      break;
256   
257    }
258 }
259
260 //____________________________________________________________________________
261 Float_t  AliPHOSEmcRecPoint::GetDispersion() 
262 {
263   // Calculates the dispersion of the shower at the origine of the RecPoint
264
265   Float_t d    = 0 ;
266   Float_t wtot = 0 ;
267
268   TVector3 locpos;
269   GetLocalPosition(locpos);
270   Float_t x = locpos.X() ;
271   Float_t z = locpos.Z() ;
272
273   AliPHOSDigit * digit ;
274   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
275   
276   Int_t iDigit;
277   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
278     digit = (AliPHOSDigit *) fDigitsList[iDigit];
279     Int_t relid[4] ;
280     Float_t xi ;
281     Float_t zi ;
282     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
283     phosgeom->RelPosInModule(relid, xi, zi);
284     Float_t w = TMath::Max(0.,fW0+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
285     d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ; 
286     wtot+=w ;
287   }
288
289   d /= wtot ;
290
291   return TMath::Sqrt(d) ;
292 }
293
294 //____________________________________________________________________________
295 void  AliPHOSEmcRecPoint::GetElipsAxis(Float_t * lambda)
296 {
297   // Calculates the axis of the shower ellipsoid
298   
299   Float_t wtot = 0. ;
300   Float_t x    = 0.;
301   Float_t z    = 0.;
302   Float_t dxx  = 0.;
303   Float_t dzz  = 0.;
304   Float_t dxz  = 0.;
305
306   AliPHOSDigit * digit ;
307   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
308   Int_t iDigit;
309
310   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
311     digit = (AliPHOSDigit *) fDigitsList[iDigit];
312     Int_t relid[4] ;
313     Float_t xi ;
314     Float_t zi ;
315     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
316     phosgeom->RelPosInModule(relid, xi, zi);
317     Float_t w = TMath::Max(0.,fW0+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
318     dxx  += w * xi * xi ;
319     x    += w * xi ;
320     dzz  += w * zi * zi ;
321     z    += w * zi ; 
322     dxz  += w * xi * zi ; 
323     wtot += w ;
324   }
325
326   dxx /= wtot ;
327   x   /= wtot ;
328   dxx -= x * x ;
329   dzz /= wtot ;
330   z   /= wtot ;
331   dzz -= z * z ;
332   dxz /= wtot ;
333   dxz -= x * z ;
334
335   lambda[0] = TMath::Sqrt( 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ) ;
336   lambda[1] = TMath::Sqrt( 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ) ;
337 }
338
339 //____________________________________________________________________________
340 void AliPHOSEmcRecPoint::GetLocalPosition(TVector3 &LPos)
341 {
342   // Calculates the center of gravity in the local PHOS-module coordinates 
343   
344   if( fLocPos.X() < 1000000.) { // already evaluated
345    LPos = fLocPos ;
346    return ;
347   }
348
349   Float_t wtot = 0. ;
350  
351   Int_t relid[4] ;
352
353   Float_t x = 0. ;
354   Float_t z = 0. ;
355   
356   AliPHOSDigit * digit ;
357
358   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
359
360   Int_t iDigit;
361
362
363   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
364     digit = (AliPHOSDigit *) fDigitsList[iDigit];
365
366     Float_t xi ;
367     Float_t zi ;
368     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
369     phosgeom->RelPosInModule(relid, xi, zi);
370     Float_t w = TMath::Max( 0., fW0 + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
371     x    += xi * w ;
372     z    += zi * w ;
373     wtot += w ;
374
375   }
376
377   x /= wtot ;
378   z /= wtot ;
379   fLocPos.SetX(x)  ;
380   fLocPos.SetY(0.) ;
381   fLocPos.SetZ(z)  ;
382
383   LPos = fLocPos ;
384 }
385
386 //____________________________________________________________________________
387 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void)
388 {
389   // Finds the maximum energy in the cluster
390   
391   Float_t menergy = 0. ;
392
393   Int_t iDigit;
394
395   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
396  
397     if(fEnergyList[iDigit] > menergy) 
398       menergy = fEnergyList[iDigit] ;
399   }
400   return menergy ;
401 }
402
403 //____________________________________________________________________________
404 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(Float_t H) 
405 {
406   // Calculates the multiplicity of digits with energy larger than H*energy 
407   
408   Int_t multipl   = 0 ;
409   Int_t iDigit ;
410   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
411
412     if(fEnergyList[iDigit] > H * fAmp) 
413       multipl++ ;
414   }
415   return multipl ;
416 }
417
418 //____________________________________________________________________________
419 Int_t  AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t *  maxAt, Float_t * maxAtEnergy) 
420
421   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
422   //  energy difference between two local maxima
423
424   AliPHOSDigit * digit ;
425   AliPHOSDigit * digitN ;
426   
427
428   Int_t iDigitN ;
429   Int_t iDigit ;
430
431   for(iDigit = 0; iDigit < fMulDigit; iDigit++){
432     maxAt[iDigit] = fDigitsList[iDigit] ;
433   }
434   
435   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
436     if(maxAt[iDigit] != -1) {
437       digit = (AliPHOSDigit *) maxAt[iDigit] ;
438          
439       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
440         digitN = (AliPHOSDigit *) fDigitsList[iDigitN] ; 
441         
442         if ( AreNeighbours(digit, digitN) ) {
443           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
444             maxAt[iDigitN] = -1 ;
445             // but may be digit too is not local max ?
446             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + fLocMaxCut) 
447               maxAt[iDigit] = -1 ;
448           }
449           else {
450             maxAt[iDigit] = -1 ;
451             // but may be digitN too is not local max ?
452             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - fLocMaxCut) 
453               maxAt[iDigitN] = -1 ; 
454           } 
455         } // if Areneighbours
456       } // while digitN
457     } // slot not empty
458   } // while digit
459   
460   iDigitN = 0 ;
461   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
462     if(maxAt[iDigit] != -1){
463       maxAt[iDigitN] = maxAt[iDigit] ;
464       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
465       iDigitN++ ; 
466     }
467   }
468   return iDigitN ;
469 }
470
471
472 // //____________________________________________________________________________
473 // AliPHOSEmcRecPoint& AliPHOSEmcRecPoint::operator = (AliPHOSEmcRecPoint Clu) 
474 // {
475 //   int * dl = Clu.GetDigitsList() ; 
476  
477 //  if(fDigitsList) 
478 //     delete fDigitsList  ;
479
480 //   AliPHOSDigit * digit ;
481  
482 //   Int_t iDigit;
483
484 //   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
485 //     digit = (AliPHOSDigit *) dl[iDigit];
486 //     AddDigit(*digit) ;
487 //   }
488
489 //   fAmp       = Clu.GetTotalEnergy() ;
490 //   fGeom      = Clu.GetGeom() ;
491 //   TVector3 locpos;
492 //   Clu.GetLocalPosition(locpos) ;
493 //   fLocPos    = locpos;
494 //   fMulDigit       = Clu.GetMultiplicity() ;
495 //   fMaxDigit  = Clu.GetMaximumMultiplicity() ;
496 //   fPHOSMod   = Clu.GetPHOSMod() ;
497 //   fW0        = Clu.GetLogWeightCut() ; 
498 //   fDelta     = Clu.GetDelta() ;
499 //   fLocMaxCut = Clu.GetLocMaxCut() ;
500   
501 //   delete dl ; 
502  
503 //   return *this ;
504 // }
505
506 //____________________________________________________________________________
507 void AliPHOSEmcRecPoint::Print(Option_t * option) 
508 {
509   // Print the list of digits belonging to the cluster
510   
511   cout << "AliPHOSEmcRecPoint: " << endl ;
512
513   AliPHOSDigit * digit ; 
514   Int_t iDigit;
515   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
516
517   Float_t xi ;
518   Float_t zi ;
519   Int_t relid[4] ; 
520  
521   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
522     digit = (AliPHOSDigit *) fDigitsList[iDigit];
523     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
524     phosgeom->RelPosInModule(relid, xi, zi);
525     cout << " Id = " << digit->GetId() ;  
526     cout << "   module  = " << relid[0] ;  
527     cout << "   x  = " << xi ;  
528     cout << "   z  = " << zi ;  
529     cout << "   Energy = " << fEnergyList[iDigit] << endl ;
530   }
531   cout << "       Multiplicity    = " << fMulDigit  << endl ;
532   cout << "       Cluster Energy  = " << fAmp << endl ;
533   
534 }
535
536 //______________________________________________________________________________
537 void AliPHOSEmcRecPoint::Streamer(TBuffer &R__b)
538 {
539   // Stream an object of class AliPHOSEmcRecPoint.
540   // Needed because of the array fEnergyList
541   
542    if (R__b.IsReading()) {
543       Version_t R__v = R__b.ReadVersion(); if (R__v) { }
544       AliPHOSRecPoint::Streamer(R__b);
545       R__b >> fDelta;
546       R__b >> fLocMaxCut;
547       R__b.ReadArray(fEnergyList);
548       R__b >> fW0;
549    } else {
550       R__b.WriteVersion(AliPHOSEmcRecPoint::IsA());
551       AliPHOSRecPoint::Streamer(R__b);
552       R__b << fDelta;
553       R__b << fLocMaxCut;
554       R__b.WriteArray(fEnergyList, GetMaximumDigitMultiplicity() );
555       R__b << fW0;
556    }
557 }
558
559