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