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