]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSgeomMatrix.cxx
Checking RecoParam before creating the ITSdebug.root file (A. Dainese)
[u/mrichter/AliRoot.git] / ITS / AliITSgeomMatrix.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 $Id$ 
18 */
19 ////////////////////////////////////////////////////////////////////////
20 // This is the implementation file for AliITSgeomMatrix class. It 
21 // contains the routines to manipulate, setup, and queary the geometry 
22 // of a given ITS module. An ITS module may be one of at least three
23 // ITS detector technologies, Silicon Pixel, Drift, or Strip Detectors,
24 // and variations of these in size and/or layout. These routines let
25 // one go between ALICE global coordiantes (cm) to a given modules 
26 // specific local coordinates (cm).
27 ////////////////////////////////////////////////////////////////////////
28
29 #include <Riostream.h>
30 #include <TClass.h>
31 #include <TMath.h>
32 #include <TBuffer.h>
33 #include <TCanvas.h>
34 #if ROOT_VERSION_CODE>= 331523
35 #include <TView3D.h>
36 #else
37 #include <TView.h>
38 #endif
39 #include <TPolyLine3D.h>
40 #include <TNode.h>
41 #include <TPCON.h>
42 #include <TBRIK.h>
43 #include <TXTRU.h>
44
45 #include "AliITSgeomMatrix.h"
46
47 ClassImp(AliITSgeomMatrix)
48 //----------------------------------------------------------------------
49 AliITSgeomMatrix::AliITSgeomMatrix():
50 TObject(),         // Base Class.
51 fDetectorIndex(0), // Detector type index (like fShapeIndex was)
52 fid(),       // layer, ladder, detector numbers.
53 frot(),      //! vector of rotations about x,y,z [radians].
54 ftran(),     // Translation vector of module x,y,z.
55 fCylR(0.0),  //! R Translation in Cylinderical coordinates
56 fCylPhi(0.0),//! Phi Translation vector in Cylindrical coord.
57 fm(),        // Rotation matrix based on frot.
58 fPath(){     // Path in geometry to this module
59     // The Default constructor for the AliITSgeomMatrix class. By Default
60     // the angles of rotations are set to zero, meaning that the rotation
61     // matrix is the unit matrix. The translation vector is also set to 
62     // zero as are the module id number. The detector type is set to -1 
63     // (an undefined value). The full rotation matrix is kept so that 
64     // the evaluation  of a coordinate transformation can be done 
65     // quickly and with a minimum of CPU overhead. The basic coordinate 
66     // systems are the ALICE global coordinate system and the detector 
67     // local coordinate system. In general this structure is not limited 
68     // to just those two coordinate systems.
69     //Begin_Html
70     /*
71       <img src="picts/ITS/AliITSgeomMatrix_L1.gif">
72     */
73     //End_Html
74     // Inputs:
75     //    none.
76     // Outputs:
77     //    none.
78     // Return:
79     //    A default constructes AliITSgeomMatrix class.
80     Int_t i,j;
81
82     fDetectorIndex = -1; // a value never defined.
83     for(i=0;i<3;i++){
84         fid[i] = 0;
85         frot[i] = ftran[i] = 0.0;
86         for(j=0;j<3;j++) fm[i][j] = 0.0;
87     }// end for i
88     fm[0][0] = fm[1][1] = fm[2][2] = 1.0;
89 }
90
91 //----------------------------------------------------------------------
92 AliITSgeomMatrix::AliITSgeomMatrix(const AliITSgeomMatrix &source) : 
93 TObject(source),         // Base Class.
94 fDetectorIndex(source.fDetectorIndex),// Detector type index (like 
95                                       // fShapeIndex was)
96 fid(),       // layer, ladder, detector numbers.
97 frot(),      //! vector of rotations about x,y,z [radians].
98 ftran(),     // Translation vector of module x,y,z.
99 fCylR(source.fCylR),  //! R Translation in Cylinderical coordinates
100 fCylPhi(source.fCylPhi),//! Phi Translation vector in Cylindrical coord.
101 fm(),        // Rotation matrix based on frot.
102 fPath(source.fPath){
103     // The standard Copy constructor. This make a full / proper copy of
104     // this class.
105     // Inputs:
106     //    AliITSgeomMatrix &source   The source of this copy
107     // Outputs:
108     //    none.
109     // Return:
110     //    A copy constructes AliITSgeomMatrix class.
111     Int_t i,j;
112
113     for(i=0;i<3;i++){
114         this->fid[i]     = source.fid[i];
115         this->frot[i]    = source.frot[i];
116         this->ftran[i]   = source.ftran[i];
117         for(j=0;j<3;j++) this->fm[i][j] = source.fm[i][j];
118     }// end for i
119 }
120 //----------------------------------------------------------------------
121 AliITSgeomMatrix& AliITSgeomMatrix::operator=(const AliITSgeomMatrix &source){
122     // The standard = operator. This make a full / proper copy of
123     // this class.
124     // The standard Copy constructor. This make a full / proper copy of
125     // this class.
126     // Inputs:
127     //    AliITSgeomMatrix &source   The source of this copy
128     // Outputs:
129     //    none.
130     // Return:
131     //    A copy of the source AliITSgeomMatrix class.
132
133     if(this == &source)return *this;
134     Int_t i,j;
135
136     this->fDetectorIndex = source.fDetectorIndex;
137     this->fCylR      = source.fCylR;
138     this->fCylPhi    = source.fCylPhi;
139     for(i=0;i<3;i++){
140         this->fid[i]     = source.fid[i];
141         this->frot[i]    = source.frot[i];
142         this->ftran[i]   = source.ftran[i];
143         for(j=0;j<3;j++) this->fm[i][j] = source.fm[i][j];
144     } // end for i
145     this->fPath   = source.fPath;
146     return *this;
147 }
148 //----------------------------------------------------------------------
149 AliITSgeomMatrix::AliITSgeomMatrix(Int_t idt,const Int_t id[3],
150                         const Double_t rot[3],const Double_t tran[3]):
151 TObject(),           // Base class
152 fDetectorIndex(idt), // Detector type index (like fShapeIndex was)
153 fid(),       // layer, ladder, detector numbers.
154 frot(),      //! vector of rotations about x,y,z [radians].
155 ftran(),     // Translation vector of module x,y,z.
156 fCylR(0.0),  //! R Translation in Cylinderical coordinates
157 fCylPhi(0.0),//! Phi Translation vector in Cylindrical coord.
158 fm(),        // Rotation matrix based on frot.
159 fPath(){     // Path in geometry to this moduel
160     // This is a constructor for the AliITSgeomMatrix class. The matrix is
161     // defined by 3 standard rotation angles [radians], and the translation
162     // vector tran [cm]. In addition the layer, ladder, and detector number
163     // for this particular module and the type of module must be given.
164     // The full rotation matrix is kept so that the evaluation 
165     // of a coordinate transformation can be done quickly and with a minimum
166     // of CPU overhead. The basic coordinate systems are the ALICE global
167     // coordinate system and the detector local coordinate system. In general
168     // this structure is not limited to just those two coordinate systems.
169     //Begin_Html
170     /*
171       <img src="picts/ITS/AliITSgeomMatrix_L1.gif">
172     */
173     //End_Html
174     // Inputs:
175     //    Int_t idt        The detector index value
176     //    Int_t id[3]      The layer, ladder, and detector numbers
177     //    Double_t rot[3]  The 3 Cartician rotaion angles [radians]
178     //    Double_t tran[3] The 3 Cartician translation distnaces
179     // Outputs:
180     //    none.
181     // Return:
182     //    A properly inilized AliITSgeomMatrix class.
183     Int_t i;
184
185     for(i=0;i<3;i++){
186         fid[i]   = id[i];
187         frot[i]  = rot[i];
188         ftran[i] = tran[i];
189     }// end for i
190     fCylR   = TMath::Sqrt(ftran[0]*ftran[0]+ftran[1]*ftran[1]);
191     fCylPhi = TMath::ATan2(ftran[1],ftran[0]);
192     if(fCylPhi<0.0) fCylPhi += 2.*TMath::Pi();
193     this->MatrixFromAngle();
194 }
195 //----------------------------------------------------------------------
196 AliITSgeomMatrix::AliITSgeomMatrix(Int_t idt, const Int_t id[3],
197                                    Double_t matrix[3][3],
198                                    const Double_t tran[3]):
199 TObject(),            // Base class
200 fDetectorIndex(idt), // Detector type index (like fShapeIndex was)
201 fid(),       // layer, ladder, detector numbers.
202 frot(),      //! vector of rotations about x,y,z [radians].
203 ftran(),     // Translation vector of module x,y,z.
204 fCylR(0.0),  //! R Translation in Cylinderical coordinates
205 fCylPhi(0.0),//! Phi Translation vector in Cylindrical coord.
206 fm(),        // Rotation matrix based on frot.
207 fPath(){     // Path in geometry to this module
208     // This is a constructor for the AliITSgeomMatrix class. The 
209     // rotation matrix is given as one of the inputs, and the 
210     // translation vector tran [cm]. In  addition the layer, ladder, 
211     // and detector number for this particular module and the type of 
212     // module must be given. The full rotation matrix is kept so that 
213     // the evaluation of a coordinate transformation can be done quickly 
214     // and with a minimum of CPU overhead. The basic coordinate systems 
215     // are the ALICE global coordinate system and the detector local
216     // coordinate system. In general this structure is not limited to just
217     // those two coordinate systems.
218     //Begin_Html
219     /*
220       <img src="picts/ITS/AliITSgeomMatrix_L1.gif">
221     */
222     //End_Html
223     // Inputs:
224     //    Int_t idt          The detector index value
225     //    Int_t id[3]        The layer, ladder, and detector numbers
226     //    Double_t rot[3][3] The 3x3 Cartician rotaion matrix
227     //    Double_t tran[3]   The 3 Cartician translation distnaces
228     // Outputs:
229     //    none.
230     // Return:
231     //    A properly inilized AliITSgeomMatrix class.
232     Int_t i,j;
233
234     for(i=0;i<3;i++){
235         fid[i]   = id[i];
236         ftran[i] = tran[i];
237         for(j=0;j<3;j++) fm[i][j] = matrix[i][j];
238     }// end for i
239     fCylR   = TMath::Sqrt(ftran[0]*ftran[0]+ftran[1]*ftran[1]);
240     fCylPhi = TMath::ATan2(ftran[1],ftran[0]);
241     if(fCylPhi<0.0) fCylPhi += 2.*TMath::Pi();
242     this->AngleFromMatrix();
243 }
244 //----------------------------------------------------------------------
245 void AliITSgeomMatrix::SixAnglesFromMatrix(Double_t *ang)const{
246     // This function returns the 6 GEANT 3.21 rotation angles [degrees] in
247     // the array ang which must be at least [6] long.
248     // Inputs:
249     //   none.
250     // Outputs:
251     //   Double_t ang[6]  The 6 Geant3.21 rotation angles. [degrees]
252     // Return:
253     //   noting
254     Double_t si,c=180./TMath::Pi();
255     const Double_t epsil=1.e-15;
256
257     ang[1] = TMath::ATan2(fm[0][1],fm[0][0]);
258     if( !(TMath::AreEqualAbs(TMath::Cos(ang[1]),0.,epsil))) si = fm[0][0]/TMath::Cos(ang[1]);
259     else si = fm[0][1]/TMath::Sin(ang[1]);
260     ang[0] = TMath::ATan2(si,fm[0][2]);
261
262     ang[3] = TMath::ATan2(fm[1][1],fm[1][0]);
263     if(!(TMath::AreEqualAbs(TMath::Cos(ang[3]),0.,epsil))) si = fm[1][0]/TMath::Cos(ang[3]);
264     else si = fm[1][1]/TMath::Sin(ang[3]);
265     ang[2] = TMath::ATan2(si,fm[1][2]);
266
267     ang[5] = TMath::ATan2(fm[2][1],fm[2][0]);
268     if(!(TMath::AreEqualAbs(TMath::Cos(ang[5]),0.,epsil))) si = fm[2][0]/TMath::Cos(ang[5]);
269     else si = fm[2][1]/TMath::Sin(ang[5]);
270     ang[4] = TMath::ATan2(si,fm[2][2]);
271
272     for(Int_t i=0;i<6;i++) {ang[i] *= c; if(ang[i]<0.0) ang[i] += 360.;}
273 }
274 //----------------------------------------------------------------------
275 void AliITSgeomMatrix::MatrixFromSixAngles(const Double_t *ang){
276     // Given the 6 GEANT 3.21 rotation angles [degree], this will compute and
277     // set the rotations matrix and 3 standard rotation angles [radians].
278     // These angles and rotation matrix are overwrite the existing values in
279     // this class.
280     // Inputs:
281     //   Double_t ang[6]  The 6 Geant3.21 rotation angles. [degrees]
282     // Outputs:
283     //   none.
284     // Return:
285     //   noting
286     Int_t    i,j;
287     Double_t si,lr[9],c=TMath::Pi()/180.;
288     const Double_t epsil = 1.e-15;
289
290     si    = TMath::Sin(c*ang[0]);
291     if(TMath::AreEqualAbs(ang[0],90.,epsil)) si = +1.0;
292     if(TMath::AreEqualAbs(ang[0],270.,epsil)) si = -1.0;
293     if(TMath::AreEqualAbs(ang[0],0.,epsil) ||TMath::AreEqualAbs(ang[0],180.,epsil)) si =  0.0;
294     lr[0] = si * TMath::Cos(c*ang[1]);
295     lr[1] = si * TMath::Sin(c*ang[1]);
296     lr[2] = TMath::Cos(c*ang[0]);
297     if(TMath::AreEqualAbs(ang[0],90.,epsil)||TMath::AreEqualAbs(ang[0],270.,epsil)) lr[2] =  0.0;
298     if(TMath::AreEqualAbs(ang[0],0.,epsil))                  lr[2] = +1.0;
299     if(TMath::AreEqualAbs(ang[0],180.,epsil))                 lr[2] = -1.0;
300 //
301     si    =  TMath::Sin(c*ang[2]);
302     if(TMath::AreEqualAbs(ang[2],90.,epsil))                 si = +1.0; 
303     if(TMath::AreEqualAbs(ang[2],270.,epsil))                 si = -1.0;
304     if(TMath::AreEqualAbs(ang[2],0.,epsil) || TMath::AreEqualAbs(ang[2],180.,epsil)) si =  0.0;
305     lr[3] = si * TMath::Cos(c*ang[3]);
306     lr[4] = si * TMath::Sin(c*ang[3]);
307     lr[5] = TMath::Cos(c*ang[2]);
308     if(TMath::AreEqualAbs(ang[2],90.,epsil) || TMath::AreEqualAbs(ang[2],270.,epsil)) lr[5] =  0.0;
309     if(TMath::AreEqualAbs(ang[2],0.,epsil))                 lr[5] = +1.0;
310     if(TMath::AreEqualAbs(ang[2],180.,epsil))                 lr[5] = -1.0;
311 //
312     si    = TMath::Sin(c*ang[4]);
313     if(TMath::AreEqualAbs(ang[4],90.,epsil))                 si = +1.0;
314     if(TMath::AreEqualAbs(ang[4],270.0,epsil))                 si = -1.0;
315     if(TMath::AreEqualAbs(ang[4],0.,epsil)|| TMath::AreEqualAbs(ang[4],180.,epsil)) si =  0.0;
316     lr[6] = si * TMath::Cos(c*ang[5]);
317     lr[7] = si * TMath::Sin(c*ang[5]);
318     lr[8] = TMath::Cos(c*ang[4]);
319     if(TMath::AreEqualAbs(ang[4],90.0,epsil) ||TMath::AreEqualAbs(ang[4],270.,epsil)) lr[8] =  0.0;
320     if(TMath::AreEqualAbs(ang[4],0.,epsil))                  lr[8] = +1.0;
321     if(TMath::AreEqualAbs(ang[4],180.0,epsil))                  lr[8] = -1.0;
322     // Normalize these elements and fill matrix fm.
323     for(i=0;i<3;i++){// reuse si.
324         si = 0.0;
325         for(j=0;j<3;j++) si += lr[3*i+j]*lr[3*i+j];
326         si = TMath::Sqrt(1./si);
327         for(j=0;j<3;j++) fm[i][j] = si*lr[3*i+j];
328     } // end for i
329     this->AngleFromMatrix();
330 }
331 //----------------------------------------------------------------------
332 AliITSgeomMatrix::AliITSgeomMatrix(const Double_t rotd[6]/*degrees*/,
333                                    Int_t idt,const Int_t id[3],
334                                    const Double_t tran[3]):
335 TObject(),            // Base class
336 fDetectorIndex(idt), // Detector type index (like fShapeIndex was)
337 fid(),       // layer, ladder, detector numbers.
338 frot(),      //! vector of rotations about x,y,z [radians].
339 ftran(),     // Translation vector of module x,y,z.
340 fCylR(0.0),  //! R Translation in Cylinderical coordinates
341 fCylPhi(0.0),//! Phi Translation vector in Cylindrical coord.
342 fm(),        // Rotation matrix based on frot.
343 fPath(){     // Path in geometry to this module
344     // This is a constructor for the AliITSgeomMatrix class. The matrix 
345     // is defined by the 6 GEANT 3.21 rotation angles [degrees], and 
346     // the translation vector tran [cm]. In addition the layer, ladder, 
347     // and detector number for this particular module and the type of 
348     // module must be given. The full rotation matrix is kept so that 
349     // the evaluation  of a coordinate transformation can be done 
350     // quickly and with a minimum of CPU overhead. The basic coordinate 
351     // systems are the ALICE global coordinate system and the detector 
352     // local coordinate system. In general this structure is not limited 
353     // to just those two coordinate systems.
354     //Begin_Html
355     /*
356       <img src="picts/ITS/AliITSgeomMatrix_L1.gif">
357     */
358     //End_Html
359     // Inputs:
360     //    Double_t rotd[6]  The 6 Geant 3.21 rotation angles [degrees]
361     //    Int_t idt         The module Id number
362     //    Int_t id[3]       The layer, ladder and detector number
363     //    Double_t tran[3]  The translation vector
364     Int_t i;
365
366     for(i=0;i<3;i++){
367         fid[i]   = id[i];
368         ftran[i] = tran[i];
369     }// end for i
370     fCylR   = TMath::Sqrt(ftran[0]*ftran[0]+ftran[1]*ftran[1]);
371     fCylPhi = TMath::ATan2(ftran[1],ftran[0]);
372     if(fCylPhi<0.0) fCylPhi += 2.*TMath::Pi();
373     this->MatrixFromSixAngles(rotd);
374 }
375 //----------------------------------------------------------------------
376 void AliITSgeomMatrix::AngleFromMatrix(){
377     // Computes the angles from the rotation matrix up to a phase of 
378     // 180 degrees. The matrix used in AliITSgeomMatrix::MatrixFromAngle()
379     // and  its inverse AliITSgeomMatrix::AngleFromMatrix() are defined in 
380     // the following ways, R = Rz*Ry*Rx (M=R*L+T) where
381     //     1   0   0       Cy  0 +Sy       Cz -Sz  0
382     // Rx= 0   Cx -Sx  Ry=  0  1   0   Rz=+Sz  Cz  0
383     //     0  +Sx  Cx     -Sy  0  Cy        0   0  1
384     // The choice of the since of S, comes from the choice between 
385     // the rotation of the object or the coordinate system (view). I think
386     // that this choice is the first, the rotation of the object.
387     // Inputs:
388     //   none
389     // Outputs:
390     //   none
391     // Return:
392     //   none
393     Double_t rx,ry,rz;
394     // get angles from matrix up to a phase of 180 degrees.
395
396     rx = TMath::ATan2(fm[2][1],fm[2][2]);if(rx<0.0) rx += 2.0*TMath::Pi();
397     ry = TMath::ASin(-fm[0][2]);         if(ry<0.0) ry += 2.0*TMath::Pi();
398     rz = TMath::ATan2(fm[1][0],fm[0][0]);if(rz<0.0) rz += 2.0*TMath::Pi();
399     frot[0] = rx;
400     frot[1] = ry;
401     frot[2] = rz;
402     return;
403 }
404 //----------------------------------------------------------------------
405 void AliITSgeomMatrix::MatrixFromAngle(){
406     // Computes the Rotation matrix from the angles [radians] kept in this
407     // class. The matrix used in AliITSgeomMatrix::MatrixFromAngle() and 
408     // its inverse AliITSgeomMatrix::AngleFromMatrix() are defined in 
409     // the following ways, R = Rz*Ry*Rx (M=R*L+T) where
410     //     1   0   0       Cy  0 +Sy       Cz -Sz  0
411     // Rx= 0   Cx -Sx  Ry=  0  1   0   Rz=+Sz  Cz  0
412     //     0  +Sx  Cx     -Sy  0  Cy        0   0  1
413     // The choice of the since of S, comes from the choice between 
414     // the rotation of the object or the coordinate system (view). I think
415     // that this choice is the first, the rotation of the object.
416     // Inputs:
417     //   none
418     // Outputs:
419     //   none
420     // Return:
421     //   none
422     Double_t sx,sy,sz,cx,cy,cz;
423
424     sx = TMath::Sin(frot[0]); cx = TMath::Cos(frot[0]);
425     sy = TMath::Sin(frot[1]); cy = TMath::Cos(frot[1]);
426     sz = TMath::Sin(frot[2]); cz = TMath::Cos(frot[2]);
427     fm[0][0] = +cz*cy;             // fr[0]
428     fm[0][1] = +cz*sy*sx - sz*cx;  // fr[1]
429     fm[0][2] = +cz*sy*cx + sz*sx;  // fr[2]
430     fm[1][0] = +sz*cy;             // fr[3]
431     fm[1][1] = +sz*sy*sx + cz*cx;  // fr[4]
432     fm[1][2] = +sz*sy*cx - cz*sx;  // fr[5]
433     fm[2][0] = -sy;                // fr[6]
434     fm[2][1] = +cy*sx;             // fr[7]
435     fm[2][2] = +cy*cx;             // fr[8]
436 }
437 //----------------------------------------------------------------------
438 void AliITSgeomMatrix::SetEulerAnglesChi(const Double_t ang[3]){
439     // Computes the Rotation matrix from the Euler angles [radians], 
440     // Chi-convention, kept in this class. The matrix used in 
441     // AliITSgeomMatrix::SetEulerAnglesChi and 
442     // its inverse AliITSgeomMatrix::GetEulerAnglesChi() are defined in 
443     // the following ways, R = Rb*Rc*Rd (M=R*L+T) where
444     //     C2 +S2  0       1  0    0       C0 +S0  0
445     // Rb=-S2  C2  0  Rc=  0  C1 +S1   Rd=-S0  C0  0
446     //     0   0   1       0 -S1  C1       0   0   1
447     // This form is taken from Wolfram Research's Geometry>
448     // Transformations>Rotations web page (also should be
449     // found in their book).
450     // Inputs:
451     //   Double_t ang[3] The three Euler Angles Phi, Theta, Psi
452     // Outputs:
453     //   none
454     // Return:
455     //   none
456     Double_t s0,s1,s2,c0,c1,c2;
457
458     s0 = TMath::Sin(ang[0]); c0 = TMath::Cos(ang[0]);
459     s1 = TMath::Sin(ang[1]); c1 = TMath::Cos(ang[1]);
460     s2 = TMath::Sin(ang[2]); c2 = TMath::Cos(ang[2]);
461     fm[0][0] = +c2*c0-c1*s0*s2;  // fr[0]
462     fm[0][1] = +c2*s0+c1*c0*s2;  // fr[1]
463     fm[0][2] = +s2*s1;           // fr[2]
464     fm[1][0] = -s2*c0-c1*s0*c2;  // fr[3]
465     fm[1][1] = -s2*s0+c1*c0*c2;  // fr[4]
466     fm[1][2] = +c2*s1;           // fr[5]
467     fm[2][0] = s1*s0;            // fr[6]
468     fm[2][1] = -s1*c0;           // fr[7]
469     fm[2][2] = +c1;              // fr[8]
470     AngleFromMatrix();
471     return ;
472 }
473 //----------------------------------------------------------------------
474 void AliITSgeomMatrix::GtoLPosition(const Double_t g0[3],Double_t l[3]) const {
475     // Returns the local coordinates given the global coordinates [cm].
476     // Inputs:
477     //   Double_t g[3]   The position represented in the ALICE 
478     //                   global coordinate system
479     // Outputs:
480     //   Double_t l[3]  The poistion represented in the local
481     //                  detector coordiante system
482     // Return:
483     //   none
484     Int_t    i,j;
485     Double_t g[3];
486
487     for(i=0;i<3;i++) g[i] = g0[i] - ftran[i];
488     for(i=0;i<3;i++){
489         l[i] = 0.0;
490         for(j=0;j<3;j++) l[i] += fm[i][j]*g[j];
491         // g = R l + translation
492     } // end for i
493     return;
494 }
495 //----------------------------------------------------------------------
496 void AliITSgeomMatrix::LtoGPosition(const Double_t l[3],Double_t g[3]) const {
497     // Returns the global coordinates given the local coordinates [cm].
498     // Inputs:
499     //   Double_t l[3]   The poistion represented in the detector 
500     //                   local coordinate system
501     // Outputs:
502     //   Double_t g[3]   The poistion represented in the ALICE
503     //                   Global coordinate system
504     // Return:
505     //   none.
506     Int_t    i,j;
507
508     for(i=0;i<3;i++){
509         g[i] = 0.0;
510         for(j=0;j<3;j++) g[i] += fm[j][i]*l[j];
511         g[i] += ftran[i];
512         // g = R^t l + translation
513     } // end for i
514     return;
515 }
516 //----------------------------------------------------------------------
517 void AliITSgeomMatrix::GtoLMomentum(const Double_t g[3],Double_t l[3]) const{
518     // Returns the local coordinates of the momentum given the global
519     // coordinates of the momentum. It transforms just like GtoLPosition
520     // except that the translation vector is zero.
521     // Inputs:
522     //   Double_t g[3] The momentum represented in the ALICE global 
523     //                 coordinate system
524     // Outputs:
525     //   Double_t l[3] the momentum represented in the detector 
526     //                 local coordinate system
527     // Return:
528     //   none.
529     Int_t    i,j;
530
531     for(i=0;i<3;i++){
532         l[i] = 0.0;
533         for(j=0;j<3;j++) l[i] += fm[i][j]*g[j];
534         // g = R l
535     } // end for i
536     return;
537 }
538 //----------------------------------------------------------------------
539 void AliITSgeomMatrix::LtoGMomentum(const Double_t l[3],Double_t g[3]) const {
540     // Returns the Global coordinates of the momentum given the local
541     // coordinates of the momentum. It transforms just like LtoGPosition
542     // except that the translation vector is zero.
543     // Inputs:
544     //   Double_t l[3] the momentum represented in the detector 
545     //                 local coordinate system
546     // Outputs:
547     //   Double_t g[3] The momentum represented in the ALICE global 
548     //                 coordinate system
549     // Return:
550     //   none.
551     Int_t    i,j;
552
553     for(i=0;i<3;i++){
554         g[i] = 0.0;
555         for(j=0;j<3;j++) g[i] += fm[j][i]*l[j];
556         // g = R^t l
557     } // end for i
558     return;
559 }
560 //----------------------------------------------------------------------
561 void AliITSgeomMatrix::GtoLPositionError(const Double_t g[3][3],
562                                          Double_t l[3][3]) const {
563     // Given an Uncertainty matrix in Global coordinates it is 
564     // rotated so that  its representation in local coordinates can 
565     // be returned. There is no effect due to the translation vector 
566     // or its uncertainty.
567     // Inputs:
568     //   Double_t g[3][3] The error matrix represented in the ALICE global 
569     //                    coordinate system
570     // Outputs:
571     //   Double_t l[3][3] the error matrix represented in the detector 
572     //                    local coordinate system
573     // Return:
574     //   none.
575     Int_t    i,j,k,m;
576
577     for(i=0;i<3;i++)for(m=0;m<3;m++){
578         l[i][m] = 0.0;
579         for(j=0;j<3;j++)for(k=0;k<3;k++)
580             l[i][m] += fm[j][i]*g[j][k]*fm[k][m];
581     } // end for i,m
582     // g = R^t l R
583     return;
584 }
585 //----------------------------------------------------------------------
586 void AliITSgeomMatrix::LtoGPositionError(const Double_t l[3][3],
587                                                Double_t g[3][3]) const {
588     // Given an Uncertainty matrix in Local coordinates it is rotated so that 
589     // its representation in global coordinates can be returned. There is no
590     // effect due to the translation vector or its uncertainty.
591     // Inputs:
592     //   Double_t l[3][3] the error matrix represented in the detector 
593     //                    local coordinate system
594     // Outputs:
595     //   Double_t g[3][3] The error matrix represented in the ALICE global 
596     //                    coordinate system
597     // Return:
598     //   none.
599     Int_t    i,j,k,m;
600
601     for(i=0;i<3;i++)for(m=0;m<3;m++){
602         g[i][m] = 0.0;
603         for(j=0;j<3;j++)for(k=0;k<3;k++)
604             g[i][m] += fm[i][j]*l[j][k]*fm[m][k];
605     } // end for i,m
606     // g = R l R^t
607     return;
608 }
609 //----------------------------------------------------------------------
610 void AliITSgeomMatrix::GtoLPositionTracking(const Double_t g[3],
611                                             Double_t l[3]) const {
612     // A slightly different coordinate system is used when tracking.
613     // This coordinate system is only relevant when the geometry represents
614     // the cylindrical ALICE ITS geometry. For tracking the Z axis is left
615     // alone but X -> -Y and Y -> X such that X always points out of the
616     // ITS Cylinder for every layer including layer 1 (where the detector 
617     // are mounted upside down).
618     //Begin_Html
619     /*
620       <img src="picts/ITS/AliITSgeomMatrix_T1.gif">
621     */
622     //End_Html
623     // Inputs:
624     //   Double_t g[3]   The position represented in the ALICE 
625     //                   global coordinate system
626     // Outputs:
627     //   Double_t l[3]  The poistion represented in the local
628     //                  detector coordiante system
629     // Return:
630     //   none
631     Double_t l0[3];
632
633     this->GtoLPosition(g,l0);
634     if(fid[0]==1){ // for layer 1 the detector are flipped upside down
635                    // with respect to the others.
636         l[0] = +l0[1];
637         l[1] = -l0[0];
638         l[2] = +l0[2];
639     }else{
640         l[0] = -l0[1];
641         l[1] = +l0[0];
642         l[2] = +l0[2];
643     } // end if
644     return;
645 }
646 //----------------------------------------------------------------------
647 void AliITSgeomMatrix::LtoGPositionTracking(const Double_t l[3],
648                                             Double_t g[3]) const {
649     // A slightly different coordinate system is used when tracking.
650     // This coordinate system is only relevant when the geometry represents
651     // the cylindrical ALICE ITS geometry. For tracking the Z axis is left
652     // alone but X -> -Y and Y -> X such that X always points out of the
653     // ITS Cylinder for every layer including layer 1 (where the detector 
654     // are mounted upside down).
655     //Begin_Html
656     /*
657       <img src="picts/ITS/AliITSgeomMatrix_T1.gif">
658     */
659     //End_Html
660     // Inputs:
661     //   Double_t l[3]   The poistion represented in the detector 
662     //                   local coordinate system
663     // Outputs:
664     //   Double_t g[3]   The poistion represented in the ALICE
665     //                   Global coordinate system
666     // Return:
667     //   none.
668     Double_t l0[3];
669
670     if(fid[0]==1){ // for layer 1 the detector are flipped upside down
671                    // with respect to the others.
672         l0[0] = -l[1];
673         l0[1] = +l[0];
674         l0[2] = +l[2];
675     }else{
676         l0[0] = +l[1];
677         l0[1] = -l[0];
678         l0[2] = +l[2];
679     } // end if
680     this->LtoGPosition(l0,g);
681     return;
682 }
683 //----------------------------------------------------------------------
684 void AliITSgeomMatrix::GtoLMomentumTracking(const Double_t g[3],
685                                             Double_t l[3]) const {
686     // A slightly different coordinate system is used when tracking.
687     // This coordinate system is only relevant when the geometry represents
688     // the cylindrical ALICE ITS geometry. For tracking the Z axis is left
689     // alone but X -> -Y and Y -> X such that X always points out of the
690     // ITS Cylinder for every layer including layer 1 (where the detector 
691     // are mounted upside down).
692     //Begin_Html
693     /*
694       <img src="picts/ITS/AliITSgeomMatrix_T1.gif">
695     */
696     //End_Html
697     // Inputs:
698     //   Double_t g[3] The momentum represented in the ALICE global 
699     //                 coordinate system
700     // Outputs:
701     //   Double_t l[3] the momentum represented in the detector 
702     //                 local coordinate system
703     // Return:
704     //   none.
705     Double_t l0[3];
706
707     this->GtoLMomentum(g,l0);
708     if(fid[0]==1){ // for layer 1 the detector are flipped upside down
709                    // with respect to the others.
710         l[0] = +l0[1];
711         l[1] = -l0[0];
712         l[2] = +l0[2];
713     }else{
714         l[0] = -l0[1];
715         l[1] = +l0[0];
716         l[2] = +l0[2];
717     } // end if
718     return;
719 }
720 //----------------------------------------------------------------------
721 void AliITSgeomMatrix::LtoGMomentumTracking(const Double_t l[3],
722                                             Double_t g[3]) const {
723     // A slightly different coordinate system is used when tracking.
724     // This coordinate system is only relevant when the geometry represents
725     // the cylindrical ALICE ITS geometry. For tracking the Z axis is left
726     // alone but X -> -Y and Y -> X such that X always points out of the
727     // ITS Cylinder for every layer including layer 1 (where the detector 
728     // are mounted upside down).
729     //Begin_Html
730     /*
731       <img src="picts/ITS/AliITSgeomMatrix_T1.gif">
732     */
733     //End_Html
734     // Inputs:
735     //   Double_t l[3] the momentum represented in the detector 
736     //                 local coordinate system
737     // Outputs:
738     //   Double_t g[3] The momentum represented in the ALICE global 
739     //                 coordinate system
740     // Return:
741     //   none.
742     Double_t l0[3];
743
744     if(fid[0]==1){ // for layer 1 the detector are flipped upside down
745                    // with respect to the others.
746         l0[0] = -l[1];
747         l0[1] = +l[0];
748         l0[2] = +l[2];
749     }else{
750         l0[0] = +l[1];
751         l0[1] = -l[0];
752         l0[2] = +l[2];
753     } // end if
754     this->LtoGMomentum(l0,g);
755         return;
756 }
757 //----------------------------------------------------------------------
758 void AliITSgeomMatrix::GtoLPositionErrorTracking(const Double_t g[3][3],
759                                                  Double_t l[3][3]) const {
760     // A slightly different coordinate system is used when tracking.
761     // This coordinate system is only relevant when the geometry represents
762     // the cylindrical ALICE ITS geometry. For tracking the Z axis is left
763     // alone but X -> -Y and Y -> X such that X always points out of the
764     // ITS Cylinder for every layer including layer 1 (where the detector 
765     // are mounted upside down).
766     //Begin_Html
767     /*
768       <img src="picts/ITS/AliITSgeomMatrix_TE1.gif">
769     */
770     //End_Html
771     // Inputs:
772     //   Double_t g[3][3] The error matrix represented in the ALICE global 
773     //                    coordinate system
774     // Outputs:
775     //   Double_t l[3][3] the error matrix represented in the detector 
776     //                    local coordinate system
777     // Return:
778     Int_t    i,j,k,m;
779     Double_t rt[3][3];
780     Double_t a0[3][3] = {{0.,+1.,0.},{-1.,0.,0.},{0.,0.,+1.}};
781     Double_t a1[3][3] = {{0.,-1.,0.},{+1.,0.,0.},{0.,0.,+1.}};
782
783     if(fid[0]==1) for(i=0;i<3;i++)for(j=0;j<3;j++)for(k=0;k<3;k++)
784         rt[i][k] = a0[i][j]*fm[j][k];
785     else for(i=0;i<3;i++)for(j=0;j<3;j++)for(k=0;k<3;k++)
786         rt[i][k] = a1[i][j]*fm[j][k];
787     for(i=0;i<3;i++)for(m=0;m<3;m++){
788         l[i][m] = 0.0;
789         for(j=0;j<3;j++)for(k=0;k<3;k++)
790             l[i][m] += rt[j][i]*g[j][k]*rt[k][m];
791     } // end for i,m
792     // g = R^t l R
793     return;
794 }
795 //----------------------------------------------------------------------
796 void AliITSgeomMatrix::LtoGPositionErrorTracking(const Double_t l[3][3],
797                                                  Double_t g[3][3]) const {
798     // A slightly different coordinate system is used when tracking.
799     // This coordinate system is only relevant when the geometry represents
800     // the cylindrical ALICE ITS geometry. For tracking the Z axis is left
801     // alone but X -> -Y and Y -> X such that X always points out of the
802     // ITS Cylinder for every layer including layer 1 (where the detector 
803     // are mounted upside down).
804     //Begin_Html
805     /*
806       <img src="picts/ITS/AliITSgeomMatrix_TE1.gif">
807     */
808     //End_Html
809     // Inputs:
810     //   Double_t l[3][3] the error matrix represented in the detector 
811     //                    local coordinate system
812     // Outputs:
813     //   Double_t g[3][3] The error matrix represented in the ALICE global 
814     //                    coordinate system
815     // Return:
816     //   none.
817     Int_t    i,j,k,m;
818     Double_t rt[3][3];
819     Double_t a0[3][3] = {{0.,+1.,0.},{-1.,0.,0.},{0.,0.,+1.}};
820     Double_t a1[3][3] = {{0.,-1.,0.},{+1.,0.,0.},{0.,0.,+1.}};
821
822     if(fid[0]==1) for(i=0;i<3;i++)for(j=0;j<3;j++)for(k=0;k<3;k++)
823         rt[i][k] = a0[i][j]*fm[j][k];
824     else for(i=0;i<3;i++)for(j=0;j<3;j++)for(k=0;k<3;k++)
825         rt[i][k] = a1[i][j]*fm[j][k];
826     for(i=0;i<3;i++)for(m=0;m<3;m++){
827         g[i][m] = 0.0;
828         for(j=0;j<3;j++)for(k=0;k<3;k++)
829             g[i][m] += rt[i][j]*l[j][k]*rt[m][k];
830     } // end for i,m
831     // g = R l R^t
832     return;
833 }
834 //----------------------------------------------------------------------
835 void AliITSgeomMatrix::PrintTitles(ostream *os) const {
836     // Standard output format for this class but it includes variable
837     // names and formatting that makes it easer to read.
838     // Inputs:
839     //    ostream *os   The output stream to print the title on
840     // Outputs:
841     //    none.
842     // Return:
843     //    none.
844     Int_t i,j;
845
846     *os << "fDetectorIndex=" << fDetectorIndex << " fid[3]={";
847     for(i=0;i<3;i++) *os << fid[i]   << " ";
848     *os << "} frot[3]={";
849     for(i=0;i<3;i++) *os << frot[i]  << " ";
850     *os << "} ftran[3]={";
851     for(i=0;i<3;i++) *os << ftran[i] << " ";
852     *os << "} fm[3][3]={";
853     for(i=0;i<3;i++){for(j=0;j<3;j++){  *os << fm[i][j] << " ";} *os <<"}{";}
854     *os << "}" << endl;
855     return;
856 }
857 //----------------------------------------------------------------------
858 void AliITSgeomMatrix::PrintComment(ostream *os) const {
859     //  output format used by Print.
860     // Inputs:
861     //    ostream *os   The output stream to print the comments on
862     // Outputs:
863     //    none.
864     // Return:
865     //    none.
866     *os << "fDetectorIndex fid[0] fid[1] fid[2] ftran[0] ftran[1] ftran[2] ";
867     *os << "fm[0][0]  fm[0][1]  fm[0][2]  fm[1][0]  fm[1][1]  fm[1][2]  ";
868     *os << "fm[2][0]  fm[2][1]  fm[2][2] ";
869     return;
870 }
871 //----------------------------------------------------------------------
872 void AliITSgeomMatrix::Print(ostream *os)const{
873     // Standard output format for this class.
874     // Inputs:
875     //    ostream *os   The output stream to print the class data on
876     // Outputs:
877     //    none.
878     // Return:
879     //    none.
880     Int_t i,j;
881 #if defined __GNUC__
882 #if __GNUC__ > 2
883     ios::fmtflags fmt;
884 #else
885     Int_t fmt;
886 #endif
887 #else
888 #if defined __ICC || defined __ECC || defined __xlC__
889     ios::fmtflags fmt;
890 #else
891     Int_t fmt;
892 #endif
893 #endif
894
895     fmt = os->setf(ios::scientific);  // set scientific floating point output
896     *os << fDetectorIndex << " ";
897     for(i=0;i<3;i++) *os << fid[i]   << " ";
898 //    for(i=0;i<3;i++) *os << frot[i]  << " ";  // Redundant with fm[][].
899     for(i=0;i<3;i++) *os << setprecision(16) << ftran[i] << " ";
900     for(i=0;i<3;i++)for(j=0;j<3;j++)  *os << setprecision(16) << 
901                                           fm[i][j] << " ";
902     *os << fPath.Length()<< " ";
903     for(i=0;i<fPath.Length();i++) *os << fPath[i];
904     *os << endl;
905     os->flags(fmt); // reset back to old formating.
906     return;
907 }
908 //----------------------------------------------------------------------
909 void AliITSgeomMatrix::Read(istream *is){
910     // Standard input format for this class.
911     // Inputs:
912     //    istream *is   The input stream to read on
913     // Outputs:
914     //    none.
915     // Return:
916     //    none.
917     Int_t i,j;
918
919     *is >> fDetectorIndex;
920     for(i=0;i<3;i++) *is >> fid[i];
921 //    for(i=0;i<3;i++) *is >> frot[i]; // Redundant with fm[][].
922     for(i=0;i<3;i++) *is >> ftran[i];
923     for(i=0;i<3;i++)for(j=0;j<3;j++)  *is >> fm[i][j];
924     while(is->peek()==' ')is->get(); // skip white spaces
925     if(isprint(is->peek())){ // old format did not have path.
926         *is >> j; // string length
927         fPath.Resize(j);
928         for(i=0;i<j;i++) {*is >> fPath[i];}
929     } // end if
930     AngleFromMatrix(); // compute angles frot[].
931     fCylR   = TMath::Sqrt(ftran[0]*ftran[0]+ftran[1]*ftran[1]);
932     fCylPhi = TMath::ATan2(ftran[1],ftran[0]);
933     if(fCylPhi<0.0) fCylPhi += 2.*TMath::Pi();
934     return;
935 }
936 //______________________________________________________________________
937 void AliITSgeomMatrix::Streamer(TBuffer &R__b){
938    // Stream an object of class AliITSgeomMatrix.
939     // Inputs:
940     //     TBuffer &R__b   The output buffer to stream data on.
941     // Outputs:
942     //    none.
943     // Return:
944     //    none.
945
946     if (R__b.IsReading()) {
947         AliITSgeomMatrix::Class()->ReadBuffer(R__b, this);
948         fCylR   = TMath::Sqrt(ftran[0]*ftran[0]+ftran[1]*ftran[1]);
949         fCylPhi = TMath::ATan2(ftran[1],ftran[0]);
950         this->AngleFromMatrix();
951         if(fCylPhi<0.0) fCylPhi += 2.*TMath::Pi();
952     } else {
953         AliITSgeomMatrix::Class()->WriteBuffer(R__b, this);
954     } // end if
955 }
956 //______________________________________________________________________
957 void AliITSgeomMatrix::SetTranslation(const Double_t tran[3]){
958     // Sets the translation vector and computes fCylR and fCylPhi.
959     // Inputs:
960     //   Double_t trans[3]   The translation vector to be used
961     // Outputs:
962     //   none.
963     // Return:
964     //   none.
965     for(Int_t i=0;i<3;i++) ftran[i] = tran[i];
966     fCylR   = TMath::Sqrt(ftran[0]*ftran[0]+ftran[1]*ftran[1]);
967     fCylPhi = TMath::ATan2(ftran[1],ftran[0]);
968     if(fCylPhi<0.0) fCylPhi += 2.*TMath::Pi();
969 }
970 //----------------------------------------------------------------------
971 TPolyLine3D* AliITSgeomMatrix::CreateLocalAxis() const {
972     // This class is used as part of the documentation of this class
973     // Inputs:
974     //   none.
975     // Outputs:
976     //   none.
977     // Return:
978     //   A pointer to a new TPolyLine3D object showing the 3 line
979     //   segments that make up the this local axis in the global
980     //   reference system.
981     Float_t  gf[15];
982     Double_t g[5][3];
983     Double_t l[5][3]={{1.0,0.0,0.0},{0.0,0.0,0.0},{0.0,1.0,0.0},{0.0,0.0,0.0},
984                       {0.0,0.0,1.0}};
985     Int_t i;
986
987     for(i=0;i<5;i++) {
988         LtoGPosition(l[i],g[i]);
989         gf[3*i]=(Float_t)g[i][0];
990         gf[3*i+1]=(Float_t)g[i][1];
991         gf[3*i+2]=(Float_t)g[i][2];
992     } // end for i
993     return new TPolyLine3D(5,gf);
994 }
995 //----------------------------------------------------------------------
996 TPolyLine3D* AliITSgeomMatrix::CreateLocalAxisTracking() const {
997     // This class is used as part of the documentation of this class
998     // Inputs:
999     //   none.
1000     // Outputs:
1001     //   none.
1002     // Return:
1003     //   A pointer to a new TPolyLine3D object showing the 3 line
1004     //   segments that make up the this local axis in the global
1005     //   reference system.
1006     Float_t gf[15];
1007     Double_t g[5][3];
1008     Double_t l[5][3]={{1.0,0.0,0.0},{0.0,0.0,0.0},{0.0,1.0,0.0},{0.0,0.0,0.0},
1009                       {0.0,0.0,1.0}};
1010     Int_t i;
1011
1012     for(i=0;i<5;i++) {
1013         LtoGPositionTracking(l[i],g[i]);
1014         gf[3*i]=(Float_t)g[i][0];
1015         gf[3*i+1]=(Float_t)g[i][1];
1016         gf[3*i+2]=(Float_t)g[i][2];
1017     } // end for i
1018     return new TPolyLine3D(5,gf);
1019 }
1020 //----------------------------------------------------------------------
1021 TNode* AliITSgeomMatrix::CreateNode(const Char_t *nodeName,
1022                                     const Char_t *nodeTitle,TNode *mother,
1023                                     TShape *shape,Bool_t axis) const {
1024     // Creates a node inside of the node mother out of the shape shape
1025     // in the position, with respect to mother, indecated by "this". If axis
1026     // is ture, it will insert an axis within this node/shape.
1027     // Inputs:
1028     //   Char_t *nodeName  This name of this node
1029     //   Char_t *nodeTitle This node title
1030     //   TNode  *mother    The node this node will be inside of/with respect to
1031     //   TShape *shape     The shape of this node
1032     //   Bool_t axis       If ture, a set of x,y,z axis will be included
1033     // Outputs:
1034     //   none.
1035     // Return:
1036     //   A pointer to "this" node.
1037     Double_t trans[3],matrix[3][3],*matr;
1038     TRotMatrix *rot = new TRotMatrix();
1039     TString name,title;
1040
1041     matr = &(matrix[0][0]);
1042     this->GetTranslation(trans);
1043     this->GetMatrix(matrix);
1044     rot->SetMatrix(matr);
1045     //
1046     name = nodeName;
1047     title = nodeTitle;
1048     //
1049     mother->cd();
1050     TNode *node1 = new TNode(name.Data(),title.Data(),shape,
1051                              trans[0],trans[1],trans[2],rot);
1052     if(axis){
1053         Int_t i,j;
1054         const Float_t kScale=0.5,kLw=0.2;
1055         Float_t xchar[13][2]={
1056             {0.5*kLw,1.},{0.,0.5*kLw},{0.5-0.5*kLw,0.5},
1057             {0.,0.5*kLw},{0.5*kLw,0.},{0.5,0.5-0.5*kLw},
1058             {1-0.5*kLw,0.},{1.,0.5*kLw},{0.5+0.5*kLw,0.5},
1059             {1.,1.-0.5*kLw},{1.-0.5*kLw,1.},{0.5,0.5+0.5*kLw},
1060             {0.5*kLw,1.}};
1061         Float_t ychar[10][2]={
1062             {.5-0.5*kLw,0.},{.5+0.5*kLw,0.},{.5+0.5*kLw,0.5-0.5*kLw},
1063             {1.,1.-0.5*kLw},{1.-0.5*kLw,1.},{0.5+0.5*kLw,0.5},
1064             {0.5*kLw,1.}   ,{0.,1-0.5*kLw} ,{0.5-0.5*kLw,0.5},
1065             {.5-0.5*kLw,0.}};
1066         Float_t zchar[11][2]={
1067             {0.,1.},{0,1.-kLw},{1.-kLw,1.-kLw},{0.,kLw}   ,{0.,0.},
1068             {1.,0.},{1.,kLw}  ,{kLw,kLw}      ,{1.,1.-kLw},{1.,1.},
1069             {0.,1.}};
1070         for(i=0;i<13;i++)for(j=0;j<2;j++){
1071             if(i<13) xchar[i][j] = kScale*xchar[i][j];
1072             if(i<10) ychar[i][j] = kScale*ychar[i][j];
1073             if(i<11) zchar[i][j] = kScale*zchar[i][j];
1074         } // end for i,j
1075         TXTRU *axisxl = new TXTRU("x","x","text",12,2);
1076         for(i=0;i<12;i++) axisxl->DefineVertex(i,xchar[i][0],xchar[i][1]);
1077         axisxl->DefineSection(0,-0.5*kLw);axisxl->DefineSection(1,0.5*kLw);
1078         TXTRU *axisyl = new TXTRU("y","y","text",9,2);
1079         for(i=0;i<9;i++) axisyl->DefineVertex(i,ychar[i][0],ychar[i][1]);
1080         axisyl->DefineSection(0,-0.5*kLw);axisyl->DefineSection(1,0.5*kLw);
1081         TXTRU *axiszl = new TXTRU("z","z","text",10,2);
1082         for(i=0;i<10;i++) axiszl->DefineVertex(i,zchar[i][0],zchar[i][1]);
1083         axiszl->DefineSection(0,-0.5*kLw);axiszl->DefineSection(1,0.5*kLw);
1084         Float_t lxy[13][2]={
1085             {-0.5*kLw,-0.5*kLw},{0.8,-0.5*kLw},{0.8,-0.1},{1.0,0.0},
1086             {0.8,0.1},{0.8,0.5*kLw},{0.5*kLw,0.5*kLw},{0.5*kLw,0.8},
1087             {0.1,0.8},{0.0,1.0},{-0.1,0.8},{-0.5*kLw,0.8},
1088             {-0.5*kLw,-0.5*kLw}};
1089         TXTRU *axisxy = new TXTRU("axisxy","axisxy","text",13,2);
1090         for(i=0;i<13;i++) axisxy->DefineVertex(i,lxy[i][0],lxy[i][1]);
1091         axisxy->DefineSection(0,-0.5*kLw);axisxy->DefineSection(1,0.5*kLw);
1092         Float_t lz[8][2]={
1093             {0.5*kLw,-0.5*kLw},{0.8,-0.5*kLw},{0.8,-0.1},{1.0,0.0},
1094             {0.8,0.1},{0.8,0.5*kLw},{0.5*kLw,0.5*kLw},
1095             {0.5*kLw,-0.5*kLw}};
1096         TXTRU *axisz = new TXTRU("axisz","axisz","text",8,2);
1097         for(i=0;i<8;i++) axisz->DefineVertex(i,lz[i][0],lz[i][1]);
1098         axisz->DefineSection(0,-0.5*kLw);axisz->DefineSection(1,0.5*kLw);
1099         //TRotMatrix *xaxis90= new TRotMatrix("xaixis90","",90.0, 0.0, 0.0);
1100         TRotMatrix *yaxis90= new TRotMatrix("yaixis90","", 0.0,90.0, 0.0);
1101         TRotMatrix *zaxis90= new TRotMatrix("zaixis90","", 0.0, 0.0,90.0);
1102         //
1103         node1->cd();
1104         title = name.Append("axisxy");
1105         TNode *nodeaxy = new TNode(title.Data(),title.Data(),axisxy);
1106         title = name.Append("axisz");
1107         TNode *nodeaz = new TNode(title.Data(),title.Data(),axisz,
1108                                   0.,0.,0.,yaxis90);
1109         TNode *textboxX0 = new TNode("textboxX0","textboxX0",axisxl,
1110                                     lxy[3][0],lxy[3][1],0.0);
1111         TNode *textboxX1 = new TNode("textboxX1","textboxX1",axisxl,
1112                                     lxy[3][0],lxy[3][1],0.0,yaxis90);
1113         TNode *textboxX2 = new TNode("textboxX2","textboxX2",axisxl,
1114                                     lxy[3][0],lxy[3][1],0.0,zaxis90);
1115         TNode *textboxY0 = new TNode("textboxY0","textboxY0",axisyl,
1116                                     lxy[9][0],lxy[9][1],0.0);
1117         TNode *textboxY1 = new TNode("textboxY1","textboxY1",axisyl,
1118                                     lxy[9][0],lxy[9][1],0.0,yaxis90);
1119         TNode *textboxY2 = new TNode("textboxY2","textboxY2",axisyl,
1120                                     lxy[9][0],lxy[9][1],0.0,zaxis90);
1121         TNode *textboxZ0 = new TNode("textboxZ0","textboxZ0",axiszl,
1122                                     0.0,0.0,lz[3][0]);
1123         TNode *textboxZ1 = new TNode("textboxZ1","textboxZ1",axiszl,
1124                                     0.0,0.0,lz[3][0],yaxis90);
1125         TNode *textboxZ2 = new TNode("textboxZ2","textboxZ2",axiszl,
1126                                     0.0,0.0,lz[3][0],zaxis90);
1127         nodeaxy->Draw();
1128         nodeaz->Draw();
1129         textboxX0->Draw();
1130         textboxX1->Draw();
1131         textboxX2->Draw();
1132         textboxY0->Draw();
1133         textboxY1->Draw();
1134         textboxY2->Draw();
1135         textboxZ0->Draw();
1136         textboxZ1->Draw();
1137         textboxZ2->Draw();
1138     } // end if
1139     mother->cd();
1140     return node1;
1141 }
1142 //----------------------------------------------------------------------
1143 void AliITSgeomMatrix::MakeFigures() const {
1144     // make figures to help document this class
1145     // Inputs:
1146     //   none.
1147     // Outputs:
1148     //   none.
1149     // Return:
1150     //   none.
1151     const Double_t kDx0=550.,kDy0=550.,kDz0=550.; // cm
1152     const Double_t kDx=1.0,kDy=0.300,kDz=3.0,kRmax=0.1; // cm
1153     Float_t l[5][3]={{1.0,0.0,0.0},{0.0,0.0,0.0},{0.0,1.0,0.0},{0.0,0.0,0.0},
1154                       {0.0,0.0,1.0}};
1155     TCanvas *c = new TCanvas(kFALSE);// create a batch mode canvas.
1156 #if ROOT_VERSION_CODE>= 331523
1157     Double_t rmin[]={-1,-1,-1};
1158     Double_t rmax[]={ 1, 1, 1};
1159     TView *view = new TView3D(1,rmin,rmax);
1160 #else
1161     TView   *view = new TView(1); // Create Cartesian coordiante view
1162 #endif
1163     TBRIK   *mother  = new TBRIK("Mother","Mother","void",kDx0,kDy0,kDz0);
1164     TBRIK   *det  = new TBRIK("Detector","","Si",kDx,kDy,kDz);
1165     TPolyLine3D *axis = new TPolyLine3D(5,&(l[0][0]));
1166     TPCON *arrow      = new TPCON("arrow","","air",0.0,360.,2);
1167     TRotMatrix *xarrow= new TRotMatrix("xarrow","",90.,0.0,0.0);
1168     TRotMatrix *yarrow= new TRotMatrix("yarrow","",0.0,90.,0.0);
1169
1170     det->SetLineColor(0); // black
1171     det->SetLineStyle(1); // solid line
1172     det->SetLineWidth(2); // pixel units
1173     det->SetFillColor(1); // black
1174     det->SetFillStyle(4010); // window is 90% transparent
1175     arrow->SetLineColor(det->GetLineColor());
1176     arrow->SetLineWidth(det->GetLineWidth());
1177     arrow->SetLineStyle(det->GetLineStyle());
1178     arrow->SetFillColor(1); // black
1179     arrow->SetFillStyle(4100); // window is 100% opaque
1180     arrow->DefineSection(0,0.0,0.0,kRmax);
1181     arrow->DefineSection(1,2.*kRmax,0.0,0.0);
1182     view->SetRange(-kDx0,-kDy0,-kDz0,kDx0,kDy0,kDz0);
1183     //
1184     TNode *node0 = new TNode("NODE0","NODE0",mother);
1185     node0->cd();
1186     TNode *node1 = new TNode("NODE1","NODE1",det);
1187     node1->cd();
1188     TNode *nodex = new TNode("NODEx","NODEx",arrow,
1189                              l[0][0],l[0][1],l[0][2],xarrow);
1190     TNode *nodey = new TNode("NODEy","NODEy",arrow,
1191                              l[2][0],l[2][1],l[2][2],yarrow);
1192     TNode *nodez = new TNode("NODEz","NODEz",arrow,l[4][0],l[4][1],l[4][2]);
1193     //
1194     axis->Draw();
1195     nodex->Draw();
1196     nodey->Draw();
1197     nodez->Draw();
1198     
1199     //
1200     node0->cd();
1201     node0->Draw();
1202     c->Update();
1203     c->SaveAs("AliITSgeomMatrix_L1.gif");
1204 }
1205 //----------------------------------------------------------------------
1206 ostream &operator<<(ostream &os,AliITSgeomMatrix &p){
1207     // Standard output streaming function.
1208     // Inputs:
1209     //    ostream &os          The output stream to print the class data on
1210     //    AliITSgeomMatrix &p  This class
1211     // Outputs:
1212     //    none.
1213     // Return:
1214     //    none.
1215
1216     p.Print(&os);
1217     return os;
1218 }
1219 //----------------------------------------------------------------------
1220 istream &operator>>(istream &is,AliITSgeomMatrix &r){
1221     // Standard input streaming function.
1222     // Inputs:
1223     //    ostream &os          The input stream to print the class data on
1224     //    AliITSgeomMatrix &p  This class
1225     // Outputs:
1226     //    none.
1227     // Return:
1228     //    none.
1229
1230     r.Read(&is);
1231     return is;
1232 }
1233 //----------------------------------------------------------------------