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