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