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