2 // Main authors: Matevz Tadel & Alja Mrak-Tadel & Bogdan Vulpescu: 2006, 2007
4 /**************************************************************************
5 * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6 * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for *
7 * full copyright notice. *
8 **************************************************************************/
10 #include "AliEveMUONTrack.h"
12 #include <EveBase/AliEveEventManager.h>
15 #include <AliMagFMaps.h>
17 #include <AliESDMuonTrack.h>
18 #include <AliESDEvent.h>
19 #include <AliESDVertex.h>
20 #include <AliRunLoader.h>
23 #include <AliMUONTrack.h>
24 #include <AliMUONTriggerTrack.h>
25 #include <AliMUONTrackParam.h>
26 #include <AliMUONConstants.h>
27 #include <AliMUONESDInterface.h>
28 #include <AliMUONVCluster.h>
30 #include <TClonesArray.h>
33 #include <TParticle.h>
34 #include <TParticlePDG.h>
36 #include <Riostream.h>
38 //==============================================================================
40 //==============================================================================
41 //==============================================================================
43 //______________________________________________________________________________
44 // Produce TEveUtil:TEveTrack from AliMUONTrack with dipole field model
46 ClassImp(AliEveMUONTrack)
48 AliMagF* AliEveMUONTrack::fgFieldMap = 0;
50 //______________________________________________________________________________
51 AliEveMUONTrack::AliEveMUONTrack(TEveRecTrack* t, TEveTrackPropagator* rs) :
57 fIsMUONTriggerTrack(kFALSE),
67 fgFieldMap = AliEveEventManager::AssertMagField();
70 //______________________________________________________________________________
71 AliEveMUONTrack::~AliEveMUONTrack()
77 if (fIsRefTrack || fIsESDTrack) delete fTrack;
78 if (fIsMCTrack) delete fPart;
82 //______________________________________________________________________________
83 void AliEveMUONTrack::PrintMCTrackInfo()
86 // information about the MC particle
92 cout << " ! no particle ..." << endl;
97 cout << " MC track parameters at vertex" << endl;
98 cout << " -------------------------------------------------------------------------------------" << endl;
99 cout << " PDG code Vx Vy Vz Px Py Pz " << endl;
102 setw(8) << setprecision(0) <<
103 fPart->GetPdgCode() << " " <<
104 setw(8) << setprecision(3) <<
105 fPart->Vx() << " " <<
106 setw(8) << setprecision(3) <<
107 fPart->Vy() << " " <<
108 setw(8) << setprecision(3) <<
109 fPart->Vz() << " " <<
110 setw(8) << setprecision(3) <<
111 fPart->Px() << " " <<
112 setw(8) << setprecision(3) <<
113 fPart->Py() << " " <<
114 setw(8) << setprecision(4) <<
115 fPart->Pz() << " " <<
119 pt = TMath::Sqrt(fPart->Px()*fPart->Px()+fPart->Py()*fPart->Py());
120 p = TMath::Sqrt(fPart->Px()*fPart->Px()+fPart->Py()*fPart->Py()+fPart->Pz()*fPart->Pz());
124 setw(8) << setprecision(3) <<
125 pt << " GeV/c" << endl;
128 setw(8) << setprecision(4) <<
129 p << " GeV/c" << endl;
133 //______________________________________________________________________________
134 void AliEveMUONTrack::PrintMUONTrackInfo()
137 // information about the reconstructed/reference track; at hits and at vertex
140 Double_t radDeg = 180.0/TMath::Pi();
143 Float_t pt, bc, nbc, zc;
144 AliMUONTrackParam *mtp;
145 TClonesArray *trackParamAtCluster;
148 cout << " ! no reconstructed track ..." << endl;
154 cout << " TEveTrack number " << fLabel << endl;
155 cout << " ---------------------------------------------------------------------------------------------------------------------------------" << endl;
157 cout << " Number of clusters " << fTrack->GetNClusters() << endl;
158 cout << " Match to trigger " << fTrack->GetMatchTrigger() << endl;
159 if (fTrack->GetMatchTrigger()) {
160 cout << " Chi2 tracking-trigger " << fTrack->GetChi2MatchTrigger() << endl;
161 cout << " Local trigger number " << fTrack->GetLoTrgNum() << endl;
167 cout << " TEveTrack reference number " << fLabel << endl;
168 cout << " ---------------------------------------------------------------------------------------------------------------------------------" << endl;
170 cout << " Number of clusters " << fTrack->GetNClusters() << endl;
173 trackParamAtCluster = fTrack->GetTrackParamAtCluster();
174 nparam = trackParamAtCluster->GetEntries();
177 cout << " trackParamAtCluster entries " << nparam << "" << endl;
178 cout << " ---------------------------------------------------------------------------------------------------------------------------------" << endl;
179 cout << " Number InvBendMom BendSlope NonBendSlope BendCoord NonBendCoord Z Px Py Pz P" << endl;
181 for (Int_t i = 0; i < nparam; i++) {
183 mtp = (AliMUONTrackParam*)trackParamAtCluster->At(i);
186 setw(9)<< setprecision(3) <<
189 setw(8) << setprecision(3) <<
190 mtp->GetInverseBendingMomentum() << " " <<
192 setw(8) << setprecision(3) <<
193 mtp->GetBendingSlope()*radDeg << " " <<
195 setw(8) << setprecision(3) <<
196 mtp->GetNonBendingSlope()*radDeg << " " <<
198 setw(8) << setprecision(4) <<
199 mtp->GetBendingCoor() << " " <<
201 setw(8) << setprecision(4) <<
202 mtp->GetNonBendingCoor() << " " <<
204 setw(10) << setprecision(6) <<
205 mtp->GetZ() << " " <<
207 setw(8) << setprecision(4) <<
210 setw(8) << setprecision(4) <<
213 setw(8) << setprecision(4) <<
216 setw(8) << setprecision(4) <<
224 cout << " TEveTrack parameters at vertex" << endl;
225 cout << " --------------------------------------------------------------------------------------------------------------------" << endl;
226 cout << " InvBendMom BendSlope NonBendSlope BendCoord NonBendCoord Z Px Py Pz P" << endl;
228 mtp = (AliMUONTrackParam*)fTrack->GetTrackParamAtVertex();
230 bc = mtp->GetBendingCoor();
231 nbc = mtp->GetNonBendingCoor();
233 if (bc < 0.001) bc = 0.0;
234 if (nbc < 0.001) nbc = 0.0;
235 if (zc < 0.001) zc = 0.0;
238 setw(8) << setprecision(3) <<
239 mtp->GetInverseBendingMomentum() << " " <<
241 setw(8) << setprecision(3) <<
242 mtp->GetBendingSlope()*radDeg << " " <<
244 setw(8) << setprecision(3) <<
245 mtp->GetNonBendingSlope()*radDeg << " " <<
247 setw(8) << setprecision(4) <<
250 setw(8) << setprecision(4) <<
253 setw(10) << setprecision(6) <<
256 setw(8) << setprecision(4) <<
259 setw(8) << setprecision(4) <<
262 setw(8) << setprecision(4) <<
265 setw(8) << setprecision(4) <<
270 pt = TMath::Sqrt(mtp->Px()*mtp->Px()+mtp->Py()*mtp->Py());
274 setw(8) << setprecision(3) <<
275 pt << " GeV/c" << endl;
279 //______________________________________________________________________________
280 void AliEveMUONTrack::PrintMUONTriggerTrackInfo()
283 // information about the trigger track
286 // Double_t radDeg = 180.0/TMath::Pi();
290 //______________________________________________________________________________
291 void AliEveMUONTrack::PrintESDTrackInfo()
294 // information about the reconstructed ESD track at vertex
297 Double_t radDeg = 180.0/TMath::Pi();
300 AliMUONTrackParam *mtp = (AliMUONTrackParam*)fTrack->GetTrackParamAtVertex();
303 cout << " ESD muon track " << endl;
304 cout << " -----------------------------------------------------------------------------------------------------------" << endl;
305 cout << " InvBendMom BendSlope NonBendSlope BendCoord NonBendCoord Z Px Py Pz" << endl;
309 setw(8) << setprecision(4) <<
310 mtp->GetInverseBendingMomentum() << " " <<
312 setw(8) << setprecision(3) <<
313 mtp->GetBendingSlope()*radDeg << " " <<
315 setw(8) << setprecision(3) <<
316 mtp->GetNonBendingSlope()*radDeg << " " <<
318 setw(8) << setprecision(4) <<
319 mtp->GetBendingCoor() << " " <<
321 setw(8) << setprecision(4) <<
322 mtp->GetNonBendingCoor() << " " <<
324 setw(10) << setprecision(6) <<
325 mtp->GetZ() << " " <<
327 setw(8) << setprecision(3) <<
330 setw(8) << setprecision(3) <<
333 setw(8) << setprecision(3) <<
338 pt = TMath::Sqrt(mtp->Px()*mtp->Px()+mtp->Py()*mtp->Py());
342 setw(8) << setprecision(3) <<
343 pt << " GeV/c" << endl;
346 setw(8) << setprecision(4) <<
347 mtp->P() << " GeV/c" << endl;
349 AliESDEvent* esd = AliEveEventManager::AssertESD();
351 Double_t spdVertexX = 0;
352 Double_t spdVertexY = 0;
353 Double_t spdVertexZ = 0;
354 Double_t esdVertexX = 0;
355 Double_t esdVertexY = 0;
356 Double_t esdVertexZ = 0;
358 AliESDVertex* spdVertex = (AliESDVertex*) esd->GetVertex();
359 if (spdVertex->GetNContributors()) {
360 spdVertexZ = spdVertex->GetZv();
361 spdVertexY = spdVertex->GetYv();
362 spdVertexX = spdVertex->GetXv();
365 AliESDVertex* esdVertex = (AliESDVertex*) esd->GetPrimaryVertex();
366 if (esdVertex->GetNContributors()) {
367 esdVertexZ = esdVertex->GetZv();
368 esdVertexY = esdVertex->GetYv();
369 esdVertexX = esdVertex->GetXv();
372 Float_t t0v = esd->GetT0zVertex();
376 cout << "External vertex SPD: " <<
378 spdVertex->GetNContributors() << " " <<
379 setw(8) << setprecision(3) <<
382 spdVertexZ << " " << endl;
383 cout << "External vertex ESD: " <<
385 esdVertex->GetNContributors() << " " <<
386 setw(8) << setprecision(3) <<
389 esdVertexZ << " " << endl;
390 cout << "External vertex T0: " <<
391 setw(8) << setprecision(3) <<
396 //______________________________________________________________________________
397 void AliEveMUONTrack::MUONTrackInfo()
407 if (fIsMUONTrack || fIsRefTrack) {
408 PrintMUONTrackInfo();
415 if (fIsMUONTriggerTrack) {
416 PrintMUONTriggerTrackInfo();
422 cout << " (slopes [deg], coord [cm], p [GeV/c])" << endl;
426 //______________________________________________________________________________
427 void AliEveMUONTrack::MUONTriggerInfo()
434 TEveUtil::TEveUtil::LoadMacro("MUON_trigger_info.C");
435 gROOT->ProcessLine(Form("MUON_trigger_info(%d);", fLabel));
438 cout << "This is a reference track!" << endl;
441 cout << "This is a Monte-Carlo track!" << endl;
445 AliESDEvent* esd = AliEveEventManager::AssertESD();
446 ULong64_t triggerMask = esd->GetTriggerMask();
449 cout << ">>>>>#########################################################################################################################" << endl;
452 cout << " ESD track trigger info" << endl;
453 cout << " -----------------------------------------------------" << endl;
456 cout << " Match to trigger " << fTrack->GetMatchTrigger() << endl;
458 cout << " ESD trigger mask = " << triggerMask << endl;
461 cout << "#########################################################################################################################<<<<<" << endl;
468 //______________________________________________________________________________
469 void AliEveMUONTrack::MakeMUONTrack(AliMUONTrack *mtrack)
472 // builds the track with dipole field
477 fIsMUONTrack = kTRUE;
480 fTrack = new AliMUONTrack(*mtrack);
485 Float_t ax, bx, ay, by;
486 Float_t xr[28], yr[28], zr[28];
487 Float_t xrc[28], yrc[28], zrc[28];
491 TMatrixD smatrix(2,2);
495 // middle z between the two detector planes of the trigger chambers
496 Float_t zg[4] = { -1603.5, -1620.5, -1703.5, -1720.5 };
499 Float_t pv[3] = {0., 0., 0.};
502 if (mtrack->GetMatchTrigger()) {
503 sprintf(form,"MUONTrack %2d (MT)", fLabel);
505 sprintf(form,"MUONTrack %2d ", fLabel);
511 AliMUONTrackParam *trackParam = 0x0;
512 if (fIsMUONTrack || fIsESDTrack) {
513 trackParam = mtrack->GetTrackParamAtVertex();
514 SetPoint(fCount,trackParam->GetNonBendingCoor(),trackParam->GetBendingCoor(),trackParam->GetZ());
518 for (Int_t i = 0; i < 28; i++) {
519 xr[i]=yr[i]=zr[i]=0.0;
523 Int_t nTrackHits = mtrack->GetNClusters();
524 TClonesArray* trackParamAtCluster = mtrack->GetTrackParamAtCluster();
525 for (Int_t iHit = 0; iHit < nTrackHits; iHit++){
527 trackParam = (AliMUONTrackParam*) trackParamAtCluster->At(iHit);
530 if (IsMUONTrack() || IsESDTrack()) {
531 pt = TMath::Sqrt(trackParam->Px()*trackParam->Px()+trackParam->Py()*trackParam->Py());
532 printf("Set line color = %d \n",ColorIndex(pt));
533 SetLineColor(ColorIndex(pt));
535 pv[0] = trackParam->Px();
536 pv[1] = trackParam->Py();
537 pv[2] = trackParam->Pz();
541 xr[iHit] = trackParam->GetNonBendingCoor();
542 yr[iHit] = trackParam->GetBendingCoor();
543 zr[iHit] = trackParam->GetZ();
544 chr[iHit] = trackParam->GetClusterPtr()->GetChamberId();
548 SetPoint(fCount,xr[0],yr[0],zr[0]);
550 for (Int_t iHit = 1; iHit < nTrackHits; iHit++) {
551 if (chr[iHit] > 3 && chr[iHit-1] < 6) Propagate(xr,yr,zr,iHit-1,iHit);
552 SetPoint(fCount,xr[iHit],yr[iHit],zr[iHit]);
556 if (!fIsMUONTrack && !fIsESDTrack) return;
559 if (mtrack->GetMatchTrigger()) {
561 for (Int_t i = 0; i < nTrackHits; i++) {
562 if (TMath::Abs(zr[i]) > 1000.0) {
563 //printf("TEveHit %d x %f y %f z %f \n",iHit,xr[i],yr[i],zr[i]);
576 for (Int_t i = 0; i < nrc; i++) {
577 xv = (Double_t)zrc[i];
578 yv = (Double_t)xrc[i];
579 //printf("x-z: xv %f yv %f \n",xv,yv);
581 smatrix(1,1) += xv*xv;
587 res = smatrix.Invert() * sums;
594 for (Int_t i = 0; i < nrc; i++) {
595 xv = (Double_t)zrc[i];
596 yv = (Double_t)yrc[i];
597 //printf("y-z: xv %f yv %f \n",xv,yv);
599 smatrix(1,1) += xv*xv;
605 res = smatrix.Invert() * sums;
609 Float_t xtc, ytc, ztc;
610 for (Int_t ii = 0; ii < 4; ii++) {
616 //printf("tc: x %f y %f z %f \n",xtc,ytc,ztc);
618 SetPoint(fCount,xtc,ytc,ztc);
623 } // end match trigger
627 //______________________________________________________________________________
628 void AliEveMUONTrack::MakeMUONTriggerTrack(AliMUONTriggerTrack *mtrack)
631 // builds the trigger track from one point and direction
634 Float_t x1 = mtrack->GetX11();
635 Float_t y1 = mtrack->GetY11();
636 Float_t thex = mtrack->GetThetax();
637 Float_t they = mtrack->GetThetay();
639 Float_t z11 = -1600.0;
640 Float_t z22 = -1724.0;
641 Float_t dz = z22-z11;
643 Float_t x2 = x1 + dz*TMath::Tan(thex);
644 Float_t y2 = y1 + dz*TMath::Tan(they);
646 SetPoint(fCount,x1,y1,z11); fCount++;
647 SetPoint(fCount,x2,y2,z22); fCount++;
651 sprintf(form,"MUONTriggerTrack %2d",mtrack->GetLoTrgNum());
657 //______________________________________________________________________________
658 void AliEveMUONTrack::MakeESDTrack(AliESDMuonTrack *mtrack)
661 // builds the track with dipole field starting from the TParticle
667 if (mtrack->GetMatchTrigger()) {
668 sprintf(form,"ESDTrack %2d (MT)", fLabel);
670 sprintf(form,"ESDTrack %2d ", fLabel);
676 fTrack = new AliMUONTrack();
678 // create a simple track from the ESD track
679 AliMUONESDInterface::ESDToMUON(*mtrack,*fTrack);
681 // reset track parameters at vertex to the ones at DCA
682 AliMUONTrackParam paramAtDCA;
683 AliMUONESDInterface::GetParamAtDCA(*mtrack, paramAtDCA);
684 fTrack->SetTrackParamAtVertex(¶mAtDCA);
686 MakeMUONTrack(fTrack);
690 //______________________________________________________________________________
691 void AliEveMUONTrack::MakeMCTrack(TParticle *part)
694 // builds the track with dipole field starting from the TParticle
699 fPart = new TParticle(*part);
702 sprintf(form,"TEveMCTrack %2d ", fLabel);
707 Double_t vect[7], vout[7];
716 vect[0] = fPart->Vx();
717 vect[1] = fPart->Vy();
718 vect[2] = fPart->Vz();
719 vect[3] = fPart->Px()/fPart->P();
720 vect[4] = fPart->Py()/fPart->P();
721 vect[5] = fPart->Pz()/fPart->P();
722 vect[6] = fPart->P();
724 TParticlePDG *ppdg = fPart->GetPDG(1);
725 Int_t charge = (Int_t)(ppdg->Charge()/3.0);
727 Double_t zMax = -1750.0;
728 Double_t rMax = 350.0;
732 while ((vect[2] > zMax) && (nSteps < 10000) && (r < rMax)) {
734 OneStepRungekutta(charge, step, vect, vout);
735 SetPoint(fCount,vout[0],vout[1],vout[2]);
737 for (Int_t i = 0; i < 7; i++) {
740 r = TMath::Sqrt(vect[0]*vect[0]+vect[1]*vect[1]);
745 //______________________________________________________________________________
746 void AliEveMUONTrack::MakeRefTrack(AliMUONTrack *mtrack)
749 // builds the track with dipole field starting from the TParticle
755 sprintf(form,"RefTrack %2d ", fLabel);
760 MakeMUONTrack(mtrack);
764 //______________________________________________________________________________
765 void AliEveMUONTrack::Propagate(Float_t *xr, Float_t *yr, Float_t *zr, Int_t i1, Int_t i2)
768 // propagate in magnetic field between hits of indices i1 and i2
771 Double_t vect[7], vout[7];
775 AliMUONTrackParam *trackParam = 0;
776 TClonesArray *trackParamAtCluster = 0;
779 zMax = zr[i1]+1.5*step;
781 zMax = zr[i2]+1.5*step;
784 trackParamAtCluster = fTrack->GetTrackParamAtCluster();
786 if (IsMUONTrack() || IsESDTrack()) {
787 trackParam = (AliMUONTrackParam*)trackParamAtCluster->At(i1);
788 charge = (Int_t)trackParam->GetCharge();
791 trackParam = fTrack->GetTrackParamAtVertex();
792 charge = (Int_t)trackParam->GetCharge();
793 trackParam = (AliMUONTrackParam*)trackParamAtCluster->At(i1);
799 vect[6] = trackParam->P();
800 vect[3] = trackParam->Px()/vect[6];
801 vect[4] = trackParam->Py()/vect[6];
802 vect[5] = trackParam->Pz()/vect[6];
805 while ((vect[2] > zMax) && (nSteps < 10000)) {
807 OneStepRungekutta(charge, step, vect, vout);
808 SetPoint(fCount,vout[0],vout[1],vout[2]);
810 for (Int_t i = 0; i < 7; i++) {
817 //______________________________________________________________________________
818 void AliEveMUONTrack::GetField(Double_t *position, Double_t *field)
821 // returns field components at position, for a give field map
824 /// interface for arguments in double precision (Why ? ChF)
827 x[0] = position[0]; x[1] = position[1]; x[2] = position[2];
830 fgFieldMap->Field(x,b);
833 AliWarning("No field map");
834 field[0] = field[1] = field[2] = 0.0;
842 field[0] = b[0]; field[1] = b[1]; field[2] = b[2];
848 //______________________________________________________________________________
849 void AliEveMUONTrack::OneStepRungekutta(Double_t charge, Double_t step,
850 Double_t* vect, Double_t* vout)
852 /// ******************************************************************
854 /// * Runge-Kutta method for tracking a particle through a magnetic *
855 /// * field. Uses Nystroem algorithm (See Handbook Nat. Bur. of *
856 /// * Standards, procedure 25.5.20) *
858 /// * Input parameters *
859 /// * CHARGE Particle charge *
860 /// * STEP Step size *
861 /// * VECT Initial co-ords,direction cosines,momentum *
862 /// * Output parameters *
863 /// * VOUT Output co-ords,direction cosines,momentum *
864 /// * User routine called *
865 /// * CALL GUFLD(X,F) *
867 /// * ==>Called by : <USER>, GUSWIM *
868 /// * Authors R.Brun, M.Hansroul ********* *
869 /// * V.Perevoztchikov (CUT STEP implementation) *
872 /// ******************************************************************
874 Double_t h2, h4, f[4];
875 Double_t xyzt[3], a, b, c, ph,ph2;
876 Double_t secxs[4],secys[4],seczs[4],hxp[3];
877 Double_t g1, g2, g3, g4, g5, g6, ang2, dxt, dyt, dzt;
878 Double_t est, at, bt, ct, cba;
879 Double_t f1, f2, f3, f4, rho, tet, hnorm, hp, rho1, sint, cost;
889 Double_t maxit = 1992;
890 Double_t maxcut = 11;
892 const Double_t kdlt = 1e-4;
893 const Double_t kdlt32 = kdlt/32.;
894 const Double_t kthird = 1./3.;
895 const Double_t khalf = 0.5;
896 const Double_t kec = 2.9979251e-4;
898 const Double_t kpisqua = 9.86960440109;
902 const Int_t kipx = 3;
903 const Int_t kipy = 4;
904 const Int_t kipz = 5;
907 // *. ------------------------------------------------------------------
909 // * this constant is for units cm,gev/c and kgauss
913 for(Int_t j = 0; j < 7; j++)
916 Double_t pinv = kec * charge / vect[6];
924 if (TMath::Abs(h) > TMath::Abs(rest)) h = rest;
925 //cmodif: call gufld(vout,f) changed into:
930 // * start of integration
943 secxs[0] = (b * f[2] - c * f[1]) * ph2;
944 secys[0] = (c * f[0] - a * f[2]) * ph2;
945 seczs[0] = (a * f[1] - b * f[0]) * ph2;
946 ang2 = (secxs[0]*secxs[0] + secys[0]*secys[0] + seczs[0]*seczs[0]);
947 if (ang2 > kpisqua) break;
949 dxt = h2 * a + h4 * secxs[0];
950 dyt = h2 * b + h4 * secys[0];
951 dzt = h2 * c + h4 * seczs[0];
956 // * second intermediate point
959 est = TMath::Abs(dxt) + TMath::Abs(dyt) + TMath::Abs(dzt);
961 if (ncut++ > maxcut) break;
970 //cmodif: call gufld(xyzt,f) changed into:
977 secxs[1] = (bt * f[2] - ct * f[1]) * ph2;
978 secys[1] = (ct * f[0] - at * f[2]) * ph2;
979 seczs[1] = (at * f[1] - bt * f[0]) * ph2;
983 secxs[2] = (bt * f[2] - ct * f[1]) * ph2;
984 secys[2] = (ct * f[0] - at * f[2]) * ph2;
985 seczs[2] = (at * f[1] - bt * f[0]) * ph2;
986 dxt = h * (a + secxs[2]);
987 dyt = h * (b + secys[2]);
988 dzt = h * (c + seczs[2]);
992 at = a + 2.*secxs[2];
993 bt = b + 2.*secys[2];
994 ct = c + 2.*seczs[2];
996 est = TMath::Abs(dxt)+TMath::Abs(dyt)+TMath::Abs(dzt);
997 if (est > 2.*TMath::Abs(h)) {
998 if (ncut++ > maxcut) break;
1007 //cmodif: call gufld(xyzt,f) changed into:
1010 z = z + (c + (seczs[0] + seczs[1] + seczs[2]) * kthird) * h;
1011 y = y + (b + (secys[0] + secys[1] + secys[2]) * kthird) * h;
1012 x = x + (a + (secxs[0] + secxs[1] + secxs[2]) * kthird) * h;
1014 secxs[3] = (bt*f[2] - ct*f[1])* ph2;
1015 secys[3] = (ct*f[0] - at*f[2])* ph2;
1016 seczs[3] = (at*f[1] - bt*f[0])* ph2;
1017 a = a+(secxs[0]+secxs[3]+2. * (secxs[1]+secxs[2])) * kthird;
1018 b = b+(secys[0]+secys[3]+2. * (secys[1]+secys[2])) * kthird;
1019 c = c+(seczs[0]+seczs[3]+2. * (seczs[1]+seczs[2])) * kthird;
1021 est = TMath::Abs(secxs[0]+secxs[3] - (secxs[1]+secxs[2]))
1022 + TMath::Abs(secys[0]+secys[3] - (secys[1]+secys[2]))
1023 + TMath::Abs(seczs[0]+seczs[3] - (seczs[1]+seczs[2]));
1025 if (est > kdlt && TMath::Abs(h) > 1.e-4) {
1026 if (ncut++ > maxcut) break;
1032 // * if too many iterations, go to helix
1033 if (iter++ > maxit) break;
1038 cba = 1./ TMath::Sqrt(a*a + b*b + c*c);
1046 if (step < 0.) rest = -rest;
1047 if (rest < 1.e-5*TMath::Abs(step)) return;
1051 // angle too big, use helix
1056 f4 = TMath::Sqrt(f1*f1+f2*f2+f3*f3);
1065 hxp[0] = f2*vect[kipz] - f3*vect[kipy];
1066 hxp[1] = f3*vect[kipx] - f1*vect[kipz];
1067 hxp[2] = f1*vect[kipy] - f2*vect[kipx];
1069 hp = f1*vect[kipx] + f2*vect[kipy] + f3*vect[kipz];
1072 sint = TMath::Sin(tet);
1073 cost = 2.*TMath::Sin(khalf*tet)*TMath::Sin(khalf*tet);
1077 g3 = (tet-sint) * hp*rho1;
1082 vout[kix] = vect[kix] + g1*vect[kipx] + g2*hxp[0] + g3*f1;
1083 vout[kiy] = vect[kiy] + g1*vect[kipy] + g2*hxp[1] + g3*f2;
1084 vout[kiz] = vect[kiz] + g1*vect[kipz] + g2*hxp[2] + g3*f3;
1086 vout[kipx] = vect[kipx] + g4*vect[kipx] + g5*hxp[0] + g6*f1;
1087 vout[kipy] = vect[kipy] + g4*vect[kipy] + g5*hxp[1] + g6*f2;
1088 vout[kipz] = vect[kipz] + g4*vect[kipz] + g5*hxp[2] + g6*f3;
1093 //______________________________________________________________________________
1094 Int_t AliEveMUONTrack::ColorIndex(Float_t val)
1097 // returns color index in the palette for a give value
1100 Float_t threshold = 0.0;
1101 Float_t maxVal = 10.0;
1103 Float_t div = TMath::Max(1, (Int_t)(maxVal - threshold));
1104 Int_t nCol = gStyle->GetNumberOfColors();
1105 Int_t cBin = (Int_t) TMath::Nint(nCol*(val - threshold)/div);
1107 return gStyle->GetColorPalette(TMath::Min(nCol - 1, cBin));