]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSCpvRecPoint.cxx
b3ee63a8601f405eb11db5a41118d248b9c1961d
[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 //*-- Author: Yuri Kharlov
22 //  (after Dmitri Peressounko (RRC KI & SUBATECH))
23 //  30 October 2000 
24
25 // --- ROOT system ---
26 #include "TPad.h"
27 #include "TH2.h"
28 #include "TMath.h" 
29 #include "TCanvas.h" 
30 #include "TClonesArray.h" 
31 #include "AliPHOSGetter.h" 
32
33 // --- Standard library ---
34
35 #include <iostream.h> 
36
37 // --- AliRoot header files ---
38
39 #include "AliPHOSCpvRecPoint.h"
40 #include "AliPHOSPpsdRecPoint.h"
41
42 ClassImp(AliPHOSCpvRecPoint)
43
44 //____________________________________________________________________________
45 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint() : AliPHOSEmcRecPoint()
46 {
47   // ctor
48
49   fLengX = -1;
50   fLengZ = -1;
51 }
52
53 //____________________________________________________________________________
54 AliPHOSCpvRecPoint::~AliPHOSCpvRecPoint()
55 {
56   // dtor
57 }
58
59
60 //____________________________________________________________________________
61 Bool_t AliPHOSCpvRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
62 {
63   // Tells if (true) or not (false) two digits are neighbors)
64   
65   Bool_t aren = kFALSE ;
66   
67   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
68   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
69
70   Int_t relid1[4] ; 
71   phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ; 
72
73   Int_t relid2[4] ; 
74   phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ; 
75   
76   Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
77   Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
78
79   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
80     aren = kTRUE ;
81   
82   return aren ;
83 }
84
85 //____________________________________________________________________________
86 Int_t AliPHOSCpvRecPoint::Compare(const TObject * obj) const
87 {
88   // Compares two RecPoints according to their position in the PHOS modules
89
90   Float_t delta = 1 ; //Width of "Sorting row". If you changibg this 
91                       //value (what is senseless) change as vell delta in
92                       //AliPHOSTrackSegmentMakerv* and other RecPoints...
93
94   Int_t rv ; 
95
96   if( (strcmp(obj->ClassName() , "AliPHOSPpsdRecPoint" )) == 0)  // PPSD Rec Point
97     {
98       AliPHOSPpsdRecPoint * clu = (AliPHOSPpsdRecPoint *)obj ; 
99       if(this->GetPHOSMod()  < clu->GetPHOSMod() ) 
100         rv = -1 ;
101       else 
102         rv = 1 ;
103       return rv ;
104     }
105   else
106     {
107       AliPHOSCpvRecPoint * clu  = (AliPHOSCpvRecPoint *) obj ; 
108       
109       Int_t phosmod1 = GetPHOSMod() ;
110       Int_t phosmod2 = clu->GetPHOSMod() ;
111       
112       TVector3 locpos1; 
113       GetLocalPosition(locpos1) ;
114       TVector3 locpos2;  
115       clu->GetLocalPosition(locpos2) ;  
116       
117       if(phosmod1 == phosmod2 ) {
118         Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
119         if (rowdif> 0) 
120           rv = 1 ;
121         else if(rowdif < 0) 
122           rv = -1 ;
123         else if(locpos1.Z()>locpos2.Z()) 
124           rv = -1 ;
125         else 
126           rv = 1 ; 
127       }
128       
129       else {
130         if(phosmod1 < phosmod2 ) 
131           rv = -1 ;
132         else 
133           rv = 1 ;
134       }
135       
136       return rv ; 
137     }
138 }
139
140 //______________________________________________________________________________
141 void AliPHOSCpvRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
142 {
143 //   // Execute action corresponding to one event
144 //   //  This member function is called when a AliPHOSRecPoint is clicked with the locator
145 //   //
146 //   //  If Left button is clicked on AliPHOSRecPoint, the digits are switched on    
147 //   //  and switched off when the mouse button is released.
148 //   //
149
150 //   //   static Int_t pxold, pyold;
151
152 //   AliPHOSGetter * gime =  AliPHOSGetter::GetInstance() ; 
153   
154 //   static TGraph *  digitgraph = 0 ;
155   
156 //   if (!gPad->IsEditable()) return;
157   
158 //   TH2F * histo = 0 ;
159 //   TCanvas * histocanvas ; 
160   
161 //   switch (event) {
162     
163 //   case kButton1Down: {
164 //     AliPHOSDigit * digit ;
165 //     AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
166 //     Int_t iDigit;
167 //     Int_t relid[4] ;
168     
169 //     const Int_t kMulDigit = AliPHOSCpvRecPoint::GetDigitsMultiplicity() ; 
170 //     Float_t * xi = new Float_t[kMulDigit] ; 
171 //     Float_t * zi = new Float_t[kMulDigit] ; 
172     
173 //     // create the histogram for the single cluster 
174 //     // 1. gets histogram boundaries
175 //     Float_t ximax = -999. ; 
176 //     Float_t zimax = -999. ; 
177 //     Float_t ximin = 999. ; 
178 //     Float_t zimin = 999. ;
179     
180 //     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
181 //       digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
182 //       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
183 //       phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
184 //       if ( xi[iDigit] > ximax )
185 //      ximax = xi[iDigit] ; 
186 //       if ( xi[iDigit] < ximin )
187 //      ximin = xi[iDigit] ; 
188 //       if ( zi[iDigit] > zimax )
189 //      zimax = zi[iDigit] ; 
190 //       if ( zi[iDigit] < zimin )
191 //      zimin = zi[iDigit] ;     
192 //     }
193 //     ximax += phosgeom->GetCrystalSize(0) / 2. ;
194 //     zimax += phosgeom->GetCrystalSize(2) / 2. ;
195 //     ximin -= phosgeom->GetCrystalSize(0) / 2. ;
196 //     zimin -= phosgeom->GetCrystalSize(2) / 2. ;
197 //     Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5  ) ; 
198 //     Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
199     
200 //     // 2. gets the histogram title
201     
202 //     Text_t title[100] ; 
203 //     sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
204     
205 //     if (!histo) {
206 //       delete histo ; 
207 //       histo = 0 ; 
208 //     }
209 //     histo = new TH2F("cluster3D", title,  xdim, ximin, ximax, zdim, zimin, zimax)  ;
210     
211 //     Float_t x, z ; 
212 //     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
213 //       digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
214 //       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
215 //       phosgeom->RelPosInModule(relid, x, z);
216 //       histo->Fill(x, z, fEnergyList[iDigit] ) ;
217 //     }
218     
219 //     if (!digitgraph) {
220 //       digitgraph = new TGraph(kMulDigit,xi,zi);
221 //       digitgraph-> SetMarkerStyle(5) ; 
222 //       digitgraph-> SetMarkerSize(1.) ;
223 //       digitgraph-> SetMarkerColor(1) ;
224 //       digitgraph-> Paint("P") ;
225 //     }
226     
227 //     Print() ;
228 //     histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ; 
229 //     histocanvas->Draw() ; 
230 //     histo->Draw("lego1") ; 
231     
232 //     delete[] xi ; 
233 //     delete[] zi ; 
234     
235 //     break;
236 //   }
237   
238 //   case kButton1Up: 
239 //     if (digitgraph) {
240 //       delete digitgraph  ;
241 //       digitgraph = 0 ;
242 //     }
243 //     break;
244   
245 //    }
246 }
247
248 //____________________________________________________________________________
249 void AliPHOSCpvRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits)
250 {
251   // wraps other methods
252   AliPHOSEmcRecPoint::EvalAll(logWeight,digits) ;
253   EvalClusterLengths(digits) ;
254 }
255 //____________________________________________________________________________
256 void AliPHOSCpvRecPoint::EvalLocalPosition(Float_t logWeight,TClonesArray * digits)
257 {
258   // Calculates the center of gravity in the local PHOS-module coordinates 
259
260   Float_t wtot = 0. ;
261  
262   Int_t relid[4] ;
263
264   Float_t x = 0. ;
265   Float_t z = 0. ;
266   
267   AliPHOSDigit * digit ;
268
269   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
270
271   Int_t iDigit;
272
273   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
274     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]);
275
276     Float_t xi ;
277     Float_t zi ;
278     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
279     phosgeom->RelPosInModule(relid, xi, zi);
280     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
281     x    += xi * w ;
282     z    += zi * w ;
283     wtot += w ;
284   }
285
286   if (wtot != 0) {
287     x /= wtot ;
288     z /= wtot ;
289   } else {
290     x = -1e6 ;
291     z = -1e6 ;
292     if (fMulDigit != 0) cout << "AliPHOSCpvRecPoint: too low log weight factor "
293                              << "to evaluate cluster's center\n";
294   }
295   fLocPos.SetX(x)  ;
296   fLocPos.SetY(0.) ;
297   fLocPos.SetZ(z)  ;
298
299 }
300
301 //____________________________________________________________________________
302 void AliPHOSCpvRecPoint::EvalClusterLengths(TClonesArray * digits)
303 {
304   //Modified 15.03.2001 by Dmitri Peressounko
305   
306   // Calculates the cluster lengths along X and Z axes
307   // These characteristics are needed for CPV to tune
308   // digitization+reconstruction to experimental data
309   // Yuri Kharlov. 24 October 2000
310
311   Int_t relid[4] ;
312
313   AliPHOSDigit * digit ;
314
315   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
316
317   const Int_t kMaxLeng=20;
318   Int_t idX[kMaxLeng], idZ[kMaxLeng];
319   fLengX = 0;
320   fLengZ = 0;
321   Bool_t dejavu;
322
323   for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
324     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
325     Int_t absId = digit->GetId();
326     phosgeom->AbsToRelNumbering(absId, relid) ;
327
328     Int_t i;
329     dejavu=kFALSE;
330     for (i=0; i<fLengX; i++) if (relid[2]==idX[i]) { dejavu=kTRUE; break; }
331     if (!dejavu) {
332       idX[fLengX]=relid[2];
333       fLengX++;
334       fLengX = TMath::Min(fLengX,kMaxLeng);
335     }
336
337     dejavu=kFALSE;
338     for (i=0; i<fLengZ; i++) if (relid[3]==idZ[i]) { dejavu=kTRUE; break; }
339     if (!dejavu) {
340       idZ[fLengZ]=relid[3];
341       fLengZ++;
342       fLengZ = TMath::Min(fLengZ,kMaxLeng);
343     }
344   }
345 }
346
347
348
349 //____________________________________________________________________________
350 void AliPHOSCpvRecPoint::Print(Option_t * option) 
351 {
352   // Print the list of digits belonging to the cluster
353   
354   cout << "AliPHOSCpvRecPoint: " << endl ;
355
356   Int_t iDigit;
357
358   cout << "Digits #   "  ;
359   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
360     cout << fDigitsList[iDigit] << "  "  ; 
361   cout << endl ;
362
363   cout << "Energies: "  ;
364   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
365     cout << fEnergyList[iDigit] << "  "  ; 
366   cout << endl ;
367   
368   cout << "       Multiplicity    = " << fMulDigit  << endl ;
369   cout << "       Cluster Energy  = " << fAmp << endl ;
370   cout << "       Stored at position " << GetIndexInList() << endl ; 
371  
372 }