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