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