1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////
20 // Representation and extrapolation of AliTracks in a magnetic field.
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.
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.
31 // To indicate the track starting point, the memberfunction SetMarker()
33 // By default no marker will be displayed.
38 // Display and extrapolation of individual tracks
39 // ----------------------------------------------
48 // r1.SetVector(vec,"car");
53 // p.SetVector(vec,"car");
56 // t.SetBeginPoint(r1);
60 // // The magnetic field vector in Tesla
65 // b.SetVector(vec,"car");
67 // AliHelix* helix=new AliHelix();
69 // helix->SetTofmax(1e-7);
71 // TCanvas* c1=new TCanvas("c1","c1");
72 // TView* view=new TView(1);
73 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
77 // Double_t range[2]={0,600};
78 // helix->Display(&t,range,3);
79 // t.SetCharge(-t.GetCharge());
80 // helix->Display(&t);
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 // ======================================================================
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
93 // cout << " ***" << endl;
94 // cout << " *** AliCollider run for " << nevents << " events." << endl;
95 // cout << " ***" << endl;
97 // AliCollider* gen=new AliCollider();
99 // gen->OpenFortranFile(6,"dump.log");
101 // gen->SetVertexMode(2);
102 // gen->SetResolution(1e-4);
104 // gen->SetRunNumber(jrun);
105 // gen->SetPrintFreq(1);
107 // gen->SetSpectatorPmin(0.01);
114 // gen->Init("fixt",zp,ap,zt,at,158);
116 // AliHelix* helix=new AliHelix();
117 // Float_t vec[3]={0,2,0};
119 // b.SetVector(vec,"car");
122 // helix->Refresh(-1); // Refresh display after each event
124 // TCanvas* c1=new TCanvas("c1","c1");
125 // TView* view=new TView(1);
126 // view->SetRange(-200,-200,-200,200,200,200);
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);
136 // for (Int_t i=0; i<nevents; i++)
139 // gen->MakeEvent(npart);
140 // AliEvent* evt=gen->GetEvent();
143 // helix->Display(evt);
145 // gSystem->Sleep(5000); // Some delay to keep the display on screen
148 // ======================================================================
150 //--- Author: Nick van Eijndhoven 17-jun-2004 Utrecht University
151 //- Modified: NvE $Date$ Utrecht University
152 ///////////////////////////////////////////////////////////////////////////
155 #include "AliHelix.h"
156 #include "Riostream.h"
158 ClassImp(AliHelix) // Class implementation to enable ROOT I/O
160 AliHelix::AliHelix() : THelix()
162 // Default constructor
175 ///////////////////////////////////////////////////////////////////////////
176 AliHelix::~AliHelix()
178 // Destructor to delete dynamically allocated memory.
190 ///////////////////////////////////////////////////////////////////////////
191 AliHelix::AliHelix(const AliHelix& h) : THelix(h)
202 ///////////////////////////////////////////////////////////////////////////
203 void AliHelix::SetB(Ali3Vector& b)
205 // Set the magnetic field vector in Tesla.
211 fB.GetVector(axis,"car");
215 ///////////////////////////////////////////////////////////////////////////
216 Ali3Vector& AliHelix::GetB()
218 // Provide the magnetic field vector in Tesla.
221 ///////////////////////////////////////////////////////////////////////////
222 void AliHelix::SetTofmax(Float_t tof)
224 // Set the maximum time of flight for straight tracks in seconds.
225 // This maximum tof will be used for drawing etc... in case no begin
226 // and endpoints can be determined from the track info.
229 // 1) In case the user specifies an explicit range, it will override
230 // the maximum tof limit.
231 // 2) By default the tofmax is set to 10 ns in the AliHelix constructor.
234 ///////////////////////////////////////////////////////////////////////////
235 Float_t AliHelix::GetTofmax() const
237 // Provide the maximum time of flight for straight tracks in seconds.
240 ///////////////////////////////////////////////////////////////////////////
241 void AliHelix::SetMarker(Int_t style,Float_t size,Int_t col)
243 // Specify the marker (style, size and colour) to indicate the starting point
244 // of a track in a display.
245 // In case col<0 the marker will have the same color as the track itself.
247 // Defaults are style=8, size=0.2 and col=-1.
253 ///////////////////////////////////////////////////////////////////////////
254 void AliHelix::UseEndPoint(Int_t mode)
256 // Select usage of track endpoint in drawing and extrapolation.
257 // This allows correct event displays even for very long tracks.
259 // mode = 0 : Do not use the track endpoint
260 // 1 : Use the track endpoint
262 // The default value is mode=1 (which is also set in the constructor).
264 if (mode==0 || mode==1) fEnduse=mode;
266 ///////////////////////////////////////////////////////////////////////////
267 void AliHelix::MakeCurve(AliTrack* t,Double_t* range,Int_t iaxis,Double_t scale)
269 // Make the helix curve for the specified AliTrack.
270 // Detailed information of all the helix points can be obtained via the
271 // GetN() and GetP() memberfunctions of TPolyLine3D.
272 // In case one wants to display or extrapolate an AliTrack it is preferable
273 // to use the Display() or Extrapolate() memberfunctions.
274 // It is assumed that the track charge is stored in elementary units
275 // (i.e. charge=1 for a proton).
276 // The input argument "scale" specifies the unit scale for the various
277 // locations where scale=0.01 indicates unit scales in cm etc...
278 // In case scale<=0, the unit scale for locations is determined from the
279 // begin, reference or endpoint of the track. If neither of these
280 // positions is present, all locations are assumed to be given in meter.
281 // The lower and upper bounds for the range are specified by range[0] and
282 // range[1] and the argument "iaxis" indicates along which axis this range
284 // The range can be specified either in the LAB frame or in the Helix frame.
285 // The latter is the frame in which the Z axis points in the B direction.
287 // The conventions for the "iaxis" argument are the following :
288 // iaxis = 1 ==> X axis in the LAB frame
289 // 2 ==> Y axis in the LAB frame
290 // 3 ==> Z axis in the LAB frame
291 // -1 ==> X axis in the Helix frame
292 // -2 ==> Y axis in the Helix frame
293 // -3 ==> Z axis in the Helix frame
295 // In case range=0 the begin/end/reference points of the AliTrack and the
296 // maximum time of flight (see the SetTofmax() memberfunction) will be used
297 // and an appropriate choice for the iaxis parameter will be made automatically
298 // based on the track kinematics.
299 // In case the reference point is not present, the begin or endpoint will be used
300 // as reference point for the 3-momentum specification. If neither of these positions
301 // is present, (0,0,0) will be taken as the reference point.
303 // The default values are range=0, iaxis=3 and scale=-1.
305 SetPolyLine(0); // Reset the polyline data points
307 if (!t || (range && !iaxis)) return;
309 Double_t energy=t->GetEnergy(1); // Track energy in GeV
310 Double_t betanorm=t->GetBeta();
312 if (energy<=0 || betanorm<=0) return;
314 AliPosition* rbeg=t->GetBeginPoint();
316 if (fEnduse) rend=t->GetEndPoint();
317 AliPosition* rref=t->GetReferencePoint();
319 // Magnetic field vector or default Z-direction
320 Double_t bvec[3]={0,0,1};
321 if (fB.GetNorm()>0) fB.GetVector(bvec,"car");
323 // The unit scale for locations if not specified by the user
326 scale=1; // Set default to meter
329 scale=rbeg->GetUnitScale();
333 scale=rend->GetUnitScale();
337 scale=rref->GetUnitScale();
341 Double_t c=2.99792458e8/scale; // Lightspeed in the selected unit scale
343 // The helix angular frequency
344 Double_t w=9e7*(t->GetCharge()*fB.GetNorm())/energy;
346 // The particle velocity in the LAB frame
347 Ali3Vector beta=t->GetBetaVector();
350 v.GetVector(vel,"car");
352 // The particle velocity in the Helix frame
353 Ali3Vector betaprim=beta.GetPrimed(fRotMat);
354 v=v.GetPrimed(fRotMat);
356 v.GetVector(velprim,"car");
358 // Check compatibility of velocity and range specification.
362 if (iaxis>0) beta.GetVector(betavec,"car");
363 if (iaxis<0) betaprim.GetVector(betavec,"car");
364 if (fabs(betavec[abs(iaxis)-1])/betanorm<1e-10) return;
367 // The LAB location in which the velocity of the particle is defined
368 Double_t loc[3]={0,0,0};
373 rx=(Ali3Vector)(*rref);
374 scalex=rref->GetUnitScale();
378 rx=(Ali3Vector)(*rbeg);
379 scalex=rbeg->GetUnitScale();
383 rx=(Ali3Vector)(*rend);
384 scalex=rend->GetUnitScale();
387 if (scalex>0 && (scalex/scale>1.1 || scale/scalex>1.1)) rx*=scalex/scale;
388 rx.GetVector(loc,"car");
390 // Initialisation of Helix kinematics
391 SetHelix(loc,vel,w,0,kUnchanged,bvec);
394 if (fabs(w)>0 && fabs(fVt)>0) bend=1;
396 // Flight time boundaries.
397 // The time origin t=0 is chosen to indicate the position in which
398 // the particle velocity was defined.
399 // The total flight time is initialised to the (user specified) tofmax.
400 Double_t tmin=0,tmax=0;
401 Double_t tof=fTofmax;
404 // The trajectory begin and end points
405 Double_t vec1[3]={0,0,0};
406 Double_t vec2[3]={0,0,0};
414 ////////////////////////////////////////
415 // Treatment of straight trajectories //
416 ////////////////////////////////////////
418 if (range) // Specified range allows for exact flight time boundaries
422 tmin=(range[0]-loc[iaxis-1])/vel[iaxis-1];
423 tmax=(range[1]-loc[iaxis-1])/vel[iaxis-1];
430 tmin=(range[0]-loc[abs(iaxis)-1])/velprim[abs(iaxis)-1];
431 tmax=(range[1]-loc[abs(iaxis)-1])/velprim[abs(iaxis)-1];
439 // Make the 'curve' in the LAB frame and exit.
440 // Use the parametrisation : r(t)=r0+t*v
441 // using the range based flight time boundaries.
442 // An additional point in the middle of the trajectory is
443 // generated in view of accuracy in the case of extrapolations.
449 r1.GetVector(vec1,"car");
450 SetNextPoint(float(vec1[0]),float(vec1[1]),float(vec1[2]));
453 r2.GetVector(vec2,"car");
454 SetNextPoint(float(vec2[0]),float(vec2[1]),float(vec2[2]));
457 r2.GetVector(vec2,"car");
458 SetNextPoint(float(vec2[0]),float(vec2[1]),float(vec2[2]));
460 else // Automatic range determination
462 // Initially the point with Z=0 in the Helix frame is taken as a starting point.
463 // In case this point can't be reached, the point in which the particle velocity
464 // was defined is taken as the starting point.
465 // The endpoint is initially obtained by applying the tofmax from the start point.
467 if (fabs(fVz)>0) tmin=-fZ0/fVz;
473 // Override the initial begin and endpoint settings by the track data
476 r1=(Ali3Vector)(*rbeg);
477 scale1=rbeg->GetUnitScale();
478 // All coordinates in the selected unit scale
479 if (scale1/scale>1.1 || scale/scale1>1.1) r1*=scale1/scale;
486 r2=(Ali3Vector)(*rend);
487 scale2=rend->GetUnitScale();
488 // All coordinates in the selected unit scale
489 if (scale2/scale>1.1 || scale/scale2>1.1) r2*=scale2/scale;
492 r1.GetVector(vec1,"car");
493 r2.GetVector(vec2,"car");
495 // Make the 'curve' in the LAB frame and exit.
496 SetNextPoint(float(vec1[0]),float(vec1[1]),float(vec1[2]));
497 SetNextPoint(float(vec2[0]),float(vec2[1]),float(vec2[2]));
502 //////////////////////////////////////
503 // Treatment of curved trajectories //
504 //////////////////////////////////////
506 // Initialisation of the flight time boundaries.
507 // Based on the constant motion of the particle along the Helix Z-axis,
508 // the parametrisation z(t)=z0+fVz*t in the Helix frame is used.
509 // If possible the point with Z=0 in the Helix frame is taken as a starting point.
510 // In case this point can't be reached, the point in which the particle velocity
511 // was defined is taken as the starting point.
513 if (fabs(fVz)>0) tmin=-fZ0/fVz;
523 // Determination of the range in the helix frame
525 if (!range) // Automatic range determination
531 r1=rbeg->GetPrimed(fRotMat);
532 scale1=rbeg->GetUnitScale();
533 // All coordinates in the selected unit scale
534 if (scale1/scale>1.1 || scale/scale1>1.1) r1*=scale1/scale;
535 // Re-calculate the tmin for this new starting point
536 r1.GetVector(vec1,"car");
537 if (fabs(fVz)>0) tmin=(vec1[2]-fZ0)/fVz;
542 r2=rend->GetPrimed(fRotMat);
543 scale2=rend->GetUnitScale();
544 // All coordinates in the selected unit scale
545 if (scale2/scale>1.1 || scale/scale2>1.1) r2*=scale2/scale;
546 r2.GetVector(vec2,"car");
547 if (fabs(fVz)>0) tmax=(vec2[2]-fZ0)/fVz;
549 // Make the curve on basis of the flight time boundaries and exit
556 SetRange(tmin,tmax,kHelixT);
558 else // User explicitly specified range
560 vec1[abs(iaxis)-1]=range[0];
561 vec2[abs(iaxis)-1]=range[1];
562 r1.SetVector(vec1,"car");
563 r2.SetVector(vec2,"car");
564 if (iaxis>0) // Range specified in LAB frame
566 r1=r1.GetPrimed(fRotMat);
567 r1.GetVector(vec1,"car");
568 r2=r2.GetPrimed(fRotMat);
569 r2.GetVector(vec2,"car");
571 // Determination of the axis component with the
572 // largest range difference
576 for (Int_t i=0; i<3; i++)
578 test=fabs(vec1[i]-vec2[i]);
586 Double_t rmin=vec1[imax];
587 Double_t rmax=vec2[imax];
595 // The kinematic range boundaries in the helix frame
596 Double_t xmin=fX0-fVt/fW;
597 Double_t xmax=fX0+fVt/fW;
598 Double_t ymin=fY0-fVt/fW;
599 Double_t ymax=fY0+fVt/fW;
614 // Set the range for the helix
615 if (imax==2 && dmax>0) SetRange(rmin,rmax,kHelixZ);
618 // Limit range to kinematic boundaries if needed
619 if (rmin<=ymin) rmin=ymin+1e-6*dmax;
620 if (rmax>=ymax) rmax=ymax-1e-6*dmax;
621 if (rmin<rmax) SetRange(rmin,rmax,kHelixY);
625 // Limit range to kinematic boundaries if needed
626 if (rmin<=xmin) rmin=xmin+1e-6*dmax;
627 if (rmax>=xmax) rmax=xmax-1e-6*dmax;
628 if (rmin<rmax) SetRange(rmin,rmax,kHelixX);
634 ///////////////////////////////////////////////////////////////////////////
635 void AliHelix::Display(AliTrack* t,Double_t* range,Int_t iaxis,Double_t scale)
637 // Display the helix curve of an AliTrack.
638 // Various curves can be displayed together or individually; please refer to
639 // the memberfunction Refresh() for further details.
640 // It is assumed that the track charge is stored in elementary units
641 // (i.e. charge=1 for a proton).
642 // The input argument "scale" specifies the unit scale for the various
643 // locations where scale=0.01 indicates unit scales in cm etc...
644 // In case scale<=0, the unit scale for locations is determined from the
645 // begin, reference or endpoint of the track. If neither of these
646 // positions is present, all locations are assumed to be given in meter.
647 // The lower and upper bounds for the range are specified by range[0] and
648 // range[1] and the argument "iaxis" indicates along which axis this range
650 // The range can be specified either in the LAB frame or in the Helix frame.
651 // The latter is the frame in which the Z axis points in the B direction.
653 // The conventions for the "iaxis" argument are the following :
654 // iaxis = 1 ==> X axis in the LAB frame
655 // 2 ==> Y axis in the LAB frame
656 // 3 ==> Z axis in the LAB frame
657 // -1 ==> X axis in the Helix frame
658 // -2 ==> Y axis in the Helix frame
659 // -3 ==> Z axis in the Helix frame
661 // In case range=0 the begin/end/reference points of the AliTrack and the
662 // maximum time of flight (see the SetTofmax() memberfunction) will be used
663 // and an appropriate choice for the iaxis parameter will be made automatically
664 // based on the track kinematics.
665 // In case the reference point is not present, the begin or endpoint will be used
666 // as reference point for the 3-momentum specification. If neither of these positions
667 // is present, (0,0,0) will be taken as the reference point.
669 // The default values are range=0, iaxis=3 and scale=-1.
673 // Before any display activity, a TCanvas and a TView have to be initiated
674 // first by the user like for instance
676 // TCanvas* c1=new TCanvas("c1","c1");
677 // TView* view=new TView(1);
678 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
681 // The user can also use the 3D viewing facilities from the TCanvas menu
682 // to open an appropriate view.
684 if (!t || (range && !iaxis)) return;
686 MakeCurve(t,range,iaxis,scale);
688 if (fRefresh>0) Refresh(fRefresh);
690 Int_t np=GetLastPoint()+1;
693 Float_t* points=GetP();
694 TPolyLine3D* curve=new TPolyLine3D(np,points);
696 curve->SetLineWidth(fLineWidth);
699 Float_t q=t->GetCharge();
700 curve->SetLineColor(kGreen);
701 if (q>0) curve->SetLineColor(kRed);
702 if (q<0) curve->SetLineColor(kBlue);
706 curve->SetLineColor(fLineColor);
712 fCurves=new TObjArray();
717 // Display the marker for the track starting point
720 TPolyMarker3D* m=new TPolyMarker3D();
721 m->SetPoint(0,points[0],points[1],points[2]);
722 m->SetMarkerStyle(fMstyle);
723 m->SetMarkerSize(fMsize);
724 Int_t col=curve->GetLineColor();
725 if (fMcol>0) col=fMcol;
726 m->SetMarkerColor(col);
731 ///////////////////////////////////////////////////////////////////////////
732 void AliHelix::Refresh(Int_t mode)
734 // Refresh the display screen before showing the next curve.
736 // mode = 0 : refreshing fully under user control.
737 // 1 : the display screen will be refreshed automatically
738 // at each individual track display.
739 // -1 : the display screen will be refreshed automatically
740 // at each event display.
742 // The default is mode=0.
744 if (abs(mode)<2) fRefresh=mode;
745 if (fCurves) fCurves->Delete();
747 ///////////////////////////////////////////////////////////////////////////
748 void AliHelix::Display(AliEvent* evt,Double_t* range,Int_t iaxis,Double_t scale)
750 // Display the helix curves of all tracks of the specified event.
751 // Various events can be displayed together or individually; please refer to
752 // the memberfunction Refresh() for further details.
753 // Please refer to the track display memberfunction for further details
754 // on the input arguments.
756 // The default values are range=0, iaxis=3 and scale=-1.
760 // Before any display activity, a TCanvas and a TView have to be initiated
761 // first by the user like for instance
763 // TCanvas* c1=new TCanvas("c1","c1");
764 // TView* view=new TView(1);
765 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
768 // The user can also use the 3D viewing facilities from the TCanvas menu
769 // to open an appropriate view.
773 if (fRefresh<0) Refresh(fRefresh);
775 Int_t ntk=evt->GetNtracks();
776 for (Int_t jtk=1; jtk<=ntk; jtk++)
778 AliTrack* tx=evt->GetTrack(jtk);
779 if (tx) Display(tx,range,iaxis,scale);
782 ///////////////////////////////////////////////////////////////////////////
783 void AliHelix::Display(TObjArray* arr,Double_t* range,Int_t iaxis,Double_t scale)
785 // Display the helix curves of all tracks in the specified array.
786 // A convenient way to obtain an array with selected tracks from e.g. an AliEvent
787 // is to make use of its GetTracks() selection facility.
788 // Various arrays can be displayed together or individually; please refer to
789 // the memberfunction Refresh() for further details.
790 // Please refer to the track display memberfunction for further details
791 // on the input arguments.
793 // The default values are range=0, iaxis=3 and scale=-1.
797 // Before any display activity, a TCanvas and a TView have to be initiated
798 // first by the user like for instance
800 // TCanvas* c1=new TCanvas("c1","c1");
801 // TView* view=new TView(1);
802 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
805 // The user can also use the 3D viewing facilities from the TCanvas menu
806 // to open an appropriate view.
810 Int_t ntk=arr->GetEntries();
811 for (Int_t jtk=0; jtk<ntk; jtk++)
813 TObject* obj=arr->At(jtk);
815 if (!(obj->InheritsFrom("AliTrack"))) continue;
816 AliTrack* tx=(AliTrack*)obj;
817 Display(tx,range,iaxis,scale);
820 ///////////////////////////////////////////////////////////////////////////
821 AliPosition* AliHelix::Extrapolate(AliTrack* t,Double_t* pars,Double_t scale)
823 // Extrapolate an AliTrack according to the corresponding helix curve
824 // and provide a pointer to the impact position w.r.t. a specified plane.
825 // In case the track can never reach the specified plane, the returned
826 // position pointer is zero.
827 // Detailed information of all the helix points used in the extrapolation
828 // can be obtained via the GetN() and GetP() memberfunctions of TPolyLine3D.
829 // It is assumed that the track charge is stored in elementary units
830 // (i.e. charge=1 for a proton).
831 // The input argument "scale" specifies the unit scale for the various
832 // locations where scale=0.01 indicates unit scales in cm etc...
833 // In case scale<=0, the unit scale for locations is determined from the
834 // begin, reference or endpoint of the track. If neither of these
835 // positions is present, all locations are assumed to be given in meter.
836 // The extrapolation parameters for the impact plane and required accuracy
837 // are specified by pars[0], pars[1] and pars[2], respectively.
838 // pars[0] = coordinate value of the plane for the impact point
839 // pars[1] = required accuracy on the specified impact plane coordinate
840 // pars[2] = the axis along which the value of par[0] is specified
842 // The parameters can be specified either w.r.t. the LAB frame or the Helix frame.
843 // The latter is the frame in which the Z axis points in the B direction.
845 // The conventions for the par[2] argument are the following :
846 // par[2] = 1 ==> X axis in the LAB frame
847 // 2 ==> Y axis in the LAB frame
848 // 3 ==> Z axis in the LAB frame
849 // -1 ==> X axis in the Helix frame
850 // -2 ==> Y axis in the Helix frame
851 // -3 ==> Z axis in the Helix frame
855 // To obtain an extrapolation to the plane Z=0 in the LAB frame
856 // with an accuracy of 0.001 cm the input arguments would be
857 // pars[0]=0 pars[1]=0.001 pars[2]=3 scale=0.01
859 // Note : The default value for the scale is -1.
867 if (!t || !pars) return fExt;
869 AliPosition* rbeg=t->GetBeginPoint();
870 AliPosition* rend=t->GetEndPoint();
871 AliPosition* rref=t->GetReferencePoint();
873 // The unit scale for locations if not specified by the user
876 scale=1; // Set default to meter
879 scale=rbeg->GetUnitScale();
883 scale=rend->GetUnitScale();
887 scale=rref->GetUnitScale();
892 range[0]=pars[0]-fabs(pars[1])/2.;
893 range[1]=pars[0]+fabs(pars[1])/2.;
895 Int_t iaxis=int(pars[2]);
897 MakeCurve(t,range,iaxis,scale);
899 Int_t np=GetLastPoint()+1;
900 if (!np) return fExt;
902 Float_t* points=GetP();
904 // First point of the curve around the impact
906 Float_t first[3]={points[3*ip],points[3*ip+1],points[3*ip+2]};
908 // Last point of the curve around the impact
910 Float_t last[3]={points[3*ip],points[3*ip+1],points[3*ip+2]};
912 // The accuracy on the impact point
914 err[0]=fabs(first[0]-last[0]);
915 err[1]=fabs(first[1]-last[1]);
916 err[2]=fabs(first[2]-last[2]);
918 // Take the middle point as impact location
920 Float_t imp[3]={points[3*ip],points[3*ip+1],points[3*ip+2]};
922 fExt=new AliPosition();
923 fExt->SetUnitScale(scale);
924 fExt->SetPosition(imp,"car");
925 fExt->SetPositionErrors(err,"car");
929 ///////////////////////////////////////////////////////////////////////////