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