]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSv11Geometry.cxx
Optionally produce a tree with cluster topologies (for the compression studies)
[u/mrichter/AliRoot.git] / ITS / AliITSv11Geometry.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17  $Id$ 
18 */
19
20
21 ////////////////////////////////////////////////////////////////////////
22 //  This class is a base class for the ITS geometry version 11. It 
23 //  contains common/standard functions used in many places in defining 
24 //  the ITS geometry, version 11. Large posions of the ITS geometry, 
25 //  version 11, should be derived from this class so as to make maximum 
26 //  use of these common functions. This class also defines the proper 
27 //  conversion valuse such, to cm and degrees, such that the most usefull 
28 //  units, those used in the Engineering drawings, can be used.
29 ////////////////////////////////////////////////////////////////////////
30
31
32 #include <Riostream.h>
33 #include <TMath.h>
34 #include <TArc.h>
35 #include <TLine.h>
36 #include <TArrow.h>
37 #include <TCanvas.h>
38 #include <TText.h>
39 #include <TGeoPcon.h>
40 #include <TGeoCone.h>
41 #include <TGeoTube.h> // contaings TGeoTubeSeg
42 #include <TGeoArb8.h>
43 #include <TGeoElement.h>
44 #include <TGeoMaterial.h>
45 #include <TPolyMarker.h>
46 #include <TPolyLine.h>
47 #include <AliMagF.h>
48 #include <AliRun.h>
49 #include "AliITSv11Geometry.h"
50
51 using std::endl;
52 using std::cout;
53 using std::cin;
54 ClassImp(AliITSv11Geometry)
55
56 const Double_t AliITSv11Geometry::fgkmicron = 1.0E-4;
57 const Double_t AliITSv11Geometry::fgkmm = 0.10;
58 const Double_t AliITSv11Geometry::fgkcm = 1.00;
59 const Double_t AliITSv11Geometry::fgkDegree = 1.0;
60 const Double_t AliITSv11Geometry::fgkRadian = 180./3.14159265358979323846;
61 const Double_t AliITSv11Geometry::fgkgcm3 = 1.0; // assume default is g/cm^3
62 const Double_t AliITSv11Geometry::fgkKgm3 = 1.0E+3;// assume Kg/m^3
63 const Double_t AliITSv11Geometry::fgkKgdm3 = 1.0; // assume Kg/dm^3
64 const Double_t AliITSv11Geometry::fgkCelsius = 1.0; // Assume default is C
65 const Double_t AliITSv11Geometry::fgkPascal  = 1.0E-3; // Assume kPascal
66 const Double_t AliITSv11Geometry::fgkKPascal = 1.0;    // Asume kPascal
67 const Double_t AliITSv11Geometry::fgkeV      = 1.0E-9; // GeV default
68 const Double_t AliITSv11Geometry::fgkKeV     = 1.0e-6; // GeV default
69 const Double_t AliITSv11Geometry::fgkMeV     = 1.0e-3; // GeV default
70 const Double_t AliITSv11Geometry::fgkGeV     = 1.0;    // GeV default
71 //______________________________________________________________________
72 void AliITSv11Geometry::IntersectLines(Double_t m, Double_t x0, Double_t y0,
73                                        Double_t n, Double_t x1, Double_t y1,
74                                        Double_t &xi, Double_t &yi)const{
75     // Given the two lines, one passing by (x0,y0) with slope m and
76     // the other passing by (x1,y1) with slope n, returns the coordinates
77     // of the intersecting point (xi,yi)
78     // Inputs:
79     //    Double_t   m     The slope of the first line
80     //    Double_t  x0,y0  The x and y coord. of the first point
81     //    Double_t   n     The slope of the second line
82     //    Double_t  x1,y1  The x and y coord. of the second point
83     // Outputs:
84     //    The coordinates xi and yi of the intersection point
85     // Return:
86     //    none.
87     // Created:      14 Dec 2009  Mario Sitta
88
89     if (TMath::Abs(m-n) < 0.000001) {
90       AliError(Form("Lines are parallel: m = %f n = %f\n",m,n));
91       return;
92     }
93
94     xi = (y1 - n*x1 - y0 + m*x0)/(m - n);
95     yi = y0 + m*(xi - x0);
96
97     return;
98 }
99 //______________________________________________________________________
100 Bool_t AliITSv11Geometry::IntersectCircle(Double_t m, Double_t x0, Double_t y0,
101                                           Double_t rr, Double_t xc, Double_t yc,
102                                           Double_t &xi1, Double_t &yi1,
103                                           Double_t &xi2, Double_t &yi2){
104     // Given a lines  passing by (x0,y0) with slope m and a circle with
105     // radius rr and center (xc,yc), returns the coordinates of the
106     // intersecting points (xi1,yi1) and (xi2,yi2) (xi1 > xi2)
107     // Inputs:
108     //    Double_t   m     The slope of the line
109     //    Double_t  x0,y0  The x and y coord. of the point
110     //    Double_t   rr     The radius of the circle
111     //    Double_t  xc,yc  The x and y coord. of the center of circle
112     // Outputs:
113     //    The coordinates xi and yi of the intersection points
114     // Return:
115     //    kFALSE if the line does not intercept the circle, otherwise kTRUE
116     // Created:      18 Dec 2009  Mario Sitta
117
118     Double_t p = m*x0 - y0;
119     Double_t q = m*m + 1;
120
121     p = p-m*xc+yc;
122
123     Double_t delta = m*m*p*p - q*(p*p - rr*rr);
124
125     if (delta < 0)
126       return kFALSE;
127     else {
128       Double_t root = TMath::Sqrt(delta);
129       xi1 = (m*p + root)/q + xc;
130       xi2 = (m*p - root)/q + xc;
131       yi1 = m*(xi1 - x0) + y0;
132       yi2 = m*(xi2 - x0) + y0;
133       return kTRUE;
134     }
135 }
136 //______________________________________________________________________
137 Double_t AliITSv11Geometry::Yfrom2Points(Double_t x0,Double_t y0,
138                                          Double_t x1,Double_t y1,
139                                          Double_t x)const{
140     // Given the two points (x0,y0) and (x1,y1) and the location x, returns
141     // the value y corresponding to that point x on the line defined by the
142     // two points.
143     // Inputs:
144     //    Double_t  x0  The first x value defining the line
145     //    Double_t  y0  The first y value defining the line
146     //    Double_t  x1  The second x value defining the line
147     //    Double_t  y1  The second y value defining the line
148     //    Double_t   x  The x value for which the y value is wanted.
149     // Outputs:
150     //    none.
151     // Return:
152     //    The value y corresponding to the point x on the line defined by
153     //    the two points (x0,y0) and (x1,y1).
154
155     if(x0==x1 && y0==y1) {
156         printf("Error: AliITSv11Geometry::Yfrom2Ponts The two points are "
157                "the same (%e,%e) and (%e,%e)",x0,y0,x1,y1);
158         return 0.0;
159     } // end if
160     if(x0==x1){
161         printf("Warning: AliITSv11Geometry::Yfrom2Points x0=%e == x1=%e. "
162                "line vertical ""returning mean y",x0,x1);
163         return 0.5*(y0+y1);
164     }// end if x0==x1
165     Double_t m = (y0-y1)/(x0-x1);
166     return m*(x-x0)+y0;
167 }
168 //______________________________________________________________________
169 Double_t AliITSv11Geometry::Xfrom2Points(Double_t x0,Double_t y0,
170                                          Double_t x1,Double_t y1,
171                                          Double_t y)const{
172     // Given the two points (x0,y0) and (x1,y1) and the location y, returns
173     // the value x corresponding to that point y on the line defined by the
174     // two points.
175     // Inputs:
176     //    Double_t  x0  The first x value defining the line
177     //    Double_t  y0  The first y value defining the line
178     //    Double_t  x1  The second x value defining the line
179     //    Double_t  y1  The second y value defining the line
180     //    Double_t   y  The y value for which the x value is wanted.
181     // Outputs:
182     //    none.
183     // Return:
184     //    The value x corresponding to the point y on the line defined by
185     //    the two points (x0,y0) and (x1,y1).
186
187     if(x0==x1 && y0==y1) {
188         printf("Error: AliITSv11Geometry::Yfrom2Ponts The two points are "
189                "the same (%e,%e) and (%e,%e)",x0,y0,x1,y1);
190         return 0.0;
191     } // end if
192     if(y0==y1){
193         printf("Warrning: AliITSv11Geometry::Yfrom2Points y0=%e == y1=%e. "
194                "line horizontal returning mean x",y0,y1);
195         return 0.5*(x0+x1);
196     }// end if y0==y1
197     Double_t m = (x0-x1)/(y0-y1);
198     return m*(y-y0)+x0;
199 }
200 //______________________________________________________________________
201 Double_t AliITSv11Geometry::RmaxFrom2Points(const TGeoPcon *p,Int_t i1,
202                                             Int_t i2,Double_t z)const{
203     // functions Require at parts of Volume A to be already defined.
204     // Retruns the value of Rmax corresponding to point z alone the line
205     // defined by the two points p.Rmax(i1),p-GetZ(i1) and p->GetRmax(i2),
206     // p->GetZ(i2).
207     // Inputs:
208     //    TGeoPcon *p  The Polycone where the two points come from
209     //    Int_t    i1  Point 1
210     //    Int_t    i2  Point 2
211     //    Double_t  z  The value of z for which Rmax is to be found
212     // Outputs:
213     //    none.
214     // Return:
215     //    Double_t Rmax the value corresponding to z
216     Double_t d0,d1,d2,r;
217
218     d0 = p->GetRmax(i1)-p->GetRmax(i2);// cout <<"L263: d0="<<d0<<endl;
219     d1 = z-p->GetZ(i2);// cout <<"L264: d1="<<d1<<endl;
220     d2 = p->GetZ(i1)-p->GetZ(i2);// cout <<"L265: d2="<<d2<<endl;
221     r  = p->GetRmax(i2) + d1*d0/d2;// cout <<"L266: r="<<r<<endl;
222     return r;
223 }
224 //______________________________________________________________________
225 Double_t AliITSv11Geometry::RminFrom2Points(const TGeoPcon *p,Int_t i1,
226                                             Int_t i2,Double_t z)const{
227     // Retruns the value of Rmin corresponding to point z alone the line
228     // defined by the two points p->GetRmin(i1),p->GetZ(i1) and 
229     // p->GetRmin(i2),  p->GetZ(i2).
230     // Inputs:
231     //    TGeoPcon *p  The Polycone where the two points come from
232     //    Int_t    i1  Point 1
233     //    Int_t    i2  Point 2
234     //    Double_t  z  The value of z for which Rmax is to be found
235     // Outputs:
236     //    none.
237     // Return:
238     //    Double_t Rmax the value corresponding to z
239
240     return p->GetRmin(i2)+(p->GetRmin(i1)-p->GetRmin(i2))*(z-p->GetZ(i2))/
241      (p->GetZ(i1)-p->GetZ(i2));
242 }
243 //______________________________________________________________________
244 Double_t AliITSv11Geometry::RFrom2Points(const Double_t *p,const Double_t *az,
245                                          Int_t i1,Int_t i2,Double_t z)const{
246     // Retruns the value of Rmin corresponding to point z alone the line
247     // defined by the two points p->GetRmin(i1),p->GetZ(i1) and 
248     // p->GetRmin(i2), p->GetZ(i2).
249     // Inputs:
250     //    Double_t az  Array of z values
251     //    Double_t  r  Array of r values
252     //    Int_t    i1  First Point in arrays
253     //    Int_t    i2  Second Point in arrays
254     //    Double_t z   Value z at which r is to be found
255     // Outputs:
256     //    none.
257     // Return:
258     //    The value r corresponding to z and the line defined by the two points
259
260     return p[i2]+(p[i1]-p[i2])*(z-az[i2])/(az[i1]-az[i2]);
261 }
262 //______________________________________________________________________
263 Double_t AliITSv11Geometry::Zfrom2MinPoints(const TGeoPcon *p,Int_t i1,
264                                             Int_t i2,Double_t r)const{
265     // Retruns the value of Z corresponding to point R alone the line
266     // defined by the two points p->GetRmin(i1),p->GetZ(i1) and 
267     // p->GetRmin(i2),p->GetZ(i2)
268     // Inputs:
269     //    TGeoPcon *p  The Poly cone where the two points come from.
270     //    Int_t    i1  First Point in arrays
271     //    Int_t    i2  Second Point in arrays
272     //    Double_t r   Value r min at which z is to be found
273     // Outputs:
274     //    none.
275     // Return:
276     //    The value z corresponding to r min and the line defined by 
277     //    the two points
278
279     return p->GetZ(i2)+(p->GetZ(i1)-p->GetZ(i2))*(r-p->GetRmin(i2))/
280      (p->GetRmin(i1)-p->GetRmin(i2));
281 }
282 //______________________________________________________________________
283 Double_t AliITSv11Geometry::Zfrom2MaxPoints(const TGeoPcon *p,Int_t i1,
284                                             Int_t i2,Double_t r)const{
285     // Retruns the value of Z corresponding to point R alone the line
286     // defined by the two points p->GetRmax(i1),p->GetZ(i1) and 
287     // p->GetRmax(i2),p->GetZ(i2)
288     // Inputs:
289     //    TGeoPcon *p  The Poly cone where the two points come from.
290     //    Int_t    i1  First Point in arrays
291     //    Int_t    i2  Second Point in arrays
292     //    Double_t r   Value r max at which z is to be found
293     // Outputs:
294     //    none.
295     // Return:
296     //    The value z corresponding to r max and the line defined by 
297     //    the two points
298
299     return p->GetZ(i2)+(p->GetZ(i1)-p->GetZ(i2))*(r-p->GetRmax(i2))/
300      (p->GetRmax(i1)-p->GetRmax(i2));
301 }
302 //______________________________________________________________________
303 Double_t AliITSv11Geometry::Zfrom2Points(const Double_t *z,const Double_t *ar,
304                                          Int_t i1,Int_t i2,Double_t r)const{
305     // Retruns the value of z corresponding to point R alone the line
306     // defined by the two points p->GetRmax(i1),p->GetZ(i1) and 
307     // p->GetRmax(i2),p->GetZ(i2)
308     // Inputs:
309     //    Double_t  z  Array of z values
310     //    Double_t ar  Array of r values
311     //    Int_t    i1  First Point in arrays
312     //    Int_t    i2  Second Point in arrays
313     //    Double_t r   Value r at which z is to be found
314     // Outputs:
315     //    none.
316     // Return:
317     //    The value z corresponding to r and the line defined by the two points
318
319     return z[i2]+(z[i1]-z[i2])*(r-ar[i2])/(ar[i1]-ar[i2]);
320 }
321 //______________________________________________________________________
322 Double_t AliITSv11Geometry::RmaxFromZpCone(const TGeoPcon *p,int ip,
323                                            Double_t tc,Double_t z,
324                                            Double_t th)const{
325     // General Outer Cone surface equation Rmax.
326     // Intputs:
327     //     TGeoPcon  *p   The poly cone where the initial point comes from
328     //     Int_t     ip   The index in p to get the point location
329     //     Double_t  tc   The angle of that part of the cone is at
330     //     Double_t   z   The value of z to compute Rmax from
331     //     Double_t  th   The perpendicular distance the parralell line is
332     //                    from the point ip.
333     // Outputs:
334     //     none.
335     // Return:
336     //     The value Rmax correstponding to the line at angle th, offeset by
337     //     th, and the point p->GetZ/Rmin[ip] at the location z.
338     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
339     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
340
341     return -tantc*(z-p->GetZ(ip))+p->GetRmax(ip)+th/costc;
342 }
343 //______________________________________________________________________
344 Double_t AliITSv11Geometry::RFromZpCone(const Double_t *ar,
345                                         const Double_t *az,int ip,
346                                         Double_t tc,Double_t z,
347                                         Double_t th)const{
348     // General Cone surface equation R(z).
349     // Intputs:
350     //     Double_t  ar   The array of R values
351     //     Double_t  az   The array of Z values
352     //     Int_t     ip   The index in p to get the point location
353     //     Double_t  tc   The angle of that part of the cone is at
354     //     Double_t   z   The value of z to compute R from
355     //     Double_t  th   The perpendicular distance the parralell line is
356     //                    from the point ip.
357     // Outputs:
358     //     none.
359     // Return:
360     //     The value R correstponding to the line at angle th, offeset by
361     //     th, and the point p->GetZ/Rmax[ip] at the locatin z.
362     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
363     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
364
365     return -tantc*(z-az[ip])+ar[ip]+th/costc;
366 }
367 //______________________________________________________________________
368 Double_t AliITSv11Geometry::RminFromZpCone(const TGeoPcon *p,Int_t ip,
369                                            Double_t tc,Double_t z,
370                                            Double_t th)const{
371     // General Inner Cone surface equation Rmin.
372     // Intputs:
373     //     TGeoPcon  *p   The poly cone where the initial point comes from
374     //     Int_t     ip   The index in p to get the point location
375     //     Double_t  tc   The angle of that part of the cone is at
376     //     Double_t   z   The value of z to compute Rmin from
377     //     Double_t  th   The perpendicular distance the parralell line is
378     //                    from the point ip.
379     // Outputs:
380     //     none.
381     // Return:
382     //     The value Rmin correstponding to the line at angle th, offeset by
383     //     th, and the point p->GetZ/Rmin[ip] at the location z.
384     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
385     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
386
387     return -tantc*(z-p->GetZ(ip))+p->GetRmin(ip)+th/costc;
388 }
389 //______________________________________________________________________
390 Double_t AliITSv11Geometry::ZFromRmaxpCone(const TGeoPcon *p,int ip,
391                                            Double_t tc,Double_t r,
392                                            Double_t th)const{
393     // General Outer cone Surface equation for z.
394     // Intputs:
395     //     TGeoPcon  *p   The poly cone where the initial point comes from
396     //     Int_t     ip   The index in p to get the point location
397     //     Double_t  tc   The angle of that part of the cone is at
398     //     Double_t   r   The value of Rmax to compute z from
399     //     Double_t  th   The perpendicular distance the parralell line is
400     //                    from the point ip.
401     // Outputs:
402     //     none.
403     // Return:
404     //     The value Z correstponding to the line at angle th, offeset by
405     //     th, and the point p->GetZ/Rmax[ip] at the location r.
406     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
407     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
408
409     return p->GetZ(ip)+(p->GetRmax(ip)+th/costc-r)/tantc;
410 }
411 //______________________________________________________________________
412 Double_t AliITSv11Geometry::ZFromRmaxpCone(const Double_t *ar,
413                                            const Double_t *az,int ip,
414                                            Double_t tc,Double_t r,
415                                            Double_t th)const{
416     // General Outer cone Surface equation for z.
417     // Intputs:
418     //     Double_t  ar   The array of R values
419     //     Double_t  az   The array of Z values
420     //     Int_t     ip   The index in p to get the point location
421     //     Double_t  tc   The angle of that part of the cone is at
422     //     Double_t   r   The value of Rmax to compute z from
423     //     Double_t  th   The perpendicular distance the parralell line is
424     //                    from the point ip.
425     // Outputs:
426     //     none.
427     // Return:
428     //     The value Z correstponding to the line at angle th, offeset by
429     //     th, and the point p->GetZ/Rmax[ip] at the locatin r.
430     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
431     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
432
433     return az[ip]+(ar[ip]+th/costc-r)/tantc;
434 }
435 //______________________________________________________________________
436 Double_t AliITSv11Geometry::ZFromRminpCone(const TGeoPcon *p,int ip,
437                                            Double_t tc,Double_t r,
438                                            Double_t th)const{
439     // General Inner cone Surface equation for z.
440     // Intputs:
441     //     TGeoPcon  *p   The poly cone where the initial point comes from
442     //     Int_t     ip   The index in p to get the point location
443     //     Double_t  tc   The angle of that part of the cone is at
444     //     Double_t   r   The value of Rmin to compute z from
445     //     Double_t  th   The perpendicular distance the parralell line is
446     //                    from the point ip.
447     // Outputs:
448     //     none.
449     // Return:
450     //     The value Z correstponding to the line at angle th, offeset by
451     //     th, and the point p->GetZ/Rmin[ip] at the location r.
452     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
453     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
454
455     return p->GetZ(ip)+(p->GetRmin(ip)+th/costc-r)/tantc;
456 }
457 //______________________________________________________________________
458 void AliITSv11Geometry::RadiusOfCurvature(Double_t rc,Double_t theta0,
459                                           Double_t z0,Double_t r0,
460                                           Double_t theta1,Double_t &z1,
461                                           Double_t &r1)const{
462     // Given a initial point z0,r0, the initial angle theta0, and the radius
463     // of curvature, returns the point z1, r1 at the angle theta1. Theta
464     // measured from the r axis in the clock wise direction [degrees].
465     // Inputs:
466     //    Double_t rc     The radius of curvature
467     //    Double_t theta0 The starting angle (degrees)
468     //    Double_t z0     The value of z at theta0
469     //    Double_t r0     The value of r at theta0
470     //    Double_t theta1 The ending angle (degrees)
471     // Outputs:
472     //    Double_t &z1  The value of z at theta1
473     //    Double_t &r1  The value of r at theta1
474     // Return:
475     //    none.
476
477     z1 = rc*(TMath::Sin(theta1*TMath::DegToRad())-TMath::Sin(theta0*TMath::DegToRad()))+z0;
478     r1 = rc*(TMath::Cos(theta1*TMath::DegToRad())-TMath::Cos(theta0*TMath::DegToRad()))+r0;
479     return;
480 }
481 //______________________________________________________________________
482 void AliITSv11Geometry::InsidePoint(const TGeoPcon *p,Int_t i1,Int_t i2,
483                                     Int_t i3,Double_t c,TGeoPcon *q,Int_t j1,
484                                     Bool_t max)const{
485     // Given two lines defined by the points i1, i2,i3 in the TGeoPcon 
486     // class p that intersect at point p->GetZ(i2) return the point z,r 
487     // that is Cthick away in the TGeoPcon class q. If points i1=i2
488     // and max == kTRUE, then p->GetRmin(i1) and p->GetRmax(i2) are used.
489     // if points i2=i3 and max=kTRUE then points p->GetRmax(i2) and
490     // p->GetRmin(i3) are used. If i2=i3 and max=kFALSE, then p->GetRmin(i2)
491     // and p->GetRmax(i3) are used.
492     // Inputs:
493     //    TGeoPcon  *p  Class where points i1, i2, and i3 are taken from
494     //    Int_t     i1  First point in class p
495     //    Int_t     i2  Second point in class p
496     //    Int_t     i3  Third point in class p
497     //    Double_t  c   Distance inside the outer surface/inner suface
498     //                  that the point j1 is to be computed for.
499     //    TGeoPcon  *q  Pointer to class for results to be put into.
500     //    Int_t     j1  Point in class q where data is to be stored.
501     //    Bool_t    max if kTRUE, then a Rmax value is computed,
502     //                  else a Rmin valule is computed.
503     // Output:
504     //    TGeoPcon  *q  Pointer to class for results to be put into.
505     // Return:
506     //    none.
507     Double_t x0,y0,x1,y1,x2,y2,x,y;
508
509     if(max){
510         c = -c; //cout <<"L394 c="<<c<<endl;
511         y0 = p->GetRmax(i1);
512         if(i1==i2) y0 = p->GetRmin(i1); //cout <<"L396 y0="<<y0<<endl;
513         y1 = p->GetRmax(i2);  //cout <<"L397 y1="<<y1<<endl;
514         y2 = p->GetRmax(i3); //cout <<"L398 y2="<<y2<<endl;
515         if(i2==i3) y2 = p->GetRmin(i3); //cout <<"L399 y2="<<y2<<endl;
516     }else{ // min
517         y0 = p->GetRmin(i1); //cout <<"L401 y0="<<y0<<endl;
518         y1 = p->GetRmin(i2); //cout <<"L402 y1="<<y1<<endl;
519         y2 = p->GetRmin(i3);
520         if(i2==i3) y2 = p->GetRmax(i3); //cout <<"L404 y2="<<y2<<endl;
521     } // end if
522     x0 = p->GetZ(i1); //cout <<"L406 x0="<<x0<<endl;
523     x1 = p->GetZ(i2); //cout <<"L407 x1="<<x1<<endl;
524     x2 = p->GetZ(i3); //cout <<"L408 x2="<<x2<<endl;
525     //
526     InsidePoint(x0,y0,x1,y1,x2,y2,c,x,y);
527     q->Z(j1) = x;
528     if(max) q->Rmax(j1) = y;
529     else    q->Rmin(j1) = y;
530     return;
531 }
532 //----------------------------------------------------------------------
533 void AliITSv11Geometry::InsidePoint(Double_t x0,Double_t y0,
534                                     Double_t x1,Double_t y1,
535                                     Double_t x2,Double_t y2,Double_t c,
536                                     Double_t &x,Double_t &y)const{
537     // Given two intersecting lines defined by the points (x0,y0), (x1,y1) and
538     // (x1,y1), (x2,y2) {intersecting at (x1,y1)} the point (x,y) a distance
539     // c away is returned such that two lines a distance c away from the
540     // lines defined above intersect at (x,y).
541     // Inputs:
542     //    Double_t  x0 X point on the first intersecting sets of lines
543     //    Double_t  y0 Y point on the first intersecting sets of lines
544     //    Double_t  x1 X point on the first/second intersecting sets of lines
545     //    Double_t  y1 Y point on the first/second intersecting sets of lines
546     //    Double_t  x2 X point on the second intersecting sets of lines
547     //    Double_t  y2 Y point on the second intersecting sets of lines
548     //    Double_t  c  Distance the two sets of lines are from each other
549     // Output:
550     //    Double_t  x  X point for the intersecting sets of parellel lines
551     //    Double_t  y  Y point for the intersecting sets of parellel lines
552     // Return:
553     //    none.
554     Double_t dx01,dx12,dy01,dy12,r01,r12,m;
555
556     //printf("InsidePoint: x0=% #12.7g y0=% #12.7g x1=% #12.7g y1=% #12.7g "
557     //       "x2=% #12.7g y2=% #12.7g c=% #12.7g ",x0,y0,x1,y2,x2,y2,c);
558     dx01 = x0-x1; //cout <<"L410 dx01="<<dx01<<endl;
559     dx12 = x1-x2; //cout <<"L411 dx12="<<dx12<<endl;
560     dy01 = y0-y1; //cout <<"L412 dy01="<<dy01<<endl;
561     dy12 = y1-y2; //cout <<"L413 dy12="<<dy12<<endl;
562     r01  = TMath::Sqrt(dy01*dy01+dx01*dx01); //cout <<"L414 r01="<<r01<<endl;
563     r12  = TMath::Sqrt(dy12*dy12+dx12*dx12); //cout <<"L415 r12="<<r12<<endl;
564     m = dx12*dy01-dy12*dx01;
565     if(m*m<DBL_EPSILON){ // m == n
566         if(dy01==0.0){ // line are =
567             x = x1+c; //cout <<"L419 x="<<x<<endl;
568             y = y1; //cout <<"L420 y="<<y<<endl;
569             //printf("dy01==0.0 x=% #12.7g y=% #12.7g\n",x,y);
570             return;
571         }else if(dx01==0.0){
572             x = x1;
573             y = y1+c;
574             //printf("dx01==0.0 x=% #12.7g y=% #12.7g\n",x,y);
575             return;
576         }else{ // dx01!=0 and dy01 !=0.
577             x = x1-0.5*c*r01/dy01; //cout <<"L434 x="<<x<<endl;
578             y = y1+0.5*c*r01/dx01; //cout <<"L435 y="<<y<<endl;
579             //printf("m*m<DBL_E x=% #12.7g y=% #12.7g\n",x,y);
580         } // end if
581         return;
582     } //
583     x = x1+c*(dx12*r01-dx01*r12)/m; //cout <<"L442 x="<<x<<endl;
584     y = y1+c*(dy12*r01-dy01*r12)/m; //cout <<"L443 y="<<y<<endl;
585     //printf("          x=% #12.7g y=% #12.7g\n",x,y);
586     //cout <<"=============================================="<<endl;
587     return;
588 }
589 //----------------------------------------------------------------------
590 void AliITSv11Geometry:: PrintArb8(const TGeoArb8 *a)const{
591     // Prints out the content of the TGeoArb8. Usefull for debugging.
592     // Inputs:
593     //   TGeoArb8 *a
594     // Outputs:
595     //   none.
596     // Return:
597     //   none.
598
599     if(!GetDebug()) return;
600     printf("%s",a->GetName());
601     a->InspectShape();
602     return;
603 }
604 //----------------------------------------------------------------------
605 void AliITSv11Geometry:: PrintPcon(const TGeoPcon *a)const{
606     // Prints out the content of the TGeoPcon. Usefull for debugging.
607     // Inputs:
608     //   TGeoPcon *a
609     // Outputs:
610     //   none.
611     // Return:
612     //   none.
613   
614     if(!GetDebug()) return;
615     cout << a->GetName() << ": N=" << a->GetNz() << " Phi1=" << a->GetPhi1()
616          << ", Dphi=" << a->GetDphi() << endl;
617     cout << "i\t   Z   \t  Rmin \t  Rmax" << endl;
618     for(Int_t iii=0;iii<a->GetNz();iii++){
619         cout << iii << "\t" << a->GetZ(iii) << "\t" << a->GetRmin(iii)
620              << "\t" << a->GetRmax(iii) << endl;
621     } // end for iii
622     return;
623 }
624 //----------------------------------------------------------------------
625 void AliITSv11Geometry::PrintTube(const TGeoTube *a)const{
626     // Prints out the content of the TGeoTube. Usefull for debugging.
627     // Inputs:
628     //   TGeoTube *a
629     // Outputs:
630     //   none.
631     // Return:
632     //   none.
633
634     if(!GetDebug()) return;
635     cout << a->GetName() <<": Rmin="<<a->GetRmin()
636          <<" Rmax=" <<a->GetRmax()<<" Dz="<<a->GetDz()<<endl;
637     return;
638 }
639 //----------------------------------------------------------------------
640 void AliITSv11Geometry::PrintTubeSeg(const TGeoTubeSeg *a)const{
641     // Prints out the content of the TGeoTubeSeg. Usefull for debugging.
642     // Inputs:
643     //   TGeoTubeSeg *a
644     // Outputs:
645     //   none.
646     // Return:
647     //   none.
648
649     if(!GetDebug()) return;
650     cout << a->GetName() <<": Phi1="<<a->GetPhi1()<<
651         " Phi2="<<a->GetPhi2()<<" Rmin="<<a->GetRmin()
652          <<" Rmax=" <<a->GetRmax()<<" Dz="<<a->GetDz()<<endl;
653     return;
654 }
655 //----------------------------------------------------------------------
656 void AliITSv11Geometry::PrintConeSeg(const TGeoConeSeg *a)const{
657     // Prints out the content of the TGeoConeSeg. Usefull for debugging.
658     // Inputs:
659     //   TGeoConeSeg *a
660     // Outputs:
661     //   none.
662     // Return:
663     //   none.
664
665     if(!GetDebug()) return;
666     cout << a->GetName() <<": Phi1="<<a->GetPhi1()<<
667         " Phi2="<<a->GetPhi2()<<" Rmin1="<<a->GetRmin1()
668          <<" Rmax1=" <<a->GetRmax1()<<" Rmin2="<<a->GetRmin2()
669          <<" Rmax2=" <<a->GetRmax2()<<" Dz="<<a->GetDz()<<endl;
670     return;
671 }
672 //----------------------------------------------------------------------
673 void AliITSv11Geometry::PrintBBox(const TGeoBBox *a)const{
674     // Prints out the content of the TGeoBBox. Usefull for debugging.
675     // Inputs:
676     //   TGeoBBox *a
677     // Outputs:
678     //   none.
679     // Return:
680     //   none.
681
682     if(!GetDebug()) return;
683     cout << a->GetName() <<": Dx="<<a->GetDX()<<
684         " Dy="<<a->GetDY()<<" Dz="<<a->GetDZ() <<endl;
685     return;
686 }
687 //---------------------------------------------------------------------
688 void AliITSv11Geometry::CreateDefaultMaterials(){
689     // Create ITS materials
690     // Defined media here should correspond to the one defined in galice.cuts
691     // File which is red in (AliMC*) fMCApp::Init() { ReadTransPar(); }
692     // Inputs:
693     //    none.
694     // Outputs:
695     //   none.
696     // Return:
697     //   none.
698     Int_t i;
699     Double_t w;
700
701     // Define some elements
702     TGeoElement *itsH  = new TGeoElement("ITS_H","Hydrogen",1,1.00794);
703     TGeoElement *itsHe = new TGeoElement("ITS_He","Helium",2,4.002602);
704     TGeoElement *itsC  = new TGeoElement("ITS_C","Carbon",6,12.0107);
705     TGeoElement *itsN  = new TGeoElement("ITS_N","Nitrogen",7,14.0067);
706     TGeoElement *itsO  = new TGeoElement("ITS_O","Oxygen",8,15.994);
707     TGeoElement *itsF  = new TGeoElement("ITS_F","Florine",9,18.9984032);
708     TGeoElement *itsNe = new TGeoElement("ITS_Ne","Neon",10,20.1797);
709     TGeoElement *itsMg = new TGeoElement("ITS_Mg","Magnesium",12,24.3050);
710     TGeoElement *itsAl = new TGeoElement("ITS_Al","Aluminum",13,26981538);
711     TGeoElement *itsSi = new TGeoElement("ITS_Si","Silicon",14,28.0855);
712     TGeoElement *itsP  = new TGeoElement("ITS_P" ,"Phosphorous",15,30.973761);
713     TGeoElement *itsS  = new TGeoElement("ITS_S" ,"Sulfur",16,32.065);
714     TGeoElement *itsAr = new TGeoElement("ITS_Ar","Argon",18,39.948);
715     TGeoElement *itsTi = new TGeoElement("ITS_Ti","Titanium",22,47.867);
716     TGeoElement *itsCr = new TGeoElement("ITS_Cr","Chromium",24,51.9961);
717     TGeoElement *itsMn = new TGeoElement("ITS_Mn","Manganese",25,54.938049);
718     TGeoElement *itsFe = new TGeoElement("ITS_Fe","Iron",26,55.845);
719     TGeoElement *itsCo = new TGeoElement("ITS_Co","Cobalt",27,58.933200);
720     TGeoElement *itsNi = new TGeoElement("ITS_Ni","Nickrl",28,56.6930);
721     TGeoElement *itsCu = new TGeoElement("ITS_Cu","Copper",29,63.546);
722     TGeoElement *itsZn = new TGeoElement("ITS_Zn","Zinc",30,65.39);
723     TGeoElement *itsKr = new TGeoElement("ITS_Kr","Krypton",36,83.80);
724     TGeoElement *itsMo = new TGeoElement("ITS_Mo","Molylibdium",42,95.94);
725     TGeoElement *itsXe = new TGeoElement("ITS_Xe","Zeon",54,131.293);
726
727     // Start with the Materials since for any one material there
728     // can be defined more than one Medium.
729     // Air, dry. at 15degree C, 101325Pa at sea-level, % by volume
730     // (% by weight). Density is 351 Kg/m^3
731     // N2 78.084% (75.47%), O2 20.9476% (23.20%), Ar 0.934 (1.28%)%,
732     // C02 0.0314% (0.0590%), Ne 0.001818% (0.0012%, CH4 0.002% (),
733     // He 0.000524% (0.00007%), Kr 0.000114% (0.0003%), H2 0.00005% (3.5E-6%), 
734     // Xe 0.0000087% (0.00004 %), H2O 0.0% (dry) + trace amounts at the ppm
735     // levels.
736     TGeoMixture *itsAir = new TGeoMixture("ITS_Air",9);
737     w = 75.47E-2;
738     itsAir->AddElement(itsN,w);// Nitorgen, atomic
739     w = 23.29E-2 + // O2
740          5.90E-4 * 2.*15.994/(12.0107+2.*15.994);// CO2.
741     itsAir->AddElement(itsO,w);// Oxygen, atomic
742     w = 1.28E-2;
743     itsAir->AddElement(itsAr,w);// Argon, atomic
744     w =  5.90E-4*12.0107/(12.0107+2.*15.994)+  // CO2
745          2.0E-5 *12.0107/(12.0107+4.* 1.00794);  // CH4
746     itsAir->AddElement(itsC,w);// Carbon, atomic
747     w = 1.818E-5;
748     itsAir->AddElement(itsNe,w);// Ne, atomic
749     w = 3.5E-8;
750     itsAir->AddElement(itsHe,w);// Helium, atomic
751     w = 7.0E-7;
752     itsAir->AddElement(itsKr,w);// Krypton, atomic
753     w = 3.0E-6;
754     itsAir->AddElement(itsH,w);// Hydrogen, atomic
755     w = 4.0E-7;
756     itsAir->AddElement(itsXe,w);// Xenon, atomic
757     itsAir->SetDensity(351.0*fgkKgm3); //
758     itsAir->SetPressure(101325*fgkPascal);
759     itsAir->SetTemperature(15.0*fgkCelsius);
760     itsAir->SetState(TGeoMaterial::kMatStateGas);
761     //
762     // Silicone
763     TGeoMaterial *itsSiDet = new TGeoMaterial("ITS_Si",itsSi,2.33*fgkgcm3);
764     itsSiDet->SetTemperature(15.0*fgkCelsius);
765     itsSiDet->SetState(TGeoMaterial::kMatStateSolid);
766     //
767     // Epoxy C18 H19 O3
768     TGeoMixture *itsEpoxy = new TGeoMixture("ITS_Epoxy",3);
769     itsEpoxy->AddElement(itsC,18);
770     itsEpoxy->AddElement(itsH,19);
771     itsEpoxy->AddElement(itsO,3);
772     itsEpoxy->SetDensity(1.8*fgkgcm3);
773     itsEpoxy->SetTemperature(15.0*fgkCelsius);
774     itsEpoxy->SetState(TGeoMaterial::kMatStateSolid);
775     //
776     // Carbon Fiber, M55J, 60% fiber by volume. Fiber density
777     // 1.91 g/cm^3. See ToryaCA M55J data sheet.
778     //Begin_Html
779     /*
780        <A HREF="http://torayusa.com/cfa/pdfs/M55JDataSheet.pdf"> Data Sheet
781        </A>
782     */
783     //End_Html
784     TGeoMixture *itsCarbonFiber = new TGeoMixture("ITS_CarbonFiber-M55J",4);
785     // Assume that the epoxy fill in the space between the fibers and so
786     // no change in the total volume. To compute w, assume 1cm^3 total 
787     // volume.
788     w = 1.91/(1.91+(1.-.60)*itsEpoxy->GetDensity());
789     itsCarbonFiber->AddElement(itsC,w);
790     w = (1.-.60)*itsEpoxy->GetDensity()/(1.91+(1.-.06)*itsEpoxy->GetDensity());
791     for(i=0;i<itsEpoxy->GetNelements();i++)
792         itsCarbonFiber->AddElement(itsEpoxy->GetElement(i),
793                                    itsEpoxy->GetWmixt()[i]*w);
794     itsCarbonFiber->SetDensity((1.91+(1.-.60)*itsEpoxy->GetDensity())*fgkgcm3);
795     itsCarbonFiber->SetTemperature(22.0*fgkCelsius);
796     itsCarbonFiber->SetState(TGeoMaterial::kMatStateSolid);
797     //
798     // 
799     //
800     // Rohacell 51A  millable foam product.
801     // C9 H13 N1 O2  52Kg/m^3
802     // Elemental composition, Private comunications with
803     // Bjorn S. Nilsen
804     //Begin_Html
805     /*
806       <A HREF="http://www.rohacell.com/en/performanceplastics8344.html">
807        Rohacell-A see Properties
808        </A>
809      */
810     //End_Html
811     TGeoMixture *itsFoam = new TGeoMixture("ITS_Foam",4);
812     itsFoam->AddElement(itsC,9);
813     itsFoam->AddElement(itsH,13);
814     itsFoam->AddElement(itsN,1);
815     itsFoam->AddElement(itsO,2);
816     itsFoam->SetTitle("Rohacell 51 A");
817     itsFoam->SetDensity(52.*fgkKgm3);
818     itsFoam->SetTemperature(22.0*fgkCelsius);
819     itsFoam->SetState(TGeoMaterial::kMatStateSolid);
820     //
821     // Kapton % by weight, H 2.6362, C69.1133, N 7.3270, O 20.0235
822     // Density 1.42 g/cm^3
823     //Begin_Html
824     /*
825         <A HREF="http://www2.dupont.com/Kapton/en_US/assets/downloads/pdf/summaryofprop.pdf">
826         Kapton. also see </A>
827         <A HREF="http://physics.nist.gov/cgi-bin/Star/compos.pl?matno=179">
828         </A>
829      */
830     //End_Html
831     TGeoMixture *itsKapton = new TGeoMixture("ITS_Kapton",4);
832     itsKapton->AddElement(itsH,0.026362);
833     itsKapton->AddElement(itsC,0.691133);
834     itsKapton->AddElement(itsN,0.073270);
835     itsKapton->AddElement(itsO,0.200235);
836     itsKapton->SetTitle("Kapton ribon and cable base");
837     itsKapton->SetDensity(1.42*fgkgcm3);
838     itsKapton->SetTemperature(22.0*fgkCelsius);
839     itsKapton->SetState(TGeoMaterial::kMatStateSolid);
840     //
841     // UPILEX-S C16 H6 O4 N2 polymer (a Kapton like material)
842     // Density 1.47 g/cm^3
843     //Begin_Html
844     /*
845         <A HREF="http://northamerica.ube.com/page.php?pageid=9">
846         UPILEX-S. also see </A>
847         <A HREF="http://northamerica.ube.com/page.php?pageid=81">
848         </A>
849      */
850     //End_Html
851     TGeoMixture *itsUpilex = new TGeoMixture("ITS_Upilex",4);
852     itsUpilex->AddElement(itsC,16);
853     itsUpilex->AddElement(itsH,6);
854     itsUpilex->AddElement(itsN,2);
855     itsUpilex->AddElement(itsO,4);
856     itsUpilex->SetTitle("Upilex ribon, cable, and pcb base");
857     itsUpilex->SetDensity(1.47*fgkgcm3);
858     itsUpilex->SetTemperature(22.0*fgkCelsius);
859     itsUpilex->SetState(TGeoMaterial::kMatStateSolid);
860     //
861     // Aluminum 6061 (Al used by US groups)
862     // % by weight, Cr 0.04-0.35 range [0.0375 nominal value used]
863     // Cu 0.15-0.4 [0.275], Fe Max 0.7 [0.35], Mg 0.8-1.2 [1.0],
864     // Mn Max 0.15 [0.075] Si 0.4-0.8 [0.6], Ti Max 0.15 [0.075],
865     // Zn Max 0.25 [0.125], Rest Al [97.4625]. Density 2.7 g/cm^3
866     //Begin_Html
867     /*
868       <A HREG="http://www.matweb.com/SpecificMaterial.asp?bassnum=MA6016&group=General">
869       Aluminum 6061 specifications
870       </A>
871      */
872     //End_Html
873     TGeoMixture *itsAl6061 = new TGeoMixture("ITS_Al6061",9);
874     itsAl6061->AddElement(itsCr,0.000375);
875     itsAl6061->AddElement(itsCu,0.00275);
876     itsAl6061->AddElement(itsFe,0.0035);
877     itsAl6061->AddElement(itsMg,0.01);
878     itsAl6061->AddElement(itsMn,0.00075);
879     itsAl6061->AddElement(itsSi,0.006);
880     itsAl6061->AddElement(itsTi,0.00075);
881     itsAl6061->AddElement(itsZn,0.00125);
882     itsAl6061->AddElement(itsAl,0.974625);
883     itsAl6061->SetTitle("Aluminum Alloy 6061");
884     itsAl6061->SetDensity(2.7*fgkgcm3);
885     itsAl6061->SetTemperature(22.0*fgkCelsius);
886     itsAl6061->SetState(TGeoMaterial::kMatStateSolid);
887     //
888     // Aluminum 7075  (Al used by Italian groups)
889     // % by weight, Cr 0.18-0.28 range [0.23 nominal value used]
890     // Cu 1.2-2.0 [1.6], Fe Max 0.5 [0.25], Mg 2.1-2.9 [2.5],
891     // Mn Max 0.3 [0.125] Si Max 0.4 [0.2], Ti Max 0.2 [0.1],
892     // Zn 5.1-6.1 [5.6], Rest Al [89.395]. Density 2.81 g/cm^3
893     //Begin_Html
894     /*
895       <A HREG="http://asm.matweb.com/search/SpecificMaterial.asp?bassnum=MA7075T6">
896       Aluminum 7075 specifications
897       </A>
898      */
899     //End_Html
900     TGeoMixture *itsAl7075 = new TGeoMixture("ITS_Al7075",9);
901     itsAl7075->AddElement(itsCr,0.0023);
902     itsAl7075->AddElement(itsCu,0.016);
903     itsAl7075->AddElement(itsFe,0.0025);
904     itsAl7075->AddElement(itsMg,0.025);
905     itsAl7075->AddElement(itsMn,0.00125);
906     itsAl7075->AddElement(itsSi,0.002);
907     itsAl7075->AddElement(itsTi,0.001);
908     itsAl7075->AddElement(itsZn,0.056);
909     itsAl7075->AddElement(itsAl,0.89395);
910     itsAl7075->SetTitle("Aluminum Alloy 7075");
911     itsAl7075->SetDensity(2.81*fgkgcm3);
912     itsAl7075->SetTemperature(22.0*fgkCelsius);
913     itsAl7075->SetState(TGeoMaterial::kMatStateSolid);
914     //
915     // "Ruby" spheres, Al2 O3
916     // "Ruby" Sphere posts, Ryton R-4 04
917     //Begin_Html
918     /*
919       <A HREF="">
920       Ruby Sphere Posts
921       </A>
922      */
923     //End_Html
924     TGeoMixture *itsRuby = new TGeoMixture("ITS_RubySphere",2);
925     itsRuby->AddElement(itsAl,2);
926     itsRuby->AddElement(itsO,3);
927     itsRuby->SetTitle("Ruby reference sphere");
928     itsRuby->SetDensity(2.81*fgkgcm3);
929     itsRuby->SetTemperature(22.0*fgkCelsius);
930     itsRuby->SetState(TGeoMaterial::kMatStateSolid);
931     //
932     //
933     // Inox, AISI 304L, compoistion % by weight (assumed) 
934     // C Max 0.03 [0.015], Mn Max 2.00 [1.00], Si Max 1.00 [0.50]
935     // P Max 0.045 [0.0225], S Max 0.03 [0.015], Ni 8.0-10.5 [9.25]
936     // Cr 18-20 [19.], Mo 2.-2.5 [2.25], rest Fe: density 7.93 Kg/dm^3
937     //Begin_Html
938     /*
939       <A HREF="http://www.cimap.fr/caracter.pdf">
940       Stainless steal (INOX) AISI 304L composition
941       </A>
942      */
943     //End_Html
944     TGeoMixture *itsInox304L = new TGeoMixture("ITS_Inox304L",9);
945     itsInox304L->AddElement(itsC,0.00015);
946     itsInox304L->AddElement(itsMn,0.010);
947     itsInox304L->AddElement(itsSi,0.005);
948     itsInox304L->AddElement(itsP,0.000225);
949     itsInox304L->AddElement(itsS,0.00015);
950     itsInox304L->AddElement(itsNi,0.0925);
951     itsInox304L->AddElement(itsCr,0.1900);
952     itsInox304L->AddElement(itsMo,0.0225);
953     itsInox304L->AddElement(itsFe,0.679475); // Rest Fe
954     itsInox304L->SetTitle("ITS Stainless Steal (Inox) type AISI 304L");
955     itsInox304L->SetDensity(7.93*fgkKgdm3);
956     itsInox304L->SetTemperature(22.0*fgkCelsius);
957     itsInox304L->SetState(TGeoMaterial::kMatStateSolid);
958     //
959     // Inox, AISI 316L, composition % by weight (assumed)
960     // C Max 0.03 [0.015], Mn Max 2.00 [1.00], Si Max 1.00 [0.50]
961     // P Max 0.045 [0.0225], S Max 0.03 [0.015], Ni 10.0-14. [12.]
962     // Cr 16-18 [17.], Mo 2-3 [2.5]: density 7.97 Kg/dm^3
963     //Begin_Html
964     /*
965       <A HREF="http://www.cimap.fr/caracter.pdf">
966       Stainless steal (INOX) AISI 316L composition
967       </A>
968      */
969     //End_Html
970     TGeoMixture *itsInox316L = new TGeoMixture("ITS_Inox316L",9);
971     itsInox316L->AddElement(itsC,0.00015);
972     itsInox316L->AddElement(itsMn,0.010);
973     itsInox316L->AddElement(itsSi,0.005);
974     itsInox316L->AddElement(itsP,0.000225);
975     itsInox316L->AddElement(itsS,0.00015);
976     itsInox316L->AddElement(itsNi,0.12);
977     itsInox316L->AddElement(itsCr,0.17);
978     itsInox316L->AddElement(itsMo,0.025);
979     itsInox316L->AddElement(itsFe,0.66945); // Rest Fe
980     itsInox316L->SetTitle("ITS Stainless Steal (Inox) type AISI 316L");
981     itsInox316L->SetDensity(7.97*fgkKgdm3);
982     itsInox316L->SetTemperature(22.0*fgkCelsius);
983     itsInox316L->SetState(TGeoMaterial::kMatStateSolid);
984     //
985     // Inox, Phynox or Elgiloy AMS 5833, composition % by weight
986     // C Max 0.15 [0.15], Mn Max 2.00 [2.00], Be max 0.0001 [none]
987     // Ni 18. [18.], Cr 21.5 [21.5], Mo 7.5 [7.5], Co 42 [42.]: 
988     // density 8.3 Kg/dm^3
989     //Begin_Html
990     /*
991       <A HREF="http://www.freepatentsonline.com/20070032816.html">
992       Compostion of Phynox or Elgiloy AMS 5833, also see
993       </A>
994       <A HREF="http://www.alloywire.com/phynox_alloy.html">
995       under corss reference number [0024].
996       </A>
997      */
998     //End_Html
999     TGeoMixture *itsPhynox = new TGeoMixture("ITS_Phynox",7);
1000     itsPhynox->AddElement(itsC,0.0015);
1001     itsPhynox->AddElement(itsMn,0.020);
1002     itsPhynox->AddElement(itsNi,0.18);
1003     itsPhynox->AddElement(itsCr,0.215);
1004     itsPhynox->AddElement(itsMo,0.075);
1005     itsPhynox->AddElement(itsCo,0.42);
1006     itsPhynox->AddElement(itsFe,0.885);
1007     itsPhynox->SetTitle("ITS Cooling tube alloy");
1008     itsPhynox->SetDensity(8.3*fgkgcm3);
1009     itsPhynox->SetTemperature(22.0*fgkCelsius);
1010     itsPhynox->SetState(TGeoMaterial::kMatStateSolid);
1011     //
1012     // G10FR4
1013     //
1014     // Demineralized Water H2O SDD & SSD Cooling liquid
1015     TGeoMixture *itsWater = new TGeoMixture("ITS_Water",2);
1016     itsWater->AddElement(itsH,2);
1017     itsWater->AddElement(itsO,1);
1018     itsWater->SetTitle("ITS Cooling Water");
1019     itsWater->SetDensity(1.0*fgkgcm3);
1020     itsWater->SetTemperature(22.0*fgkCelsius);
1021     itsWater->SetState(TGeoMaterial::kMatStateLiquid);
1022     //
1023     // Freon SPD Cooling liquid PerFluorobuthane C4F10
1024     //Begin_Html
1025     /*
1026       <A HREF=" http://st-support-cooling-electronics.web.cern.ch/st-support-cooling-electronics/default.htm">
1027       SPD 2 phase cooling using PerFluorobuthane
1028       </A>
1029      */
1030     //End_Html
1031     TGeoMixture *itsFreon = new TGeoMixture("ITS_SPD_Freon",2);
1032     itsFreon->AddElement(itsC,4);
1033     itsFreon->AddElement(itsF,10);
1034     itsFreon->SetTitle("ITS SPD 2 phase Cooling freon");
1035     itsFreon->SetDensity(1.52*fgkgcm3);
1036     itsFreon->SetTemperature(22.0*fgkCelsius);
1037     itsFreon->SetState(TGeoMaterial::kMatStateLiquid);
1038     //
1039     //    Int_t   ifield = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ();
1040     //    Float_t fieldm = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();
1041
1042     //    Float_t tmaxfd = 0.1;//  1.0;//  Degree
1043     //    Float_t stemax = 1.0;//  cm
1044     //   Float_t deemax = 0.1;// 30.0;// Fraction of particle's energy 0<deemax<=1
1045     //    Float_t epsil  = 1.0E-4;//  1.0;  cm
1046     //    Float_t stmin  = 0.0; // cm "Default value used"
1047
1048     //    Float_t tmaxfdSi = 0.1; // .10000E+01; // Degree
1049     //   Float_t stemaxSi = 0.0075; //  .10000E+01; // cm
1050     //   Float_t deemaxSi = 0.1; // Fraction of particle's energy 0<deemax<=1
1051     //    Float_t epsilSi  = 1.0E-4;// .10000E+01;
1052     /*
1053     Float_t stminSi  = 0.0; // cm "Default value used"
1054
1055     Float_t tmaxfdAir = 0.1; // .10000E+01; // Degree
1056     Float_t stemaxAir = .10000E+01; // cm
1057     Float_t deemaxAir = 0.1; // 0.30000E-02; // Fraction of particle's energy 0<deemax<=1
1058     Float_t epsilAir  = 1.0E-4;// .10000E+01;
1059     Float_t stminAir  = 0.0; // cm "Default value used"
1060
1061     Float_t tmaxfdServ = 1.0; // 10.0; // Degree
1062     Float_t stemaxServ = 1.0; // 0.01; // cm
1063     Float_t deemaxServ = 0.5; // 0.1; // Fraction of particle's energy 0<deemax<=1
1064     Float_t epsilServ  = 1.0E-3; // 0.003; // cm
1065     Float_t stminServ  = 0.0; //0.003; // cm "Default value used"
1066
1067     // Freon PerFluorobuthane C4F10 see 
1068     // http://st-support-cooling-electronics.web.cern.ch/
1069     //        st-support-cooling-electronics/default.htm
1070     Float_t afre[2]  = { 12.011,18.9984032 };
1071     Float_t zfre[2]  = { 6., 9. };
1072     Float_t wfre[2]  = { 4.,10. };
1073     Float_t densfre  = 1.52;
1074
1075     //CM55J
1076     Float_t aCM55J[4]={12.0107,14.0067,15.9994,1.00794};
1077     Float_t zCM55J[4]={6.,7.,8.,1.};
1078     Float_t wCM55J[4]={0.908508078,0.010387573,0.055957585,0.025146765};
1079     Float_t dCM55J = 1.63;
1080
1081     //ALCM55J
1082     Float_t aALCM55J[5]={12.0107,14.0067,15.9994,1.00794,26.981538};
1083     Float_t zALCM55J[5]={6.,7.,8.,1.,13.};
1084     Float_t wALCM55J[5]={0.817657902,0.0093488157,0.0503618265,0.0226320885,0.1};
1085     Float_t dALCM55J = 1.9866;
1086
1087     //Si Chips
1088     Float_t aSICHIP[6]={12.0107,14.0067,15.9994,1.00794,28.0855,107.8682};
1089     Float_t zSICHIP[6]={6.,7.,8.,1.,14., 47.};
1090     Float_t wSICHIP[6]={0.039730642,0.001396798,0.01169634,
1091                         0.004367771,0.844665,0.09814344903};
1092     Float_t dSICHIP = 2.36436;
1093
1094     //Inox
1095     Float_t aINOX[9]={12.0107,54.9380, 28.0855,30.9738,32.066,
1096                       58.6928,55.9961,95.94,55.845};
1097     Float_t zINOX[9]={6.,25.,14.,15.,16., 28.,24.,42.,26.};
1098     Float_t wINOX[9]={0.0003,0.02,0.01,0.00045,0.0003,0.12,0.17,0.025,0.654};
1099     Float_t dINOX = 8.03;
1100
1101     //SDD HV microcable
1102     Float_t aHVm[5]={12.0107,1.00794,14.0067,15.9994,26.981538};
1103     Float_t zHVm[5]={6.,1.,7.,8.,13.};
1104     Float_t wHVm[5]={0.520088819984,0.01983871336,0.0551367996,0.157399667056, 0.247536};
1105     Float_t dHVm = 1.6087;
1106
1107     //SDD LV+signal cable
1108     Float_t aLVm[5]={12.0107,1.00794,14.0067,15.9994,26.981538};
1109     Float_t zLVm[5]={6.,1.,7.,8.,13.};
1110     Float_t wLVm[5]={0.21722436468,0.0082859922,0.023028867,0.06574077612, 0.68572};
1111     Float_t dLVm = 2.1035;
1112
1113     //SDD hybrid microcab
1114     Float_t aHLVm[5]={12.0107,1.00794,14.0067,15.9994,26.981538};
1115     Float_t zHLVm[5]={6.,1.,7.,8.,13.};
1116     Float_t wHLVm[5]={0.24281879711,0.00926228815,0.02574224025,0.07348667449, 0.64869};
1117     Float_t dHLVm = 2.0502;
1118
1119     //SDD anode microcab
1120     Float_t aALVm[5]={12.0107,1.00794,14.0067,15.9994,26.981538};
1121     Float_t zALVm[5]={6.,1.,7.,8.,13.};
1122     Float_t wALVm[5]={0.392653705471,0.0128595919215,
1123                       0.041626868025,0.118832707289, 0.431909};
1124     Float_t dALVm = 2.0502;
1125
1126     //X7R capacitors
1127     Float_t aX7R[7]={137.327,47.867,15.9994,58.6928,63.5460,118.710,207.2};
1128     Float_t zX7R[7]={56.,22.,8.,28.,29.,50.,82.};
1129     Float_t wX7R[7]={0.251639432,0.084755042,0.085975822,
1130                      0.038244751,0.009471271,0.321736471,0.2081768};
1131     Float_t dX7R = 7.14567;
1132
1133     // AIR
1134     Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
1135     Float_t zAir[4]={6.,7.,8.,18.};
1136     Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
1137     Float_t dAir = 1.20479E-3;
1138
1139     // Water
1140     Float_t aWater[2]={1.00794,15.9994};
1141     Float_t zWater[2]={1.,8.};
1142     Float_t wWater[2]={0.111894,0.888106};
1143     Float_t dWater   = 1.0;
1144
1145     // CERAMICS
1146     //     94.4% Al2O3 , 2.8% SiO2 , 2.3% MnO , 0.5% Cr2O3
1147     Float_t acer[5]  = { 26.981539,15.9994,28.0855,54.93805,51.9961 };
1148     Float_t zcer[5]  = {       13.,     8.,    14.,     25.,    24. };
1149     Float_t wcer[5]  = {.4443408,.5213375,.0130872,.0178135,.003421};
1150     Float_t denscer  = 3.6;
1151
1152     // G10FR4
1153     Float_t zG10FR4[14] = {14.00,       20.00,  13.00,  12.00,  5.00,
1154                            22.00,       11.00,  19.00,  26.00,  9.00,
1155                            8.00,        6.00,   7.00,   1.00};
1156     Float_t aG10FR4[14] = {28.0855000,40.0780000,26.9815380,24.3050000,
1157                            10.8110000,47.8670000,22.9897700,39.0983000,
1158                            55.8450000,18.9984000,15.9994000,12.0107000,
1159                            14.0067000,1.0079400};
1160     Float_t wG10FR4[14] = {0.15144894,0.08147477,0.04128158,0.00904554,
1161                            0.01397570,0.00287685,0.00445114,0.00498089,
1162                            0.00209828,0.00420000,0.36043788,0.27529426,
1163                            0.01415852,0.03427566};
1164     Float_t densG10FR4= 1.8;
1165
1166     //--- EPOXY  --- C18 H19 O3
1167     Float_t aEpoxy[3] = {15.9994, 1.00794, 12.0107} ; 
1168     Float_t zEpoxy[3] = {     8.,      1.,      6.} ; 
1169     Float_t wEpoxy[3] = {     3.,     19.,     18.} ; 
1170     Float_t dEpoxy = 1.8 ;
1171
1172     // rohacell: C9 H13 N1 O2
1173     Float_t arohac[4] = {12.01,  1.01, 14.010, 16.};
1174     Float_t zrohac[4] = { 6.,    1.,    7.,     8.};
1175     Float_t wrohac[4] = { 9.,   13.,    1.,     2.};
1176     Float_t drohac    = 0.05;
1177
1178     // If he/she means stainless steel (inox) + Aluminium and Zeff=15.3383 then
1179     // %Al=81.6164 %inox=100-%Al
1180     Float_t aInAl[5] = {27., 55.847,51.9961,58.6934,28.0855 };
1181     Float_t zInAl[5] = {13., 26.,24.,28.,14. };
1182     Float_t wInAl[5] = {.816164, .131443,.0330906,.0183836,.000919182};
1183     Float_t dInAl    = 3.075;
1184
1185     // Kapton
1186     Float_t aKapton[4]={1.00794,12.0107, 14.010,15.9994};
1187     Float_t zKapton[4]={1.,6.,7.,8.};
1188     Float_t wKapton[4]={0.026362,0.69113,0.07327,0.209235};
1189     Float_t dKapton   = 1.42;
1190
1191     //SDD ruby sph.
1192     Float_t aAlOxide[2]  = { 26.981539,15.9994};
1193     Float_t zAlOxide[2]  = {       13.,     8.};
1194     Float_t wAlOxide[2]  = {0.4707, 0.5293};
1195     Float_t dAlOxide     = 3.97;
1196     */
1197 }
1198 //---------------------------------------------------------------------
1199 void AliITSv11Geometry::DrawCrossSection(const TGeoPcon *p,
1200                             Int_t fillc,Int_t fills,
1201                             Int_t linec,Int_t lines,Int_t linew,
1202                             Int_t markc,Int_t marks,Float_t marksize)const{
1203     // Draws a cross sectional view of the TGeoPcon, Primarily for debugging.
1204     // A TCanvas should exist first.
1205     //  Inputs:
1206     //    TGeoPcon  *p  The TGeoPcon to be "drawn"
1207     //    Int_t  fillc  The fill color to be used
1208     //    Int_t  fills  The fill style to be used
1209     //    Int_t  linec  The line color to be used
1210     //    Int_t  lines  The line style to be used
1211     //    Int_t  linew  The line width to be used
1212     //    Int_t  markc  The markder color to be used
1213     //    Int_t  marks  The markder style to be used
1214     //    Float_t marksize The marker size
1215     // Outputs:
1216     //   none.
1217     // Return:
1218     //   none.
1219     Int_t n=0,m=0,i=0;
1220     Double_t *z=0,*r=0;
1221     TPolyMarker *pts=0;
1222     TPolyLine   *line=0;
1223
1224     n = p->GetNz();
1225     if(n<=0) return;
1226     m = 2*n+1;
1227     z = new Double_t[m];
1228     r = new Double_t[m];
1229
1230     for(i=0;i<n;i++){
1231         z[i] = p->GetZ(i);
1232         r[i] = p->GetRmax(i);
1233         z[i+n] = p->GetZ(n-1-i);
1234         r[i+n] = p->GetRmin(n-1-i);
1235     } //  end for i
1236     z[n-1] = z[0];
1237     r[n-1] = r[0];
1238
1239     line = new TPolyLine(n,z,r);
1240     pts  = new TPolyMarker(n,z,r);
1241
1242     line->SetFillColor(fillc);
1243     line->SetFillStyle(fills);
1244     line->SetLineColor(linec);
1245     line->SetLineStyle(lines);
1246     line->SetLineWidth(linew);
1247     pts->SetMarkerColor(markc);
1248     pts->SetMarkerStyle(marks);
1249     pts->SetMarkerSize(marksize);
1250
1251     line->Draw("f");
1252     line->Draw();
1253     pts->Draw();
1254
1255     delete[] z;
1256     delete[] r;
1257
1258     cout<<"Hit Return to continue"<<endl;
1259     cin >> n;
1260     delete line;
1261     delete pts;
1262     return;
1263 }
1264 //______________________________________________________________________
1265 Bool_t AliITSv11Geometry::AngleOfIntersectionWithLine(Double_t x0,Double_t y0,
1266                                                       Double_t x1,Double_t y1,
1267                                                       Double_t xc,Double_t yc,
1268                                                       Double_t rc,Double_t &t0,
1269                                                       Double_t &t1)const{
1270     // Computes the angles, t0 and t1 corresponding to the intersection of
1271     // the line, defined by {x0,y0} {x1,y1}, and the circle, defined by
1272     // its center {xc,yc} and radius r. If the line does not intersect the
1273     // line, function returns kFALSE, otherwise it returns kTRUE. If the
1274     // line is tangent to the circle, the angles t0 and t1 will be the same.
1275     // Inputs:
1276     //   Double_t x0   X of first point defining the line
1277     //   Double_t y0   Y of first point defining the line
1278     //   Double_t x1   X of Second point defining the line
1279     //   Double_t y1   Y of Second point defining the line
1280     //   Double_t xc   X of Circle center point defining the line
1281     //   Double_t yc   Y of Circle center point defining the line
1282     //   Double_t r    radius of circle
1283     // Outputs:
1284     //   Double_t &t0  First angle where line intersects circle
1285     //   Double_t &t1  Second angle where line intersects circle
1286     // Return:
1287     //    kTRUE, line intersects circle, kFALSE line does not intersect circle
1288     //           or the line is not properly defined point {x0,y0} and {x1,y1}
1289     //           are the same point.
1290     Double_t dx,dy,cx,cy,s2,t[4];
1291     Double_t a0,b0,c0,a1,b1,c1,sinthp,sinthm,costhp,costhm;
1292     Int_t i,j;
1293
1294     t0 = 400.0;
1295     t1 = 400.0;
1296     dx = x1-x0;
1297     dy = y1-y0;
1298     cx = xc-x0;
1299     cy = yc-y0;
1300     s2 = dx*dx+dy*dy;
1301     if(s2==0.0) return kFALSE;
1302
1303     a0 = rc*rc*s2;
1304     if(a0==0.0) return kFALSE;
1305     b0 = 2.0*rc*dx*(dx*cy-cx*dy);
1306     c0 = dx*dx*cy*cy-2.0*dy*dx*cy*cx+cx*cx*dy*dy-rc*rc*dy*dy;
1307     c0 = 0.25*b0*b0/(a0*a0)-c0/a0;
1308     if(c0<0.0) return kFALSE;
1309     sinthp = -0.5*b0/a0+TMath::Sqrt(c0);
1310     sinthm = -0.5*b0/a0-TMath::Sqrt(c0);
1311
1312     a1 = rc*rc*s2;
1313     if(a1==0.0) return kFALSE;
1314     b1 = 2.0*rc*dy*(dy*cx-dx*cy);
1315     c1 = dy*dy*cx*cx-2.0*dy*dx*cy*cx+dx*dx*cy*cy-rc*rc*dx*dx;
1316     c1 = 0.25*b1*b1/(a1*a1)-c1/a1;
1317     if(c1<0.0) return kFALSE;
1318     costhp = -0.5*b1/a1+TMath::Sqrt(c1);
1319     costhm = -0.5*b1/a1-TMath::Sqrt(c1);
1320
1321     t[0] = t[1] = t[2] = t[3] = 400.;
1322     a0 = TMath::ATan2(sinthp,costhp); if(a0<0.0) a0 += 2.0*TMath::Pi();
1323     a1 = TMath::ATan2(sinthp,costhm); if(a1<0.0) a1 += 2.0*TMath::Pi();
1324     b0 = TMath::ATan2(sinthm,costhp); if(b0<0.0) b0 += 2.0*TMath::Pi();
1325     b1 = TMath::ATan2(sinthm,costhm); if(b1<0.0) b1 += 2.0*TMath::Pi();
1326     x1 = xc+rc*TMath::Cos(a0);
1327     y1 = yc+rc*TMath::Sin(a0);
1328     s2 = dx*(y1-y0)-dy*(x1-x0);
1329     if(s2*s2<DBL_EPSILON) t[0] = a0*TMath::RadToDeg();
1330     x1 = xc+rc*TMath::Cos(a1);
1331     y1 = yc+rc*TMath::Sin(a1);
1332     s2 = dx*(y1-y0)-dy*(x1-x0);
1333     if(s2*s2<DBL_EPSILON) t[1] = a1*TMath::RadToDeg();
1334     x1 = xc+rc*TMath::Cos(b0);
1335     y1 = yc+rc*TMath::Sin(b0);
1336     s2 = dx*(y1-y0)-dy*(x1-x0);
1337     if(s2*s2<DBL_EPSILON) t[2] = b0*TMath::RadToDeg();
1338     x1 = xc+rc*TMath::Cos(b1);
1339     y1 = yc+rc*TMath::Sin(b1);
1340     s2 = dx*(y1-y0)-dy*(x1-x0);
1341     if(s2*s2<DBL_EPSILON) t[3] = b1*TMath::RadToDeg();
1342     for(i=0;i<4;i++)for(j=i+1;j<4;j++){
1343         if(t[i]>t[j]) {t0 = t[i];t[i] = t[j];t[j] = t0;}
1344     } // end for i,j
1345     t0 = t[0];
1346     t1 = t[1];
1347     //
1348     return kTRUE;
1349 }
1350 //______________________________________________________________________
1351 Double_t AliITSv11Geometry::AngleForRoundedCorners0(Double_t dx,Double_t dy,
1352                                                     Double_t sdr)const{
1353     // Basic function used to determine the ending angle and starting angles
1354     // for rounded corners given the relative distance between the centers
1355     // of the circles and the difference/sum of their radii. Case 0.
1356     // Inputs:
1357     //   Double_t dx    difference in x locations of the circle centers
1358     //   Double_t dy    difference in y locations of the circle centers
1359     //   Double_t sdr   difference or sum of the circle radii
1360     // Outputs:
1361     //   none.
1362     // Return:
1363     //   the angle in Degrees
1364     Double_t a,b;
1365
1366     b = dy*dy+dx*dx-sdr*sdr;
1367     if(b<0.0) Error("AngleForRoundedCorners0",
1368                     "dx^2(%e)+dy^2(%e)-sdr^2(%e)=b=%e<0",dx,dy,sdr,b);
1369     b = TMath::Sqrt(b);
1370     a = -sdr*dy+dx*b;
1371     b = -sdr*dx-dy*b;
1372     return TMath::ATan2(a,b)*TMath::RadToDeg();
1373     
1374 }
1375 //______________________________________________________________________
1376 Double_t AliITSv11Geometry::AngleForRoundedCorners1(Double_t dx,Double_t dy,
1377                                                     Double_t sdr)const{
1378     // Basic function used to determine the ending angle and starting angles
1379     // for rounded corners given the relative distance between the centers
1380     // of the circles and the difference/sum of their radii. Case 1.
1381     // Inputs:
1382     //   Double_t dx    difference in x locations of the circle centers
1383     //   Double_t dy    difference in y locations of the circle centers
1384     //   Double_t sdr   difference or sum of the circle radii
1385     // Outputs:
1386     //   none.
1387     // Return:
1388     //   the angle in Degrees
1389     Double_t a,b;
1390
1391     b = dy*dy+dx*dx-sdr*sdr;
1392     if(b<0.0) Error("AngleForRoundedCorners1",
1393                     "dx^2(%e)+dy^2(%e)-sdr^2(%e)=b=%e<0",dx,dy,sdr,b);
1394     b = TMath::Sqrt(b);
1395     a = -sdr*dy-dx*b;
1396     b = -sdr*dx+dy*b;
1397     return TMath::ATan2(a,b)*TMath::RadToDeg();
1398     
1399 }
1400 //----------------------------------------------------------------------
1401 void AliITSv11Geometry::AnglesForRoundedCorners(Double_t x0,Double_t y0,
1402                                                 Double_t r0,Double_t x1,
1403                                                 Double_t y1,Double_t r1,
1404                                                 Double_t &t0,Double_t &t1)
1405     const{
1406     // Function to compute the ending angle, for arc 0, and starting angle,
1407     // for arc 1, such that a straight line will connect them with no
1408     // discontinuities.
1409     //Begin_Html
1410     /*
1411       <img src="picts/ITS/AliITSv11Geometry_AnglesForRoundedCorners.gif">
1412      */
1413     //End_Html
1414     // Inputs:
1415     //    Double_t x0  X Coordinate of arc 0 center.
1416     //    Double_t y0  Y Coordinate of arc 0 center.
1417     //    Double_t r0  Radius of curvature of arc 0. For signe see figure.
1418     //    Double_t x1  X Coordinate of arc 1 center.
1419     //    Double_t y1  Y Coordinate of arc 1 center.
1420     //    Double_t r1  Radius of curvature of arc 1. For signe see figure.
1421     // Outputs:
1422     //    Double_t t0  Ending angle of arch 0, with respect to x axis, Degrees.
1423     //    Double_t t1  Starting angle of arch 1, with respect to x axis, 
1424     //                 Degrees.
1425     // Return:
1426     //    none.
1427     Double_t t;
1428
1429     if(r0>=0.0&&r1>=0.0) { // Inside to inside    ++
1430         t = AngleForRoundedCorners1(x1-x0,y1-y0,r1-r0);
1431         t0 = t1 = t;
1432         return;
1433     }else if(r0>=0.0&&r1<=0.0){ // Inside to Outside  +-
1434         r1 = -r1; // make positive
1435         t = AngleForRoundedCorners0(x1-x0,y1-y0,r1+r0);
1436         t0 = 180.0 + t;
1437         if(t0<0.0) t += 360.;
1438         if(t<0.0) t += 360.;
1439         t1 = t;
1440         return;
1441     }else if(r0<=0.0&&r1>=0.0){ // Outside to Inside  -+
1442         r0 = - r0; // make positive
1443         t = AngleForRoundedCorners1(x1-x0,y1-y0,r1+r0);
1444         t0 = 180.0 + t;
1445         if(t0>180.) t0 -= 360.;
1446         if(t >180.) t  -= 360.;
1447         t1 = t;
1448         return;
1449     }else if(r0<=0.0&&r1<=0.0) { // Outside to outside --
1450         r0 = -r0; // make positive
1451         r1 = -r1; // make positive
1452         t = AngleForRoundedCorners0(x1-x0,y1-y0,r1-r0);
1453         t0 = t1 = t;
1454         return;
1455     } // end if
1456     return;
1457 }
1458 //----------------------------------------------------------------------
1459 void AliITSv11Geometry::MakeFigure1(Double_t x0,Double_t y0,Double_t r0,
1460                                     Double_t x1,Double_t y1,Double_t r1){
1461     // Function to create the figure discribing how the function 
1462     // AnglesForRoundedCorners works.
1463     //
1464     // Inputs:
1465     //    Double_t x0  X Coordinate of arc 0 center.
1466     //    Double_t y0  Y Coordinate of arc 0 center.
1467     //    Double_t r0  Radius of curvature of arc 0. For signe see figure.
1468     //    Double_t x1  X Coordinate of arc 1 center.
1469     //    Double_t y1  Y Coordinate of arc 1 center.
1470     //    Double_t r1  Radius of curvature of arc 1. For signe see figure.
1471     // Outputs:
1472     //    none.
1473     // Return:
1474     //    none.
1475     Double_t t0[4],t1[4],xa0[4],ya0[4],xa1[4],ya1[4],ra0[4],ra1[4];
1476     Double_t xmin,ymin,xmax,ymax,h;
1477     Int_t j;
1478
1479     for(j=0;j<4;j++) {
1480         ra0[j] = r0; if(j%2) ra0[j] = -r0;
1481         ra1[j] = r1; if(j>1) ra1[j] = -r1;
1482         AnglesForRoundedCorners(x0,y0,ra0[j],x1,y1,ra1[j],t0[j],t1[j]);
1483         xa0[j] = TMath::Abs(r0)*CosD(t0[j])+x0;
1484         ya0[j] = TMath::Abs(r0)*SinD(t0[j])+y0;
1485         xa1[j] = TMath::Abs(r1)*CosD(t1[j])+x1;
1486         ya1[j] = TMath::Abs(r1)*SinD(t1[j])+y1;
1487     } // end for j
1488     if(r0<0.0) r0 = -r0;
1489     if(r1<0.0) r1 = -r1;
1490     xmin = TMath::Min(x0 - r0,x1-r1);
1491     ymin = TMath::Min(y0 - r0,y1-r1);
1492     xmax = TMath::Max(x0 + r0,x1+r1);
1493     ymax = TMath::Max(y0 + r0,y1+r1);
1494     for(j=1;j<4;j++) {
1495         xmin = TMath::Min(xmin,xa0[j]);
1496         xmin = TMath::Min(xmin,xa1[j]);
1497         ymin = TMath::Min(ymin,ya0[j]);
1498         ymin = TMath::Min(ymin,ya1[j]);
1499
1500         xmax = TMath::Max(xmax,xa0[j]);
1501         xmax = TMath::Max(xmax,xa1[j]);
1502         ymax = TMath::Max(ymax,ya0[j]);
1503         ymax = TMath::Max(ymax,ya1[j]);
1504     } // end for j
1505     if(xmin<0.0) xmin *= 1.1; else xmin *= 0.9;
1506     if(ymin<0.0) ymin *= 1.1; else ymin *= 0.9;
1507     if(xmax<0.0) xmax *= 0.9; else xmax *= 1.1;
1508     if(ymax<0.0) ymax *= 0.9; else ymax *= 1.1;
1509     j = (Int_t)(500.0*(ymax-ymin)/(xmax-xmin));
1510     TCanvas *can = new TCanvas("AliITSv11Geometry_AnglesForRoundedCorners",
1511                                "Figure for AliITSv11Geometry",500,j);
1512     h = ymax-ymin; if(h<0) h = -h;
1513     can->Range(xmin,ymin,xmax,ymax);
1514     TArc *c0 = new TArc(x0,y0,r0);
1515     TArc *c1 = new TArc(x1,y1,r1);
1516     TLine *line[4];
1517     TArrow *ar0[4];
1518     TArrow *ar1[4];
1519     for(j=0;j<4;j++){
1520         ar0[j] = new TArrow(x0,y0,xa0[j],ya0[j]);
1521         ar1[j] = new TArrow(x1,y1,xa1[j],ya1[j]);
1522         line[j] = new TLine(xa0[j],ya0[j],xa1[j],ya1[j]);
1523         ar0[j]->SetLineColor(j+1);
1524         ar0[j]->SetArrowSize(0.1*r0/h);
1525         ar1[j]->SetLineColor(j+1);
1526         ar1[j]->SetArrowSize(0.1*r1/h);
1527         line[j]->SetLineColor(j+1);
1528     } // end for j
1529     c0->Draw();
1530     c1->Draw();
1531     for(j=0;j<4;j++){
1532         ar0[j]->Draw();
1533         ar1[j]->Draw();
1534         line[j]->Draw();
1535     } // end for j
1536     TText *t = new TText();
1537     t->SetTextSize(0.02);
1538     Char_t txt[100];
1539     snprintf(txt,99,"(x0=%5.2f,y0=%5.2f)",x0,y0);
1540     t->DrawText(x0,y0,txt);
1541     snprintf(txt,99,"(x1=%5.2f,y1=%5.2f)",x1,y1);
1542     for(j=0;j<4;j++) {
1543         t->SetTextColor(j+1);
1544         t->DrawText(x1,y1,txt);
1545         snprintf(txt,99,"r0=%5.2f",ra0[j]);
1546         t->DrawText(0.5*(x0+xa0[j]),0.5*(y0+ya0[j]),txt);
1547         snprintf(txt,99,"r1=%5.2f",ra1[j]);
1548         t->DrawText(0.5*(x1+xa1[j]),0.5*(y1+ya1[j]),txt);
1549     } // end for j
1550 }