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