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