]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliHelix.cxx
Improved version of online detector-algorithm makefile
[u/mrichter/AliRoot.git] / RALICE / AliHelix.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 // $Id$
17
18 ///////////////////////////////////////////////////////////////////////////
19 // Class AliHelix
20 // Representation and extrapolation of AliTracks in a magnetic field.
21 //
22 // This class is meant to provide a means to display and extrapolate
23 // AliTrack objects in the presence of a constant homogeneous magnetic field. 
24 //
25 // For track/event displays the line width, colour etc... can be set using the
26 // standard facilities (see TAttLine).
27 // By default the linewith is set to 2 and the colour set to -1 in the constructor.
28 // The latter results in an automatic colour coding according to the track charge
29 // with the convention positive=red neutral=green negative=blue.
30 //
31 // To indicate the track starting point, the memberfunction SetMarker()
32 // may be used.
33 // By default no marker will be displayed. 
34 //
35 // Examples :
36 // ==========
37 //
38 // Display and extrapolation of individual tracks 
39 // ----------------------------------------------
40 // Float_t vec[3];
41 // AliPosition r1;
42 // Ali3Vector p;
43 // AliTrack t;
44 //
45 // vec[0]=0;
46 // vec[1]=0;
47 // vec[2]=0;
48 // r1.SetVector(vec,"car");
49 //
50 // vec[0]=1;
51 // vec[1]=0;
52 // vec[2]=0.3;
53 // p.SetVector(vec,"car");
54 //
55 // t.Set3Momentum(p);
56 // t.SetBeginPoint(r1);
57 // t.SetCharge(-1);
58 // t.SetMass(0.139);
59 //
60 // // The magnetic field vector in Tesla
61 // Ali3Vector b;
62 // vec[0]=0;
63 // vec[1]=0;
64 // vec[2]=1;
65 // b.SetVector(vec,"car");
66 //
67 // AliHelix* helix=new AliHelix();
68 // helix->SetB(b);
69 // helix->SetTofmax(1e-7);
70 //
71 // TCanvas* c1=new TCanvas("c1","c1");
72 // TView* view=new TView(1);
73 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
74 // view->ShowAxis();
75 //
76 // // Track displays 
77 // Double_t range[2]={0,600};
78 // helix->Display(&t,range,3);
79 // t.SetCharge(-t.GetCharge());
80 // helix->Display(&t);
81 //
82 // // Track extrapolation
83 // Double_t pars[3]={550,0.001,3};
84 // AliPosition* rext=helix->Extrapolate(&t,pars);
85 // if (rext) rext->Data();
86 // ======================================================================
87 //
88 // Online display of events generated via AliCollider
89 // -------------------------------------------------- 
90 // Int_t nevents=5;   // Number of events to be generated
91 // Int_t jrun=1;      // The run number of this batch of generated events
92 //
93 // cout << " ***" << endl;
94 // cout << " *** AliCollider run for " << nevents << " events." << endl; 
95 // cout << " ***" << endl;
96 //
97 // AliCollider* gen=new AliCollider();
98 //
99 // gen->OpenFortranFile(6,"dump.log");
100 //
101 // gen->SetVertexMode(2);
102 // gen->SetResolution(1e-4);
103 //
104 // gen->SetRunNumber(jrun);
105 // gen->SetPrintFreq(1);
106 //
107 // gen->SetSpectatorPmin(0.01);
108 //
109 // Int_t zp=1;
110 // Int_t ap=1;
111 // Int_t zt=2;
112 // Int_t at=4;
113 //
114 // gen->Init("fixt",zp,ap,zt,at,158);
115 //
116 // AliHelix* helix=new AliHelix();
117 // Float_t vec[3]={0,2,0};
118 // Ali3Vector b;
119 // b.SetVector(vec,"car");
120 // helix->SetB(b);
121 //
122 // helix->Refresh(-1); // Refresh display after each event
123 //
124 // TCanvas* c1=new TCanvas("c1","c1");
125 // TView* view=new TView(1);
126 // view->SetRange(-200,-200,-200,200,200,200);
127 // view->ShowAxis();
128 //
129 // // Prepare random number sequence for this run
130 // // to obtain the number of participants for each event
131 // AliRandom rndm(abs(jrun));
132 // Float_t* rans=new Float_t[nevents];
133 // rndm.Uniform(rans,nevents,2,ap+at);
134 // Int_t npart=0;
135 // Int_t ntk=0;
136 // for (Int_t i=0; i<nevents; i++)
137 // {
138 //  npart=rans[i];
139 //  gen->MakeEvent(npart);
140 //  AliEvent* evt=gen->GetEvent();
141 //  if (evt)
142 //  {
143 //   helix->Display(evt);
144 //   c1->Update();
145 //   gSystem->Sleep(5000); // Some delay to keep the display on screen
146 //  }
147 // }
148 // ======================================================================
149 //
150 //--- Author: Nick van Eijndhoven 17-jun-2004 Utrecht University
151 //- Modified: NvE $Date$ Utrecht University
152 ///////////////////////////////////////////////////////////////////////////
153
154 #include "AliHelix.h"
155 #include "Riostream.h"
156  
157 ClassImp(AliHelix) // Class implementation to enable ROOT I/O
158  
159 AliHelix::AliHelix() : THelix()
160 {
161 // Default constructor
162  fRefresh=0;
163  fCurves=0;
164  fExt=0;
165  fTofmax=1e-8;
166  fMstyle=-1;
167  fMsize=0;
168  fMcol=0;
169  fEnduse=1;
170
171  fLineWidth=2;
172  fLineColor=-1;
173 }
174 ///////////////////////////////////////////////////////////////////////////
175 AliHelix::~AliHelix()
176 {
177 // Destructor to delete dynamically allocated memory.
178  if (fCurves)
179  {
180   delete fCurves;
181   fCurves=0;
182  }
183  if (fExt)
184  {
185   delete fExt;
186   fExt=0;
187  }
188 }
189 ///////////////////////////////////////////////////////////////////////////
190 AliHelix::AliHelix(const AliHelix& h) : THelix(h)
191 {
192 // Copy constructor
193  fB=h.fB;
194  fRefresh=h.fRefresh;
195  fTofmax=h.fTofmax;
196  fMstyle=h.fMstyle;
197  fMsize=h.fMsize;
198  fMcol=h.fMcol;
199  fEnduse=h.fEnduse;
200 }
201 ///////////////////////////////////////////////////////////////////////////
202 void AliHelix::SetB(Ali3Vector& b)
203 {
204 // Set the magnetic field vector in Tesla.
205  fB=b;
206
207  if (fB.GetNorm()>0)
208  {
209   Double_t axis[3];
210   fB.GetVector(axis,"car");
211   SetAxis(axis);
212  }
213 }
214 ///////////////////////////////////////////////////////////////////////////
215 Ali3Vector& AliHelix::GetB()
216 {
217 // Provide the magnetic field vector in Tesla.
218  return fB;
219 }
220 ///////////////////////////////////////////////////////////////////////////
221 void AliHelix::SetTofmax(Float_t tof)
222 {
223 // Set the maximum time of flight for straight tracks in seconds.
224 // This maximum tof will be used for drawing etc... in case no begin
225 // and endpoints can be determined from the track info.
226 // Notes :
227 // -------
228 // 1) In case the user specifies an explicit range, it will override
229 //    the maximum tof limit.
230 // 2) By default the tofmax is set to 10 ns in the AliHelix constructor.
231  fTofmax=tof;
232 }
233 ///////////////////////////////////////////////////////////////////////////
234 Float_t AliHelix::GetTofmax() const
235 {
236 // Provide the maximum time of flight for straight tracks in seconds.
237  return fTofmax;
238 }
239 ///////////////////////////////////////////////////////////////////////////
240 void AliHelix::SetMarker(Int_t style,Float_t size,Int_t col)
241 {
242 // Specify the marker (style, size and colour) to indicate the starting point
243 // of a track in a display.
244 // In case col<0 the marker will have the same color as the track itself.
245 // 
246 // Defaults are style=8, size=0.2 and col=-1.
247  
248  fMstyle=style;
249  fMsize=size;
250  fMcol=col;
251 }
252 ///////////////////////////////////////////////////////////////////////////
253 void AliHelix::UseEndPoint(Int_t mode)
254 {
255 // Select usage of track endpoint in drawing and extrapolation.
256 // This allows correct event displays even for very long tracks.
257 //
258 // mode = 0 : Do not use the track endpoint
259 //        1 : Use the track endpoint
260 // 
261 // The default value is mode=1 (which is also set in the constructor).
262
263  if (mode==0 || mode==1) fEnduse=mode; 
264 }
265 ///////////////////////////////////////////////////////////////////////////
266 void AliHelix::MakeCurve(AliTrack* t,Double_t* range,Int_t iaxis,Double_t scale)
267 {
268 // Make the helix curve for the specified AliTrack.
269 // Detailed information of all the helix points can be obtained via the
270 // GetN() and GetP() memberfunctions of TPolyLine3D.
271 // In case one wants to display or extrapolate an AliTrack it is preferable
272 // to use the Display() or Extrapolate() memberfunctions.
273 // It is assumed that the track charge is stored in elementary units
274 // (i.e. charge=1 for a proton).
275 // The input argument "scale" specifies the unit scale for the various
276 // locations where scale=0.01 indicates unit scales in cm etc...
277 // In case scale<=0, the unit scale for locations is determined from the
278 // begin, reference or endpoint of the track. If neither of these
279 // positions is present, all locations are assumed to be given in meter.
280 // The lower and upper bounds for the range are specified by range[0] and
281 // range[1] and the argument "iaxis" indicates along which axis this range
282 // is specified.
283 // The range can be specified either in the LAB frame or in the Helix frame.
284 // The latter is the frame in which the Z axis points in the B direction.
285 //
286 // The conventions for the "iaxis" argument are the following :
287 // iaxis = 1 ==> X axis in the LAB frame
288 //         2 ==> Y axis in the LAB frame
289 //         3 ==> Z axis in the LAB frame
290 //        -1 ==> X axis in the Helix frame
291 //        -2 ==> Y axis in the Helix frame
292 //        -3 ==> Z axis in the Helix frame
293 //
294 // In case range=0 the begin/end/reference points of the AliTrack and the
295 // maximum time of flight (see the SetTofmax() memberfunction) will be used
296 // and an appropriate choice for the iaxis parameter will be made automatically
297 // based on the track kinematics.
298 // In case the reference point is not present, the begin or endpoint will be used
299 // as reference point for the 3-momentum specification. If neither of these positions
300 // is present, (0,0,0) will be taken as the reference point.
301 // 
302 // The default values are range=0, iaxis=3 and scale=-1.
303
304  SetPolyLine(0); // Reset the polyline data points
305
306  if (!t || (range && !iaxis)) return;
307
308  Double_t energy=t->GetEnergy(1); // Track energy in GeV
309  Double_t betanorm=t->GetBeta();
310
311  if (energy<=0 || betanorm<=0) return;
312
313  AliPosition* rbeg=t->GetBeginPoint();
314  AliPosition* rend=0;
315  if (fEnduse) rend=t->GetEndPoint();
316  AliPosition* rref=t->GetReferencePoint();
317
318  // Magnetic field vector or default Z-direction
319  Double_t bvec[3]={0,0,1};
320  if (fB.GetNorm()>0) fB.GetVector(bvec,"car");
321
322  // The unit scale for locations if not specified by the user
323  if (scale<=0)
324  {
325   scale=1; // Set default to meter
326   if (rbeg)
327   {
328    scale=rbeg->GetUnitScale();
329   }
330   else if (rend)
331   {
332    scale=rend->GetUnitScale();
333   }
334   else if (rref)
335   {
336    scale=rref->GetUnitScale();
337   }
338  }
339
340  Double_t c=2.99792458e8/scale; // Lightspeed in the selected unit scale
341
342  // The helix angular frequency
343  Double_t w=9e7*(t->GetCharge()*fB.GetNorm())/energy;
344
345  // The particle velocity in the LAB frame
346  Ali3Vector beta=t->GetBetaVector();
347  Ali3Vector v=beta*c;
348  Double_t vel[3];
349  v.GetVector(vel,"car");
350
351  // The particle velocity in the Helix frame
352  Ali3Vector betaprim=beta.GetPrimed(fRotMat);
353  v=v.GetPrimed(fRotMat);
354  Double_t velprim[3];
355  v.GetVector(velprim,"car");
356
357  // Check compatibility of velocity and range specification.
358  if (range)
359  {
360   Double_t betavec[3];
361   if (iaxis>0) beta.GetVector(betavec,"car");
362   if (iaxis<0) betaprim.GetVector(betavec,"car");
363   if (fabs(betavec[abs(iaxis)-1])/betanorm<1e-10) return;
364  }
365
366  // The LAB location in which the velocity of the particle is defined
367  Double_t loc[3]={0,0,0};
368  Ali3Vector rx;
369  Double_t scalex=-1;
370  if (rref)
371  {
372   rx=(Ali3Vector)(*rref);
373   scalex=rref->GetUnitScale();
374  }
375  else if (rbeg)
376  {
377   rx=(Ali3Vector)(*rbeg);
378   scalex=rbeg->GetUnitScale();
379  }
380  else if (rend)
381  {
382   rx=(Ali3Vector)(*rend);
383   scalex=rend->GetUnitScale();
384  }
385
386  if (scalex>0 && (scalex/scale>1.1 || scale/scalex>1.1)) rx*=scalex/scale;
387  rx.GetVector(loc,"car");
388
389  // Initialisation of Helix kinematics
390  SetHelix(loc,vel,w,0,kUnchanged,bvec);
391
392  Int_t bend=0;
393  if (fabs(w)>0 && fabs(fVt)>0) bend=1;
394
395  // Flight time boundaries.
396  // The time origin t=0 is chosen to indicate the position in which
397  // the particle velocity was defined.
398  // The total flight time is initialised to the (user specified) tofmax.
399  Double_t tmin=0,tmax=0;
400  Double_t tof=fTofmax;
401  Double_t dum=0;
402
403  // The trajectory begin and end points
404  Double_t vec1[3]={0,0,0};
405  Double_t vec2[3]={0,0,0};
406  Ali3Vector r1;
407  Ali3Vector r2;
408  Double_t scale1=1;
409  Double_t scale2=1;
410
411  if (!bend)
412  {
413   ////////////////////////////////////////
414   // Treatment of straight trajectories //
415   ////////////////////////////////////////
416   Ali3Vector r;
417   if (range) // Specified range allows for exact flight time boundaries
418   {
419    if (iaxis>0)
420    {
421     tmin=(range[0]-loc[iaxis-1])/vel[iaxis-1];
422     tmax=(range[1]-loc[iaxis-1])/vel[iaxis-1];
423    }
424    else
425    {
426     loc[0]=fX0;
427     loc[1]=fY0;
428     loc[2]=fZ0;
429     tmin=(range[0]-loc[abs(iaxis)-1])/velprim[abs(iaxis)-1];
430     tmax=(range[1]-loc[abs(iaxis)-1])/velprim[abs(iaxis)-1];
431    }
432    if (tmax<tmin)
433    {
434     dum=tmin;
435     tmin=tmax;
436     tmax=dum;
437    }
438    // Make the 'curve' in the LAB frame and exit.
439    // Use the parametrisation : r(t)=r0+t*v
440    // using the range based flight time boundaries.
441    // An additional point in the middle of the trajectory is
442    // generated in view of accuracy in the case of extrapolations.
443    tof=tmax-tmin;
444    v=beta*c;
445    r1=rx;
446    r=v*tmin;
447    r1=r1+r;
448    r1.GetVector(vec1,"car");
449    SetNextPoint(float(vec1[0]),float(vec1[1]),float(vec1[2]));
450    r=v*(tof/2.);
451    r2=r1+r;
452    r2.GetVector(vec2,"car");
453    SetNextPoint(float(vec2[0]),float(vec2[1]),float(vec2[2]));
454    r=v*tof;
455    r2=r1+r;
456    r2.GetVector(vec2,"car");
457    SetNextPoint(float(vec2[0]),float(vec2[1]),float(vec2[2]));
458   }
459   else // Automatic range determination
460   {
461    // Initially the point with Z=0 in the Helix frame is taken as a starting point.
462    // In case this point can't be reached, the point in which the particle velocity
463    // was defined is taken as the starting point.
464    // The endpoint is initially obtained by applying the tofmax from the start point.
465    tmin=0;
466    if (fabs(fVz)>0) tmin=-fZ0/fVz;
467    v=beta*c;
468    r1=rx;
469    r=v*tmin;
470    r1=r1+r;
471
472    // Override the initial begin and endpoint settings by the track data
473    if (rbeg)
474    {
475     r1=(Ali3Vector)(*rbeg); 
476     scale1=rbeg->GetUnitScale();
477     // All coordinates in the selected unit scale
478     if (scale1/scale>1.1 || scale/scale1>1.1) r1*=scale1/scale;
479    }
480
481    r=v*fTofmax;
482    r2=r1+r;
483    if (rend)
484    {
485     r2=(Ali3Vector)(*rend); 
486     scale2=rend->GetUnitScale();
487     // All coordinates in the selected unit scale
488     if (scale2/scale>1.1 || scale/scale2>1.1) r2*=scale2/scale;
489    }
490    
491    r1.GetVector(vec1,"car");
492    r2.GetVector(vec2,"car");
493
494    // Make the 'curve' in the LAB frame and exit.
495    SetNextPoint(float(vec1[0]),float(vec1[1]),float(vec1[2]));
496    SetNextPoint(float(vec2[0]),float(vec2[1]),float(vec2[2]));
497   }
498  }
499  else
500  {
501   //////////////////////////////////////
502   // Treatment of curved trajectories //
503   //////////////////////////////////////
504
505   // Initialisation of the flight time boundaries.
506   // Based on the constant motion of the particle along the Helix Z-axis,
507   // the parametrisation z(t)=z0+fVz*t in the Helix frame is used. 
508   // If possible the point with Z=0 in the Helix frame is taken as a starting point.
509   // In case this point can't be reached, the point in which the particle velocity
510   // was defined is taken as the starting point.
511   tmin=0;
512   if (fabs(fVz)>0) tmin=-fZ0/fVz;
513   tmax=tmin+fTofmax;
514
515   if (tmax<tmin)
516   {
517    dum=tmin;
518    tmin=tmax;
519    tmax=dum;
520   }
521
522   // Determination of the range in the helix frame
523
524   if (!range) // Automatic range determination
525   {
526    scale1=1;
527    scale2=1;
528    if (rbeg)
529    {
530     r1=rbeg->GetPrimed(fRotMat);
531     scale1=rbeg->GetUnitScale();
532     // All coordinates in the selected unit scale
533     if (scale1/scale>1.1 || scale/scale1>1.1) r1*=scale1/scale;
534     // Re-calculate the tmin for this new starting point
535     r1.GetVector(vec1,"car");
536     if (fabs(fVz)>0) tmin=(vec1[2]-fZ0)/fVz;
537     tmax=tmin+fTofmax;
538    }
539    if (rend)
540    {
541     r2=rend->GetPrimed(fRotMat);
542     scale2=rend->GetUnitScale();
543     // All coordinates in the selected unit scale
544     if (scale2/scale>1.1 || scale/scale2>1.1) r2*=scale2/scale;
545     r2.GetVector(vec2,"car");
546     if (fabs(fVz)>0) tmax=(vec2[2]-fZ0)/fVz;
547    }
548    // Make the curve on basis of the flight time boundaries and exit
549    if (tmax<tmin)
550    {
551     dum=tmin;
552     tmin=tmax;
553     tmax=dum;
554    }
555    SetRange(tmin,tmax,kHelixT);
556   }
557   else // User explicitly specified range
558   {
559    vec1[abs(iaxis)-1]=range[0];
560    vec2[abs(iaxis)-1]=range[1];
561    r1.SetVector(vec1,"car");
562    r2.SetVector(vec2,"car");
563    if (iaxis>0) // Range specified in LAB frame
564    {
565     r1=r1.GetPrimed(fRotMat);
566     r1.GetVector(vec1,"car");
567     r2=r2.GetPrimed(fRotMat);
568     r2.GetVector(vec2,"car");
569    } 
570    // Determination of the axis component with the
571    // largest range difference
572    Double_t dmax=0;
573    Int_t imax=0;
574    Double_t test=0;
575    for (Int_t i=0; i<3; i++)
576    {
577     test=fabs(vec1[i]-vec2[i]);
578     if (test>dmax)
579     {
580      dmax=test;
581      imax=i;
582     }
583    }
584
585    Double_t rmin=vec1[imax];
586    Double_t rmax=vec2[imax];
587    if (rmax<rmin)
588    {
589     dum=rmin;
590     rmin=rmax;
591     rmax=dum;
592    }
593
594    // The kinematic range boundaries in the helix frame
595    Double_t xmin=fX0-fVt/fW;
596    Double_t xmax=fX0+fVt/fW;
597    Double_t ymin=fY0-fVt/fW;
598    Double_t ymax=fY0+fVt/fW;
599
600    if (xmax<xmin)
601    {
602     dum=xmin;
603     xmin=xmax;
604     xmax=dum;
605    }
606    if (ymax<ymin)
607    {
608     dum=ymin;
609     ymin=ymax;
610     ymax=dum;
611    }
612
613    // Set the range for the helix
614    if (imax==2 && dmax>0) SetRange(rmin,rmax,kHelixZ);
615    if (imax==1)
616    {
617     // Limit range to kinematic boundaries if needed
618     if (rmin<=ymin) rmin=ymin+1e-6*dmax;
619     if (rmax>=ymax) rmax=ymax-1e-6*dmax;
620     if (rmin<rmax) SetRange(rmin,rmax,kHelixY);
621    }
622    if (imax==0)
623    {
624     // Limit range to kinematic boundaries if needed
625     if (rmin<=xmin) rmin=xmin+1e-6*dmax;
626     if (rmax>=xmax) rmax=xmax-1e-6*dmax;
627     if (rmin<rmax) SetRange(rmin,rmax,kHelixX);
628    }
629   }
630  }
631  return;
632 }
633 ///////////////////////////////////////////////////////////////////////////
634 void AliHelix::Display(AliTrack* t,Double_t* range,Int_t iaxis,Double_t scale)
635 {
636 // Display the helix curve of an AliTrack.
637 // Various curves can be displayed together or individually; please refer to
638 // the memberfunction Refresh() for further details.
639 // It is assumed that the track charge is stored in elementary units
640 // (i.e. charge=1 for a proton).
641 // The input argument "scale" specifies the unit scale for the various
642 // locations where scale=0.01 indicates unit scales in cm etc...
643 // In case scale<=0, the unit scale for locations is determined from the
644 // begin, reference or endpoint of the track. If neither of these
645 // positions is present, all locations are assumed to be given in meter.
646 // The lower and upper bounds for the range are specified by range[0] and
647 // range[1] and the argument "iaxis" indicates along which axis this range
648 // is specified.
649 // The range can be specified either in the LAB frame or in the Helix frame.
650 // The latter is the frame in which the Z axis points in the B direction.
651 //
652 // The conventions for the "iaxis" argument are the following :
653 // iaxis = 1 ==> X axis in the LAB frame
654 //         2 ==> Y axis in the LAB frame
655 //         3 ==> Z axis in the LAB frame
656 //        -1 ==> X axis in the Helix frame
657 //        -2 ==> Y axis in the Helix frame
658 //        -3 ==> Z axis in the Helix frame
659 //
660 // In case range=0 the begin/end/reference points of the AliTrack and the
661 // maximum time of flight (see the SetTofmax() memberfunction) will be used
662 // and an appropriate choice for the iaxis parameter will be made automatically
663 // based on the track kinematics.
664 // In case the reference point is not present, the begin or endpoint will be used
665 // as reference point for the 3-momentum specification. If neither of these positions
666 // is present, (0,0,0) will be taken as the reference point.
667 // 
668 // The default values are range=0, iaxis=3 and scale=-1.
669 //
670 // Note :
671 // ------
672 // Before any display activity, a TCanvas and a TView have to be initiated
673 // first by the user like for instance
674 // 
675 // TCanvas* c1=new TCanvas("c1","c1");
676 // TView* view=new TView(1);
677 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
678 // view->ShowAxis();
679 //
680 // The user can also use the 3D viewing facilities from the TCanvas menu
681 // to open an appropriate view. 
682
683  if (!t || (range && !iaxis)) return;
684
685  MakeCurve(t,range,iaxis,scale);
686
687  if (fRefresh>0) Refresh(fRefresh);
688
689  Int_t np=GetLastPoint()+1;
690  if (!np) return;
691
692  Float_t* points=GetP();
693  TPolyLine3D* curve=new TPolyLine3D(np,points);
694
695  curve->SetLineWidth(fLineWidth);
696  if (fLineColor<0)
697  {
698   Float_t q=t->GetCharge();
699   curve->SetLineColor(kGreen);
700   if (q>0) curve->SetLineColor(kRed);
701   if (q<0) curve->SetLineColor(kBlue);
702  }
703  else
704  {
705   curve->SetLineColor(fLineColor);
706  }
707  curve->Draw();
708
709  if (!fCurves)
710  {
711   fCurves=new TObjArray();
712   fCurves->SetOwner();
713  }
714  fCurves->Add(curve);
715
716  // Display the marker for the track starting point
717  if (fMstyle>0)
718  {
719   TPolyMarker3D* m=new TPolyMarker3D();
720   m->SetPoint(0,points[0],points[1],points[2]);
721   m->SetMarkerStyle(fMstyle);
722   m->SetMarkerSize(fMsize);
723   Int_t col=curve->GetLineColor();
724   if (fMcol>0) col=fMcol;
725   m->SetMarkerColor(col);
726   m->Draw();
727   fCurves->Add(m);
728  }
729 }
730 ///////////////////////////////////////////////////////////////////////////
731 void AliHelix::Refresh(Int_t mode)
732 {
733 // Refresh the display screen before showing the next curve.
734 //
735 // mode = 0 : refreshing fully under user control.
736 //        1 : the display screen will be refreshed automatically
737 //            at each individual track display.
738 //       -1 : the display screen will be refreshed automatically
739 //            at each event display.
740 //
741 // The default is mode=0.
742
743  if (abs(mode)<2) fRefresh=mode;
744  if (fCurves) fCurves->Delete();
745 }
746 ///////////////////////////////////////////////////////////////////////////
747 void AliHelix::Display(AliEvent* evt,Double_t* range,Int_t iaxis,Double_t scale)
748 {
749 // Display the helix curves of all tracks of the specified event.
750 // Various events can be displayed together or individually; please refer to
751 // the memberfunction Refresh() for further details.
752 // Please refer to the track display memberfunction for further details
753 // on the input arguments.
754 // 
755 // The default values are range=0, iaxis=3 and scale=-1.
756 //
757 // Note :
758 // ------
759 // Before any display activity, a TCanvas and a TView have to be initiated
760 // first by the user like for instance
761 // 
762 // TCanvas* c1=new TCanvas("c1","c1");
763 // TView* view=new TView(1);
764 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
765 // view->ShowAxis();
766 //
767 // The user can also use the 3D viewing facilities from the TCanvas menu
768 // to open an appropriate view. 
769
770  if (!evt) return;
771
772  if (fRefresh<0) Refresh(fRefresh);
773
774  Int_t ntk=evt->GetNtracks();
775  for (Int_t jtk=1; jtk<=ntk; jtk++)
776  {
777   AliTrack* tx=evt->GetTrack(jtk);
778   if (tx) Display(tx,range,iaxis,scale);
779  }
780 }
781 ///////////////////////////////////////////////////////////////////////////
782 void AliHelix::Display(TObjArray* arr,Double_t* range,Int_t iaxis,Double_t scale)
783 {
784 // Display the helix curves of all tracks in the specified array.
785 // A convenient way to obtain an array with selected tracks from e.g. an AliEvent
786 // is to make use of its GetTracks() selection facility.
787 // Various arrays can be displayed together or individually; please refer to
788 // the memberfunction Refresh() for further details.
789 // Please refer to the track display memberfunction for further details
790 // on the input arguments.
791 // 
792 // The default values are range=0, iaxis=3 and scale=-1.
793 //
794 // Note :
795 // ------
796 // Before any display activity, a TCanvas and a TView have to be initiated
797 // first by the user like for instance
798 // 
799 // TCanvas* c1=new TCanvas("c1","c1");
800 // TView* view=new TView(1);
801 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
802 // view->ShowAxis();
803 //
804 // The user can also use the 3D viewing facilities from the TCanvas menu
805 // to open an appropriate view. 
806
807  if (!arr) return;
808
809  Int_t ntk=arr->GetEntries();
810  for (Int_t jtk=0; jtk<ntk; jtk++)
811  {
812   TObject* obj=arr->At(jtk);
813   if (!obj) continue;
814   if (!(obj->InheritsFrom("AliTrack"))) continue;
815   AliTrack* tx=(AliTrack*)obj;
816   Display(tx,range,iaxis,scale);
817  }
818 }
819 ///////////////////////////////////////////////////////////////////////////
820 AliPosition* AliHelix::Extrapolate(AliTrack* t,Double_t* pars,Double_t scale)
821 {
822 // Extrapolate an AliTrack according to the corresponding helix curve
823 // and provide a pointer to the impact position w.r.t. a specified plane.
824 // In case the track can never reach the specified plane, the returned
825 // position pointer is zero.
826 // Detailed information of all the helix points used in the extrapolation
827 // can be obtained via the GetN() and GetP() memberfunctions of TPolyLine3D.
828 // It is assumed that the track charge is stored in elementary units
829 // (i.e. charge=1 for a proton).
830 // The input argument "scale" specifies the unit scale for the various
831 // locations where scale=0.01 indicates unit scales in cm etc...
832 // In case scale<=0, the unit scale for locations is determined from the
833 // begin, reference or endpoint of the track. If neither of these
834 // positions is present, all locations are assumed to be given in meter.
835 // The extrapolation parameters for the impact plane and required accuracy
836 // are specified by pars[0], pars[1] and pars[2], respectively.
837 // pars[0] = coordinate value of the plane for the impact point
838 // pars[1] = required accuracy on the specified impact plane coordinate
839 // pars[2] = the axis along which the value of par[0] is specified
840 //
841 // The parameters can be specified either w.r.t. the LAB frame or the Helix frame.
842 // The latter is the frame in which the Z axis points in the B direction.
843 //
844 // The conventions for the par[2] argument are the following :
845 // par[2] = 1 ==> X axis in the LAB frame
846 //          2 ==> Y axis in the LAB frame
847 //          3 ==> Z axis in the LAB frame
848 //         -1 ==> X axis in the Helix frame
849 //         -2 ==> Y axis in the Helix frame
850 //         -3 ==> Z axis in the Helix frame
851 //
852 // Example :
853 // ---------
854 // To obtain an extrapolation to the plane Z=0 in the LAB frame
855 // with an accuracy of 0.001 cm the input arguments would be
856 // pars[0]=0  pars[1]=0.001  pars[2]=3  scale=0.01
857 //
858 // Note : The default value for the scale is -1.
859
860  if (fExt)
861  {
862   delete fExt;
863   fExt=0;
864  }
865
866  if (!t || !pars) return fExt;
867
868  AliPosition* rbeg=t->GetBeginPoint();
869  AliPosition* rend=t->GetEndPoint();
870  AliPosition* rref=t->GetReferencePoint();
871
872  // The unit scale for locations if not specified by the user
873  if (scale<=0)
874  {
875   scale=1; // Set default to meter
876   if (rbeg)
877   {
878    scale=rbeg->GetUnitScale();
879   }
880   else if (rend)
881   {
882    scale=rend->GetUnitScale();
883   }
884   else if (rref)
885   {
886    scale=rref->GetUnitScale();
887   }
888  }
889
890  Double_t range[2];
891  range[0]=pars[0]-fabs(pars[1])/2.;
892  range[1]=pars[0]+fabs(pars[1])/2.;
893
894  Int_t iaxis=int(pars[2]);
895
896  MakeCurve(t,range,iaxis,scale);
897
898  Int_t np=GetLastPoint()+1;
899  if (!np) return fExt;
900
901  Float_t* points=GetP();
902
903  // First point of the curve around the impact
904  Int_t ip=0;
905  Float_t first[3]={points[3*ip],points[3*ip+1],points[3*ip+2]};
906
907  // Last point of the curve around the impact
908  ip=GetLastPoint();
909  Float_t last[3]={points[3*ip],points[3*ip+1],points[3*ip+2]};
910
911  // The accuracy on the impact point 
912  Float_t err[3];
913  err[0]=fabs(first[0]-last[0]);
914  err[1]=fabs(first[1]-last[1]);
915  err[2]=fabs(first[2]-last[2]);
916
917  // Take the middle point as impact location
918  ip=np/2;
919  Float_t imp[3]={points[3*ip],points[3*ip+1],points[3*ip+2]};
920
921  fExt=new AliPosition();
922  fExt->SetUnitScale(scale);
923  fExt->SetPosition(imp,"car");
924  fExt->SetPositionErrors(err,"car");
925
926  return fExt;
927 }
928 ///////////////////////////////////////////////////////////////////////////