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