]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSBaseGeometry.cxx
Removed warning from part of code not properly implimneted yet.
[u/mrichter/AliRoot.git] / ITS / AliITSBaseGeometry.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 $Log$
18
19 $Id$
20 */
21
22 /*
23   A base geometry class defining all of the ITS volumes that make up an ITS
24 geometry.
25 Auhors: B. S. Nilsen
26 Version 0
27 Created February 2003.
28 */
29
30 #include <Riostream.h>
31 #include <stdio.h>
32 #include <stdlib.h>
33 #include <string.h>
34 #include <TMath.h>
35 #include <TGeometry.h>
36 #include <TNode.h>
37 #include <TTUBE.h>
38 #include <TTUBS.h>
39 #include <TPCON.h>
40 #include <TFile.h>    // only required for Tracking function?
41 #include <TCanvas.h>
42 #include <TObjArray.h>
43 #include <TLorentzVector.h>
44 #include <TObjString.h>
45 #include <TClonesArray.h>
46 #include <TBRIK.h>
47 #include <TSystem.h>
48 #include <AliRun.h>
49 #include <AliMagF.h>
50 #include <AliConst.h>
51 #include "AliITSBaseGeometry.h"
52
53 ClassImp(AliITSBaseGeometry)
54 //______________________________________________________________________
55 AliITSBaseGeometry::AliITSBaseGeometry(){
56     // Default construtor for the ITS Base Geometry class.
57     // Inputs:
58     //    none.
59     // Outputs:
60     //    none.
61     // Return:
62     //    none.
63
64     fScale = 1.0; // Default value.
65     fits = 0; // zero pointers.
66     fNCreates++; // incrament this creation counter.
67 }
68 //______________________________________________________________________
69 AliITSBaseGeometry::AliITSBaseGeometry(AliModule *its,Int_t iflag){
70     // Standard construtor for the ITS Base Geometry class.
71     // Inputs:
72     //    Int_t iflag  flag to indecate specific swiches in the geometry
73     // Outputs:
74     //    none.
75     // Return:
76     //    none.
77
78     fScale = 1.0; // Default value.
79     fits = its; // get a copy of the pointer to the ITS.
80     fNCreates++; // incrament this creation counter.
81 }
82 //______________________________________________________________________
83 AliITSBaseGeometry::~AliITSBaseGeometry(){
84     // Standeard destructor for the ITS Base Geometry class.
85     // Inputs:
86     //    Int_t iflag  flag to indecate specific swiches in the geometry
87     // Outputs:
88     //    none.
89     // Return:
90     //    none.
91
92     fits = 0; // This class does not own this class. It contaitns a pointer
93     // to it for conveniance.
94     fidmed = 0; // This class does not own this array of media indexs. It
95     fNCreates--;
96     if(fNCreates==0){ // Now delete the static members
97         Int_t i;
98         if(fVolName!=0){
99             for(i=0;i<fVolNameLast;i++) delete fVolName[i];
100             fVolNameSize = 0;
101             fVolNameLast = 0;
102             delete[] fVolName;
103         }// end if
104         delete[] fidrot;
105         fidrotsize = fidrotlast = 0;
106     }// end if
107 }
108 //______________________________________________________________________
109 Int_t AliITSBaseGeometry::AddVolName(const TString name){
110     // Checks if the volume name already exist, if not it adds it to
111     // the list of volume names and returns an index to that volume name.
112     // it will create and expand the array of volume names as needed.
113     // If the volume name already exists, it will give an error message and
114     // return an index <0.
115     // Inputs:
116     //    const TString name  Volume name to be added to the list.
117     // Outputs:
118     //    none.
119     // Return:
120     //    The index where this volume name is stored.
121     Int_t i;
122
123     if(fVolName==0){ // must create array.
124         fVolNameSize = 1000;
125         fVolName = new TString[fVolNameSize];
126         fVolNameLast = 0;
127     } // end if
128     for(i=0;i<fVolNameLast;i++) if(fVolName[i].CompareTo(name)==0){ // Error
129         Error("AddVolName","Volume name already exists for volume %d",i);
130         return -1;
131     } // end for i
132     if(fVolNameSize==fVolNameLast-1){ // Array is full must expand.
133         Int_t size = fVolNameSize*2;
134         TString *old = fVolName;
135         fVolName = new TString[fVolNameSize];
136         for(i=0;i<fVolNameLast;i++) fVolName[i] = old[i];
137         delete[] old;
138         fVolNameSize = size;
139     } // end if
140     if(strcmp(ITSIndexToITSG3name(fVolNameLast),"ITSV")==0){
141         // Special Reserved Geant 3 volumen name. Skip it
142         // fill it with explination for conveniance.
143         fVolName[fVolNameLast] = "ITS Master Mother Volume";
144         fVolNameLast++;
145     } // end if
146     fVolName[fVolNameLast] = name;
147     fVolNameLast++;
148     return fVolNameLast-1; // return the index
149 }
150 //______________________________________________________________________
151 char* AliITSBaseGeometry::ITSIndexToITSG3name(const Int_t i){
152     // Given the ITS volume index i, it returns the Geant3 ITS volume
153     // name. The valid characters must be in the range
154     // '0' through 'Z'. This will include all upper case letter and the
155     // numbers 0-9. In addition it does not will include the following simbols
156     // ":;<=>?@"
157     // Inputs:
158     //    const Int_t i  the ITS volume index
159     // Output:
160     //    none.
161     // Return:
162     //    char[4] with the ITS volume name starting from "I000" to "IZZZ"
163     const Int_t rangen=(Int_t)('9'-'0'+1); // range of numbers
164     const Int_t rangel=(Int_t)('Z'-'A'+1); // range of letters
165     const Int_t range = rangen+rangel; // the number of characters between 
166                                        // 0-9 and A-Z.
167     char a[4];
168     Int_t j = i;
169
170     a[0] = (char)('I');
171     a[1] = (char)('0'+j/(range*range));
172     if(a[1]>'9') a[1] += 'A'-'0'; // if it is a letter add in gap for simples.
173     j -= range*range*(a[1]-'0');
174     a[2] = (char)('0'+j/range);
175     if(a[2]>'9') a[2] += 'A'-'0'; // if it is a letter add in gap for simples.
176     j -= range*(a[2]-'0');
177     a[3] = (char)('0'+j);
178     if(a[3]>'9') a[3] += 'A'-'0'; // if it is a letter add in gap for simples.
179     return a;
180 }
181 //______________________________________________________________________
182 Int_t AliITSBaseGeometry::ITSG3VnameToIndex(const char name[3])const{
183     // Given the last three characters of the ITS Geant3 volume name,
184     // this returns the index. The valid characters must be in the range
185     // '0' through 'Z'. This will include all upper case letter and the
186     // numbers 0-9. In addition it will include the following simbles
187     // ":;<=>?@"
188     // Inputs:
189     //    const char name[3]  The last three characters of the ITS Geant3
190     //                        volume name
191     // Output:
192     //    none.
193     // Return:
194     //    Int_t the index.
195     const Int_t rangen=(Int_t)('9'-'0'+1); // range of numbers
196     const Int_t rangel=(Int_t)('Z'-'A'+1); // range of letters
197     const Int_t range = rangen+rangel; // the number of characters between 
198                                        // 0-9 and A-Z.
199     Int_t i,j;
200
201     i = 0;
202     for(j=3;j>-1;j--){
203         if(isdigit(name[j])){ // number
204             i += (Int_t)(name[j]-'0')*TMath::Power(range,(Double_t)j);
205         }else{ // Letter
206             i += (Int_t)(name[j]-'A'+rangen)*TMath::Power(range,(Double_t)j);
207         } // end if
208     } // end for j
209     return i;
210 }
211 //______________________________________________________________________
212 TString AliITSBaseGeometry::GetVolName(const Int_t i)const{
213     // Returns the volume name at a given index i. Index must be in
214     // range and the array of volume names must exist. If there is an
215     // error, a message is written and 0 is returned.
216     // Inputs:
217     //   const Int_t i Index
218     // Output:
219     //   none.
220     // Return:
221     //   A TString contianing the ITS volume name.
222
223     if(i<0||i>=fVolNameLast){
224         Error("GetVolName","Index=%d out of range but be witin 0<%d",i,
225               fVolName-1);
226         return 0;
227     } // end if Error
228     return fVolName[i];
229 }
230 //______________________________________________________________________
231 Int_t AliITSBaseGeometry::GetVolumeIndex(const TString &a){
232     // Return the index corresponding the the volume name a. If the
233     // Volumen name is not found, return -1, and a warning message given.
234     // Inputs:
235     //   const TString &a  Name of volume for which index is wanted.
236     // Output:
237     //   none.
238     // Return:
239     //   Int_t Index corresponding the volume a. If not found -1 is returned.
240     Int_t i;
241
242     for(i=0;i<fVolNameLast;i++) if(fVolName[i].CompareTo(a)==0) return i;
243     Info("GetVolumeIndex","Volume name %s not found",a.Data());
244     return -1;
245 }
246 //______________________________________________________________________
247 void AliITSBaseGeometry::Box(const char gnam[3],const TString &dis,
248                              Double_t dx,Double_t dy,Double_t dz,Int_t med){
249     // Interface to TMC->Gsvolu() for ITS bos geometries. Box with faces
250     // perpendicular to the axes. It has 3 paramters. See SetScale() for
251     // units. Default units are geant 3 [cm].
252     // Inputs:
253     //    const char gnam[3]  3 character geant volume name. The letter "I"
254     //                        is appended to the front to indecate that this
255     //                        is an ITS volume.
256     //    TString &dis        String containging part discription.
257     //    Double_t dx         half-length of box in x-axis
258     //    Double_t dy         half-length of box in y-axis
259     //    Double_t dz         half-length of box in z-axis
260     //    Int_t    med        media index number.
261     // Output:
262     //    none.
263     // Return.
264     //    none.
265     char name[4];
266     Float_t param[3];
267
268     if(fidmed==0) SetMedArray();
269     param[0] = fScale*dx;
270     param[1] = fScale*dy;
271     param[2] = fScale*dz;
272     name[3] = 'I';
273     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
274     gMC->Gsvolu(name,"BOX ",fidmed[med],param,3);
275 }
276 //______________________________________________________________________
277 void AliITSBaseGeometry::Trapezoid1(const char gnam[3],const TString &dis,
278                                     Double_t dxn,Double_t dxp,Double_t dy,
279                                     Double_t dz,Int_t med){
280     // Interface to TMC->Gsvolu() for ITS TRD1 geometries. Trapezoid with the 
281     // x dimension varing along z. It has 4 parameters. See SetScale() for
282     // units. Default units are geant 3 [cm].
283     // Inputs:
284     //    const char gnam[3]  3 character geant volume name. The letter "I"
285     //                        is appended to the front to indecate that this
286     //                        is an ITS volume.
287     //    TString &dis        String containging part discription.
288     //    Double_t dxn        half-length along x at the z surface positioned 
289     //                        at -DZ
290     //    Double_t dxp        half-length along x at the z surface positioned 
291     //                        at +DZ
292     //    Double_t dy         half-length along the y-axis
293     //    Double_t dz         half-length along the z-axis
294     //    Int_t    med        media index number.
295     // Output:
296     //    none.
297     // Return.
298     //    none.
299     char name[4];
300     Float_t param[4];
301
302     if(fidmed==0) SetMedArray();
303     param[0] = fScale*dxn;
304     param[1] = fScale*dxp;
305     param[2] = fScale*dy;
306     param[3] = fScale*dz;
307     name[3] = 'I';
308     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
309     gMC->Gsvolu(name,"TRD1",fidmed[med],param,4);
310 }
311 //______________________________________________________________________
312 void AliITSBaseGeometry::Trapezoid2(const char gnam[3],const TString &dis,
313                                     Double_t dxn,Double_t dxp,Double_t dyn,
314                                     Double_t dyp,Double_t dz,Int_t med){
315     // Interface to TMC->Gsvolu() for ITS TRD2 geometries. Trapezoid with the 
316     // x and y dimension varing along z. It has 5 parameters. See SetScale() 
317     // for units. Default units are geant 3 [cm].
318     // Inputs:
319     //    const char gnam[3]  3 character geant volume name. The letter "I"
320     //                        is appended to the front to indecate that this
321     //                        is an ITS volume.
322     //    TString &dis        String containging part discription.
323     //    Double_t dxn        half-length along x at the z surface positioned 
324     //                        at -DZ
325     //    Double_t dxp        half-length along x at the z surface positioned 
326     //                        at +DZ
327     //    Double_t dyn        half-length along x at the z surface positioned 
328     //                        at -DZ
329     //    Double_t dyp        half-length along x at the z surface positioned 
330     //                        at +DZ
331     //    Double_t dz         half-length along the z-axis
332     //    Int_t    med        media index number.
333     // Output:
334     //    none.
335     // Return.
336     //    none.
337     char name[4];
338     Float_t param[5];
339
340     if(fidmed==0) SetMedArray();
341     param[0] = fScale*dxn;
342     param[1] = fScale*dxp;
343     param[2] = fScale*dyn;
344     param[3] = fScale*dyp;
345     param[4] = fScale*dz;
346     name[3] = 'I';
347     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
348     gMC->Gsvolu(name,"TRD2",fidmed[med],param,5);
349 }
350 //______________________________________________________________________
351 void AliITSBaseGeometry::Trapezoid(const char gnam[3],const TString &dis,
352                                    Double_t dz,Double_t thet,Double_t phi,
353                                    Double_t h1,Double_t bl1,Double_t tl1,
354                                    Double_t alp1,Double_t h2,Double_t bl2,
355                                    Double_t tl2,Double_t alp2,Int_t med){
356     // Interface to TMC->Gsvolu() for ITS TRAP geometries. General Trapezoid, 
357     // The faces perpendicular to z are trapezia and their centers are not 
358     // necessarily on a line parallel to the z axis. This shape has 11 
359     // parameters, but only cosidering that the faces should be planar, only 9 
360     // are really independent. A check is performed on the user parameters and 
361     // a message is printed in case of non-planar faces. Ignoring this warning 
362     // may cause unpredictable effects at tracking time. See SetScale() 
363     // for units. Default units are geant 3 [cm].
364     // Inputs:
365     //    const char gnam[3]  3 character geant volume name. The letter "I"
366     //                        is appended to the front to indecate that this
367     //                        is an ITS volume.
368     //    TString &dis        String containging part discription.
369     //    Double_t dz         Half-length along the z-asix
370     //    Double_t thet       Polar angle of the line joing the center of the 
371     //                        face at -dz to the center of the one at dz 
372     //                        [degree].
373     //    Double_t phi        aximuthal angle of the line joing the center of 
374     //                        the face at -dz to the center of the one at +dz 
375     //                        [degree].
376     //    Double_t h1         half-length along y of the face at -dz.
377     //    Double_t bl1        half-length along x of the side at -h1 in y of 
378     //                        the face at -dz in z.
379     //    Double_t tl1        half-length along x of teh side at +h1 in y of 
380     //                        the face at -dz in z.
381     //    Double_t alp1       angle with respect to the y axis from the center 
382     //                        of the side at -h1 in y to the cetner of the 
383     //                        side at +h1 in y of the face at -dz in z 
384     //                        [degree].
385     //    Double_t h2         half-length along y of the face at +dz
386     //    Double_t bl2        half-length along x of the side at -h2 in y of
387     //                        the face at +dz in z.
388     //    Double_t tl2        half-length along x of the side at _h2 in y of 
389     //                        the face at +dz in z.
390     //    Double_t alp2       angle with respect to the y axis from the center 
391     //                        of the side at -h2 in y to the center of the 
392     //                        side at +h2 in y of the face at +dz in z 
393     //                        [degree].
394     //    Int_t    med        media index number.
395     // Output:
396     //    none.
397     // Return.
398     //    none.
399     char name[4];
400     Float_t param[11];
401
402     if(fidmed==0) SetMedArray();
403     param[0] = fScale*dz;
404     param[1] = thet;
405     param[2] = phi;
406     param[3] = fScale*h1;
407     param[4] = fScale*bl1;
408     param[5] = fScale*tl1;
409     param[6] = alp1;
410     param[7] = fScale*h2;
411     param[8] = fScale*bl2;
412     param[9] = fScale*tl2;
413     param[10] = alp2;
414     name[3] = 'I';
415     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
416     gMC->Gsvolu(name,"TRAP",fidmed[med],param,11);
417 }
418 //______________________________________________________________________
419 void AliITSBaseGeometry::Tube(const char gnam[3],const TString &dis,
420                               Double_t rmin,Double_t rmax,Double_t dz,
421                               Int_t med){
422     // Interface to TMC->Gsvolu() for ITS TUBE geometries. Simple Tube. It has
423     // 3 parameters. See SetScale() 
424     // for units. Default units are geant 3 [cm].
425     // Inputs:
426     //    const char gnam[3]  3 character geant volume name. The letter "I"
427     //                        is appended to the front to indecate that this
428     //                        is an ITS volume.
429     //    TString &dis        String containging part discription.
430     //    Double_t rmin       Inside Radius.
431     //    Double_t rmax       Outside Radius.
432     //    Double_t dz         half-length along the z-axis
433     //    Int_t    med        media index number.
434     // Output:
435     //    none.
436     // Return.
437     //    none.
438     char name[4];
439     Float_t param[3];
440
441     if(fidmed==0) SetMedArray();
442     param[0] = fScale*rmin;
443     param[1] = fScale*rmax;
444     param[2] = fScale*dz;
445     name[3] = 'I';
446     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
447     gMC->Gsvolu(name,"TUBE",fidmed[med],param,3);
448 }
449 //______________________________________________________________________
450 void AliITSBaseGeometry::TubeSegment(const char gnam[3],const TString &dis,
451                                      Double_t rmin,Double_t rmax,Double_t dz,
452                                      Double_t phi1,Double_t phi2,Int_t med){
453     // Interface to TMC->Gsvolu() for ITS TUBE geometries. Phi segment of a 
454     // tube. It has 5  parameters. Phi1 should be smaller than phi2. If this is
455     // not the case, the system adds 360 degrees to phi2. See SetScale() 
456     // for units. Default units are geant 3 [cm].
457     // Inputs:
458     //    const char gnam[3]  3 character geant volume name. The letter "I"
459     //                        is appended to the front to indecate that this
460     //                        is an ITS volume.
461     //    TString &dis        String containging part discription.
462     //    Double_t rmin       Inside Radius.
463     //    Double_t rmax       Outside Radius.
464     //    Double_t dz         half-length along the z-axis
465     //    Double_t phi1       Starting angle of the segment [degree].
466     //    Double_t phi2       Ending angle of the segment [degree].
467     //    Int_t    med        media index number.
468     // Output:
469     //    none.
470     // Return.
471     //    none.
472     char name[4];
473     Float_t param[5];
474
475     if(fidmed==0) SetMedArray();
476     param[0] = fScale*rmin;
477     param[1] = fScale*rmax;
478     param[2] = fScale*dz;
479     param[3] = phi1;
480     param[4] = phi2;
481     name[3] = 'I';
482     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
483     gMC->Gsvolu(name,"TUBS",fidmed[med],param,5);
484 }
485 //______________________________________________________________________
486 void AliITSBaseGeometry::Cone(const char gnam[3],const TString &dis,
487                               Double_t dz,Double_t rmin1,Double_t rmax1,
488                               Double_t rmin2,Double_t rmax2,Int_t med){
489     // Interface to TMC->Gsvolu() for ITS Cone geometries. Conical tube. It 
490     // has 5 parameters. See SetScale() 
491     // for units. Default units are geant 3 [cm].
492     // Inputs:
493     //    const char gnam[3]  3 character geant volume name. The letter "I"
494     //                        is appended to the front to indecate that this
495     //                        is an ITS volume.
496     //    TString &dis        String containging part discription.
497     //    Double_t dz         half-length along the z-axis
498     //    Double_t rmin1      Inside Radius at -dz.
499     //    Double_t rmax1      Outside Radius at -dz.
500     //    Double_t rmin2      inside radius at +dz.
501     //    Double_t rmax2      outside radius at +dz.
502     //    Int_t    med        media index number.
503     // Output:
504     //    none.
505     // Return.
506     //    none.
507     char name[4];
508     Float_t param[5];
509
510     if(fidmed==0) SetMedArray();
511     param[0] = fScale*dz;
512     param[1] = fScale*rmin1;
513     param[2] = fScale*rmax1;
514     param[3] = fScale*rmin2;
515     param[4] = fScale*rmax2;
516     name[3] = 'I';
517     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
518     gMC->Gsvolu(name,"CONS",fidmed[med],param,5);
519 }
520 //______________________________________________________________________
521 void AliITSBaseGeometry::ConeSegment(const char gnam[3],const TString &dis,
522                                      Double_t dz,Double_t rmin1,Double_t rmax1,
523                                      Double_t rmin2,Double_t rmax2,
524                                      Double_t phi1,Double_t phi2,Int_t med){
525     // Interface to TMC->Gsvolu() for ITS ConS geometries. One segment of a 
526     // conical tube. It has 7 parameters. Phi1 should be smaller than phi2. If 
527     // this is not the case, the system adds 360 degrees to phi2. See 
528     // SetScale() for units. Default units are geant 3 [cm].
529     // Inputs:
530     //    const char gnam[3]  3 character geant volume name. The letter "I"
531     //                        is appended to the front to indecate that this
532     //                        is an ITS volume.
533     //    TString &dis        String containging part discription.
534     //    Double_t dz         half-length along the z-axis
535     //    Double_t rmin1      Inside Radius at -dz.
536     //    Double_t rmax1      Outside Radius at -dz.
537     //    Double_t rmin2      inside radius at +dz.
538     //    Double_t rmax2      outside radius at +dz.
539     //    Double_t phi1       Starting angle of the segment [degree].
540     //    Double_t phi2       Ending angle of the segment [degree].
541     //    Int_t    med        media index number.
542     // Output:
543     //    none.
544     // Return.
545     //    none.
546     char name[4];
547     Float_t param[7];
548
549     if(fidmed==0) SetMedArray();
550     param[0] = fScale*dz;
551     param[1] = fScale*rmin1;
552     param[2] = fScale*rmax1;
553     param[3] = fScale*rmin2;
554     param[4] = fScale*rmax2;
555     param[5] = phi1;
556     param[6] = phi2;
557     name[3] = 'I';
558     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
559     gMC->Gsvolu(name,"CONS",fidmed[med],param,7);
560 }
561 //______________________________________________________________________
562 void AliITSBaseGeometry::Sphere(const char gnam[3],const TString &dis,
563                                 Double_t rmin,Double_t rmax,Double_t the1,
564                                 Double_t the2,Double_t phi1,Double_t phi2,
565                                 Int_t med){
566     // Interface to TMC->Gsvolu() for ITS SPHE geometries. Segment of a 
567     // sphereical shell. It has 6 parameters. See SetScale() 
568     // for units. Default units are geant 3 [cm].
569     // Inputs:
570     //    const char gnam[3]  3 character geant volume name. The letter "I"
571     //                        is appended to the front to indecate that this
572     //                        is an ITS volume.
573     //    TString &dis        String containging part discription.
574     //    Double_t rmin       Inside Radius.
575     //    Double_t rmax       Outside Radius.
576     //    Double_t the1       staring polar angle of the shell [degree].
577     //    Double_t the2       ending polar angle of the shell [degree].
578     //    Double_t phui       staring asimuthal angle of the shell [degree].
579     //    Double_t phi2       ending asimuthal angle of the shell [degree].
580     //    Int_t    med        media index number.
581     // Output:
582     //    none.
583     // Return.
584     //    none.
585     char name[4];
586     Float_t param[6];
587
588     if(fidmed==0) SetMedArray();
589     param[0] = fScale*rmin;
590     param[1] = fScale*rmax;
591     param[2] = the1;
592     param[3] = the2;
593     param[4] = phi1;
594     param[5] = phi2;
595     name[3] = 'I';
596     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
597     gMC->Gsvolu(name,"SPHE",fidmed[med],param,6);
598 }
599 //______________________________________________________________________
600 void AliITSBaseGeometry::Parallelepiped(const char gnam[3],const TString &dis,
601                                         Double_t dx,Double_t dy,Double_t dz,
602                                         Double_t alpha,Double_t thet,
603                                         Double_t phi,Int_t med){
604     // Interface to TMC->Gsvolu() for ITS PARA geometries. Parallelepiped. It 
605     // has 6 parameters. See SetScale() for units. Default units are geant 3 
606     // [cm].
607     // Inputs:
608     //    const char gnam[3]  3 character geant volume name. The letter "I"
609     //                        is appended to the front to indecate that this
610     //                        is an ITS volume.
611     //    TString &dis        String containging part discription.
612     //    Double_t dx         half-length allong x-axis
613     //    Double_t dy         half-length allong y-axis
614     //    Double_t dz         half-length allong z-axis
615     //    Double_t alpha      angle formed by the y axis and by the plane 
616     //                        joining the center of teh faces parallel to the 
617     //                        z-x plane at -dY and +dy [degree].
618     //    Double_t thet       polar angle of the line joining the centers of 
619     //                        the faces at -dz and +dz in z [degree].
620     //    Double_t phi        azimuthal angle of teh line joing the centers of 
621     //                        the faaces at -dz and +dz in z [degree].
622     //    Int_t    med        media index number.
623     // Output:
624     //    none.
625     // Return.
626     //    none.
627     char name[4];
628     Float_t param[6];
629
630     if(fidmed==0) SetMedArray();
631     param[0] = fScale*dx;
632     param[1] = fScale*dy;
633     param[2] = fScale*dz;
634     param[3] = alpha;
635     param[4] = thet;
636     param[5] = phi;
637     name[3] = 'I';
638     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
639     gMC->Gsvolu(name,"PARA",fidmed[med],param,6);
640 }
641 //______________________________________________________________________
642 void AliITSBaseGeometry::Polygon(const char gnam[3],const TString &dis,
643                                  Double_t phi1,Double_t dphi,Int_t npdv,
644                                  Int_t nz,Double_t *z,Double_t *rmin,
645                                  Double_t *rmax,Int_t med){
646     // Interface to TMC->Gsvolu() for ITS PGON geometry. Polygon It has 10 
647     // parameters or more. See SetScale() for units. Default units are geant 3 
648     // [cm].
649     // Inputs:
650     //    const char gnam[3]  3 character geant volume name. The letter "I"
651     //                        is appended to the front to indecate that this
652     //                        is an ITS volume.
653     //    TString &dis        String containging part discription.
654     //    Double_t phi1       the azimuthal angle at which the volume begins 
655     //                        (angles are counted clouterclockwise) [degrees].
656     //    Double_t dphi       opening angle of the volume, which extends from 
657     //                        phi1 to phi1+dphi [degree].
658     //    Int_t npdv          the number of sides of teh cross section between 
659     //                        the given phi limits.
660     //    Int_t nz            number of planes perpendicular to the z axis 
661     //                        where the dimension of the section is given - 
662     //                        this number should be at least 2 and NP triples 
663     //                        of number must follow.
664     //    Double_t *z         array [nz] of z coordiates of the sections..
665     //    Double_t *rmin      array [nz] of radius of teh circle tangent to 
666     //                        the sides of the inner polygon in teh 
667     //                        cross-section.
668     //    Double_t *rmax      array [nz] of radius of the circle tangent to 
669     //                        the sides of the outer polygon in the 
670     //                       cross-section.
671     //    Int_t    med        media index number.
672     // Output:
673     //    none.
674     // Return.
675     //    none.
676     char name[4];
677     Float_t *param;
678     Int_t n,i;
679
680     if(fidmed==0) SetMedArray();
681     n = 4+3*nz;
682     param = new Float_t[n];
683     param[0] = phi1;
684     param[1] = dphi;
685     param[2] = (Float_t)npdv;
686     param[3] = (Float_t)nz;
687     for(i=0;i<nz;i++){
688         param[4+3*i] = z[i];
689         param[5+3*i] = rmin[i];
690         param[6+3*i] = rmax[i];
691     } // end for i
692     name[3] = 'I';
693     for(i=0;i<3;i++) name[i+1] = gnam[i];
694     gMC->Gsvolu(name,"PGON",fidmed[med],param,n);
695
696     delete[] param;
697 }
698 //______________________________________________________________________
699 void AliITSBaseGeometry::PolyCone(const char gnam[3],const TString &dis,
700                                   Double_t phi1,Double_t dphi,Int_t nz,
701                                   Double_t *z,Double_t *rmin,Double_t *rmax,
702                                   Int_t med){
703     // Interface to TMC->Gsvolu() for ITS PCON geometry. Poly-cone It has 9 
704     // parameters or more. See SetScale() for units. Default units are geant 3 
705     // [cm].
706     // Inputs:
707     //    const char gnam[3]  3 character geant volume name. The letter "I"
708     //                        is appended to the front to indecate that this
709     //                        is an ITS volume.
710     //    TString &dis        String containging part discription.
711     //    Double_t phi1       the azimuthal angle at which the volume begins 
712     //                        (angles are counted clouterclockwise) [degrees].
713     //    Double_t dphi       opening angle of the volume, which extends from 
714     //                        phi1 to phi1+dphi [degree].
715     //    Int_t nz            number of planes perpendicular to the z axis 
716     //                        where the dimension of the section is given - 
717     //                        this number should be at least 2 and NP triples 
718     //                        of number must follow.
719     //    Double_t *z         Array [nz] of z coordinate of the section.
720     //    Double_t *rmin      Array [nz] of radius of teh inner circle in the 
721     //                        cross-section.
722     //    Double_t *rmax      Array [nz] of radius of the outer circle in the 
723     //                        cross-section.
724     //    Int_t    med        media index number.
725     // Output:
726     //    none.
727     // Return.
728     //    none.
729     char name[4];
730     Float_t *param;
731     Int_t n,i;
732
733     if(fidmed==0) SetMedArray();
734     n = 3+3*nz;
735     param = new Float_t[n];
736     param[0] = phi1;
737     param[1] = dphi;
738     param[2] = (Float_t) nz;
739     for(i=0;i<nz;i++){
740         param[3+3*i] = z[i];
741         param[4+3*i] = rmin[i];
742         param[5+3*i] = rmax[i];
743     } // end for i
744     name[3] = 'I';
745     for(i=0;i<3;i++) name[i+1] = gnam[i];
746     gMC->Gsvolu(name,"PCON",fidmed[med],param,n);
747
748     delete[] param;
749 }
750 //______________________________________________________________________
751 void AliITSBaseGeometry::TubeElliptical(const char gnam[3],const TString &dis,
752                                Double_t p1,Double_t p2,Double_t dz,Int_t med){
753     // Interface to TMC->Gsvolu() for ITS ELTU geometries. Elliptical 
754     // cross-section Tube. It has 3 parameters. See SetScale() 
755     // for units. Default units are geant 3 [cm]. The equation of the surface 
756     // is x^2 * p1^-2 + y^2 * p2^-2 = 1.
757     // Inputs:
758     //    const char gnam[3]  3 character geant volume name. The letter "I"
759     //                        is appended to the front to indecate that this
760     //                        is an ITS volume.
761     //    TString &dis        String containging part discription.
762     //    Double_t p1         semi-axis of the elipse along x.
763     //    Double_t p2         semi-axis of the elipse along y.
764     //    Double_t dz         half-length along the z-axis
765     //    Int_t    med        media index number.
766     // Output:
767     //    none.
768     // Return.
769     //    none.
770     char name[4];
771     Float_t param[3];
772
773     if(fidmed==0) SetMedArray();
774     param[0] = fScale*p1;
775     param[1] = fScale*p2;
776     param[2] = fScale*dz;
777     name[3] = 'I';
778     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
779     gMC->Gsvolu(name,"ELTU",fidmed[med],param,3);
780 }
781 //______________________________________________________________________
782 void AliITSBaseGeometry::HyperbolicTube(const char gnam[3],const TString &dis,
783                                Double_t rmin,Double_t rmax,Double_t dz,
784                                Double_t thet,Int_t med){
785     // Interface to TMC->Gsvolu() for ITS HYPE geometries. Hyperbolic tube. 
786     // Fore example the inner and outer surfaces are hyperboloids, as would be 
787     // foumed by a system of cylinderical wires which were then rotated 
788     // tangentially about their centers. It has 4 parameters. See SetScale() 
789     // for units. Default units are geant 3 [cm]. The hyperbolic surfaces are 
790     // given by r^2 = (ztan(thet)^2 + r(z=0)^2.
791     // Inputs:
792     //    const char gnam[3]  3 character geant volume name. The letter "I"
793     //                        is appended to the front to indecate that this
794     //                        is an ITS volume.
795     //    TString &dis        String containging part discription.
796     //    Double_t rmin       Inner radius at z=0 where tube is narrowest.
797     //    Double_t rmax       Outer radius at z=0 where tube is narrowest.
798     //    Double_t dz         half-length along the z-axis
799     //    Double_t thet       stero angel of rotation of the two faces 
800     //                       [degrees].
801     //    Int_t    med        media index number.
802     // Output:
803     //    none.
804     // Return.
805     //    none.
806     char name[4];
807     Float_t param[4];
808
809     if(fidmed==0) SetMedArray();
810     param[0] = fScale*rmin;
811     param[1] = fScale*rmax;
812     param[2] = fScale*dz;
813     param[3] = thet;
814     name[3] = 'I';
815     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
816     gMC->Gsvolu(name,"HYPE",fidmed[med],param,4);
817 }
818 //______________________________________________________________________
819 void AliITSBaseGeometry::TwistedTrapezoid(const char gnam[3],
820                                           const TString &dis,
821                                  Double_t dz,Double_t thet,Double_t phi,
822                                  Double_t twist,Double_t h1,Double_t bl1,
823                                  Double_t tl1,Double_t apl1,Double_t h2,
824                                  Double_t bl2,Double_t tl2,Double_t apl2,
825                                  Int_t med){
826     // Interface to TMC->Gsvolu() for ITS GTRA geometries. General twisted 
827     // trapazoid. The faces perpendicular to z are trapazia and their centers 
828     // are not necessarily on a line parallel to the z axis as the TRAP. 
829     // Additionally, the faces may be twisted so that none of their edges are 
830     // parallel. It is a TRAP shape, exept that it is twisted in the x-y plane 
831     // as a function of z. The parallel sides perpendicular to the x axis are 
832     // rotated with respect to the x axis by an angle TWIST, which is one of 
833     // the parameters. The shape is defined by the eight corners and is assumed
834     // to be constructed of straight lines joingin points on the boundry of the
835     // trapezoidal face at Z=-dz to the coresponding points on the face at 
836     // z=+dz. Divisions are not allowed. It has 12 parameters. See SetScale() 
837     // for units. Default units are geant 3 [cm]. Note: This shape suffers from
838     // the same limitations than the TRAP. The tracking routines assume that 
839     // the faces are planar, but htis constraint is not easily expressed in 
840     // terms of the 12 parameters. Additionally, no check on th efaces is 
841     // performed in this case. Users should avoid to use this shape as much as 
842     // possible, and if they have to do so, they should make sure that the 
843     // faces are really planes. If this is not the case, the result of the 
844     // trasport is unpredictable. To accelerat ethe computations necessary for 
845     // trasport, 18 additioanl parameters are calculated for this shape are
846     // 1 DXODZ dx/dz of the line joing the centers of the faces at z=+_dz.
847     // 2 DYODZ dy/dz of the line joing the centers of the faces at z=+_dz.
848     // 3 XO1    x at z=0 for line joing the + on parallel side, perpendicular 
849     //          corners at z=+_dz.
850     // 4 YO1    y at z=0 for line joing the + on parallel side, + on 
851     //          perpendicular corners at z=+-dz.
852     // 5 DXDZ1  dx/dz for line joing the + on parallel side, + on 
853     //          perpendicular corners at z=+-dz.
854     // 6 DYDZ1  dy/dz for line joing the + on parallel side, + on 
855     //          perpendicular corners at z=+-dz.
856     // 7 X02    x at z=0 for line joing the - on parallel side, + on 
857     //          perpendicular corners at z=+-dz.
858     // 8 YO2    y at z=0 for line joing the - on parallel side, + on 
859     //          perpendicular corners at z=+-dz.
860     // 9 DXDZ2  dx/dz for line joing the - on parallel side, + on 
861     //          perpendicular corners at z=+-dz.
862     // 10 DYDZ2dy/dz for line joing the - on parallel side, + on 
863     //          perpendicular corners at z=+-dz.
864     // 11 XO3   x at z=0 for line joing the - on parallel side, - on 
865     //          perpendicular corners at z=+-dz.
866     // 12 YO3   y at z=0 for line joing the - on parallel side, - on 
867     //          perpendicular corners at z=+-dz.
868     // 13 DXDZ3 dx/dzfor line joing the - on parallel side, - on 
869     //          perpendicular corners at z=+-dz.
870     // 14 DYDZ3 dydz for line joing the - on parallel side, - on 
871     //          perpendicular corners at z=+-dz.
872     // 15 XO4   x at z=0 for line joing the + on parallel side, - on 
873     //          perpendicular corners at z=+-dz.
874     // 16 YO4   y at z=0 for line joing the + on parallel side, - on 
875     //          perpendicular corners at z=+-dz.
876     // 17 DXDZ4 dx/dz for line joing the + on parallel side, - on 
877     //          perpendicular corners at z=+-dz.
878     // 18 DYDZ4 dydz for line joing the + on parallel side, - on 
879     //          perpendicular corners at z=+-dz.
880     // Inputs:
881     //    const char gnam[3]  3 character geant volume name. The letter "I"
882     //                        is appended to the front to indecate that this
883     //                        is an ITS volume.
884     //    TString &dis        String containging part discription.
885     //    Double_t dz         half-length along the z axis.
886     //    Double_t thet       polar angle of the line joing the center of the 
887     //                        face at -dz to the center of the one at +dz 
888     //                        [degrees].
889     //    Double_t phi        Azymuthal angle of teh line joing the centre of 
890     //                        the face at -dz to the center of the one at +dz 
891     //                        [degrees].
892     //    Double_t twist      Twist angle of the faces parallel to the x-y 
893     //                        plane at z=+-dz around an axis parallel to z 
894     //                        passing through their centre [degrees].
895     //    Double_t h1         Half-length along y of the face at -dz.
896     //    Double_t bl1        half-length along x of the side -h1 in y of the 
897     //                        face at -dz in z.
898     //    Double_t tl1        half-length along x of the side at +h1 in y of 
899     //                        the face at -dz in z.
900     //    Double_t apl1       Angle with respect to the y ais from the center 
901     //                        of the side at -h1 in y to the centere of the 
902     //                        side at +h1 in y of the face at -dz in z 
903     //                        [degrees].
904     //    Double_t h2         half-length along the face at +dz.
905     //    Double_t bl2        half-length along x of the side at -h2 in y of 
906     //                        the face at -dz in z.
907     //    Double_t tl2        half-length along x of the side at +h2 in y of 
908     //                        the face at +dz in z.
909     //    Double_t apl2       angle with respect to the y axis from the center 
910     //                        of the side at -h2 in y to the center of the side
911     //                        at +h2 in y of the face at +dz in z [degrees].
912     //    Int_t    med        media index number.
913     // Output:
914     //    none.
915     // Return.
916     //    none.
917     char name[4];
918     Float_t param[12];
919
920     if(fidmed==0) SetMedArray();
921     param[0] = fScale*dz;
922     param[1] = thet;
923     param[2] = phi;
924     param[3] = twist;
925     param[4] = fScale*h1;
926     param[5] = fScale*bl1;
927     param[6] = fScale*tl1;
928     param[7] = apl1;
929     param[8] = fScale*h2;
930     param[9] = fScale*bl2;
931     param[10] = fScale*tl2;
932     param[11] = apl2;
933     name[3] = 'I';
934     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
935     gMC->Gsvolu(name,"GTRA",fidmed[med],param,12);
936 }
937 //______________________________________________________________________
938 void AliITSBaseGeometry::CutTube(const char gnam[3],const TString &dis,
939                                  Double_t rmin,Double_t rmax,Double_t dz,
940                                  Double_t phi1,Double_t phi2,Double_t lx,
941                                  Double_t ly,Double_t lz,Double_t hx,
942                                  Double_t hy,Double_t hz,Int_t med){
943     // Interface to TMC->Gsvolu() for ITS CTUB geometries. Cut tube. A tube cut
944     // at the extremities with planes not necessarily perpendicular tot he z 
945     // axis. It has 11 parameters. See SetScale() for units. Default units are 
946     // geant 3 [cm]. phi1 should be smaller than phi2. If this is not the case,
947     // the system adds 360 degrees to phi2.
948     // Inputs:
949     //    const char gnam[3]  3 character geant volume name. The letter "I"
950     //                        is appended to the front to indecate that this
951     //                        is an ITS volume.
952     //    TString &dis        String containging part discription.
953     //    Double_t rmin       Inner radius at z=0 where tube is narrowest.
954     //    Double_t rmax       Outer radius at z=0 where tube is narrowest.
955     //    Double_t dz         half-length along the z-axis
956     //    Double_t dz         half-length along the z-axis
957     //    Double_t phi1       Starting angle of the segment [degree].
958     //    Double_t phi2       Ending angle of the segment [degree].
959     //    Double_t lx         x component of a unit vector perpendicular to 
960     //                        the face at -dz.
961     //    Double_t ly         y component of a unit vector perpendicular to 
962     //                        the face at -dz.
963     //    Double_t lz         z component of a unit vector perpendicular to 
964     //                        the face at -dz.
965     //    Double_t hx         x component of a unit vector perpendicular to 
966     //                        the face at +dz.
967     //    Double_t hy         y component of a unit vector perpendicular to 
968     //                        the face at +dz.
969     //    Double_t hz         z component of a unit vector perpendicular to 
970     //                        the face at +dz.
971     //    Int_t    med        media index number.
972     // Output:
973     //    none.
974     // Return.
975     //    none.
976     char name[4];
977     Float_t param[11];
978
979     if(fidmed==0) SetMedArray();
980     param[0] = fScale*rmin;
981     param[1] = fScale*rmax;
982     param[2] = fScale*dz;
983     param[3] = phi1;
984     param[4] = phi2;
985     param[5] = lx;
986     param[6] = ly;
987     param[7] = lz;
988     param[8] = hx;
989     param[9] = hy;
990     param[10] = hz;
991     name[3] = 'I';
992     for(Int_t i=0;i<3;i++) name[i+1] = gnam[i];
993     gMC->Gsvolu(name,"CTUB",fidmed[med],param,11);
994 }
995 //______________________________________________________________________
996 void AliITSBaseGeometry::Pos(const char vol[3],Int_t cn,const char moth[3],
997                              Double_t x,Double_t y,Double_t z,Int_t irot){
998     // Place a copy of a volume previously defined by a call to GSVOLU inside 
999     // its mother volulme moth.
1000     // Inputs:
1001     //   const char vol[3]  3 character geant volume name. The letter "I"
1002     //                      is appended to the front to indecate that this
1003     //                      is an ITS volume.
1004     //   const char moth[3] 3 character geant volume name of the mother volume 
1005     //                      in which vol will be placed. The letter "I" is 
1006     //                      appended to the front to indecate that this is an 
1007     //                      ITS volume.
1008     //   Double_t x         The x positon of the volume in the mother's 
1009     //                      reference system
1010     //   Double_t y         The y positon of the volume in the mother's 
1011     //                      reference system
1012     //   Double_t z         The z positon of the volume in the mother's 
1013     //                      reference system
1014     //   Int_t irot         the index for the rotation matrix to be used.
1015     //                      irot=-1 => unit rotation.
1016     // Outputs:
1017     //    none.
1018     // Return:
1019     //    none.
1020     char name[4],mother[4];
1021     Float_t param[3];
1022     Int_t r=0,i;
1023
1024     param[0] = x;
1025     param[1] = y;
1026     param[2] = z;
1027     name[3] = 'I';
1028     for(i=0;i<3;i++) name[i+1] = vol[i];
1029     mother[3] = 'I';
1030     for(i=0;i<3;i++) mother[i+1] = moth[i];
1031     if(irot>=0) r=fidrot[irot];
1032     gMC->Gspos(name,1,mother,param[0],param[1],param[2],r,"ONLY");
1033 }
1034 //______________________________________________________________________
1035 void AliITSBaseGeometry::Matrix(Int_t irot,Double_t thet1,Double_t phi1,
1036                        Double_t thet2,Double_t phi2,
1037                        Double_t thet3,Double_t phi3){
1038     // Defines a Geant rotation matrix. checks to see if it is the unit
1039     // matrix. If so, then no additonal matrix is defined. Stores rotation 
1040     // matrix irot in the data structure JROTM. If the matrix is not 
1041     // orthonormal, it will be corrected by setting y' perpendicular to x' 
1042     // and z' = x' X y'. A warning message is printed in this case.
1043     // Inputs:
1044     //   Int_t irot     Intex specifing which rotation matrix.
1045     //   Double_t thet1 Polar angle for axisw x [degrees].
1046     //   Double_t phi1  azimuthal angle for axis x [degrees].
1047     //   Double_t thet12Polar angle for axisw y [degrees].
1048     //   Double_t phi2  azimuthal angle for axis y [degrees].
1049     //   Double_t thet3 Polar angle for axisw z [degrees].
1050     //   Double_t phi3  azimuthal angle for axis z [degrees].
1051     // Outputs:
1052     //    none.
1053     // Return:
1054     //    none.
1055     Float_t t1,p1,t2,p2,t3,p3;
1056
1057     if(thet1==90.0&&phi1==0.0&&thet2==90.0&&phi2==90.0&&thet3==0.0&&phi3==0.0){
1058         fidrot[irot] = 0; // Unit matrix
1059     }else{
1060         t1 = thet1;
1061         p1 = phi1;
1062         t2 = thet2;
1063         p2 = phi2;
1064         t3 = thet3;
1065         p3 = phi3;
1066         fits->AliMatrix(fidrot[irot],t1,p1,t2,p2,t3,p3);
1067     } // end if
1068 }
1069 //______________________________________________________________________
1070 void AliITSBaseGeometry::Matrix(Int_t irot,Int_t axis,Double_t thet){
1071     // Defines a Geant rotation matrix. checks to see if it is the unit
1072     // matrix. If so, then no additonal matrix is defined. Stores rotation 
1073     // matrix irot in the data structure JROTM. If the matrix is not 
1074     // orthonormal, it will be corrected by setting y' perpendicular to x' 
1075     // and z' = x' X y'. A warning message is printed in this case.
1076     // Inputs:
1077     //   Int_t irot         Intex specifing which rotation matrix.
1078     //   Int_t axis         Axis about which rotation is to be done.
1079     //   Double_t thet      Angle to rotate by [degrees].
1080     // Outputs:
1081     //    none.
1082     // Return:
1083     //    none.
1084
1085     if(thet==0.0){
1086         fidrot[irot] = 0; // Unit matrix
1087     }else{
1088         switch (irot) {
1089         case 0: //Rotate about x-axis, x-axis does not change.
1090             fits->AliMatrix(fidrot[irot],90.0,0.0,90.0+thet,90.0,thet,90.0);
1091             break;
1092         case 1: //Rotate about y-axis, y-axis does not change.
1093             fits->AliMatrix(fidrot[irot],-90.0-thet,0.0,90.0,90.0,thet,90.0);
1094             break;
1095         case 2: //Rotate about z-axis, z-axis does not change.
1096             fits->AliMatrix(fidrot[irot],90.0,thet,90.0,-thet-90.0,0.0,0.0);
1097             break;
1098         default:
1099             Error("Matrix","axis must be either 0, 1, or 2. for matrix=%d",
1100                   irot);
1101             break;
1102         } // end switch
1103     } // end if
1104 }
1105 //______________________________________________________________________
1106 void AliITSBaseGeometry::Matrix(Int_t irot,Double_t rot[3][3]){
1107     // Defines a Geant rotation matrix. checks to see if it is the unit
1108     // matrix. If so, then no additonal matrix is defined. Stores rotation 
1109     // matrix irot in the data structure JROTM. If the matrix is not 
1110     // orthonormal, it will be corrected by setting y' perpendicular to x' 
1111     // and z' = x' X y'. A warning message is printed in this case.
1112     // Inputs:
1113     //   Int_t irot         Intex specifing which rotation matrix.
1114     //   Double_t rot[3][3] The 3 by 3 rotation matrix.
1115     // Outputs:
1116     //    none.
1117     // Return:
1118     //    none.
1119
1120     if(rot[0][0]==1.0&&rot[1][1]==1.0&&rot[2][2]==1.0&&
1121        rot[0][1]==0.0&&rot[0][2]==0.0&&rot[1][0]==0.0&&
1122        rot[1][2]==0.0&&rot[2][0]==0.0&&rot[2][1]==0.0){
1123         fidrot[irot] = 0; // Unit matrix
1124     }else{
1125         Double_t si,c=180./TMath::Pi();
1126         Double_t ang[6];
1127
1128         ang[1] = TMath::ATan2(rot[0][1],rot[0][0]);
1129         if(TMath::Cos(ang[1])!=0.0) si = rot[0][0]/TMath::Cos(ang[1]);
1130         else si = rot[0][1]/TMath::Sin(ang[1]);
1131         ang[0] = TMath::ATan2(si,rot[0][2]);
1132
1133         ang[3] = TMath::ATan2(rot[1][1],rot[1][0]);
1134         if(TMath::Cos(ang[3])!=0.0) si = rot[1][0]/TMath::Cos(ang[3]);
1135         else si = rot[1][1]/TMath::Sin(ang[3]);
1136         ang[2] = TMath::ATan2(si,rot[1][2]);
1137
1138         ang[5] = TMath::ATan2(rot[2][1],rot[2][0]);
1139         if(TMath::Cos(ang[5])!=0.0) si = rot[2][0]/TMath::Cos(ang[5]);
1140         else si = rot[2][1]/TMath::Sin(ang[5]);
1141         ang[4] = TMath::ATan2(si,rot[2][2]);
1142
1143         for(Int_t i=0;i<6;i++) {ang[i] *= c; if(ang[i]<0.0) ang[i] += 360.;}
1144         fits->AliMatrix(fidrot[irot],ang[0],ang[1],ang[2],ang[3],
1145                         ang[4],ang[5]);
1146     } // end if
1147 }
1148 //______________________________________________________________________
1149 Float_t AliITSBaseGeometry::GetA(Int_t z){
1150     // Returns the isotopicaly averaged atomic number.
1151     // Inputs:
1152     //    Int_t z  Elemental number
1153     // Outputs:
1154     //    none.
1155     // Return:
1156     //    The atomic mass number.
1157     const Float_t A[]={ 1.00794 ,  4.0026902,  6.941   ,  9.012182 , 10.811   ,
1158                        12.01007 , 14.00674  , 15.9994  , 18.9984032, 20.1797  ,
1159                        22.98970 , 24.3050   , 26.981538, 28.0855   , 30.973761,
1160                        32.066   , 35.4527   , 39.948   , 39.0983   , 40.078   ,
1161                        44.95591 , 47.867    , 50.9415  , 51.9961   , 54.938049,
1162                        55.845   , 58.933200 , 58.6934  , 63.546    , 65.39    ,
1163                        69.723   , 72.61     , 74.92160 , 78.96     , 79.904   ,
1164                        83.80    , 85.4678   , 87.62    , 88.9085   , 91.224   ,
1165                        92.90638 , 95.94     , 97.907215, 101.07    ,102.90550 ,
1166                       106.42    ,107.8682   ,112.411   ,114.818    ,118.710   ,
1167                       121.760   ,127.60     ,126.90447 ,131.29     ,132.90545 ,
1168                       137.327   ,138.9055   ,140.116   ,140.90765  ,144.24    ,
1169                       144.912746,150.36     ,151.964   ,157.25     ,158.92534 ,
1170                       162.50    ,164.93032  ,167.26    ,168.93421  ,173.04    ,
1171                       174.967   ,178.49     ,180.9479 ,183.84      ,186.207   ,
1172                       190.23    ,192.217    ,195.078  ,196.96655   ,200.59    ,
1173                       204.3833  ,207.2      ,208.98038,208.982415  ,209.987131,
1174                       222.017570,223.019731 ,226.025402,227.027747 ,232.0381  ,
1175                       231.03588 ,238.0289};
1176
1177     if(z<1||z>92){
1178         Error("GetA","z must be 0<z<93. z=%d",z);
1179         return 0.0;
1180     } // end if
1181     return A[z-1];
1182 }
1183 //______________________________________________________________________
1184 Float_t AliITSBaseGeometry::GetStandardMaxStepSize(Int_t istd){
1185     // Returns one of a set of standard Maximum Step Size values.
1186     // Inputs:
1187     //   Int_t istd  Index to indecate which standard.
1188     // Outputs:
1189     //    none.
1190     // Return:
1191     //    The appropreate standard Maximum Step Size value [cm].
1192     Float_t t[]={1.0, // default
1193                  0.0075, // Silicon detectors...
1194                  1.0, // Air in central detectors region
1195                  1.0  // Material in non-centeral region
1196     };
1197     return t[istd];
1198 }
1199 //______________________________________________________________________
1200 Float_t AliITSBaseGeometry::GetStandardThetaMax(Int_t istd){
1201     // Returns one of a set of standard Theata Max values.
1202     // Inputs:
1203     //   Int_t istd  Index to indecate which standard.
1204     // Outputs:
1205     //    none.
1206     // Return:
1207     //    The appropreate standard Theta max value [degrees].
1208     Float_t t[]={0.1, // default
1209                  0.1, // Silicon detectors...
1210                  0.1, // Air in central detectors region
1211                  1.0  // Material in non-centeral region
1212     };
1213     return t[istd];
1214 }
1215 //______________________________________________________________________
1216 Float_t AliITSBaseGeometry::GetStandardEfraction(Int_t istd){
1217     // Returns one of a set of standard E fraction values.
1218     // Inputs:
1219     //   Int_t istd  Index to indecate which standard.
1220     // Outputs:
1221     //    none.
1222     // Return:
1223     //    The appropreate standard E fraction value [#].
1224     Float_t t[]={0.1, // default
1225                  0.1, // Silicon detectors...
1226                  0.1, // Air in central detectors region
1227                  0.5  // Material in non-centeral region
1228     };
1229     return t[istd];
1230 }
1231 Float_t AliITSBaseGeometry::GetStandardEpsilon(Int_t istd){
1232     // Returns one of the standard Epsilon valuse
1233     // Inputs:
1234     //    Int_t istd  index of standard cuts to get
1235     // Output:
1236     //    none.
1237     // Return:
1238     //    Float_t the standard Epsilon cut value.
1239     Float_t t[]={1.0E-4, // default
1240                  1.0E-4, // Silicon detectors...
1241                  1.0E-4, // Air in central detector region
1242                  1.0E-3, // Material in non-cneteral regions
1243     };
1244
1245     return t[istd];
1246 }
1247 //______________________________________________________________________
1248 void AliITSBaseGeometry::Element(Int_t imat,const char* name,Int_t z,
1249                                  Double_t dens,Int_t istd){
1250     // Defines a Geant single element material and sets its Geant medium
1251     // proporties. The average atomic A is assumed to be given by their
1252     // natural abundances. Things like the radiation length are calculated
1253     // for you.
1254     // Inputs:
1255     //    Int_t imat       Material number.
1256     //    const char* name Material name. No need to add a $ at the end.
1257     //    Int_t z          The elemental number.
1258     //    Double_t dens    The density of the material [g/cm^3].
1259     //    Int_t istd       Defines which standard set of transport parameters
1260     //                     which should be used.
1261     // Output:
1262     //     none.
1263     // Return:
1264     //     none.
1265     Float_t rad,Z,A=GetA(z),tmax,stemax,deemax,epsilon;
1266     char *name2;
1267     Int_t len;
1268
1269     len = strlen(name)+1;
1270     name2 = new char[len];
1271     strncpy(name2,name,len-1);
1272     name2[len-1] = '\0';
1273     name2[len-2] = '$';
1274     Z = (Float_t)z;
1275     rad = GetRadLength(z)/dens;
1276     fits->AliMaterial(imat,name2,A,Z,dens,rad,0.0,0,0);
1277     tmax    = GetStandardThetaMax(istd);    // degree
1278     stemax  = GetStandardMaxStepSize(istd);  // cm
1279     deemax  = GetStandardEfraction(istd);     // ratio
1280     epsilon = GetStandardEpsilon(istd);       //
1281     fits->AliMedium(imat,name2,imat,0,gAlice->Field()->Integ(),
1282                     gAlice->Field()->Max(),tmax,stemax,deemax,epsilon,0.0);
1283     delete[] name2;
1284 }
1285 //______________________________________________________________________
1286 void AliITSBaseGeometry::MixtureByWeight(Int_t imat,const char* name,Int_t *z,
1287                                 Double_t *w,Double_t dens,Int_t n,Int_t istd){
1288     // Defines a Geant material by a set of elements and weights, and sets 
1289     // its Geant medium proporties. The average atomic A is assumed to be 
1290     // given by their natural abundances. Things like the radiation length 
1291     // are calculated for you.
1292     // Inputs:
1293     //    Int_t imat       Material number.
1294     //    const char* name Material name. No need to add a $ at the end.
1295     //    Int_t *z         Array of The elemental numbers.
1296     //    Double_t *w      Array of relative weights.
1297     //    Double_t dens    The density of the material [g/cm^3].
1298     //    Int_t n          the number of elements making up the mixture.
1299     //    Int_t istd       Defines which standard set of transport parameters
1300     //                     which should be used.   
1301     // Output:
1302     //     none.
1303     // Return:
1304     //     none.
1305     Float_t *Z,*A,*W,tmax,stemax,deemax,epsilon;
1306     char *name2;
1307     Int_t len,i;
1308     Z = new Float_t[n];
1309     A = new Float_t[n];
1310     W = new Float_t[n];
1311
1312     len = strlen(name)+1;
1313     name2 = new char[len];
1314     strncpy(name2,name,len-1);
1315     name2[len-1] = '\0';
1316     name2[len-2] = '$';
1317     for(i=0;i<n;i++){Z[i] = (Float_t)z[i];A[i] = (Float_t)GetA(z[i]);
1318                      W[i] = (Float_t)w[i];}
1319     fits->AliMixture(imat,name2,A,Z,dens,n,W);
1320     tmax    = GetStandardThetaMax(istd);    // degree
1321     stemax  = GetStandardMaxStepSize(istd);  // cm
1322     deemax  = GetStandardEfraction(istd);     // #
1323     epsilon = GetStandardEpsilon(istd);
1324     fits->AliMedium(imat,name2,imat,0,gAlice->Field()->Integ(),
1325               gAlice->Field()->Max(),tmax,stemax,deemax,epsilon,0.0);
1326     delete[] name2;
1327     delete[] Z;
1328     delete[] A;
1329     delete[] W;
1330 }
1331 //______________________________________________________________________
1332 void AliITSBaseGeometry::MixtureByNumber(Int_t imat,const char* name,Int_t *z,
1333                                 Int_t *w,Double_t dens,Int_t n,Int_t istd){
1334     // Defines a Geant material by a set of elements and number, and sets 
1335     // its Geant medium proporties. The average atomic A is assumed to be 
1336     // given by their natural abundances. Things like the radiation length 
1337     // are calculated for you.
1338     // Inputs:
1339     //    Int_t imat       Material number.
1340     //    const char* name Material name. No need to add a $ at the end.
1341     //    Int_t *z         Array of The elemental numbers.
1342     //    Int_t_t *w       Array of relative number.
1343     //    Double_t dens    The density of the material [g/cm^3].
1344     //    Int_t n          the number of elements making up the mixture.
1345     //    Int_t istd       Defines which standard set of transport parameters
1346     //                     which should be used.   
1347     // Output:
1348     //     none.
1349     // Return:
1350     //     none.
1351     Float_t *Z,*A,*W,tmax,stemax,deemax,epsilon;
1352     char *name2;
1353     Int_t len,i;
1354     Z = new Float_t[n];
1355     A = new Float_t[n];
1356     W = new Float_t[n];
1357
1358     len = strlen(name)+1;
1359     name2 = new char[len];
1360     strncpy(name2,name,len-1);
1361     name2[len-1] = '\0';
1362     name2[len-2] = '$';
1363     for(i=0;i<n;i++){Z[i] = (Float_t)z[i];A[i] = (Float_t)GetA(z[i]);
1364                      W[i] = (Float_t)w[i];}
1365     fits->AliMixture(imat,name2,A,Z,dens,-n,W);
1366     tmax    = GetStandardThetaMax(istd);    // degree
1367     stemax  = GetStandardMaxStepSize(istd);  // cm
1368     deemax  = GetStandardEfraction(istd);     // #
1369     epsilon = GetStandardEpsilon(istd);
1370     fits->AliMedium(imat,name2,imat,0,gAlice->Field()->Integ(),
1371                     gAlice->Field()->Max(),tmax,stemax,deemax,epsilon,0.0);
1372     delete[] name2;
1373     delete[] Z;
1374     delete[] A;
1375     delete[] W;
1376 }
1377 //______________________________________________________________________
1378 Double_t AliITSBaseGeometry::RadLength(Int_t iz,Double_t a){
1379     // Computes the radiation length in accordance to the PDG 2000 Section
1380     // 23.4.1 p. 166. Transladed from the c code of Flavio Tosello.
1381     // Inputs:
1382     //    Int_t iz    The elemental number
1383     //    Dougle_t    The elemental average atomic mass number
1384     // Outputs:
1385     // Return:
1386     //    Double_t returns the radiation length of the element iz in
1387     //             [gm/cm^2].
1388     Double_t z = (Double_t)iz;
1389     Double_t alphaz = fAlpha*z;
1390     Double_t alphaz2 = alphaz*alphaz;
1391     Double_t c0 = +0.20206,c1 = -0.0369,c2 = +0.0083,c3 = -0.0020;
1392     Double_t z12,z23,l,lp,c;
1393
1394     c = alphaz2*(1./(1.+alphaz2) + c0 + c1*alphaz2 + c2*alphaz2*alphaz2
1395                   +c3*alphaz2*alphaz2*alphaz2);
1396     z12 = TMath::Exp(TMath::Log(z)/3.0);
1397     z23 = z12*z12;
1398     switch (iz){
1399     case 1: //Hydrogen
1400         l  = 5.31;
1401         lp = 6.144;
1402         break;
1403     case 2: //Helium
1404         l  = 4.79;
1405         lp = 5,621;
1406         break;
1407     case 3: //Lithium
1408         l  = 4.74;
1409         lp = 5.805;
1410         break;
1411     case 4: //Berilium
1412         l  = 4.71;
1413         lp = 5.924;
1414         break;
1415     default: //Others
1416         l  = TMath::Log(184.15/z12);
1417         lp = TMath::Log(1194.0/z23);
1418         break;
1419     } // end switch
1420     Double_t re2,b,r,xz;
1421
1422     re2 = fRe*fRe;
1423     b = 4.0*fAlpha*re2*fNa/a;
1424     r = b*z*(z*(l-c)+lp);
1425     xz = 1.0/r;
1426     return xz; // [gm/cm^2]
1427 }