]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/Alieve/AliEveTrackFitter.cxx
Comment code, hide public data-members.
[u/mrichter/AliRoot.git] / EVE / Alieve / AliEveTrackFitter.cxx
1 // $Id$
2 // Main authors: Matevz Tadel & Alja Mrak-Tadel: 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 "AliEveTrackFitter.h"
11
12 #include "TCanvas.h"
13 #include "TGraph.h"
14 #include "TGraphErrors.h"
15
16 #include "AliCluster.h"
17 #include "AliRieman.h"
18 #include "AliExternalTrackParam.h"
19
20 #include <TEveTrack.h>
21 #include <TEveTrackPropagator.h>
22 #include <TEveVSDStructs.h>
23 #include <TEveManager.h>
24
25
26 //______________________________________________________________________
27 // AliEveTrackFitter
28 //
29 // AliEveTrackFitter is an interface to helix fit. It creates a set of
30 // points, listening to signal PointCtrlClicked() of any
31 // TEvePointSet. Via editor it fits selected points and creates a
32 // reconstructed track.
33 // 
34
35 ClassImp(AliEveTrackFitter)
36
37 AliEveTrackFitter::AliEveTrackFitter(const Text_t* name, Int_t n_points) :
38     TEvePointSet(name, n_points),
39
40     fGraphSelected(0),
41     fGraphFitted(0),
42
43     fRieman(0),
44     fConnected(kFALSE),
45     fTrackList(0)
46 {
47   // Constructor.
48
49   SetMarkerColor(3);
50   SetOwnIds(kFALSE);
51
52   fGraphSelected = new TGraph();
53   fGraphSelected->SetName("Selected points");
54   fGraphSelected->SetMarkerColor(4);
55   fGraphSelected->SetMarkerStyle(4);  
56   fGraphSelected->SetMarkerSize(2);
57
58   fGraphFitted = new TGraphErrors();
59   fGraphFitted->SetName("Fitted points");
60   fGraphFitted->SetMarkerColor(2);
61
62   fTrackList = new TEveTrackList("Tracks");
63   fTrackList->SetLineWidth(2);
64   fTrackList->SetLineColor(8);
65   fTrackList->IncDenyDestroy();
66   fTrackList->GetPropagator()->SetEditPathMarks(kTRUE);
67   gEve->AddElement(fTrackList, this);
68   UpdateItems();
69 }
70
71 AliEveTrackFitter::~AliEveTrackFitter()
72 {
73   // Destructor.
74
75   if(fRieman) delete fRieman;
76   fTrackList->DecDenyDestroy();
77   delete fTrackList;
78 }
79
80 /**************************************************************************/
81 void AliEveTrackFitter::DestroyElements()
82 {
83   // Virtual method of base class TEveElement.
84   // It preserves track list to have coomon track propagator attributes.
85
86   TEveElement::DestroyElements();
87   gEve->AddElement(fTrackList, this);
88   fTrackList->DestroyElements();
89   UpdateItems();
90 }
91
92 /**************************************************************************/
93 void AliEveTrackFitter::Start()
94 {
95   // Start selection of points.
96
97   Reset();
98   if(fConnected == kFALSE)
99   {
100     TQObject::Connect("TEvePointSet", "PointCtrlClicked(TEvePointSet*,Int_t)",
101                       "AliEveTrackFitter", this, "AddFitPoint(TEvePointSet*,Int_t)");
102
103     fConnected = kTRUE;
104   }
105 }
106
107 void AliEveTrackFitter::Stop()
108 {
109   // Stop selection of points.
110
111   if(fConnected)
112   {
113     TQObject::Disconnect("TEvePointSet", "AddFitPoint(TEvePointSet*,Int_t)");
114     fConnected = kFALSE;
115   }
116 }
117
118 /**************************************************************************/
119
120 void AliEveTrackFitter::AddFitPoint(TEvePointSet* ps, Int_t n)
121
122   // Add/remove given point depending if exists in the fMapPS.
123  
124   Float_t x, y, z;
125
126   std::map<Point_t, Int_t>::iterator g = fMapPS.find(Point_t(ps, n));
127   if(g != fMapPS.end())
128   {
129     Int_t idx = g->second;
130     if(idx != fLastPoint)
131     {
132       GetPoint(fLastPoint, x, y, z);
133       SetPoint(idx, x, y, z);
134     }
135     fMapPS.erase(g);
136     fLastPoint--;
137   }
138   else 
139   {
140     fMapPS[Point_t(ps, n)] = Size();
141     ps->GetPoint(n, x, y, z);
142     SetNextPoint(x, y, z); 
143     SetPointId(ps->GetPointId(n));
144   }
145   ResetBBox();
146   ElementChanged(kTRUE, kTRUE);
147 }
148
149 /**************************************************************************/
150
151 void AliEveTrackFitter::FitTrack()
152 {
153   // Fit selected points with AliRieman fitter.
154
155   using namespace TMath;
156
157   if(fRieman) delete fRieman;
158   fRieman = new AliRieman(Size());
159
160   Float_t x, y, z;
161   Int_t alphaIdx = 0;
162   GetPoint(alphaIdx, x, y, z);
163   Float_t minR2=x*x + y*y;
164   for (Int_t i=1; i<=fLastPoint; i++)
165   {
166     GetPoint(i, x, y, z);
167     Float_t cR2 = x*x + y*y;
168     if ( minR2 > cR2 )
169     {
170       minR2 = cR2;
171       alphaIdx = i;
172     }
173   }
174   GetPoint(alphaIdx, x, y, z);
175   fAlpha = ATan2(y, x);
176   Float_t sin = Sin(-fAlpha);
177   Float_t cos = Cos(-fAlpha);  
178   for (Int_t i=0; i<=fLastPoint; i++) { 
179     GetPoint(i, x, y, z);
180     fRieman->AddPoint(cos*x - sin*y, cos*y + sin*x, z, 1, 1);
181   }
182   fRieman->Update();
183
184   Double_t r = Sqrt(minR2);
185   Double_t param[5];
186   Double_t cov[15];
187   fRieman->GetExternalParameters(r, param, cov);
188   // curvature to pt
189   param[4] /= TEveTrackPropagator::fgDefMagField*TEveTrackPropagator::fgkB2C;
190   // sign in tang
191   if(param[4] < 0) param[3] *= -1;
192   AliExternalTrackParam trackParam(r, fAlpha, param, cov);
193   trackParam.Print();
194
195   // make track
196   Double_t AliEveV0[3];
197   trackParam.GetXYZAt(r, TEveTrackPropagator::fgDefMagField, AliEveV0);
198   Double_t P0[3];
199   trackParam.GetPxPyPzAt(r, TEveTrackPropagator::fgDefMagField, P0);
200   TEveRecTrack rc;
201   rc.fV.Set(AliEveV0); 
202   rc.fP.Set(P0);
203   rc.fSign = trackParam.Charge();
204
205   TEveTrack* track = new TEveTrack(&rc, fTrackList->GetPropagator());
206   track->SetName(Form("track %f", fAlpha));
207   TEvePathMark* pm = new TEvePathMark(TEvePathMark::kDaughter);
208   for(Int_t i=0; i==fLastPoint; i++)
209   {
210     GetPoint(i, x, y, z);
211     pm->fV.Set(x, y, z);
212     pm->fP.Set(P0); 
213     track->AddPathMark(pm);
214   }
215   track->MakeTrack();
216   track->SetAttLineAttMarker(fTrackList); 
217   gEve->AddElement(track, fTrackList);
218 }
219
220
221 void AliEveTrackFitter::Reset(Int_t n, Int_t ids)
222 {
223   // Reset selection.
224
225   if(fRieman) fRieman->Reset();
226   TEvePointSet::Reset(n, ids);
227   fMapPS.clear();
228 }
229
230 /**************************************************************************/
231 void AliEveTrackFitter::DrawRiemanGraph()
232 {
233   // Draw graph of rieman fit.
234
235    static const TEveException eH("AliEveTrackFitter::DrawRiemanGraph ");
236
237   if(fRieman == 0)
238     throw(eH + "fitter not set.");
239
240   Int_t nR = fRieman->GetN();
241   fGraphSelected->Set(nR);
242   fGraphFitted->Set(nR);
243
244   Double_t* x =  fRieman->GetX();  
245   Double_t* y =  fRieman->GetY();
246   Double_t* sy =  fRieman->GetSy();
247   for (Int_t i=0; i<nR; i++)
248   {
249     fGraphSelected->SetPoint(i, x[i], y[i]);
250     fGraphFitted->SetPoint(i, x[i], fRieman->GetYat(x[i]));
251     fGraphFitted->SetPointError(i, 0.1, sy[i]);
252   }
253   
254   if (gPad) gPad->Clear();
255   fGraphSelected->Draw("AP");
256   fGraphFitted->Draw("SAME P");
257   gPad->GetCanvas()->SetTitle(Form("AliRieman alpha: %f", fAlpha));
258   gPad->Modified();
259   gPad->Update();
260 }