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