]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/AliEveMUONTrack.cxx
New code for full visualization of the vzero raw data. ctrl-alt-letf-button shows...
[u/mrichter/AliRoot.git] / EVE / EveDet / AliEveMUONTrack.cxx
1
2 // Main authors: Matevz Tadel & Alja Mrak-Tadel & Bogdan Vulpescu: 2006, 2007
3
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  **************************************************************************/
9
10 #include "AliEveMUONTrack.h"
11
12 #include <EveBase/AliEveEventManager.h>
13
14 #include <AliMagF.h>
15 #include <AliMagFMaps.h>
16 #include <AliLog.h>
17 #include <AliESDMuonTrack.h>
18 #include <AliESDEvent.h>
19 #include <AliESDVertex.h>
20 #include <AliRunLoader.h>
21 #include <AliRun.h>
22
23 #include <AliMUONTrack.h>
24 #include <AliMUONTriggerTrack.h>
25 #include <AliMUONTrackParam.h>
26 #include <AliMUONConstants.h>
27 #include <AliMUONESDInterface.h>
28 #include <AliMUONVCluster.h>
29
30 #include <TClonesArray.h>
31 #include <TStyle.h>
32 #include <TROOT.h>
33 #include <TParticle.h>
34 #include <TParticlePDG.h>
35
36 #include <Riostream.h>
37
38 //==============================================================================
39 // AliEveMUONTrack
40 //==============================================================================
41 //==============================================================================
42
43 //______________________________________________________________________________
44 // Produce TEveUtil:TEveTrack from AliMUONTrack with dipole field model
45
46 ClassImp(AliEveMUONTrack)
47
48 AliMagF* AliEveMUONTrack::fgFieldMap = 0;
49
50 //______________________________________________________________________________
51 AliEveMUONTrack::AliEveMUONTrack(TEveRecTrack* t, TEveTrackPropagator* rs) :
52   TEveTrack(t,rs),
53   fTrack(0),
54   fPart(0),
55   fCount(0),
56   fIsMUONTrack(kFALSE),
57   fIsMUONTriggerTrack(kFALSE),
58   fIsESDTrack(kFALSE),
59   fIsMCTrack(kFALSE),
60   fIsRefTrack(kFALSE)
61 {
62   //
63   // constructor
64   //
65
66   if (fgFieldMap == 0)
67     fgFieldMap = AliEveEventManager::AssertMagField();
68 }
69
70 //______________________________________________________________________________
71 AliEveMUONTrack::~AliEveMUONTrack()
72 {
73   //
74   // destructor
75   //
76
77   if (fIsRefTrack || fIsESDTrack) delete fTrack;
78   if (fIsMCTrack) delete fPart;
79
80 }
81
82 //______________________________________________________________________________
83 void AliEveMUONTrack::PrintMCTrackInfo()
84 {
85   //
86   // information about the MC particle
87   //
88
89   Float_t pt, p;
90
91   if (!fPart) {
92     cout << "   ! no particle ..." << endl;
93     return;
94   }
95
96   cout << endl;
97   cout << "   MC track parameters at vertex" << endl;
98   cout << "   -------------------------------------------------------------------------------------" << endl;
99   cout << "   PDG code          Vx           Vy           Vz           Px           Py           Pz   " << endl;
100
101   cout << "   " <<
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() << "     " <<
116
117     endl;
118
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());
121
122   cout << endl;
123   cout << "   Pt = " <<
124     setw(8) << setprecision(3) <<
125     pt << "  GeV/c" << endl;
126
127   cout << "   P  = " <<
128     setw(8) << setprecision(4) <<
129     p  << "  GeV/c" << endl;
130
131 }
132
133 //______________________________________________________________________________
134 void AliEveMUONTrack::PrintMUONTrackInfo()
135 {
136   //
137   // information about the reconstructed/reference track; at hits and at vertex
138   //
139
140   Double_t radDeg = 180.0/TMath::Pi();
141
142   Int_t nparam;
143   Float_t pt, bc, nbc, zc;
144   AliMUONTrackParam *mtp;
145   TClonesArray *trackParamAtCluster;
146
147   if (!fTrack) {
148     cout << "   ! no reconstructed track ..." << endl;
149     return;
150   }
151
152   if (fIsMUONTrack) {
153     cout << endl;
154     cout << "   TEveTrack number " << fLabel << endl;
155     cout << "   ---------------------------------------------------------------------------------------------------------------------------------" << endl;
156     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;
162     }
163   }
164
165   if (fIsRefTrack) {
166     cout << endl;
167     cout << "   TEveTrack reference number " << fLabel << endl;
168     cout << "   ---------------------------------------------------------------------------------------------------------------------------------" << endl;
169     cout << endl;
170     cout << "   Number of clusters       " << fTrack->GetNClusters() << endl;
171   }
172
173   trackParamAtCluster = fTrack->GetTrackParamAtCluster();
174   nparam = trackParamAtCluster->GetEntries();
175
176   cout << endl;
177   cout << "   trackParamAtCluster entries  " << nparam << "" << endl;
178   cout << "   ---------------------------------------------------------------------------------------------------------------------------------" << endl;
179   cout << "   Number   InvBendMom   BendSlope   NonBendSlope   BendCoord   NonBendCoord               Z        Px        Py        Pz         P" << endl;
180
181   for (Int_t i = 0; i < nparam; i++) {
182
183     mtp = (AliMUONTrackParam*)trackParamAtCluster->At(i);
184
185     cout <<
186       setw(9)<< setprecision(3) <<
187       i << "     " <<
188
189       setw(8) << setprecision(3) <<
190       mtp->GetInverseBendingMomentum() << "    " <<
191
192       setw(8) << setprecision(3) <<
193       mtp->GetBendingSlope()*radDeg << "       " <<
194
195       setw(8) << setprecision(3) <<
196       mtp->GetNonBendingSlope()*radDeg << "    " <<
197
198       setw(8) << setprecision(4) <<
199       mtp->GetBendingCoor() << "       " <<
200
201       setw(8) << setprecision(4) <<
202       mtp->GetNonBendingCoor() << "      " <<
203
204       setw(10) << setprecision(6) <<
205       mtp->GetZ() << "  " <<
206
207       setw(8) << setprecision(4) <<
208       mtp->Px() << "  " <<
209
210       setw(8) << setprecision(4) <<
211       mtp->Py() << "  " <<
212
213       setw(8) << setprecision(4) <<
214       mtp->Pz() << "  " <<
215
216       setw(8) << setprecision(4) <<
217       mtp->P() << "  " <<
218
219       endl;
220
221   }
222
223   cout << endl;
224   cout << "   TEveTrack parameters at vertex" << endl;
225   cout << "   --------------------------------------------------------------------------------------------------------------------" << endl;
226   cout << "   InvBendMom   BendSlope   NonBendSlope   BendCoord   NonBendCoord           Z        Px        Py        Pz         P" << endl;
227
228   mtp = (AliMUONTrackParam*)fTrack->GetTrackParamAtVertex();
229
230   bc  = mtp->GetBendingCoor();
231   nbc = mtp->GetNonBendingCoor();
232   zc  = mtp->GetZ();
233   if (bc  < 0.001) bc  = 0.0;
234   if (nbc < 0.001) nbc = 0.0;
235   if (zc  < 0.001) zc  = 0.0;
236
237   cout << "     " <<
238     setw(8) << setprecision(3) <<
239     mtp->GetInverseBendingMomentum() << "    " <<
240
241     setw(8) << setprecision(3) <<
242     mtp->GetBendingSlope()*radDeg << "       " <<
243
244     setw(8) << setprecision(3) <<
245     mtp->GetNonBendingSlope()*radDeg << "    " <<
246
247     setw(8) << setprecision(4) <<
248     bc << "       " <<
249
250     setw(8) << setprecision(4) <<
251     nbc << "  " <<
252
253     setw(10) << setprecision(6) <<
254     zc << "  " <<
255
256     setw(8) << setprecision(4) <<
257     mtp->Px() << "  " <<
258
259     setw(8) << setprecision(4) <<
260     mtp->Py() << "  " <<
261
262     setw(8) << setprecision(4) <<
263     mtp->Pz() << "  " <<
264
265     setw(8) << setprecision(4) <<
266     mtp->P() << "  " <<
267
268     endl;
269
270   pt = TMath::Sqrt(mtp->Px()*mtp->Px()+mtp->Py()*mtp->Py());
271
272   cout << endl;
273   cout << "   Pt = " <<
274     setw(8) << setprecision(3) <<
275     pt << "  GeV/c" << endl;
276
277 }
278
279 //______________________________________________________________________________
280 void AliEveMUONTrack::PrintMUONTriggerTrackInfo()
281 {
282   //
283   // information about the trigger track
284   //
285
286   // Double_t radDeg = 180.0/TMath::Pi();
287
288 }
289
290 //______________________________________________________________________________
291 void AliEveMUONTrack::PrintESDTrackInfo()
292 {
293   //
294   // information about the reconstructed ESD track at vertex
295   //
296
297   Double_t radDeg = 180.0/TMath::Pi();
298   Float_t pt;
299
300   AliMUONTrackParam *mtp = (AliMUONTrackParam*)fTrack->GetTrackParamAtVertex();
301
302   cout << endl;
303   cout << "   ESD muon track " << endl;
304   cout << "   -----------------------------------------------------------------------------------------------------------" << endl;
305   cout << "   InvBendMom   BendSlope   NonBendSlope    BendCoord   NonBendCoord           Z        Px        Py        Pz" << endl;
306
307   cout << "     " <<
308
309     setw(8) << setprecision(4) <<
310     mtp->GetInverseBendingMomentum() << "    " <<
311
312     setw(8) << setprecision(3) <<
313     mtp->GetBendingSlope()*radDeg << "       " <<
314
315     setw(8) << setprecision(3) <<
316     mtp->GetNonBendingSlope()*radDeg << "     " <<
317
318     setw(8) << setprecision(4) <<
319     mtp->GetBendingCoor() << "       " <<
320
321     setw(8) << setprecision(4) <<
322     mtp->GetNonBendingCoor() << "  " <<
323
324     setw(10) << setprecision(6) <<
325     mtp->GetZ() << "  " <<
326
327     setw(8) << setprecision(3) <<
328     mtp->Px() << "  " <<
329
330     setw(8) << setprecision(3) <<
331     mtp->Py() << "  " <<
332
333     setw(8) << setprecision(3) <<
334     mtp->Pz() << "  " <<
335
336     endl;
337
338   pt = TMath::Sqrt(mtp->Px()*mtp->Px()+mtp->Py()*mtp->Py());
339
340   cout << endl;
341   cout << "   Pt = " <<
342     setw(8) << setprecision(3) <<
343     pt << "  GeV/c" << endl;
344
345   cout << "   P  = " <<
346     setw(8) << setprecision(4) <<
347     mtp->P()  << "  GeV/c" << endl;
348
349   AliESDEvent* esd = AliEveEventManager::AssertESD();
350
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;
357
358   AliESDVertex* spdVertex = (AliESDVertex*) esd->GetVertex();
359   if (spdVertex->GetNContributors()) {
360     spdVertexZ = spdVertex->GetZv();
361     spdVertexY = spdVertex->GetYv();
362     spdVertexX = spdVertex->GetXv();
363   }
364
365   AliESDVertex* esdVertex = (AliESDVertex*) esd->GetPrimaryVertex();
366   if (esdVertex->GetNContributors()) {
367     esdVertexZ = esdVertex->GetZv();
368     esdVertexY = esdVertex->GetYv();
369     esdVertexX = esdVertex->GetXv();
370   }
371
372   Float_t t0v = esd->GetT0zVertex();
373
374   cout << endl;
375   cout << endl;
376   cout << "External vertex SPD: " <<
377     setw(3) <<
378     spdVertex->GetNContributors() << "   " <<
379     setw(8) << setprecision(3) <<
380     spdVertexX << "   " <<
381     spdVertexY << "   " <<
382     spdVertexZ << "   " << endl;
383   cout << "External vertex ESD: " <<
384     setw(3) <<
385     esdVertex->GetNContributors() << "   " <<
386     setw(8) << setprecision(3) <<
387     esdVertexX << "   " <<
388     esdVertexY << "   " <<
389     esdVertexZ << "   " << endl;
390   cout << "External vertex T0: " <<
391     setw(8) << setprecision(3) <<
392     t0v << "   " << endl;
393
394 }
395
396 //______________________________________________________________________________
397 void AliEveMUONTrack::MUONTrackInfo()
398 {
399   //
400   // MENU function
401   //
402
403   if (fIsMCTrack) {
404     PrintMCTrackInfo();
405   }
406
407   if (fIsMUONTrack || fIsRefTrack) {
408     PrintMUONTrackInfo();
409   }
410
411   if (fIsESDTrack) {
412     PrintESDTrackInfo();
413   }
414
415   if (fIsMUONTriggerTrack) {
416     PrintMUONTriggerTrackInfo();
417   }
418
419   cout << endl;
420   cout << endl;
421   cout << endl;
422   cout << "   (slopes [deg], coord [cm], p [GeV/c])" << endl;
423
424 }
425
426 //______________________________________________________________________________
427 void AliEveMUONTrack::MUONTriggerInfo()
428 {
429   //
430   // MENU function
431   //
432
433   if (fIsMUONTrack) {
434     TEveUtil::TEveUtil::LoadMacro("MUON_trigger_info.C");
435     gROOT->ProcessLine(Form("MUON_trigger_info(%d);", fLabel));
436   }
437   if (fIsRefTrack) {
438     cout << "This is a reference track!" << endl;
439   }
440   if (fIsMCTrack) {
441     cout << "This is a Monte-Carlo track!" << endl;
442   }
443   if (fIsESDTrack) {
444
445     AliESDEvent* esd = AliEveEventManager::AssertESD();
446     ULong64_t triggerMask = esd->GetTriggerMask();
447
448     cout << endl;
449     cout << ">>>>>#########################################################################################################################" << endl;
450     cout << endl;
451
452     cout << "   ESD track trigger info" << endl;
453     cout << "   -----------------------------------------------------" << endl;
454     cout << endl;
455
456     cout << "   Match to trigger         " << fTrack->GetMatchTrigger() << endl;
457     cout << endl;
458     cout << "   ESD trigger mask = " << triggerMask << endl;
459
460     cout << endl;
461     cout << "#########################################################################################################################<<<<<" << endl;
462     cout << endl;
463
464   }
465
466 }
467
468 //______________________________________________________________________________
469 void AliEveMUONTrack::MakeMUONTrack(AliMUONTrack *mtrack)
470 {
471   //
472   // builds the track with dipole field
473   //
474
475   if (!fIsESDTrack) {
476     if (!fIsRefTrack) {
477       fIsMUONTrack = kTRUE;
478       fTrack = mtrack;
479     } else {
480       fTrack = new AliMUONTrack(*mtrack);
481     }
482   }
483
484   Double_t xv, yv;
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];
488   Int_t chr[28];
489   char form[1000];
490
491   TMatrixD smatrix(2,2);
492   TMatrixD sums(2,1);
493   TMatrixD res(2,1);
494
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 };
497
498   Float_t pt    = 0.0;
499   Float_t pv[3] = {0., 0., 0.};
500
501   if (fIsMUONTrack) {
502     if (mtrack->GetMatchTrigger()) {
503       sprintf(form,"MUONTrack %2d (MT)", fLabel);
504     } else {
505       sprintf(form,"MUONTrack %2d     ", fLabel);
506     }
507     SetName(form);
508     SetLineStyle(1);
509   }
510
511   AliMUONTrackParam *trackParam = 0x0;
512   if (fIsMUONTrack || fIsESDTrack) {
513     trackParam = mtrack->GetTrackParamAtVertex();
514     SetPoint(fCount,trackParam->GetNonBendingCoor(),trackParam->GetBendingCoor(),trackParam->GetZ());
515     fCount++;
516   }
517
518   for (Int_t i = 0; i < 28; i++) {
519     xr[i]=yr[i]=zr[i]=0.0;
520     chr[i]=-1;
521   }
522   
523   Int_t nTrackHits = mtrack->GetNClusters();
524   TClonesArray* trackParamAtCluster = mtrack->GetTrackParamAtCluster();
525   for (Int_t iHit = 0; iHit < nTrackHits; iHit++){
526
527     trackParam = (AliMUONTrackParam*) trackParamAtCluster->At(iHit);
528
529     if (iHit == 0) {
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));
534       }
535       pv[0] = trackParam->Px();
536       pv[1] = trackParam->Py();
537       pv[2] = trackParam->Pz();
538       fP.Set(pv);
539     }
540
541     xr[iHit] = trackParam->GetNonBendingCoor();
542     yr[iHit] = trackParam->GetBendingCoor();
543     zr[iHit] = trackParam->GetZ();
544     chr[iHit] = trackParam->GetClusterPtr()->GetChamberId();
545
546   }
547   
548   SetPoint(fCount,xr[0],yr[0],zr[0]);
549   fCount++;
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]);
553     fCount++;
554   }
555   
556   if (!fIsMUONTrack && !fIsESDTrack) return;
557
558   Int_t nrc = 0;
559   if (mtrack->GetMatchTrigger()) {
560
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]);
564         xrc[nrc] = xr[i];
565         yrc[nrc] = yr[i];
566         zrc[nrc] = zr[i];
567         nrc++;
568       }
569     }
570
571     if (nrc < 2) return;
572
573     // fit x-z
574     smatrix.Zero();
575     sums.Zero();
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);
580       smatrix(0,0) += 1.0;
581       smatrix(1,1) += xv*xv;
582       smatrix(0,1) += xv;
583       smatrix(1,0) += xv;
584       sums(0,0)    += yv;
585       sums(1,0)    += xv*yv;
586     }
587     res = smatrix.Invert() * sums;
588     ax = res(0,0);
589     bx = res(1,0);
590
591     // fit y-z
592     smatrix.Zero();
593     sums.Zero();
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);
598       smatrix(0,0) += 1.0;
599       smatrix(1,1) += xv*xv;
600       smatrix(0,1) += xv;
601       smatrix(1,0) += xv;
602       sums(0,0)    += yv;
603       sums(1,0)    += xv*yv;
604     }
605     res = smatrix.Invert() * sums;
606     ay = res(0,0);
607     by = res(1,0);
608
609     Float_t xtc, ytc, ztc;
610     for (Int_t ii = 0; ii < 4; ii++) {
611
612       ztc = zg[ii];
613       ytc = ay+by*zg[ii];
614       xtc = ax+bx*zg[ii];
615
616       //printf("tc: x %f y %f z %f \n",xtc,ytc,ztc);
617
618       SetPoint(fCount,xtc,ytc,ztc);
619       fCount++;
620
621     }
622
623   }  // end match trigger
624
625 }
626
627 //______________________________________________________________________________
628 void AliEveMUONTrack::MakeMUONTriggerTrack(AliMUONTriggerTrack *mtrack)
629 {
630   //
631   // builds the trigger track from one point and direction
632   //
633
634   Float_t x1   = mtrack->GetX11();
635   Float_t y1   = mtrack->GetY11();
636   Float_t thex = mtrack->GetThetax();
637   Float_t they = mtrack->GetThetay();
638
639   Float_t z11 = -1600.0;
640   Float_t z22 = -1724.0;
641   Float_t dz  = z22-z11;
642
643   Float_t x2 = x1 + dz*TMath::Tan(thex);
644   Float_t y2 = y1 + dz*TMath::Tan(they);
645
646   SetPoint(fCount,x1,y1,z11); fCount++;
647   SetPoint(fCount,x2,y2,z22); fCount++;
648
649   char form[1000];
650
651   sprintf(form,"MUONTriggerTrack %2d",mtrack->GetLoTrgNum());
652   SetName(form);
653   SetLineStyle(1);
654
655 }
656
657 //______________________________________________________________________________
658 void AliEveMUONTrack::MakeESDTrack(AliESDMuonTrack *mtrack)
659 {
660   //
661   // builds the track with dipole field starting from the TParticle
662   //
663
664   fIsESDTrack = kTRUE;
665
666   char form[1000];
667   if (mtrack->GetMatchTrigger()) {
668     sprintf(form,"ESDTrack %2d (MT)", fLabel);
669   } else {
670     sprintf(form,"ESDTrack %2d     ", fLabel);
671   }
672   SetName(form);
673   SetLineStyle(3);
674   SetLineColor(0);
675
676   fTrack = new AliMUONTrack();
677
678   // create a simple track from the ESD track
679   AliMUONESDInterface::ESDToMUON(*mtrack,*fTrack);
680   
681   // reset track parameters at vertex to the ones at DCA
682   AliMUONTrackParam paramAtDCA;
683   AliMUONESDInterface::GetParamAtDCA(*mtrack, paramAtDCA);
684   fTrack->SetTrackParamAtVertex(&paramAtDCA);
685   
686   MakeMUONTrack(fTrack);
687
688 }
689
690 //______________________________________________________________________________
691 void AliEveMUONTrack::MakeMCTrack(TParticle *part)
692 {
693   //
694   // builds the track with dipole field starting from the TParticle
695   //
696
697   fIsMCTrack = kTRUE;
698
699   fPart     = new TParticle(*part);
700
701   char form[1000];
702   sprintf(form,"TEveMCTrack %2d ", fLabel);
703   SetName(form);
704   SetLineStyle(2);
705   SetLineColor(8);
706
707   Double_t vect[7], vout[7];
708   Double_t step = 1.0;
709
710   Float_t pv[3];
711   pv[0] = fPart->Px();
712   pv[1] = fPart->Py();
713   pv[2] = fPart->Pz();
714   fP.Set(pv);
715
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();
723
724   TParticlePDG *ppdg = fPart->GetPDG(1);
725   Int_t charge = (Int_t)(ppdg->Charge()/3.0);
726
727   Double_t zMax = -1750.0;
728   Double_t rMax =   350.0;
729   Double_t r    =     0.0;
730
731   Int_t nSteps = 0;
732   while ((vect[2] > zMax) && (nSteps < 10000) && (r < rMax)) {
733     nSteps++;
734     OneStepRungekutta(charge, step, vect, vout);
735     SetPoint(fCount,vout[0],vout[1],vout[2]);
736     fCount++;
737     for (Int_t i = 0; i < 7; i++) {
738       vect[i] = vout[i];
739     }
740     r = TMath::Sqrt(vect[0]*vect[0]+vect[1]*vect[1]);
741   }
742
743 }
744
745 //______________________________________________________________________________
746 void AliEveMUONTrack::MakeRefTrack(AliMUONTrack *mtrack)
747 {
748   //
749   // builds the track with dipole field starting from the TParticle
750   //
751
752   fIsRefTrack = kTRUE;
753
754   char form[1000];
755   sprintf(form,"RefTrack %2d ", fLabel);
756   SetName(form);
757   SetLineStyle(2);
758   SetLineColor(0);
759
760   MakeMUONTrack(mtrack);
761
762 }
763
764 //______________________________________________________________________________
765 void AliEveMUONTrack::Propagate(Float_t *xr, Float_t *yr, Float_t *zr, Int_t i1, Int_t i2)
766 {
767   //
768   // propagate in magnetic field between hits of indices i1 and i2
769   //
770
771   Double_t vect[7], vout[7];
772   Double_t step = 1.0;
773   Double_t zMax = 0.0;
774   Int_t  charge =   0;
775   AliMUONTrackParam *trackParam = 0;
776   TClonesArray *trackParamAtCluster = 0;
777
778   if (i2 == 9999) {
779     zMax = zr[i1]+1.5*step;
780   } else {
781     zMax = zr[i2]+1.5*step;
782   }
783
784   trackParamAtCluster = fTrack->GetTrackParamAtCluster();
785
786   if (IsMUONTrack() || IsESDTrack()) {
787     trackParam = (AliMUONTrackParam*)trackParamAtCluster->At(i1);
788     charge = (Int_t)trackParam->GetCharge();
789   }
790   if (IsRefTrack()) {
791     trackParam = fTrack->GetTrackParamAtVertex();
792     charge = (Int_t)trackParam->GetCharge();
793     trackParam = (AliMUONTrackParam*)trackParamAtCluster->At(i1);
794   }
795
796   vect[0] = xr[i1];
797   vect[1] = yr[i1];
798   vect[2] = zr[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];
803
804   Int_t nSteps = 0;
805   while ((vect[2] > zMax) && (nSteps < 10000)) {
806     nSteps++;
807     OneStepRungekutta(charge, step, vect, vout);
808     SetPoint(fCount,vout[0],vout[1],vout[2]);
809     fCount++;
810     for (Int_t i = 0; i < 7; i++) {
811       vect[i] = vout[i];
812     }
813   }
814
815 }
816
817 //______________________________________________________________________________
818 void AliEveMUONTrack::GetField(Double_t *position, Double_t *field)
819 {
820   //
821   // returns field components at position, for a give field map
822   //
823
824   /// interface for arguments in double precision (Why ? ChF)
825   Float_t x[3], b[3];
826
827   x[0] = position[0]; x[1] = position[1]; x[2] = position[2];
828
829   if (fgFieldMap) {
830     fgFieldMap->Field(x,b);
831   }
832   else {
833     AliWarning("No field map");
834     field[0] = field[1] = field[2] = 0.0;
835     return;
836   }
837
838   // force components
839   //b[1] = 0.0;
840   //b[2] = 0.0;
841
842   field[0] = b[0]; field[1] = b[1]; field[2] = b[2];
843
844   return;
845
846 }
847
848 //______________________________________________________________________________
849 void AliEveMUONTrack::OneStepRungekutta(Double_t charge, Double_t step,
850                                   Double_t* vect, Double_t* vout)
851 {
852 ///     ******************************************************************
853 ///     *                                                                *
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)                                 *
857 ///     *                                                                *
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)                                          *
866 ///     *                                                                *
867 ///     *    ==>Called by : <USER>, GUSWIM                               *
868 ///     *       Authors    R.Brun, M.Hansroul  *********                 *
869 ///     *                  V.Perevoztchikov (CUT STEP implementation)    *
870 ///     *                                                                *
871 ///     *                                                                *
872 ///     ******************************************************************
873
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;
880
881     Double_t x;
882     Double_t y;
883     Double_t z;
884
885     Double_t xt;
886     Double_t yt;
887     Double_t zt;
888
889     Double_t maxit = 1992;
890     Double_t maxcut = 11;
891
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;
897
898     const Double_t kpisqua = 9.86960440109;
899     const Int_t kix  = 0;
900     const Int_t kiy  = 1;
901     const Int_t kiz  = 2;
902     const Int_t kipx = 3;
903     const Int_t kipy = 4;
904     const Int_t kipz = 5;
905
906     // *.
907     // *.    ------------------------------------------------------------------
908     // *.
909     // *             this constant is for units cm,gev/c and kgauss
910     // *
911     Int_t iter = 0;
912     Int_t ncut = 0;
913     for(Int_t j = 0; j < 7; j++)
914       vout[j] = vect[j];
915
916     Double_t  pinv   = kec * charge / vect[6];
917     Double_t tl = 0.;
918     Double_t h = step;
919     Double_t rest;
920
921
922     do {
923       rest  = step - tl;
924       if (TMath::Abs(h) > TMath::Abs(rest)) h = rest;
925       //cmodif: call gufld(vout,f) changed into:
926
927       GetField(vout,f);
928
929       // *
930       // *             start of integration
931       // *
932       x      = vout[0];
933       y      = vout[1];
934       z      = vout[2];
935       a      = vout[3];
936       b      = vout[4];
937       c      = vout[5];
938
939       h2     = khalf * h;
940       h4     = khalf * h2;
941       ph     = pinv * h;
942       ph2    = khalf * ph;
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;
948
949       dxt    = h2 * a + h4 * secxs[0];
950       dyt    = h2 * b + h4 * secys[0];
951       dzt    = h2 * c + h4 * seczs[0];
952       xt     = x + dxt;
953       yt     = y + dyt;
954       zt     = z + dzt;
955       // *
956       // *              second intermediate point
957       // *
958
959       est = TMath::Abs(dxt) + TMath::Abs(dyt) + TMath::Abs(dzt);
960       if (est > h) {
961         if (ncut++ > maxcut) break;
962         h *= khalf;
963         continue;
964       }
965
966       xyzt[0] = xt;
967       xyzt[1] = yt;
968       xyzt[2] = zt;
969
970       //cmodif: call gufld(xyzt,f) changed into:
971       GetField(xyzt,f);
972
973       at     = a + secxs[0];
974       bt     = b + secys[0];
975       ct     = c + seczs[0];
976
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;
980       at     = a + secxs[1];
981       bt     = b + secys[1];
982       ct     = c + seczs[1];
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]);
989       xt     = x + dxt;
990       yt     = y + dyt;
991       zt     = z + dzt;
992       at     = a + 2.*secxs[2];
993       bt     = b + 2.*secys[2];
994       ct     = c + 2.*seczs[2];
995
996       est = TMath::Abs(dxt)+TMath::Abs(dyt)+TMath::Abs(dzt);
997       if (est > 2.*TMath::Abs(h)) {
998         if (ncut++ > maxcut) break;
999         h *= khalf;
1000         continue;
1001       }
1002
1003       xyzt[0] = xt;
1004       xyzt[1] = yt;
1005       xyzt[2] = zt;
1006
1007       //cmodif: call gufld(xyzt,f) changed into:
1008       GetField(xyzt,f);
1009
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;
1013
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;
1020
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]));
1024
1025       if (est > kdlt && TMath::Abs(h) > 1.e-4) {
1026         if (ncut++ > maxcut) break;
1027         h *= khalf;
1028         continue;
1029       }
1030
1031       ncut = 0;
1032       // *               if too many iterations, go to helix
1033       if (iter++ > maxit) break;
1034
1035       tl += h;
1036       if (est < kdlt32)
1037         h *= 2.;
1038       cba    = 1./ TMath::Sqrt(a*a + b*b + c*c);
1039       vout[0] = x;
1040       vout[1] = y;
1041       vout[2] = z;
1042       vout[3] = cba*a;
1043       vout[4] = cba*b;
1044       vout[5] = cba*c;
1045       rest = step - tl;
1046       if (step < 0.) rest = -rest;
1047       if (rest < 1.e-5*TMath::Abs(step)) return;
1048
1049     } while(1);
1050
1051     // angle too big, use helix
1052
1053     f1  = f[0];
1054     f2  = f[1];
1055     f3  = f[2];
1056     f4  = TMath::Sqrt(f1*f1+f2*f2+f3*f3);
1057     rho = -f4*pinv;
1058     tet = rho * step;
1059
1060     hnorm = 1./f4;
1061     f1 = f1*hnorm;
1062     f2 = f2*hnorm;
1063     f3 = f3*hnorm;
1064
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];
1068
1069     hp = f1*vect[kipx] + f2*vect[kipy] + f3*vect[kipz];
1070
1071     rho1 = 1./rho;
1072     sint = TMath::Sin(tet);
1073     cost = 2.*TMath::Sin(khalf*tet)*TMath::Sin(khalf*tet);
1074
1075     g1 = sint*rho1;
1076     g2 = cost*rho1;
1077     g3 = (tet-sint) * hp*rho1;
1078     g4 = -cost;
1079     g5 = sint;
1080     g6 = cost * hp;
1081
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;
1085
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;
1089
1090     return;
1091 }
1092
1093 //______________________________________________________________________________
1094 Int_t AliEveMUONTrack::ColorIndex(Float_t val)
1095 {
1096   //
1097   // returns color index in the palette for a give value
1098   //
1099
1100   Float_t threshold =  0.0;
1101   Float_t maxVal    = 10.0;
1102
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);
1106
1107   return gStyle->GetColorPalette(TMath::Min(nCol - 1, cBin));
1108
1109 }