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