]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSv11Geometry.cxx
DIPO added
[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 <TPolyMarker.h>
44 #include <TPolyLine.h>
45 #include "AliITSv11Geometry.h"
46
47 ClassImp(AliITSv11Geometry)
48
49 const Double_t AliITSv11Geometry::fgkmicron = 1.0E-4;
50 const Double_t AliITSv11Geometry::fgkmm = 0.10;
51 const Double_t AliITSv11Geometry::fgkcm = 1.00;
52 const Double_t AliITSv11Geometry::fgkDegree = 1.0;
53 const Double_t AliITSv11Geometry::fgkRadian = 180./3.14159265358979323846;
54 const Double_t AliITSv11Geometry::fgkgcm3 = 1.0; // assume default is g/cm^3
55 const Double_t AliITSv11Geometry::fgkCelsius = 1.0; // Assume default is C
56 const Double_t AliITSv11Geometry::fgkPascal  = 1.0E-3; // Assume kPascal
57 const Double_t AliITSv11Geometry::fgkKPascal = 1.0;    // Asume kPascal
58 const Double_t AliITSv11Geometry::fgkeV      = 1.0E-9; // GeV default
59 const Double_t AliITSv11Geometry::fgkKeV     = 1.0e-6; // GeV default
60 const Double_t AliITSv11Geometry::fgkMeV     = 1.0e-3; // GeV default
61 const Double_t AliITSv11Geometry::fgkGeV     = 1.0;    // GeV default
62 //______________________________________________________________________
63 Double_t AliITSv11Geometry::Yfrom2Points(Double_t x0,Double_t y0,
64                                          Double_t x1,Double_t y1,
65                                          Double_t x)const{
66     // Given the two points (x0,y0) and (x1,y1) and the location x, returns
67     // the value y corresponding to that point x on the line defined by the
68     // two points.
69     // Inputs:
70     //    Double_t  x0  The first x value defining the line
71     //    Double_t  y0  The first y value defining the line
72     //    Double_t  x1  The second x value defining the line
73     //    Double_t  y1  The second y value defining the line
74     //    Double_t   x  The x value for which the y value is wanted.
75     // Outputs:
76     //    none.
77     // Return:
78     //    The value y corresponding to the point x on the line defined by
79     //    the two points (x0,y0) and (x1,y1).
80
81     if(x0==x1 && y0==y1) {
82         printf("Error: AliITSv11Geometry::Yfrom2Ponts The two points are "
83                "the same (%e,%e) and (%e,%e)",x0,y0,x1,y1);
84         return 0.0;
85     } // end if
86     if(x0==x1){
87         printf("Warning: AliITSv11Geometry::Yfrom2Points x0=%e == x1=%e. "
88                "line vertical ""returning mean y",x0,x1);
89         return 0.5*(y0+y1);
90     }// end if x0==x1
91     Double_t m = (y0-y1)/(x0-x1);
92     return m*(x-x0)+y0;
93 }
94 //______________________________________________________________________
95 Double_t AliITSv11Geometry::Xfrom2Points(Double_t x0,Double_t y0,
96                                          Double_t x1,Double_t y1,
97                                          Double_t y)const{
98     // Given the two points (x0,y0) and (x1,y1) and the location y, returns
99     // the value x corresponding to that point y on the line defined by the
100     // two points.
101     // Inputs:
102     //    Double_t  x0  The first x value defining the line
103     //    Double_t  y0  The first y value defining the line
104     //    Double_t  x1  The second x value defining the line
105     //    Double_t  y1  The second y value defining the line
106     //    Double_t   y  The y value for which the x value is wanted.
107     // Outputs:
108     //    none.
109     // Return:
110     //    The value x corresponding to the point y on the line defined by
111     //    the two points (x0,y0) and (x1,y1).
112
113     if(x0==x1 && y0==y1) {
114         printf("Error: AliITSv11Geometry::Yfrom2Ponts The two points are "
115                "the same (%e,%e) and (%e,%e)",x0,y0,x1,y1);
116         return 0.0;
117     } // end if
118     if(y0==y1){
119         printf("Warrning: AliITSv11Geometry::Yfrom2Points y0=%e == y1=%e. "
120                "line horizontal returning mean x",y0,y1);
121         return 0.5*(x0+x1);
122     }// end if y0==y1
123     Double_t m = (x0-x1)/(y0-y1);
124     return m*(y-y0)+x0;
125 }
126 //______________________________________________________________________
127 Double_t AliITSv11Geometry::RmaxFrom2Points(const TGeoPcon *p,Int_t i1,
128                                             Int_t i2,Double_t z)const{
129     // functions Require at parts of Volume A to be already defined.
130     // Retruns the value of Rmax corresponding to point z alone the line
131     // defined by the two points p.Rmax(i1),p-GetZ(i1) and p->GetRmax(i2),
132     // p->GetZ(i2).
133     // Inputs:
134     //    TGeoPcon *p  The Polycone where the two points come from
135     //    Int_t    i1  Point 1
136     //    Int_t    i2  Point 2
137     //    Double_t  z  The value of z for which Rmax is to be found
138     // Outputs:
139     //    none.
140     // Return:
141     //    Double_t Rmax the value corresponding to z
142     Double_t d0,d1,d2,r;
143
144     d0 = p->GetRmax(i1)-p->GetRmax(i2);// cout <<"L263: d0="<<d0<<endl;
145     d1 = z-p->GetZ(i2);// cout <<"L264: d1="<<d1<<endl;
146     d2 = p->GetZ(i1)-p->GetZ(i2);// cout <<"L265: d2="<<d2<<endl;
147     r  = p->GetRmax(i2) + d1*d0/d2;// cout <<"L266: r="<<r<<endl;
148     return r;
149 }
150 //______________________________________________________________________
151 Double_t AliITSv11Geometry::RminFrom2Points(const TGeoPcon *p,Int_t i1,
152                                             Int_t i2,Double_t z)const{
153     // Retruns the value of Rmin corresponding to point z alone the line
154     // defined by the two points p->GetRmin(i1),p->GetZ(i1) and 
155     // p->GetRmin(i2),  p->GetZ(i2).
156     // Inputs:
157     //    TGeoPcon *p  The Polycone where the two points come from
158     //    Int_t    i1  Point 1
159     //    Int_t    i2  Point 2
160     //    Double_t  z  The value of z for which Rmax is to be found
161     // Outputs:
162     //    none.
163     // Return:
164     //    Double_t Rmax the value corresponding to z
165
166     return p->GetRmin(i2)+(p->GetRmin(i1)-p->GetRmin(i2))*(z-p->GetZ(i2))/
167      (p->GetZ(i1)-p->GetZ(i2));
168 }
169 //______________________________________________________________________
170 Double_t AliITSv11Geometry::RFrom2Points(const Double_t *p,const Double_t *az,
171                                          Int_t i1,Int_t i2,Double_t z)const{
172     // Retruns the value of Rmin corresponding to point z alone the line
173     // defined by the two points p->GetRmin(i1),p->GetZ(i1) and 
174     // p->GetRmin(i2), p->GetZ(i2).
175     // Inputs:
176     //    Double_t az  Array of z values
177     //    Double_t  r  Array of r values
178     //    Int_t    i1  First Point in arrays
179     //    Int_t    i2  Second Point in arrays
180     //    Double_t z   Value z at which r is to be found
181     // Outputs:
182     //    none.
183     // Return:
184     //    The value r corresponding to z and the line defined by the two points
185
186     return p[i2]+(p[i1]-p[i2])*(z-az[i2])/(az[i1]-az[i2]);
187 }
188 //______________________________________________________________________
189 Double_t AliITSv11Geometry::Zfrom2MinPoints(const TGeoPcon *p,Int_t i1,
190                                             Int_t i2,Double_t r)const{
191     // Retruns the value of Z corresponding to point R alone the line
192     // defined by the two points p->GetRmin(i1),p->GetZ(i1) and 
193     // p->GetRmin(i2),p->GetZ(i2)
194     // Inputs:
195     //    TGeoPcon *p  The Poly cone where the two points come from.
196     //    Int_t    i1  First Point in arrays
197     //    Int_t    i2  Second Point in arrays
198     //    Double_t r   Value r min at which z is to be found
199     // Outputs:
200     //    none.
201     // Return:
202     //    The value z corresponding to r min and the line defined by 
203     //    the two points
204
205     return p->GetZ(i2)+(p->GetZ(i1)-p->GetZ(i2))*(r-p->GetRmin(i2))/
206      (p->GetRmin(i1)-p->GetRmin(i2));
207 }
208 //______________________________________________________________________
209 Double_t AliITSv11Geometry::Zfrom2MaxPoints(const TGeoPcon *p,Int_t i1,
210                                             Int_t i2,Double_t r)const{
211     // Retruns the value of Z corresponding to point R alone the line
212     // defined by the two points p->GetRmax(i1),p->GetZ(i1) and 
213     // p->GetRmax(i2),p->GetZ(i2)
214     // Inputs:
215     //    TGeoPcon *p  The Poly cone where the two points come from.
216     //    Int_t    i1  First Point in arrays
217     //    Int_t    i2  Second Point in arrays
218     //    Double_t r   Value r max at which z is to be found
219     // Outputs:
220     //    none.
221     // Return:
222     //    The value z corresponding to r max and the line defined by 
223     //    the two points
224
225     return p->GetZ(i2)+(p->GetZ(i1)-p->GetZ(i2))*(r-p->GetRmax(i2))/
226      (p->GetRmax(i1)-p->GetRmax(i2));
227 }
228 //______________________________________________________________________
229 Double_t AliITSv11Geometry::Zfrom2Points(const Double_t *z,const Double_t *ar,
230                                          Int_t i1,Int_t i2,Double_t r)const{
231     // Retruns the value of z corresponding to point R alone the line
232     // defined by the two points p->GetRmax(i1),p->GetZ(i1) and 
233     // p->GetRmax(i2),p->GetZ(i2)
234     // Inputs:
235     //    Double_t  z  Array of z values
236     //    Double_t ar  Array of r values
237     //    Int_t    i1  First Point in arrays
238     //    Int_t    i2  Second Point in arrays
239     //    Double_t r   Value r at which z is to be found
240     // Outputs:
241     //    none.
242     // Return:
243     //    The value z corresponding to r and the line defined by the two points
244
245     return z[i2]+(z[i1]-z[i2])*(r-ar[i2])/(ar[i1]-ar[i2]);
246 }
247 //______________________________________________________________________
248 Double_t AliITSv11Geometry::RmaxFromZpCone(const TGeoPcon *p,int ip,
249                                            Double_t tc,Double_t z,
250                                            Double_t th)const{
251     // General Outer Cone surface equation Rmax.
252     // Intputs:
253     //     TGeoPcon  *p   The poly cone where the initial point comes from
254     //     Int_t     ip   The index in p to get the point location
255     //     Double_t  tc   The angle of that part of the cone is at
256     //     Double_t   z   The value of z to compute Rmax from
257     //     Double_t  th   The perpendicular distance the parralell line is
258     //                    from the point ip.
259     // Outputs:
260     //     none.
261     // Return:
262     //     The value Rmax correstponding to the line at angle th, offeset by
263     //     th, and the point p->GetZ/Rmin[ip] at the location z.
264     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
265     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
266
267     return -tantc*(z-p->GetZ(ip))+p->GetRmax(ip)+th/costc;
268 }
269 //______________________________________________________________________
270 Double_t AliITSv11Geometry::RFromZpCone(const Double_t *ar,
271                                         const Double_t *az,int ip,
272                                         Double_t tc,Double_t z,
273                                         Double_t th)const{
274     // General Cone surface equation R(z).
275     // Intputs:
276     //     Double_t  ar   The array of R values
277     //     Double_t  az   The array of Z values
278     //     Int_t     ip   The index in p to get the point location
279     //     Double_t  tc   The angle of that part of the cone is at
280     //     Double_t   z   The value of z to compute R from
281     //     Double_t  th   The perpendicular distance the parralell line is
282     //                    from the point ip.
283     // Outputs:
284     //     none.
285     // Return:
286     //     The value R correstponding to the line at angle th, offeset by
287     //     th, and the point p->GetZ/Rmax[ip] at the locatin z.
288     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
289     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
290
291     return -tantc*(z-az[ip])+ar[ip]+th/costc;
292 }
293 //______________________________________________________________________
294 Double_t AliITSv11Geometry::RminFromZpCone(const TGeoPcon *p,Int_t ip,
295                                            Double_t tc,Double_t z,
296                                            Double_t th)const{
297     // General Inner Cone surface equation Rmin.
298     // Intputs:
299     //     TGeoPcon  *p   The poly cone where the initial point comes from
300     //     Int_t     ip   The index in p to get the point location
301     //     Double_t  tc   The angle of that part of the cone is at
302     //     Double_t   z   The value of z to compute Rmin from
303     //     Double_t  th   The perpendicular distance the parralell line is
304     //                    from the point ip.
305     // Outputs:
306     //     none.
307     // Return:
308     //     The value Rmin correstponding to the line at angle th, offeset by
309     //     th, and the point p->GetZ/Rmin[ip] at the location z.
310     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
311     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
312
313     return -tantc*(z-p->GetZ(ip))+p->GetRmin(ip)+th/costc;
314 }
315 //______________________________________________________________________
316 Double_t AliITSv11Geometry::ZFromRmaxpCone(const TGeoPcon *p,int ip,
317                                            Double_t tc,Double_t r,
318                                            Double_t th)const{
319     // General Outer cone Surface equation for z.
320     // Intputs:
321     //     TGeoPcon  *p   The poly cone where the initial point comes from
322     //     Int_t     ip   The index in p to get the point location
323     //     Double_t  tc   The angle of that part of the cone is at
324     //     Double_t   r   The value of Rmax to compute z from
325     //     Double_t  th   The perpendicular distance the parralell line is
326     //                    from the point ip.
327     // Outputs:
328     //     none.
329     // Return:
330     //     The value Z correstponding to the line at angle th, offeset by
331     //     th, and the point p->GetZ/Rmax[ip] at the location r.
332     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
333     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
334
335     return p->GetZ(ip)+(p->GetRmax(ip)+th/costc-r)/tantc;
336 }
337 //______________________________________________________________________
338 Double_t AliITSv11Geometry::ZFromRmaxpCone(const Double_t *ar,
339                                            const Double_t *az,int ip,
340                                            Double_t tc,Double_t r,
341                                            Double_t th)const{
342     // General Outer cone Surface equation for z.
343     // Intputs:
344     //     Double_t  ar   The array of R values
345     //     Double_t  az   The array of Z values
346     //     Int_t     ip   The index in p to get the point location
347     //     Double_t  tc   The angle of that part of the cone is at
348     //     Double_t   r   The value of Rmax to compute z from
349     //     Double_t  th   The perpendicular distance the parralell line is
350     //                    from the point ip.
351     // Outputs:
352     //     none.
353     // Return:
354     //     The value Z correstponding to the line at angle th, offeset by
355     //     th, and the point p->GetZ/Rmax[ip] at the locatin r.
356     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
357     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
358
359     return az[ip]+(ar[ip]+th/costc-r)/tantc;
360 }
361 //______________________________________________________________________
362 Double_t AliITSv11Geometry::ZFromRminpCone(const TGeoPcon *p,int ip,
363                                            Double_t tc,Double_t r,
364                                            Double_t th)const{
365     // General Inner cone Surface equation for z.
366     // Intputs:
367     //     TGeoPcon  *p   The poly cone where the initial point comes from
368     //     Int_t     ip   The index in p to get the point location
369     //     Double_t  tc   The angle of that part of the cone is at
370     //     Double_t   r   The value of Rmin to compute z from
371     //     Double_t  th   The perpendicular distance the parralell line is
372     //                    from the point ip.
373     // Outputs:
374     //     none.
375     // Return:
376     //     The value Z correstponding to the line at angle th, offeset by
377     //     th, and the point p->GetZ/Rmin[ip] at the location r.
378     Double_t tantc = TMath::Tan(tc*TMath::DegToRad());
379     Double_t costc = TMath::Cos(tc*TMath::DegToRad());
380
381     return p->GetZ(ip)+(p->GetRmin(ip)+th/costc-r)/tantc;
382 }
383 //______________________________________________________________________
384 void AliITSv11Geometry::RadiusOfCurvature(Double_t rc,Double_t theta0,
385                                           Double_t z0,Double_t r0,
386                                           Double_t theta1,Double_t &z1,
387                                           Double_t &r1)const{
388     // Given a initial point z0,r0, the initial angle theta0, and the radius
389     // of curvature, returns the point z1, r1 at the angle theta1. Theta
390     // measured from the r axis in the clock wise direction [degrees].
391     // Inputs:
392     //    Double_t rc     The radius of curvature
393     //    Double_t theta0 The starting angle (degrees)
394     //    Double_t z0     The value of z at theta0
395     //    Double_t r0     The value of r at theta0
396     //    Double_t theta1 The ending angle (degrees)
397     // Outputs:
398     //    Double_t &z1  The value of z at theta1
399     //    Double_t &r1  The value of r at theta1
400     // Return:
401     //    none.
402
403     z1 = rc*(TMath::Sin(theta1*TMath::DegToRad())-TMath::Sin(theta0*TMath::DegToRad()))+z0;
404     r1 = rc*(TMath::Cos(theta1*TMath::DegToRad())-TMath::Cos(theta0*TMath::DegToRad()))+r0;
405     return;
406 }
407 //______________________________________________________________________
408 void AliITSv11Geometry::InsidePoint(const TGeoPcon *p,Int_t i1,Int_t i2,
409                                     Int_t i3,Double_t c,TGeoPcon *q,Int_t j1,
410                                     Bool_t max)const{
411     // Given two lines defined by the points i1, i2,i3 in the TGeoPcon 
412     // class p that intersect at point p->GetZ(i2) return the point z,r 
413     // that is Cthick away in the TGeoPcon class q. If points i1=i2
414     // and max == kTRUE, then p->GetRmin(i1) and p->GetRmax(i2) are used.
415     // if points i2=i3 and max=kTRUE then points p->GetRmax(i2) and
416     // p->GetRmin(i3) are used. If i2=i3 and max=kFALSE, then p->GetRmin(i2)
417     // and p->GetRmax(i3) are used.
418     // Inputs:
419     //    TGeoPcon  *p  Class where points i1, i2, and i3 are taken from
420     //    Int_t     i1  First point in class p
421     //    Int_t     i2  Second point in class p
422     //    Int_t     i3  Third point in class p
423     //    Double_t  c   Distance inside the outer surface/inner suface
424     //                  that the point j1 is to be computed for.
425     //    TGeoPcon  *q  Pointer to class for results to be put into.
426     //    Int_t     j1  Point in class q where data is to be stored.
427     //    Bool_t    max if kTRUE, then a Rmax value is computed,
428     //                  else a Rmin valule is computed.
429     // Output:
430     //    TGeoPcon  *q  Pointer to class for results to be put into.
431     // Return:
432     //    none.
433     Double_t x0,y0,x1,y1,x2,y2,x,y;
434
435     if(max){
436         c = -c; //cout <<"L394 c="<<c<<endl;
437         y0 = p->GetRmax(i1);
438         if(i1==i2) y0 = p->GetRmin(i1); //cout <<"L396 y0="<<y0<<endl;
439         y1 = p->GetRmax(i2);  //cout <<"L397 y1="<<y1<<endl;
440         y2 = p->GetRmax(i3); //cout <<"L398 y2="<<y2<<endl;
441         if(i2==i3) y2 = p->GetRmin(i3); //cout <<"L399 y2="<<y2<<endl;
442     }else{ // min
443         y0 = p->GetRmin(i1); //cout <<"L401 y0="<<y0<<endl;
444         y1 = p->GetRmin(i2); //cout <<"L402 y1="<<y1<<endl;
445         y2 = p->GetRmin(i3);
446         if(i2==i3) y2 = p->GetRmax(i3); //cout <<"L404 y2="<<y2<<endl;
447     } // end if
448     x0 = p->GetZ(i1); //cout <<"L406 x0="<<x0<<endl;
449     x1 = p->GetZ(i2); //cout <<"L407 x1="<<x1<<endl;
450     x2 = p->GetZ(i3); //cout <<"L408 x2="<<x2<<endl;
451     //
452     InsidePoint(x0,y0,x1,y1,x2,y2,c,x,y);
453     q->Z(j1) = x;
454     if(max) q->Rmax(j1) = y;
455     else    q->Rmin(j1) = y;
456     return;
457 }
458 //----------------------------------------------------------------------
459 void AliITSv11Geometry::InsidePoint(Double_t x0,Double_t y0,
460                                     Double_t x1,Double_t y1,
461                                     Double_t x2,Double_t y2,Double_t c,
462                                     Double_t &x,Double_t &y)const{
463     // Given two intersecting lines defined by the points (x0,y0), (x1,y1) and
464     // (x1,y1), (x1,y2) {intersecting at (x1,y1)} the point (x,y) a distance
465     // c away is returned such that two lines a distance c away from the
466     // lines defined above intersect at (x,y).
467     // Inputs:
468     //    Double_t  x0 X point on the first intersecting sets of lines
469     //    Double_t  y0 Y point on the first intersecting sets of lines
470     //    Double_t  x1 X point on the first/second intersecting sets of lines
471     //    Double_t  y1 Y point on the first/second intersecting sets of lines
472     //    Double_t  x2 X point on the second intersecting sets of lines
473     //    Double_t  y2 Y point on the second intersecting sets of lines
474     //    Double_t  c  Distance the two sets of lines are from each other
475     // Output:
476     //    Double_t  x  X point for the intersecting sets of parellel lines
477     //    Double_t  y  Y point for the intersecting sets of parellel lines
478     // Return:
479     //    none.
480     Double_t dx01,dx12,dy01,dy12,r01,r12,m;
481     dx01 = x0-x1; //cout <<"L410 dx01="<<dx01<<endl;
482     dx12 = x1-x2; //cout <<"L411 dx12="<<dx12<<endl;
483     dy01 = y0-y1; //cout <<"L412 dy01="<<dy01<<endl;
484     dy12 = y1-y2; //cout <<"L413 dy12="<<dy12<<endl;
485     r01  = TMath::Sqrt(dy01*dy01+dx01*dx01); //cout <<"L414 r01="<<r01<<endl;
486     r12  = TMath::Sqrt(dy12*dy12+dx12*dx12); //cout <<"L415 r12="<<r12<<endl;
487     m = dx12*dy01-dy12*dx01;
488     if(m*m<DBL_EPSILON){ // m == n
489         if(dy01==0.0){ // line are =
490             x = x1+c; //cout <<"L419 x="<<x<<endl;
491             y = y1; //cout <<"L420 y="<<y<<endl;
492             return;
493         }else if(dx01==0.0){
494             x = x1;
495             y = y1+c;
496             return;
497         }else{ // dx01!=0 and dy01 !=0.
498             x = x1-0.5*c*r01/dy01; //cout <<"L434 x="<<x<<endl;
499             y = y1+0.5*c*r01/dx01; //cout <<"L435 y="<<y<<endl;
500         } // end if
501         return;
502     } //
503     x = x1+c*(dx12*r01-dx01*r12)/m; //cout <<"L442 x="<<x<<endl;
504     y = y1+c*(dy12*r01-dy01*r12)/m; //cout <<"L443 y="<<y<<endl;
505     //cout <<"=============================================="<<endl;
506     return;
507 }
508 //----------------------------------------------------------------------
509 void AliITSv11Geometry:: PrintArb8(const TGeoArb8 *a)const{
510     // Prints out the content of the TGeoArb8. Usefull for debugging.
511     // Inputs:
512     //   TGeoArb8 *a
513     // Outputs:
514     //   none.
515     // Return:
516     //   none.
517
518     if(!GetDebug()) return;
519     printf("%s",a->GetName());
520     a->InspectShape();
521     return;
522 }
523 //----------------------------------------------------------------------
524 void AliITSv11Geometry:: PrintPcon(const TGeoPcon *a)const{
525     // Prints out the content of the TGeoPcon. Usefull for debugging.
526     // Inputs:
527     //   TGeoPcon *a
528     // Outputs:
529     //   none.
530     // Return:
531     //   none.
532   
533     if(!GetDebug()) return;
534     cout << a->GetName() << ": N=" << a->GetNz() << " Phi1=" << a->GetPhi1()
535          << ", Dphi=" << a->GetDphi() << endl;
536     cout << "i\t   Z   \t  Rmin \t  Rmax" << endl;
537     for(Int_t iii=0;iii<a->GetNz();iii++){
538         cout << iii << "\t" << a->GetZ(iii) << "\t" << a->GetRmin(iii)
539              << "\t" << a->GetRmax(iii) << endl;
540     } // end for iii
541     return;
542 }
543 //----------------------------------------------------------------------
544 void AliITSv11Geometry::PrintTube(const TGeoTube *a)const{
545     // Prints out the content of the TGeoTube. Usefull for debugging.
546     // Inputs:
547     //   TGeoTube *a
548     // Outputs:
549     //   none.
550     // Return:
551     //   none.
552
553     if(!GetDebug()) return;
554     cout << a->GetName() <<": Rmin="<<a->GetRmin()
555          <<" Rmax=" <<a->GetRmax()<<" Dz="<<a->GetDz()<<endl;
556     return;
557 }
558 //----------------------------------------------------------------------
559 void AliITSv11Geometry::PrintTubeSeg(const TGeoTubeSeg *a)const{
560     // Prints out the content of the TGeoTubeSeg. Usefull for debugging.
561     // Inputs:
562     //   TGeoTubeSeg *a
563     // Outputs:
564     //   none.
565     // Return:
566     //   none.
567
568     if(!GetDebug()) return;
569     cout << a->GetName() <<": Phi1="<<a->GetPhi1()<<
570         " Phi2="<<a->GetPhi2()<<" Rmin="<<a->GetRmin()
571          <<" Rmax=" <<a->GetRmax()<<" Dz="<<a->GetDz()<<endl;
572     return;
573 }
574 //----------------------------------------------------------------------
575 void AliITSv11Geometry::PrintConeSeg(const TGeoConeSeg *a)const{
576     // Prints out the content of the TGeoConeSeg. Usefull for debugging.
577     // Inputs:
578     //   TGeoConeSeg *a
579     // Outputs:
580     //   none.
581     // Return:
582     //   none.
583
584     if(!GetDebug()) return;
585     cout << a->GetName() <<": Phi1="<<a->GetPhi1()<<
586         " Phi2="<<a->GetPhi2()<<" Rmin1="<<a->GetRmin1()
587          <<" Rmax1=" <<a->GetRmax1()<<" Rmin2="<<a->GetRmin2()
588          <<" Rmax2=" <<a->GetRmax2()<<" Dz="<<a->GetDz()<<endl;
589     return;
590 }
591 //----------------------------------------------------------------------
592 void AliITSv11Geometry::PrintBBox(const TGeoBBox *a)const{
593     // Prints out the content of the TGeoBBox. Usefull for debugging.
594     // Inputs:
595     //   TGeoBBox *a
596     // Outputs:
597     //   none.
598     // Return:
599     //   none.
600
601     if(!GetDebug()) return;
602     cout << a->GetName() <<": Dx="<<a->GetDX()<<
603         " Dy="<<a->GetDY()<<" Dz="<<a->GetDZ() <<endl;
604     return;
605 }
606 //---------------------------------------------------------------------
607 void AliITSv11Geometry::DrawCrossSection(const TGeoPcon *p,
608                             Int_t fillc,Int_t fills,
609                             Int_t linec,Int_t lines,Int_t linew,
610                             Int_t markc,Int_t marks,Float_t marksize)const{
611     // Draws a cross sectional view of the TGeoPcon, Primarily for debugging.
612     // A TCanvas should exist first.
613     //  Inputs:
614     //    TGeoPcon  *p  The TGeoPcon to be "drawn"
615     //    Int_t  fillc  The fill color to be used
616     //    Int_t  fills  The fill style to be used
617     //    Int_t  linec  The line color to be used
618     //    Int_t  lines  The line style to be used
619     //    Int_t  linew  The line width to be used
620     //    Int_t  markc  The markder color to be used
621     //    Int_t  marks  The markder style to be used
622     //    Float_t marksize The marker size
623     // Outputs:
624     //   none.
625     // Return:
626     //   none.
627     Int_t n=0,m=0,i=0;
628     Double_t *z=0,*r=0;
629     TPolyMarker *pts=0;
630     TPolyLine   *line=0;
631
632     n = p->GetNz();
633     if(n<=0) return;
634     m = 2*n+1;
635     z = new Double_t[m];
636     r = new Double_t[m];
637
638     for(i=0;i<n;i++){
639         z[i] = p->GetZ(i);
640         r[i] = p->GetRmax(i);
641         z[i+n] = p->GetZ(n-1-i);
642         r[i+n] = p->GetRmin(n-1-i);
643     } //  end for i
644     z[n-1] = z[0];
645     r[n-1] = r[0];
646
647     line = new TPolyLine(n,z,r);
648     pts  = new TPolyMarker(n,z,r);
649
650     line->SetFillColor(fillc);
651     line->SetFillStyle(fills);
652     line->SetLineColor(linec);
653     line->SetLineStyle(lines);
654     line->SetLineWidth(linew);
655     pts->SetMarkerColor(markc);
656     pts->SetMarkerStyle(marks);
657     pts->SetMarkerSize(marksize);
658
659     line->Draw("f");
660     line->Draw();
661     pts->Draw();
662
663     delete[] z;
664     delete[] r;
665
666     cout<<"Hit Return to continue"<<endl;
667     cin >> n;
668     delete line;
669     delete pts;
670     return;
671 }
672 //______________________________________________________________________
673 Bool_t AliITSv11Geometry::AngleOfIntersectionWithLine(Double_t x0,Double_t y0,
674                                                       Double_t x1,Double_t y1,
675                                                       Double_t xc,Double_t yc,
676                                                       Double_t rc,Double_t &t0,
677                                                       Double_t &t1)const{
678     // Computes the angles, t0 and t1 corresponding to the intersection of
679     // the line, defined by {x0,y0} {x1,y1}, and the circle, defined by
680     // its center {xc,yc} and radius r. If the line does not intersect the
681     // line, function returns kFALSE, otherwise it returns kTRUE. If the
682     // line is tangent to the circle, the angles t0 and t1 will be the same.
683     // Inputs:
684     //   Double_t x0   X of first point defining the line
685     //   Double_t y0   Y of first point defining the line
686     //   Double_t x1   X of Second point defining the line
687     //   Double_t y1   Y of Second point defining the line
688     //   Double_t xc   X of Circle center point defining the line
689     //   Double_t yc   Y of Circle center point defining the line
690     //   Double_t r    radius of circle
691     // Outputs:
692     //   Double_t &t0  First angle where line intersects circle
693     //   Double_t &t1  Second angle where line intersects circle
694     // Return:
695     //    kTRUE, line intersects circle, kFALSE line does not intersect circle
696     //           or the line is not properly defined point {x0,y0} and {x1,y1}
697     //           are the same point.
698     Double_t dx,dy,cx,cy,s2,t[4];
699     Double_t a0,b0,c0,a1,b1,c1,sinthp,sinthm,costhp,costhm;
700     Int_t i,j;
701
702     t0 = 400.0;
703     t1 = 400.0;
704     dx = x1-x0;
705     dy = y1-y0;
706     cx = xc-x0;
707     cy = yc-y0;
708     s2 = dx*dx+dy*dy;
709     if(s2==0.0) return kFALSE;
710
711     a0 = rc*rc*s2;
712     if(a0==0.0) return kFALSE;
713     b0 = 2.0*rc*dx*(dx*cy-cx*dy);
714     c0 = dx*dx*cy*cy-2.0*dy*dx*cy*cx+cx*cx*dy*dy-rc*rc*dy*dy;
715     c0 = 0.25*b0*b0/(a0*a0)-c0/a0;
716     if(c0<0.0) return kFALSE;
717     sinthp = -0.5*b0/a0+TMath::Sqrt(c0);
718     sinthm = -0.5*b0/a0-TMath::Sqrt(c0);
719
720     a1 = rc*rc*s2;
721     if(a1==0.0) return kFALSE;
722     b1 = 2.0*rc*dy*(dy*cx-dx*cy);
723     c1 = dy*dy*cx*cx-2.0*dy*dx*cy*cx+dx*dx*cy*cy-rc*rc*dx*dx;
724     c1 = 0.25*b1*b1/(a1*a1)-c1/a1;
725     if(c1<0.0) return kFALSE;
726     costhp = -0.5*b1/a1+TMath::Sqrt(c1);
727     costhm = -0.5*b1/a1-TMath::Sqrt(c1);
728
729     t[0] = t[1] = t[2] = t[3] = 400.;
730     a0 = TMath::ATan2(sinthp,costhp); if(a0<0.0) a0 += 2.0*TMath::Pi();
731     a1 = TMath::ATan2(sinthp,costhm); if(a1<0.0) a1 += 2.0*TMath::Pi();
732     b0 = TMath::ATan2(sinthm,costhp); if(b0<0.0) b0 += 2.0*TMath::Pi();
733     b1 = TMath::ATan2(sinthm,costhm); if(b1<0.0) b1 += 2.0*TMath::Pi();
734     x1 = xc+rc*TMath::Cos(a0);
735     y1 = yc+rc*TMath::Sin(a0);
736     s2 = dx*(y1-y0)-dy*(x1-x0);
737     if(s2*s2<DBL_EPSILON) t[0] = a0*TMath::RadToDeg();
738     x1 = xc+rc*TMath::Cos(a1);
739     y1 = yc+rc*TMath::Sin(a1);
740     s2 = dx*(y1-y0)-dy*(x1-x0);
741     if(s2*s2<DBL_EPSILON) t[1] = a1*TMath::RadToDeg();
742     x1 = xc+rc*TMath::Cos(b0);
743     y1 = yc+rc*TMath::Sin(b0);
744     s2 = dx*(y1-y0)-dy*(x1-x0);
745     if(s2*s2<DBL_EPSILON) t[2] = b0*TMath::RadToDeg();
746     x1 = xc+rc*TMath::Cos(b1);
747     y1 = yc+rc*TMath::Sin(b1);
748     s2 = dx*(y1-y0)-dy*(x1-x0);
749     if(s2*s2<DBL_EPSILON) t[3] = b1*TMath::RadToDeg();
750     for(i=0;i<4;i++)for(j=i+1;j<4;j++){
751         if(t[i]>t[j]) {t0 = t[i];t[i] = t[j];t[j] = t0;}
752     } // end for i,j
753     t0 = t[0];
754     t1 = t[1];
755     //
756     return kTRUE;
757 }
758 //______________________________________________________________________
759 Double_t AliITSv11Geometry::AngleForRoundedCorners0(Double_t dx,Double_t dy,
760                                                     Double_t sdr)const{
761     // Basic function used to determine the ending angle and starting angles
762     // for rounded corners given the relative distance between the centers
763     // of the circles and the difference/sum of their radii. Case 0.
764     // Inputs:
765     //   Double_t dx    difference in x locations of the circle centers
766     //   Double_t dy    difference in y locations of the circle centers
767     //   Double_t sdr   difference or sum of the circle radii
768     // Outputs:
769     //   none.
770     // Return:
771     //   the angle in Degrees
772     Double_t a,b;
773
774     b = dy*dy+dx*dx-sdr*sdr;
775     if(b<0.0) Error("AngleForRoundedCorners0",
776                     "dx^2(%e)+dy^2(%e)-sdr^2(%e)=b=%e<0",dx,dy,sdr,b);
777     b = TMath::Sqrt(b);
778     a = -sdr*dy+dx*b;
779     b = -sdr*dx-dy*b;
780     return TMath::ATan2(a,b)*TMath::RadToDeg();
781     
782 }
783 //______________________________________________________________________
784 Double_t AliITSv11Geometry::AngleForRoundedCorners1(Double_t dx,Double_t dy,
785                                                     Double_t sdr)const{
786     // Basic function used to determine the ending angle and starting angles
787     // for rounded corners given the relative distance between the centers
788     // of the circles and the difference/sum of their radii. Case 1.
789     // Inputs:
790     //   Double_t dx    difference in x locations of the circle centers
791     //   Double_t dy    difference in y locations of the circle centers
792     //   Double_t sdr   difference or sum of the circle radii
793     // Outputs:
794     //   none.
795     // Return:
796     //   the angle in Degrees
797     Double_t a,b;
798
799     b = dy*dy+dx*dx-sdr*sdr;
800     if(b<0.0) Error("AngleForRoundedCorners1",
801                     "dx^2(%e)+dy^2(%e)-sdr^2(%e)=b=%e<0",dx,dy,sdr,b);
802     b = TMath::Sqrt(b);
803     a = -sdr*dy-dx*b;
804     b = -sdr*dx+dy*b;
805     return TMath::ATan2(a,b)*TMath::RadToDeg();
806     
807 }
808 //----------------------------------------------------------------------
809 void AliITSv11Geometry::AnglesForRoundedCorners(Double_t x0,Double_t y0,
810                                                 Double_t r0,Double_t x1,
811                                                 Double_t y1,Double_t r1,
812                                                 Double_t &t0,Double_t &t1)
813     const{
814     // Function to compute the ending angle, for arc 0, and starting angle,
815     // for arc 1, such that a straight line will connect them with no
816     // discontinuities.
817     //Begin_Html
818     /*
819       <img src="picts/ITS/AliITSv11Geometry_AnglesForRoundedCorners.gif">
820      */
821     //End_Html
822     // Inputs:
823     //    Double_t x0  X Coordinate of arc 0 center.
824     //    Double_t y0  Y Coordinate of arc 0 center.
825     //    Double_t r0  Radius of curvature of arc 0. For signe see figure.
826     //    Double_t x1  X Coordinate of arc 1 center.
827     //    Double_t y1  Y Coordinate of arc 1 center.
828     //    Double_t r1  Radius of curvature of arc 1. For signe see figure.
829     // Outputs:
830     //    Double_t t0  Ending angle of arch 0, with respect to x axis, Degrees.
831     //    Double_t t1  Starting angle of arch 1, with respect to x axis, 
832     //                 Degrees.
833     // Return:
834     //    none.
835     Double_t t;
836
837     if(r0>=0.0&&r1>=0.0) { // Inside to inside    ++
838         t = AngleForRoundedCorners1(x1-x0,y1-y0,r1-r0);
839         t0 = t1 = t;
840         return;
841     }else if(r0>=0.0&&r1<=0.0){ // Inside to Outside  +-
842         r1 = -r1; // make positive
843         t = AngleForRoundedCorners0(x1-x0,y1-y0,r1+r0);
844         t0 = 180.0 + t;
845         if(t0<0.0) t += 360.;
846         if(t<0.0) t += 360.;
847         t1 = t;
848         return;
849     }else if(r0<=0.0&&r1>=0.0){ // Outside to Inside  -+
850         r0 = - r0; // make positive
851         t = AngleForRoundedCorners1(x1-x0,y1-y0,r1+r0);
852         t0 = 180.0 + t;
853         if(t0>180.) t0 -= 360.;
854         if(t >180.) t  -= 360.;
855         t1 = t;
856         return;
857     }else if(r0<=0.0&&r1<=0.0) { // Outside to outside --
858         r0 = -r0; // make positive
859         r1 = -r1; // make positive
860         t = AngleForRoundedCorners0(x1-x0,y1-y0,r1-r0);
861         t0 = t1 = t;
862         return;
863     } // end if
864     return;
865 }
866 //----------------------------------------------------------------------
867 void AliITSv11Geometry::MakeFigure1(Double_t x0,Double_t y0,Double_t r0,
868                                     Double_t x1,Double_t y1,Double_t r1){
869     // Function to create the figure discribing how the function 
870     // AnglesForRoundedCorners works.
871     //
872     // Inputs:
873     //    Double_t x0  X Coordinate of arc 0 center.
874     //    Double_t y0  Y Coordinate of arc 0 center.
875     //    Double_t r0  Radius of curvature of arc 0. For signe see figure.
876     //    Double_t x1  X Coordinate of arc 1 center.
877     //    Double_t y1  Y Coordinate of arc 1 center.
878     //    Double_t r1  Radius of curvature of arc 1. For signe see figure.
879     // Outputs:
880     //    none.
881     // Return:
882     //    none.
883     Double_t t0[4],t1[4],xa0[4],ya0[4],xa1[4],ya1[4],ra0[4],ra1[4];
884     Double_t xmin,ymin,xmax,ymax,h;
885     Int_t j;
886
887     for(j=0;j<4;j++) {
888         ra0[j] = r0; if(j%2) ra0[j] = -r0;
889         ra1[j] = r1; if(j>1) ra1[j] = -r1;
890         AnglesForRoundedCorners(x0,y0,ra0[j],x1,y1,ra1[j],t0[j],t1[j]);
891         xa0[j] = TMath::Abs(r0)*CosD(t0[j])+x0;
892         ya0[j] = TMath::Abs(r0)*SinD(t0[j])+y0;
893         xa1[j] = TMath::Abs(r1)*CosD(t1[j])+x1;
894         ya1[j] = TMath::Abs(r1)*SinD(t1[j])+y1;
895     } // end for j
896     if(r0<0.0) r0 = -r0;
897     if(r1<0.0) r1 = -r1;
898     xmin = TMath::Min(x0 - r0,x1-r1);
899     ymin = TMath::Min(y0 - r0,y1-r1);
900     xmax = TMath::Max(x0 + r0,x1+r1);
901     ymax = TMath::Max(y0 + r0,y1+r1);
902     for(j=1;j<4;j++) {
903         xmin = TMath::Min(xmin,xa0[j]);
904         xmin = TMath::Min(xmin,xa1[j]);
905         ymin = TMath::Min(ymin,ya0[j]);
906         ymin = TMath::Min(ymin,ya1[j]);
907
908         xmax = TMath::Max(xmax,xa0[j]);
909         xmax = TMath::Max(xmax,xa1[j]);
910         ymax = TMath::Max(ymax,ya0[j]);
911         ymax = TMath::Max(ymax,ya1[j]);
912     } // end for j
913     if(xmin<0.0) xmin *= 1.1; else xmin *= 0.9;
914     if(ymin<0.0) ymin *= 1.1; else ymin *= 0.9;
915     if(xmax<0.0) xmax *= 0.9; else xmax *= 1.1;
916     if(ymax<0.0) ymax *= 0.9; else ymax *= 1.1;
917     j = (Int_t)(500.0*(ymax-ymin)/(xmax-xmin));
918     TCanvas *can = new TCanvas("AliITSv11Geometry_AnglesForRoundedCorners",
919                                "Figure for AliITSv11Geometry",500,j);
920     h = ymax-ymin; if(h<0) h = -h;
921     can->Range(xmin,ymin,xmax,ymax);
922     TArc *c0 = new TArc(x0,y0,r0);
923     TArc *c1 = new TArc(x1,y1,r1);
924     TLine *line[4];
925     TArrow *ar0[4];
926     TArrow *ar1[4];
927     for(j=0;j<4;j++){
928         ar0[j] = new TArrow(x0,y0,xa0[j],ya0[j]);
929         ar1[j] = new TArrow(x1,y1,xa1[j],ya1[j]);
930         line[j] = new TLine(xa0[j],ya0[j],xa1[j],ya1[j]);
931         ar0[j]->SetLineColor(j+1);
932         ar0[j]->SetArrowSize(0.1*r0/h);
933         ar1[j]->SetLineColor(j+1);
934         ar1[j]->SetArrowSize(0.1*r1/h);
935         line[j]->SetLineColor(j+1);
936     } // end for j
937     c0->Draw();
938     c1->Draw();
939     for(j=0;j<4;j++){
940         ar0[j]->Draw();
941         ar1[j]->Draw();
942         line[j]->Draw();
943     } // end for j
944     TText *t = new TText();
945     t->SetTextSize(0.02);
946     Char_t txt[100];
947     sprintf(txt,"(x0=%5.2f,y0=%5.2f)",x0,y0);
948     t->DrawText(x0,y0,txt);
949     sprintf(txt,"(x1=%5.2f,y1=%5.2f)",x1,y1);
950     for(j=0;j<4;j++) {
951         t->SetTextColor(j+1);
952         t->DrawText(x1,y1,txt);
953         sprintf(txt,"r0=%5.2f",ra0[j]);
954         t->DrawText(0.5*(x0+xa0[j]),0.5*(y0+ya0[j]),txt);
955         sprintf(txt,"r1=%5.2f",ra1[j]);
956         t->DrawText(0.5*(x1+xa1[j]),0.5*(y1+ya1[j]),txt);
957     } // end for j
958 }