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