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