CbmRoot
CbmFitGlobalTracksQa.cxx
Go to the documentation of this file.
1 //-------------------------------------------------------------------------
2 //----- CbmFitGlobalTracksQa -----
3 //----- Created 07/03/06 by D. Kresan -----
4 //-------------------------------------------------------------------------
5 #include "CbmFitGlobalTracksQa.h"
6 
7 #include "CbmGlobalTrack.h"
8 
9 #include "CbmStsHit.h"
10 #include "CbmStsPoint.h"
11 #include "CbmStsTrack.h"
12 #include "CbmTrdHit.h"
13 #include "CbmTrdPoint.h"
14 #include "CbmTrdTrack.h"
15 #include "FairRootManager.h"
16 
17 #include "TClonesArray.h"
18 #include "TH1.h"
19 #include "TMath.h"
20 #include "TVector3.h"
21 
22 #include <iostream>
23 
24 using std::cout;
25 using std::endl;
26 
27 //________________________________________________________________________
28 //
29 // CbmFitGlobalTracksQa
30 //
31 // Task for performance analysis of the global track fit
32 //
33 
34 
35 // -----------------------------------------------------------------------
37  : FairTask()
38  , fArrayStsPoint(NULL)
39  , fArrayTrdPoint(NULL)
40  , fArrayStsHit(NULL)
41  , fArrayTrdHit(NULL)
42  , fArrayStsTrack(NULL)
43  , fArrayTrdTrack(NULL)
44  , fArrayGlbTrack(NULL)
45  , fEvents(0)
46  , fh_first_resx(NULL)
47  , fh_first_resy(NULL)
48  , fh_first_restx(NULL)
49  , fh_first_resty(NULL)
50  , fh_first_resz(NULL)
51  , fh_last_resx(NULL)
52  , fh_last_resy(NULL)
53  , fh_last_restx(NULL)
54  , fh_last_resty(NULL)
55  , fh_last_resz(NULL)
56  , fh_first_pullx(NULL)
57  , fh_first_pully(NULL)
58  , fh_first_pulltx(NULL)
59  , fh_first_pullty(NULL)
60  , fh_last_pullx(NULL)
61  , fh_last_pully(NULL)
62  , fh_last_pulltx(NULL)
63  , fh_last_pullty(NULL)
64  , fh_chi2ndf(NULL) {
66 }
67 // -----------------------------------------------------------------------
68 
69 
70 // -----------------------------------------------------------------------
71 CbmFitGlobalTracksQa::CbmFitGlobalTracksQa(const char* name, Int_t verbose)
72  : FairTask(name, verbose)
73  , fArrayStsPoint(NULL)
74  , fArrayTrdPoint(NULL)
75  , fArrayStsHit(NULL)
76  , fArrayTrdHit(NULL)
77  , fArrayStsTrack(NULL)
78  , fArrayTrdTrack(NULL)
79  , fArrayGlbTrack(NULL)
80  , fEvents(0)
81  , fh_first_resx(NULL)
82  , fh_first_resy(NULL)
83  , fh_first_restx(NULL)
84  , fh_first_resty(NULL)
85  , fh_first_resz(NULL)
86  , fh_last_resx(NULL)
87  , fh_last_resy(NULL)
88  , fh_last_restx(NULL)
89  , fh_last_resty(NULL)
90  , fh_last_resz(NULL)
91  , fh_first_pullx(NULL)
92  , fh_first_pully(NULL)
93  , fh_first_pulltx(NULL)
94  , fh_first_pullty(NULL)
95  , fh_last_pullx(NULL)
96  , fh_last_pully(NULL)
97  , fh_last_pulltx(NULL)
98  , fh_last_pullty(NULL)
99  , fh_chi2ndf(NULL) {
101 }
102 // -----------------------------------------------------------------------
103 
104 
105 // -----------------------------------------------------------------------
107  // Destructor
108 }
109 // -----------------------------------------------------------------------
110 
111 
112 // -----------------------------------------------------------------------
114  // Initialisation of the task
115 
116  // Get pointer to the ROOT I/O manager
117  FairRootManager* rootMgr = FairRootManager::Instance();
118  if (NULL == rootMgr) {
119  cout << "-E- CbmFitGlobalTracksQa::Init :"
120  << " ROOT manager is not instantiated" << endl;
121  return kFATAL;
122  }
123  // Get MC points
124  fArrayStsPoint = (TClonesArray*) rootMgr->GetObject("StsPoint");
125  if (NULL == fArrayStsPoint) {
126  cout << "-E- CbmFitGlobalTracksQa::Init : "
127  << "no STS point array!" << endl;
128  }
129  fArrayTrdPoint = (TClonesArray*) rootMgr->GetObject("TrdPoint");
130  if (NULL == fArrayTrdPoint) {
131  cout << "-E- CbmFitGlobalTracksQa::Init : "
132  << "no TRD point array!" << endl;
133  }
134  // Get hit arrays
135  fArrayStsHit = (TClonesArray*) rootMgr->GetObject("StsHit");
136  if (NULL == fArrayStsHit) {
137  cout << "-W- CbmFitGlobalTracksQa::Init :"
138  << " no STS hit array" << endl;
139  }
140  fArrayTrdHit = (TClonesArray*) rootMgr->GetObject("TrdHit");
141  if (NULL == fArrayTrdHit) {
142  cout << "-W- CbmFitGlobalTracksQa::Init :"
143  << " no TRD hit array" << endl;
144  }
145  // Get track arrays
146  fArrayStsTrack = (TClonesArray*) rootMgr->GetObject("StsTrack");
147  if (NULL == fArrayStsTrack) {
148  cout << "-W- CbmFitGlobalTracksQa::Init : "
149  << "no STS track array!" << endl;
150  }
151  fArrayTrdTrack = (TClonesArray*) rootMgr->GetObject("TrdTrack");
152  if (NULL == fArrayTrdTrack) {
153  cout << "-W- CbmFitGlobalTracksQa::Init : "
154  << "no TRD track array!" << endl;
155  }
156  fArrayGlbTrack = (TClonesArray*) rootMgr->GetObject("GlobalTrack");
157  if (NULL == fArrayGlbTrack) {
158  cout << "-W- CbmFitGlobalTracksQa::Init : "
159  << "no global track array!" << endl;
160  }
161  return kSUCCESS;
162 }
163 // -----------------------------------------------------------------------
164 
165 
166 // -----------------------------------------------------------------------
167 void CbmFitGlobalTracksQa::Exec(Option_t*) {
168  // Task execution
169  if (NULL == fArrayStsPoint || NULL == fArrayTrdPoint || NULL == fArrayStsHit
170  || NULL == fArrayTrdHit || NULL == fArrayStsTrack
171  || NULL == fArrayTrdTrack || NULL == fArrayGlbTrack)
172  return;
173 
174  // Loop over the global tracks
175  CbmStsPoint* stsPoint;
176  CbmTrdPoint* trdPoint;
177  CbmStsHit* stsHit;
178  CbmTrdHit* trdHit;
179  CbmStsTrack* stsTrack;
180  CbmTrdTrack* trdTrack;
181  CbmGlobalTrack* glbTrack;
182  Int_t stsTrackIndex;
183  Int_t trdTrackIndex;
184  Int_t stsHitIndex;
185  Int_t trdHitIndex;
186  Int_t stsPointIndex;
187  Int_t trdPointIndex;
188  TVector3 pos;
189  TVector3 mom;
190 
191  for (Int_t iGlbTrack = 0; iGlbTrack < fArrayGlbTrack->GetEntriesFast();
192  iGlbTrack++) {
193 
194  // Get pointer to the global track
195  glbTrack = (CbmGlobalTrack*) fArrayGlbTrack->At(iGlbTrack);
196  if (NULL == glbTrack) continue;
197 
198 
199  // Chi2/NDF
200  fh_chi2ndf->Fill(glbTrack->GetChi2() / (Double_t) glbTrack->GetNDF());
201 
202 
203  // Performance at the first plane
204  stsTrackIndex = glbTrack->GetStsTrackIndex();
205  if (stsTrackIndex >= 0) {
206  stsTrack = (CbmStsTrack*) fArrayStsTrack->At(stsTrackIndex);
207  if (NULL != stsTrack) {
208  stsPointIndex = -1;
209  if (stsTrack->GetNofHits()) {
210  stsHitIndex = stsTrack->GetHitIndex(0);
211  stsHit = (CbmStsHit*) fArrayStsHit->At(stsHitIndex);
212  if (NULL != stsHit) { stsPointIndex = stsHit->GetRefId(); }
213  } else {
214  cout << "-W- CbmFitGlobalTracksQa::Exec : "
215  << "STS track " << stsTrackIndex << " is empty!" << endl;
216  }
217 
218  if (stsPointIndex >= 0) {
219  stsPoint = (CbmStsPoint*) fArrayStsPoint->At(stsPointIndex);
220  if (NULL != stsPoint) {
221  stsPoint->Position(pos);
222  stsPoint->Momentum(mom);
223  fh_first_resx->Fill(glbTrack->GetParamFirst()->GetX() - pos.X());
224  fh_first_resy->Fill(glbTrack->GetParamFirst()->GetY() - pos.Y());
225  fh_first_resz->Fill(glbTrack->GetParamFirst()->GetZ() - pos.Z());
226  fh_first_restx->Fill(glbTrack->GetParamFirst()->GetTx()
227  - mom.X() / mom.Z());
228  fh_first_resty->Fill(glbTrack->GetParamFirst()->GetTy()
229  - mom.Y() / mom.Z());
230  fh_first_pullx->Fill(
231  (glbTrack->GetParamFirst()->GetX() - pos.X())
232  / TMath::Sqrt(glbTrack->GetParamFirst()->GetCovariance(0, 0)));
233  fh_first_pully->Fill(
234  (glbTrack->GetParamFirst()->GetY() - pos.Y())
235  / TMath::Sqrt(glbTrack->GetParamFirst()->GetCovariance(1, 1)));
236  fh_first_pulltx->Fill(
237  (glbTrack->GetParamFirst()->GetTx() - mom.X() / mom.Z())
238  / TMath::Sqrt(glbTrack->GetParamFirst()->GetCovariance(2, 2)));
239  fh_first_pullty->Fill(
240  (glbTrack->GetParamFirst()->GetTy() - mom.Y() / mom.Z())
241  / TMath::Sqrt(glbTrack->GetParamFirst()->GetCovariance(3, 3)));
242  }
243  }
244  }
245  } else {
246  trdTrackIndex = glbTrack->GetTrdTrackIndex();
247  trdTrack = (CbmTrdTrack*) fArrayTrdTrack->At(trdTrackIndex);
248  if (NULL != trdTrack) {
249  trdPointIndex = -1;
250  if (trdTrack->GetNofHits()) {
251  trdHitIndex = trdTrack->GetHitIndex(0);
252  trdHit = (CbmTrdHit*) fArrayTrdHit->At(trdHitIndex);
253  if (NULL != trdHit) { trdPointIndex = trdHit->GetRefId(); }
254  } else {
255  cout << "-W- CbmFitGlobalTracksQa::Exec : "
256  << "TRD track " << trdTrackIndex << " is empty!" << endl;
257  }
258 
259  if (trdPointIndex >= 0) {
260  trdPoint = (CbmTrdPoint*) fArrayTrdPoint->At(trdPointIndex);
261  if (NULL != trdPoint) {
262  trdPoint->Position(pos);
263  trdPoint->Momentum(mom);
264  fh_first_resx->Fill(glbTrack->GetParamFirst()->GetX() - pos.X());
265  fh_first_resy->Fill(glbTrack->GetParamFirst()->GetY() - pos.Y());
266  fh_first_resz->Fill(glbTrack->GetParamFirst()->GetZ() - pos.Z());
267  fh_first_restx->Fill(glbTrack->GetParamFirst()->GetTx()
268  - mom.X() / mom.Z());
269  fh_first_resty->Fill(glbTrack->GetParamFirst()->GetTy()
270  - mom.Y() / mom.Z());
271  fh_first_pullx->Fill(
272  (glbTrack->GetParamFirst()->GetX() - pos.X())
273  / TMath::Sqrt(glbTrack->GetParamFirst()->GetCovariance(0, 0)));
274  fh_first_pully->Fill(
275  (glbTrack->GetParamFirst()->GetY() - pos.Y())
276  / TMath::Sqrt(glbTrack->GetParamFirst()->GetCovariance(1, 1)));
277  fh_first_pulltx->Fill(
278  (glbTrack->GetParamFirst()->GetTx() - mom.X() / mom.Z())
279  / TMath::Sqrt(glbTrack->GetParamFirst()->GetCovariance(2, 2)));
280  fh_first_pullty->Fill(
281  (glbTrack->GetParamFirst()->GetTy() - mom.Y() / mom.Z())
282  / TMath::Sqrt(glbTrack->GetParamFirst()->GetCovariance(3, 3)));
283  }
284  }
285  }
286  }
287 
288 
289  // Performance at the last plane
290  trdTrackIndex = glbTrack->GetTrdTrackIndex();
291  if (trdTrackIndex >= 0) {
292  trdTrack = (CbmTrdTrack*) fArrayTrdTrack->At(trdTrackIndex);
293  if (NULL != trdTrack) {
294  trdPointIndex = -1;
295  if (trdTrack->GetNofHits()) {
296  trdHitIndex = trdTrack->GetHitIndex(trdTrack->GetNofHits() - 1);
297  trdHit = (CbmTrdHit*) fArrayTrdHit->At(trdHitIndex);
298  if (NULL != trdHit) { trdPointIndex = trdHit->GetRefId(); }
299  } else {
300  cout << "-W- CbmFitGlobalTracksQa::Exec : "
301  << "TRD track " << trdTrackIndex << " is empty!" << endl;
302  }
303 
304  if (trdPointIndex >= 0) {
305  trdPoint = (CbmTrdPoint*) fArrayTrdPoint->At(trdPointIndex);
306  if (NULL != trdPoint) {
307  trdPoint->Position(pos);
308  trdPoint->Momentum(mom);
309  fh_last_resx->Fill(glbTrack->GetParamLast()->GetX() - pos.X());
310  fh_last_resy->Fill(glbTrack->GetParamLast()->GetY() - pos.Y());
311  fh_last_resz->Fill(glbTrack->GetParamLast()->GetZ() - pos.Z());
312  fh_last_restx->Fill(glbTrack->GetParamLast()->GetTx()
313  - mom.X() / mom.Z());
314  fh_last_resty->Fill(glbTrack->GetParamLast()->GetTy()
315  - mom.Y() / mom.Z());
316  fh_last_pullx->Fill(
317  (glbTrack->GetParamLast()->GetX() - pos.X())
318  / TMath::Sqrt(glbTrack->GetParamLast()->GetCovariance(0, 0)));
319  fh_last_pully->Fill(
320  (glbTrack->GetParamLast()->GetY() - pos.Y())
321  / TMath::Sqrt(glbTrack->GetParamLast()->GetCovariance(1, 1)));
322  fh_last_pulltx->Fill(
323  (glbTrack->GetParamLast()->GetTx() - mom.X() / mom.Z())
324  / TMath::Sqrt(glbTrack->GetParamLast()->GetCovariance(2, 2)));
325  fh_last_pullty->Fill(
326  (glbTrack->GetParamLast()->GetTy() - mom.Y() / mom.Z())
327  / TMath::Sqrt(glbTrack->GetParamLast()->GetCovariance(3, 3)));
328  }
329  }
330  }
331  } else {
332  stsTrackIndex = glbTrack->GetStsTrackIndex();
333  stsTrack = (CbmStsTrack*) fArrayStsTrack->At(stsTrackIndex);
334  if (NULL != stsTrack) {
335  stsPointIndex = -1;
336  if (stsTrack->GetNofHits()) {
337  stsHitIndex = stsTrack->GetHitIndex(stsTrack->GetNofHits() - 1);
338  stsHit = (CbmStsHit*) fArrayStsHit->At(stsHitIndex);
339  if (NULL != stsHit) { stsPointIndex = stsHit->GetRefId(); }
340  } else {
341  cout << "-W- CbmFitGlobalTracksQa::Exec : "
342  << "STS track " << stsTrackIndex << " is empty!" << endl;
343  }
344 
345  if (stsPointIndex >= 0) {
346  stsPoint = (CbmStsPoint*) fArrayStsPoint->At(stsPointIndex);
347  if (NULL != stsPoint) {
348  stsPoint->Position(pos);
349  stsPoint->Momentum(mom);
350  fh_last_resx->Fill(glbTrack->GetParamLast()->GetX() - pos.X());
351  fh_last_resy->Fill(glbTrack->GetParamLast()->GetY() - pos.Y());
352  fh_last_resz->Fill(glbTrack->GetParamLast()->GetZ() - pos.Z());
353  fh_last_restx->Fill(glbTrack->GetParamLast()->GetTx()
354  - mom.X() / mom.Z());
355  fh_last_resty->Fill(glbTrack->GetParamLast()->GetTy()
356  - mom.Y() / mom.Z());
357  fh_last_pullx->Fill(
358  (glbTrack->GetParamLast()->GetX() - pos.X())
359  / TMath::Sqrt(glbTrack->GetParamLast()->GetCovariance(0, 0)));
360  fh_last_pully->Fill(
361  (glbTrack->GetParamLast()->GetY() - pos.Y())
362  / TMath::Sqrt(glbTrack->GetParamLast()->GetCovariance(1, 1)));
363  fh_last_pulltx->Fill(
364  (glbTrack->GetParamLast()->GetTx() - mom.X() / mom.Z())
365  / TMath::Sqrt(glbTrack->GetParamLast()->GetCovariance(2, 2)));
366  fh_last_pullty->Fill(
367  (glbTrack->GetParamLast()->GetTy() - mom.Y() / mom.Z())
368  / TMath::Sqrt(glbTrack->GetParamLast()->GetCovariance(3, 3)));
369  }
370  }
371  }
372  }
373  }
374 
375 
376  fEvents += 1;
377 }
378 
379 
380 // -----------------------------------------------------------------------
382  // Finish of the task execution
384 }
385 // -----------------------------------------------------------------------
386 
387 
388 // -----------------------------------------------------------------------
390  // Create control histogramms
391 
392  // Residuals at first TRD layer
393  fh_first_resx =
394  new TH1F("h_glb_first_resx", "x residual at first TRD layer", 200, -1., 1.);
395  fh_first_resy =
396  new TH1F("h_glb_first_resy", "y residual at first TRD layer", 200, -1., 1.);
397  fh_first_restx = new TH1F(
398  "h_glb_first_restx", "t_{x} residual at first TRD layer", 200, -0.01, 0.01);
399  fh_first_resty = new TH1F(
400  "h_glb_first_resty", "t_{y} residual at first TRD layer", 200, -0.01, 0.01);
401  fh_first_resz = new TH1F(
402  "h_glb_first_resz", "z residual at first TRD layer", 100, -0.5, 0.5);
403 
404  // Residuals at last TRD layer
405  fh_last_resx =
406  new TH1F("h_glb_last_resx", "x residual at last TRD layer", 200, -1., 1.);
407  fh_last_resy =
408  new TH1F("h_glb_last_resy", "y residual at last TRD layer", 200, -1., 1.);
409  fh_last_restx = new TH1F(
410  "h_glb_last_restx", "t_{x} residual at last TRD layer", 200, -0.01, 0.01);
411  fh_last_resty = new TH1F(
412  "h_glb_last_resty", "t_{y} residual at last TRD layer", 200, -0.01, 0.01);
413  fh_last_resz =
414  new TH1F("h_glb_last_resz", "z residual at last TRD layer", 100, -0.5, 0.5);
415 
416  // Pulls at first TRD layer
418  new TH1F("h_glb_first_pullx", "x pull at first TRD layer", 200, -10., 10.);
420  new TH1F("h_glb_first_pully", "y pull at first TRD layer", 200, -10., 10.);
421  fh_first_pulltx = new TH1F(
422  "h_glb_first_pulltx", "t_{x} pull at first TRD layer", 200, -10., 10.);
423  fh_first_pullty = new TH1F(
424  "h_glb_first_pullty", "t_{y} pull at first TRD layer", 200, -10., 10.);
425 
426  // Pulls at last TRD layer
427  fh_last_pullx =
428  new TH1F("h_glb_last_pullx", "x pull at last TRD layer", 200, -10., 10.);
429  fh_last_pully =
430  new TH1F("h_glb_last_pully", "y pull at last TRD layer", 200, -10., 10.);
431  fh_last_pulltx = new TH1F(
432  "h_glb_last_pulltx", "t_{x} pull at last TRD layer", 200, -10., 10.);
433  fh_last_pullty = new TH1F(
434  "h_glb_last_pullty", "t_{y} pull at last TRD layer", 200, -10., 10.);
435 
436  // Chi2/NDF of track fit
437  fh_chi2ndf =
438  new TH1F("h_glb_chi2ndf", "#chi^{2}/NDF of track fit", 500, 0., 50.);
439 }
440 // -----------------------------------------------------------------------
441 
442 
443 // -----------------------------------------------------------------------
445  // Write control histogramms to file
446 
447  // Residuals at first and last TRD layer
448  fh_first_resx->Scale(1. / fEvents);
449  fh_first_resy->Scale(1. / fEvents);
450  fh_first_restx->Scale(1. / fEvents);
451  fh_first_resty->Scale(1. / fEvents);
452  fh_last_resx->Scale(1. / fEvents);
453  fh_last_resy->Scale(1. / fEvents);
454  fh_last_restx->Scale(1. / fEvents);
455  fh_last_resty->Scale(1. / fEvents);
456  fh_first_resx->Write();
457  fh_first_resy->Write();
458  fh_first_restx->Write();
459  fh_first_resty->Write();
460  fh_first_resz->Write();
461  fh_last_resx->Write();
462  fh_last_resy->Write();
463  fh_last_restx->Write();
464  fh_last_resty->Write();
465  fh_last_resz->Write();
466 
467  // Pulls at first and last TRD layer
468  fh_first_pullx->Scale(1. / fEvents);
469  fh_first_pully->Scale(1. / fEvents);
470  fh_first_pulltx->Scale(1. / fEvents);
471  fh_first_pullty->Scale(1. / fEvents);
472  fh_last_pullx->Scale(1. / fEvents);
473  fh_last_pully->Scale(1. / fEvents);
474  fh_last_pulltx->Scale(1. / fEvents);
475  fh_last_pullty->Scale(1. / fEvents);
476  fh_first_pullx->Write();
477  fh_first_pully->Write();
478  fh_first_pulltx->Write();
479  fh_first_pullty->Write();
480  fh_last_pullx->Write();
481  fh_last_pully->Write();
482  fh_last_pulltx->Write();
483  fh_last_pullty->Write();
484 
485  // Chi2/NDF of fit
486  fh_chi2ndf->Scale(1. / fEvents);
487  fh_chi2ndf->Write();
488 }
489 // -----------------------------------------------------------------------
490 
491 
ClassImp
ClassImp(CbmFitGlobalTracksQa)
CbmFitGlobalTracksQa::fArrayTrdPoint
TClonesArray * fArrayTrdPoint
Definition: CbmFitGlobalTracksQa.h:18
CbmFitGlobalTracksQa
Definition: CbmFitGlobalTracksQa.h:14
CbmFitGlobalTracksQa::fh_last_resty
TH1F * fh_last_resty
Definition: CbmFitGlobalTracksQa.h:33
CbmFitGlobalTracksQa::fh_first_resz
TH1F * fh_first_resz
Definition: CbmFitGlobalTracksQa.h:29
CbmTrack::GetNofHits
virtual Int_t GetNofHits() const
Definition: CbmTrack.h:53
CbmFitGlobalTracksQa::fArrayGlbTrack
TClonesArray * fArrayGlbTrack
Definition: CbmFitGlobalTracksQa.h:23
CbmFitGlobalTracksQa::Exec
virtual void Exec(Option_t *option="")
Definition: CbmFitGlobalTracksQa.cxx:167
CbmFitGlobalTracksQa::fh_first_resy
TH1F * fh_first_resy
Definition: CbmFitGlobalTracksQa.h:26
CbmGlobalTrack::GetParamLast
const FairTrackParam * GetParamLast() const
Definition: CbmGlobalTrack.h:44
CbmFitGlobalTracksQa::fh_first_restx
TH1F * fh_first_restx
Definition: CbmFitGlobalTracksQa.h:27
CbmFitGlobalTracksQa::~CbmFitGlobalTracksQa
virtual ~CbmFitGlobalTracksQa()
Definition: CbmFitGlobalTracksQa.cxx:106
rootMgr
static FairRootManager * rootMgr
Definition: CbmDeviceHitBuilderTof.cxx:72
CbmFitGlobalTracksQa::fh_last_resz
TH1F * fh_last_resz
Definition: CbmFitGlobalTracksQa.h:34
CbmTrdHit
data class for a reconstructed Energy-4D measurement in the TRD
Definition: CbmTrdHit.h:35
CbmFitGlobalTracksQa::CbmFitGlobalTracksQa
CbmFitGlobalTracksQa()
Definition: CbmFitGlobalTracksQa.cxx:36
CbmGlobalTrack.h
CbmFitGlobalTracksQa::fh_last_pulltx
TH1F * fh_last_pulltx
Definition: CbmFitGlobalTracksQa.h:41
CbmFitGlobalTracksQa::fArrayStsHit
TClonesArray * fArrayStsHit
Definition: CbmFitGlobalTracksQa.h:19
CbmFitGlobalTracksQa::fh_last_restx
TH1F * fh_last_restx
Definition: CbmFitGlobalTracksQa.h:32
CbmFitGlobalTracksQa::CreateHistogramms
void CreateHistogramms()
Definition: CbmFitGlobalTracksQa.cxx:389
CbmHit::GetRefId
Int_t GetRefId() const
Definition: CbmHit.h:72
CbmGlobalTrack::GetChi2
Double_t GetChi2() const
Definition: CbmGlobalTrack.h:47
CbmStsPoint
Definition: CbmStsPoint.h:27
CbmFitGlobalTracksQa::fh_first_pulltx
TH1F * fh_first_pulltx
Definition: CbmFitGlobalTracksQa.h:37
CbmStsTrack::GetNofHits
virtual Int_t GetNofHits() const
Definition: CbmStsTrack.h:76
CbmFitGlobalTracksQa::fh_first_resty
TH1F * fh_first_resty
Definition: CbmFitGlobalTracksQa.h:28
CbmStsHit
data class for a reconstructed 3-d hit in the STS
Definition: CbmStsHit.h:31
CbmFitGlobalTracksQa::fh_first_pully
TH1F * fh_first_pully
Definition: CbmFitGlobalTracksQa.h:36
CbmStsTrack.h
Data class for STS tracks.
CbmGlobalTrack::GetStsTrackIndex
Int_t GetStsTrackIndex() const
Definition: CbmGlobalTrack.h:38
CbmFitGlobalTracksQa::WriteHistogramms
void WriteHistogramms()
Definition: CbmFitGlobalTracksQa.cxx:444
CbmFitGlobalTracksQa::fh_last_pullx
TH1F * fh_last_pullx
Definition: CbmFitGlobalTracksQa.h:39
CbmFitGlobalTracksQa::fh_last_resy
TH1F * fh_last_resy
Definition: CbmFitGlobalTracksQa.h:31
CbmTrack::GetHitIndex
Int_t GetHitIndex(Int_t iHit) const
Definition: CbmTrack.h:54
CbmGlobalTrack::GetTrdTrackIndex
Int_t GetTrdTrackIndex() const
Definition: CbmGlobalTrack.h:39
CbmTrdHit.h
Class for hits in TRD detector.
CbmTrdPoint
Definition: CbmTrdPoint.h:23
CbmTrdTrack
Definition: CbmTrdTrack.h:22
CbmFitGlobalTracksQa::fh_first_pullx
TH1F * fh_first_pullx
Definition: CbmFitGlobalTracksQa.h:35
CbmFitGlobalTracksQa::Init
virtual InitStatus Init()
Definition: CbmFitGlobalTracksQa.cxx:113
CbmGlobalTrack
Definition: CbmGlobalTrack.h:26
CbmFitGlobalTracksQa::fEvents
Int_t fEvents
Definition: CbmFitGlobalTracksQa.h:24
CbmFitGlobalTracksQa::fh_last_resx
TH1F * fh_last_resx
Definition: CbmFitGlobalTracksQa.h:30
CbmStsPoint.h
CbmFitGlobalTracksQa::fArrayTrdTrack
TClonesArray * fArrayTrdTrack
Definition: CbmFitGlobalTracksQa.h:22
CbmGlobalTrack::GetNDF
Int_t GetNDF() const
Definition: CbmGlobalTrack.h:48
CbmTrdPoint.h
CbmFitGlobalTracksQa::fh_last_pully
TH1F * fh_last_pully
Definition: CbmFitGlobalTracksQa.h:40
pos
TVector3 pos
Definition: CbmMvdSensorDigiToHitTask.cxx:60
CbmFitGlobalTracksQa::Finish
virtual void Finish()
Definition: CbmFitGlobalTracksQa.cxx:381
CbmFitGlobalTracksQa::fArrayStsPoint
TClonesArray * fArrayStsPoint
Definition: CbmFitGlobalTracksQa.h:17
CbmFitGlobalTracksQa::fh_chi2ndf
TH1F * fh_chi2ndf
Definition: CbmFitGlobalTracksQa.h:43
CbmFitGlobalTracksQa::fArrayStsTrack
TClonesArray * fArrayStsTrack
Definition: CbmFitGlobalTracksQa.h:21
CbmGlobalTrack::GetParamFirst
const FairTrackParam * GetParamFirst() const
Definition: CbmGlobalTrack.h:43
CbmFitGlobalTracksQa::fh_first_resx
TH1F * fh_first_resx
Definition: CbmFitGlobalTracksQa.h:25
CbmTrdTrack.h
CbmFitGlobalTracksQa::fh_first_pullty
TH1F * fh_first_pullty
Definition: CbmFitGlobalTracksQa.h:38
CbmFitGlobalTracksQa::fArrayTrdHit
TClonesArray * fArrayTrdHit
Definition: CbmFitGlobalTracksQa.h:20
CbmStsTrack
Definition: CbmStsTrack.h:37
CbmFitGlobalTracksQa.h
CbmFitGlobalTracksQa::fh_last_pullty
TH1F * fh_last_pullty
Definition: CbmFitGlobalTracksQa.h:42
CbmStsHit.h
Data class for a reconstructed hit in the STS.