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