]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDUtility.cxx
A new method DrawPMDModule is added
[u/mrichter/AliRoot.git] / PMD / AliPMDUtility.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 //                                                     //
17 //                                                     //
18 //  Date   : August 05 2003                            //
19 //                                                     //
20 //  Utility code for ALICE-PMD                         //
21 //                                                     //
22 //-----------------------------------------------------//
23
24 #include "Riostream.h"
25 #include "TMath.h"
26 #include "TText.h"
27 #include "TLine.h"
28
29 #include <stdio.h>
30 #include <math.h>
31
32 #include "AliPMDUtility.h"
33
34
35 ClassImp(AliPMDUtility)
36
37 AliPMDUtility::AliPMDUtility():
38   fPx(0.),
39   fPy(0.),
40   fPz(0.),
41   fTheta(0.),
42   fEta(0.),
43   fPhi(0.),
44   fWriteModule(1)
45 {
46   // Default constructor
47 }
48
49 AliPMDUtility::AliPMDUtility(Float_t px, Float_t py, Float_t pz):
50   fPx(px),
51   fPy(py),
52   fPz(pz),
53   fTheta(0.),
54   fEta(0.),
55   fPhi(0.),
56   fWriteModule(1)
57 {
58   // Constructor
59 }
60 AliPMDUtility::AliPMDUtility(const AliPMDUtility &pmdutil):
61   fPx(pmdutil.fPx),
62   fPy(pmdutil.fPy),
63   fPz(pmdutil.fPz),
64   fTheta(pmdutil.fTheta),
65   fEta(pmdutil.fEta),
66   fPhi(pmdutil.fPhi),
67   fWriteModule(pmdutil.fWriteModule)
68 {
69   // copy constructor
70 }
71 AliPMDUtility & AliPMDUtility::operator=(const AliPMDUtility &pmdutil)
72 {
73   // assignment operator
74   if(this != &pmdutil)
75     {
76       fPx = pmdutil.fPx;
77       fPy = pmdutil.fPy;
78       fPz = pmdutil.fPz;
79       fTheta = pmdutil.fTheta;
80       fEta = pmdutil.fEta;
81       fPhi = pmdutil.fPhi;
82       fWriteModule = pmdutil.fWriteModule;
83     }
84   return *this;
85 }
86 AliPMDUtility::~AliPMDUtility()
87 {
88   // Default destructor
89 }
90
91 void AliPMDUtility::RectGeomCellPos(Int_t ism, Int_t xpad, Int_t ypad, Float_t &xpos, Float_t &ypos)
92 {
93   // This routine finds the cell eta,phi for the new PMD rectangular 
94   // geometry in ALICE
95   // Authors : Bedanga Mohanty and Dipak Mishra - 29.4.2003
96   // modified by B. K. Nandi for change of coordinate sys
97   //
98   // SMA  ---> Supermodule Type A           ( SM - 0)
99   // SMAR ---> Supermodule Type A ROTATED   ( SM - 1)
100   // SMB  ---> Supermodule Type B           ( SM - 2)
101   // SMBR ---> Supermodule Type B ROTATED   ( SM - 3)
102   //
103   // ism   : Serial module number from 0 to 23 for each plane
104
105  
106   // Corner positions (x,y) of the 24 unit moudles in ALICE PMD
107
108   double xcorner[24] =
109     {
110       74.8833,  53.0045, 31.1255,    //Type-A
111       74.8833,  53.0045, 31.1255,    //Type-A
112       -74.8833, -53.0044, -31.1255,  //Type-AR
113       -74.8833, -53.0044, -31.1255,  //Type-AR
114       8.9165, -33.7471,            //Type-B
115       8.9165, -33.7471,            //Type-B
116       8.9165, -33.7471,            //Type-B
117       -8.9165, 33.7471,            //Type-BR
118       -8.9165, 33.7471,            //Type-BR
119       -8.9165, 33.7471,            //Type-BR
120     };
121
122   
123   double ycorner[24] =
124     {
125       86.225,  86.225,  86.225,      //Type-A
126       37.075,  37.075,  37.075,      //Type-A
127       -86.225, -86.225, -86.225,     //Type-AR
128       -37.075, -37.075, -37.075,     //Type-AR
129       86.225,  86.225,               //Type-B
130       61.075,  61.075,               //Type-B
131       35.925,  35.925,               //Type-B
132       -86.225, -86.225,              //Type-BR
133       -61.075, -61.075,              //Type-BR
134       -35.925, -35.925               //Type-BR
135     };
136
137   
138   const Float_t kSqroot3      = 1.73205;  // sqrt(3.);
139   const Float_t kCellRadius   = 0.25;
140   
141   //
142   //Every even row of cells is shifted and placed
143   //in geant so this condition
144   //
145   Float_t cellRadius = 0.25;
146   Float_t shift = 0.0;
147   if(xpad%2 == 0)
148     {
149       shift = -cellRadius/2.0;
150     }
151   else
152     {
153       shift = 0.0;
154     }
155
156
157   if(ism < 6)
158     {
159       ypos = ycorner[ism] - (Float_t) xpad*kCellRadius*2.0 + shift;
160       xpos = xcorner[ism] - (Float_t) ypad*kSqroot3*kCellRadius;
161     }
162   else if(ism >=6 && ism < 12)
163     {
164       ypos = ycorner[ism] + (Float_t) xpad*kCellRadius*2.0 + shift;
165       xpos = xcorner[ism] + (Float_t) ypad*kSqroot3*kCellRadius;
166     }
167   else if(ism >= 12 && ism < 18)
168     {
169       ypos = ycorner[ism] - (Float_t) xpad*kCellRadius*2.0 + shift;
170       xpos = xcorner[ism] - (Float_t) ypad*kSqroot3*kCellRadius;
171     }
172   else if(ism >= 18 && ism < 24)
173     {
174       ypos = ycorner[ism] + (Float_t) xpad*kCellRadius*2.0 + shift;
175       xpos = xcorner[ism] + (Float_t) ypad*kSqroot3*kCellRadius;
176     }
177
178 }
179
180 void AliPMDUtility::RectGeomCellPos(Int_t ism, Float_t xpad, Float_t ypad, Float_t &xpos, Float_t &ypos)
181 {
182   // If the xpad and ypad inputs are float, then 0.5 is added to it
183   // to find the layer which is shifted.
184   // This routine finds the cell eta,phi for the new PMD rectangular 
185   // geometry in ALICE
186   // Authors : Bedanga Mohanty and Dipak Mishra - 29.4.2003
187   // modified by B. K. Nnadi for change of coordinate sys
188   //
189   // SMA  ---> Supermodule Type A           ( SM - 0)
190   // SMAR ---> Supermodule Type A ROTATED   ( SM - 1)
191   // SMB  ---> Supermodule Type B           ( SM - 2)
192   // SMBR ---> Supermodule Type B ROTATED   ( SM - 3)
193   //
194   // ism   : Serial Module number from 0 to 23 for each plane
195
196   // Corner positions (x,y) of the 24 unit moudles in ALICE PMD
197
198   double xcorner[24] =
199     {
200       74.8833,  53.0045, 31.1255,    //Type-A
201       74.8833,  53.0045, 31.1255,    //Type-A
202       -74.8833, -53.0044, -31.1255,  //Type-AR
203       -74.8833, -53.0044, -31.1255,  //Type-AR
204       8.9165, -33.7471,            //Type-B
205       8.9165, -33.7471,            //Type-B
206       8.9165, -33.7471,            //Type-B
207       -8.9165, 33.7471,            //Type-BR
208       -8.9165, 33.7471,            //Type-BR
209       -8.9165, 33.7471,            //Type-BR
210     };
211
212   
213
214   double ycorner[24] =
215     {
216       86.225,  86.225,  86.225,      //Type-A
217       37.075,  37.075,  37.075,      //Type-A
218       -86.225, -86.225, -86.225,     //Type-AR
219       -37.075, -37.075, -37.075,     //Type-AR
220       86.225,  86.225,               //Type-B
221       61.075,  61.075,               //Type-B
222       35.925,  35.925,               //Type-B
223       -86.225, -86.225,              //Type-BR
224       -61.075, -61.075,              //Type-BR
225       -35.925, -35.925               //Type-BR
226     };
227
228
229   const Float_t kSqroot3    = 1.73205;  // sqrt(3.);
230   const Float_t kCellRadius = 0.25;
231   
232   //
233   //Every even row of cells is shifted and placed
234   //in geant so this condition
235   //
236   Float_t cellRadius = 0.25;
237   Float_t shift = 0.0;
238   Int_t iirow = (Int_t) (xpad+0.5);
239   if(iirow%2 == 0)
240     {
241       shift = -cellRadius/2.0;
242     }
243   else
244     {
245       shift = 0.0;
246     }
247
248   if(ism < 6)
249     {
250       ypos = ycorner[ism] - xpad*kCellRadius*2.0 + shift;
251       xpos = xcorner[ism] - ypad*kSqroot3*kCellRadius;
252     }
253   else if(ism >=6 && ism < 12)
254     {
255       ypos = ycorner[ism] + xpad*kCellRadius*2.0 + shift;
256       xpos = xcorner[ism] + ypad*kSqroot3*kCellRadius;
257     }
258   else if(ism >= 12 && ism < 18)
259     {
260       ypos = ycorner[ism] - xpad*kCellRadius*2.0 + shift;
261       xpos = xcorner[ism] - ypad*kSqroot3*kCellRadius;
262     }
263   else if(ism >= 18 && ism < 24)
264     {
265       ypos = ycorner[ism] + xpad*kCellRadius*2.0 + shift;
266       xpos = xcorner[ism] + ypad*kSqroot3*kCellRadius;
267     }
268
269 }
270 // -------------------------------------------------------- //
271
272 void AliPMDUtility::GenerateBoundaryPoints(Int_t ism, Float_t &x1ism, 
273                                            Float_t &y1ism, Float_t &x2ism,
274                                            Float_t &y2ism)
275 {
276   // Generate bounding-box.
277
278
279     Float_t xism = 0, yism = 0;
280     Float_t dxism, dyism;
281
282     const Float_t kRad     = 0.25;
283     const Float_t kSqRoot3 = 1.732050808;
284     const Float_t kDia     = 0.50;
285
286
287   const Double_t kXcorner[24] =
288     {
289       74.8833,  53.0045, 31.1255,    //Type-A
290       74.8833,  53.0045, 31.1255,    //Type-A
291       -74.8833, -53.0044, -31.1255,  //Type-AR
292       -74.8833, -53.0044, -31.1255,  //Type-AR
293       8.9165, -33.7471,            //Type-B
294       8.9165, -33.7471,            //Type-B
295       8.9165, -33.7471,            //Type-B
296       -8.9165, 33.7471,            //Type-BR
297       -8.9165, 33.7471,            //Type-BR
298       -8.9165, 33.7471,            //Type-BR
299     };
300
301
302   const Double_t kYcorner[24] =
303     {
304       86.225,  86.225,  86.225,      //Type-A
305       37.075,  37.075,  37.075,      //Type-A
306       -86.225, -86.225, -86.225,     //Type-AR
307       -37.075, -37.075, -37.075,     //Type-AR
308       86.225,  86.225,               //Type-B
309       61.075,  61.075,               //Type-B
310       35.925,  35.925,               //Type-B
311       -86.225, -86.225,              //Type-BR
312       -61.075, -61.075,              //Type-BR
313       -35.925, -35.925               //Type-BR
314     };
315
316
317   if (ism > 23) ism -= 24;
318
319
320   if (ism < 6)
321     {
322       xism  = kXcorner[ism] + kRad;
323       yism  = kYcorner[ism] + kRad;
324       dxism = -kRad*kSqRoot3*48.;
325       dyism = -kDia*96. - kRad;
326   }
327   if (ism >= 6 && ism < 12)
328     {
329       xism  = kXcorner[ism] - kRad;
330       yism  = kYcorner[ism] - kRad;
331       dxism = kRad*kSqRoot3*48.;
332       dyism = kDia*96. + kRad;
333   }
334   if (ism >= 12 && ism < 18)
335     {
336       xism  = kXcorner[ism] + kRad;
337       yism  = kYcorner[ism] + kRad;
338       dxism = -kRad*kSqRoot3*96.;
339       dyism = -kDia*48. - kRad;
340   }
341   if (ism >= 18 && ism < 24)
342     {
343       xism  = kXcorner[ism] - kRad;
344       yism  = kYcorner[ism] - kRad;
345       dxism = kRad*kSqRoot3*96.;
346       dyism = kDia*48. + kRad;
347   }
348
349   x1ism = xism;
350   x2ism = xism + dxism;
351   y1ism = yism;
352   y2ism = yism + dyism;
353
354 }
355 // ------------------------------------------------------------------- //
356
357 void AliPMDUtility::DrawPMDModule()
358 {
359
360     Float_t x1ism, x2ism, y1ism, y2ism;
361     Float_t deltaX, deltaY;
362     
363     //TH2F *h2 = new TH2F("h2","Y vs. X",200,-100.,100.,200,-100.,100.);
364     //h2->Draw();
365
366     TLine t;
367     t.SetLineColor(2);
368
369     TText tt;
370     tt.SetTextColor(4);
371
372     Char_t smnumber[10];
373
374     for(Int_t ism=0; ism < 24; ism++)
375     {
376         GenerateBoundaryPoints(ism, x1ism, y1ism, x2ism, y2ism);
377         deltaX = (x2ism - x1ism)/2.;
378         deltaY = (y2ism - y1ism)/2.;
379         if (fWriteModule == 1)
380         {
381             sprintf(smnumber,"%d",ism);
382             tt.DrawText(x1ism+deltaX,y1ism+deltaY,smnumber);
383         }
384         t.DrawLine(x1ism, y1ism, x1ism, y2ism);
385         t.DrawLine(x1ism, y1ism, x2ism, y1ism);
386         t.DrawLine(x2ism, y1ism, x2ism, y2ism);
387         t.DrawLine(x1ism, y2ism, x2ism, y2ism);
388     }
389
390 }
391
392 // ------------------------------------------------------------------- //
393
394
395 void AliPMDUtility::ApplyVertexCorrection(Float_t vertex[], Float_t xpos,
396                                           Float_t ypos, Float_t zpos)
397 {
398   // Not implemented
399   fPx = xpos - vertex[0];
400   fPy = ypos - vertex[1];
401   fPz = zpos - vertex[2];
402 }
403 void AliPMDUtility::ApplyAlignment()
404 {
405   // Not implemented
406 }
407
408 void AliPMDUtility::SetPxPyPz(Float_t px, Float_t py, Float_t pz)
409 {
410   fPx = px;
411   fPy = py;
412   fPz = pz;
413 }
414
415 void AliPMDUtility::SetXYZ(Float_t xpos, Float_t ypos, Float_t zpos)
416 {
417   fPx = xpos;
418   fPy = ypos;
419   fPz = zpos;
420 }
421 void AliPMDUtility::SetWriteModule(Int_t wrmod)
422 {
423     fWriteModule = wrmod;
424 }
425 void AliPMDUtility::CalculateEta()
426 {
427   Float_t rpxpy, theta, eta;
428
429   rpxpy  = TMath::Sqrt(fPx*fPx + fPy*fPy);
430   theta  = TMath::ATan2(rpxpy,fPz);
431   eta    = -TMath::Log(TMath::Tan(0.5*theta));
432   fTheta = theta;
433   fEta   = eta;
434 }
435 void AliPMDUtility::CalculatePhi()
436 {
437   Float_t pybypx, phi = 0., phi1;
438
439   if(fPx==0)
440     {
441       if(fPy>0) phi = 90.;
442       if(fPy<0) phi = 270.;
443     }
444   if(fPx != 0)
445     {
446       pybypx = fPy/fPx;
447       if(pybypx < 0) pybypx = - pybypx;
448       phi1 = TMath::ATan(pybypx)*180./3.14159;
449
450       if(fPx > 0 && fPy > 0) phi = phi1;        // 1st Quadrant
451       if(fPx < 0 && fPy > 0) phi = 180 - phi1;  // 2nd Quadrant
452       if(fPx < 0 && fPy < 0) phi = 180 + phi1;  // 3rd Quadrant
453       if(fPx > 0 && fPy < 0) phi = 360 - phi1;  // 4th Quadrant
454
455     }
456   phi = phi*3.14159/180.;
457
458   fPhi = phi;
459
460 }
461 void AliPMDUtility::CalculateEtaPhi()
462 {
463   Float_t rpxpy, theta, eta;
464   Float_t pybypx, phi = 0., phi1;
465
466   rpxpy = TMath::Sqrt(fPx*fPx + fPy*fPy);
467   theta = TMath::ATan2(rpxpy,fPz);
468   eta   = -TMath::Log(TMath::Tan(0.5*theta));
469   
470   if(fPx == 0)
471     {
472       if(fPy>0) phi = 90.;
473       if(fPy<0) phi = 270.;
474     }
475   if(fPx != 0)
476     {
477       pybypx = fPy/fPx;
478       if(pybypx < 0) pybypx = - pybypx;
479       phi1 = TMath::ATan(pybypx)*180./3.14159;
480       if(fPx > 0 && fPy > 0) phi = phi1;        // 1st Quadrant
481       if(fPx < 0 && fPy > 0) phi = 180 - phi1;  // 2nd Quadrant
482       if(fPx < 0 && fPy < 0) phi = 180 + phi1;  // 3rd Quadrant
483       if(fPx > 0 && fPy < 0) phi = 360 - phi1;  // 4th Quadrant
484
485     }
486   phi = phi*3.14159/180.;
487
488   fTheta = theta;
489   fEta   = eta;
490   fPhi   = phi;
491 }
492 void AliPMDUtility::CalculateXY(Float_t eta, Float_t phi, Float_t zpos)
493 {
494   // Not implemented
495
496   //  eta   = -TMath::Log(TMath::Tan(0.5*theta));
497
498   Float_t xpos = 0., ypos = 0.;
499
500   //  Float_t theta = 2.0*TMath::ATan(TMath::Log(-eta));
501
502   fEta = eta;
503   fPhi = phi;
504   fPx  = xpos;
505   fPy  = ypos;
506   fPz  = zpos;
507 }
508 Float_t AliPMDUtility::GetTheta() const
509 {
510   return fTheta;
511 }
512 Float_t AliPMDUtility::GetEta() const
513 {
514   return fEta;
515 }
516 Float_t AliPMDUtility::GetPhi() const
517 {
518   return fPhi;
519 }
520 Float_t AliPMDUtility::GetX() const
521 {
522   return fPx;
523 }
524 Float_t AliPMDUtility::GetY() const
525 {
526   return fPy;
527 }
528 Float_t AliPMDUtility::GetZ() const
529 {
530   return fPz;
531 }
532
533