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