CbmRoot
CbmLitFitQa.cxx
Go to the documentation of this file.
1 
6 #include "CbmLitFitQa.h"
7 #include "../mc/CbmLitMCTrackCreator.h"
8 #include "CbmEvent.h"
9 #include "CbmGlobalTrack.h"
10 #include "CbmHistManager.h"
11 #include "CbmLitFitQaReport.h"
12 #include "CbmMCDataArray.h"
13 #include "CbmMCDataManager.h"
14 #include "CbmMCTrack.h"
15 #include "CbmMuchGeoScheme.h"
16 #include "CbmMvdHit.h"
17 #include "CbmStsAddress.h"
18 #include "CbmStsHit.h"
19 #include "CbmStsSetup.h"
20 #include "CbmStsTrack.h"
21 #include "CbmTrack.h"
22 #include "CbmTrackMatchNew.h"
23 #include "CbmTrdAddress.h"
24 #include "TClonesArray.h"
25 #include "TH1.h"
26 #include "TH1F.h"
27 #include "TH2F.h"
28 #include <TFile.h>
29 #include <boost/assign/list_of.hpp>
30 #include <vector>
31 using std::make_pair;
32 using std::pair;
33 using std::vector;
34 
36  : fIsFixedBounds(true)
37  , fMvdMinNofHits(0)
38  , fStsMinNofHits(0)
39  , fTrdMinNofHits(0)
40  , fMuchMinNofHits(0)
41  , fOutputDir("./test/")
42  , fHM(NULL)
43  , fPRangeMin(0.)
44  , fPRangeMax(10.)
45  , fPRangeBins(20)
46  , fGlobalTracks(NULL)
47  , fStsTracks(NULL)
48  , fStsTrackMatches(NULL)
49  , fStsHits(NULL)
50  , fMvdHits(NULL)
51  , fTrdTracks(NULL)
52  , fTrdTrackMatches(NULL)
53  , fTrdHits(NULL)
54  , fMuchTracks(NULL)
55  , fMuchTrackMatches(NULL)
56  , fMuchPixelHits(NULL)
57  , fMuchStripHits(NULL)
58  , fMCTracks(NULL)
59  , fEvents(NULL)
60  , fQuota(0.7)
61  , fPrimVertex(NULL)
62  , fKFFitter()
63  , fMCTrackCreator(NULL)
64  , fDet() {}
65 
67  if (fHM) delete fHM;
68 }
69 
70 InitStatus CbmLitFitQa::Init() {
71  fHM = new CbmHistManager();
76  fKFFitter.Init();
77  return kSUCCESS;
78 }
79 
80 void CbmLitFitQa::Exec(Option_t* opt) {
81  static Int_t nofEvents = 0;
82  nofEvents++;
83  std::cout << "CbmLitFitQa::Exec: event=" << nofEvents << std::endl;
84  fMCTrackCreator->Create(nofEvents - 1);
88 }
89 
91  TDirectory* oldir = gDirectory;
92  TFile* outFile = FairRootManager::Instance()->GetOutFile();
93  if (outFile != NULL) {
94  outFile->cd();
95  fHM->WriteToFile();
96  }
97  gDirectory->cd(oldir->GetPath());
98 
99  fHM->ShrinkEmptyBinsH1ByPattern("htf_.+_WrongCov_.+");
100  CbmSimulationReport* report = new CbmLitFitQaReport();
101  report->Create(fHM, fOutputDir);
102  delete report;
103 }
104 
106  FairRootManager* ioman = FairRootManager::Instance();
107  assert(ioman != NULL);
108 
109  fGlobalTracks = (TClonesArray*) ioman->GetObject("GlobalTrack");
110  fStsTracks = (TClonesArray*) ioman->GetObject("StsTrack");
111  fStsTrackMatches = (TClonesArray*) ioman->GetObject("StsTrackMatch");
112  fStsHits = (TClonesArray*) ioman->GetObject("StsHit");
113  fMvdHits = (TClonesArray*) ioman->GetObject("MvdHit");
114  fMuchTracks = (TClonesArray*) ioman->GetObject("MuchTrack");
115  fMuchTrackMatches = (TClonesArray*) ioman->GetObject("MuchTrackMatch");
116  fMuchPixelHits = (TClonesArray*) ioman->GetObject("MuchPixelHit");
117  fMuchStripHits = (TClonesArray*) ioman->GetObject("MuchStripHit");
118  fTrdTracks = (TClonesArray*) ioman->GetObject("TrdTrack");
119  fTrdTrackMatches = (TClonesArray*) ioman->GetObject("TrdTrackMatch");
120  fTrdHits = (TClonesArray*) ioman->GetObject("TrdHit");
121 
122  CbmMCDataManager* mcManager =
123  (CbmMCDataManager*) ioman->GetObject("MCDataManager");
124 
125  if (0 == mcManager)
126  LOG(fatal)
127  << "CbmMatchRecoToMC::ReadAndCreateDataBranches() NULL MCDataManager.";
128 
129  fMCTracks = mcManager->InitBranch("MCTrack");
130  fEvents = dynamic_cast<TClonesArray*>(ioman->GetObject("Event"));
131  // fPrimVertex = (CbmVertex*) ioman->GetObject("PrimaryVertex");
132  // Get pointer to PrimaryVertex object from IOManager if it exists
133  // The old name for the object is "PrimaryVertex" the new one
134  // "PrimaryVertex." Check first for the new name
135  fPrimVertex = dynamic_cast<CbmVertex*>(ioman->GetObject("PrimaryVertex."));
136  if (nullptr == fPrimVertex) {
137  fPrimVertex = dynamic_cast<CbmVertex*>(ioman->GetObject("PrimaryVertex"));
138  }
139  if (nullptr == fPrimVertex) {
140  // LOG(fatal) << "No primary vertex";
141  }
142 }
143 
145  Int_t nofGlobalTracks = fGlobalTracks->GetEntriesFast();
146  for (Int_t iTrack = 0; iTrack < nofGlobalTracks; iTrack++) {
147  const CbmGlobalTrack* globalTrack =
148  static_cast<const CbmGlobalTrack*>(fGlobalTracks->At(iTrack));
149  ProcessStsTrack(globalTrack->GetStsTrackIndex());
150  ProcessTrdTrack(globalTrack->GetTrdTrackIndex());
151  ProcessMuchTrack(globalTrack->GetMuchTrackIndex());
152  }
153 }
154 
155 void CbmLitFitQa::ProcessStsTrack(Int_t trackId) {
156  if (NULL == fStsTracks || NULL == fStsTrackMatches || trackId < 0) return;
157 
158  CbmTrackMatchNew* match =
159  static_cast<CbmTrackMatchNew*>(fStsTrackMatches->At(trackId));
160  if (match->GetNofLinks() < 1) return;
161  Int_t mcId = match->GetMatchedLink().GetIndex();
162  Int_t mcEventId = match->GetMatchedLink().GetEntry();
163  if (mcId < 0) return; // Ghost or fake track
164 
165  // Check correctness of reconstructed track
166  if (match->GetTrueOverAllHitsRatio() < fQuota) return;
167 
168  CbmStsTrack* track = static_cast<CbmStsTrack*>(fStsTracks->At(trackId));
169  Int_t nofStsHits = track->GetNofStsHits();
170  Int_t nofMvdHits = track->GetNofMvdHits();
171  if (nofStsHits < 1) return; // No hits in STS
172  if (nofMvdHits < fMvdMinNofHits || nofStsHits + nofMvdHits < fStsMinNofHits)
173  return; // cut on number of hits in track
174 
175  const CbmLitMCTrack& mcTrack = fMCTrackCreator->GetTrack(mcEventId, mcId);
176 
177  const FairTrackParam* firstParam = track->GetParamFirst();
178  const FairTrackParam* lastParam = track->GetParamLast();
179 
180  FillTrackParamHistogramm("htp_Sts_FirstParam", firstParam);
181  FillTrackParamHistogramm("htp_Sts_LastParam", lastParam);
182 
183  // Fill histograms for first track parameters
184  if (nofMvdHits > 0) { // first track parameters in MVD
185  const CbmMvdHit* firstHit =
186  static_cast<const CbmMvdHit*>(fMvdHits->At(track->GetMvdHitIndex(0)));
187  Int_t firstStation = firstHit->GetStationNr() - 1; // to start with 0
188  if (mcTrack.GetNofPointsAtStation(ECbmModuleId::kMvd, firstStation) > 0) {
189  const CbmLitMCPoint& firstPoint =
190  mcTrack.GetPointAtStation(ECbmModuleId::kMvd, firstStation, 0);
191  FillResidualsAndPulls(firstParam,
192  &firstPoint,
193  "htf_Sts_FirstParam_",
194  nofMvdHits + nofStsHits,
196  }
197  } else { // first track parameters in STS
198  const CbmStsHit* firstHit =
199  static_cast<const CbmStsHit*>(fStsHits->At(track->GetHitIndex(0)));
200  Int_t firstStation =
202  - 1; //firstHit->GetStationNr() - 1; // to start with 0
203  if (mcTrack.GetNofPointsAtStation(ECbmModuleId::kSts, firstStation) > 0) {
204  const CbmLitMCPoint& firstPoint =
205  mcTrack.GetPointAtStation(ECbmModuleId::kSts, firstStation, 0);
206  FillResidualsAndPulls(firstParam,
207  &firstPoint,
208  "htf_Sts_FirstParam_",
209  nofMvdHits + nofStsHits,
211  }
212  }
213 
214  // Fill histograms for last track parameters
215  const CbmStsHit* lastHit = static_cast<const CbmStsHit*>(
216  fStsHits->At(track->GetHitIndex(nofStsHits - 1)));
217  Int_t lastStation =
219  - 1; //lastHit->GetStationNr() - 1; // to start with 0
220  if (mcTrack.GetNofPointsAtStation(ECbmModuleId::kSts, lastStation) > 0) {
221  const CbmLitMCPoint& lastPoint =
222  mcTrack.GetPointAtStation(ECbmModuleId::kSts, lastStation, 0);
223  FillResidualsAndPulls(lastParam,
224  &lastPoint,
225  "htf_Sts_LastParam_",
226  nofMvdHits + nofStsHits,
228  }
229 }
230 
231 void CbmLitFitQa::ProcessTrdTrack(Int_t trackId) {
232  if (NULL == fTrdTracks || NULL == fTrdTrackMatches || trackId < 0) return;
233 
234  CbmTrackMatchNew* match =
235  static_cast<CbmTrackMatchNew*>(fTrdTrackMatches->At(trackId));
236  Int_t mcId = match->GetMatchedLink().GetIndex();
237  Int_t mcEntryId = match->GetMatchedLink().GetEntry();
238  if (mcId < 0) return; // Ghost or fake track
239 
240  // Check correctness of reconstructed track
241  if (match->GetTrueOverAllHitsRatio() < fQuota) return;
242 
243  CbmTrack* track = static_cast<CbmTrack*>(fTrdTracks->At(trackId));
244  Int_t nofHits = track->GetNofHits();
245  if (nofHits < 1) return; // No hits
246  if (nofHits < fTrdMinNofHits) return; // cut on number of hits in track
247 
248  const CbmLitMCTrack& mcTrack = fMCTrackCreator->GetTrack(mcEntryId, mcId);
249 
250  const FairTrackParam* firstParam = track->GetParamFirst();
251  const FairTrackParam* lastParam = track->GetParamLast();
252 
253  FillTrackParamHistogramm("htp_Trd_FirstParam", firstParam);
254  FillTrackParamHistogramm("htp_Trd_LastParam", lastParam);
255 
256  // Fill histograms for first track parameters
257  const CbmHit* firstHit =
258  static_cast<const CbmHit*>(fTrdHits->At(track->GetHitIndex(0)));
259  //Int_t firstStation = 10 * CbmTrdAddress::GetStationNr(firstHit->GetAddress()) + CbmTrdAddress::GetLayerNr(firstHit->GetAddress());
260  Int_t firstStation = CbmTrdAddress::GetLayerId(firstHit->GetAddress());
261  if (mcTrack.GetNofPointsAtStation(ECbmModuleId::kTrd, firstStation) > 0) {
262  const CbmLitMCPoint& firstPoint =
263  mcTrack.GetPointAtStation(ECbmModuleId::kTrd, firstStation, 0);
264  FillResidualsAndPulls(firstParam,
265  &firstPoint,
266  "htf_Trd_FirstParam_",
267  nofHits,
269  }
270 
271  // Fill histograms for last track parameters
272  const CbmHit* lastHit =
273  static_cast<const CbmHit*>(fTrdHits->At(track->GetHitIndex(nofHits - 1)));
274  //Int_t lastStation = 10 * CbmTrdAddress::GetStationNr(lastHit->GetAddress()) + CbmTrdAddress::GetLayerNr(lastHit->GetAddress());
275  Int_t lastStation = CbmTrdAddress::GetLayerId(lastHit->GetAddress());
276  if (mcTrack.GetNofPointsAtStation(ECbmModuleId::kTrd, lastStation) > 0) {
277  const CbmLitMCPoint& lastPoint =
278  mcTrack.GetPointAtStation(ECbmModuleId::kTrd, lastStation, 0);
280  lastParam, &lastPoint, "htf_Trd_LastParam_", nofHits, ECbmModuleId::kTrd);
281  }
282 }
283 
284 void CbmLitFitQa::ProcessMuchTrack(Int_t trackId) {
285  if (NULL == fMuchTracks || NULL == fMuchTrackMatches || trackId < 0) return;
286 
287  const CbmTrackMatchNew* match =
288  static_cast<const CbmTrackMatchNew*>(fMuchTrackMatches->At(trackId));
289  Int_t mcId = match->GetMatchedLink().GetIndex();
290  Int_t mcEventId = match->GetMatchedLink().GetEntry();
291  if (mcId < 0) return; // Ghost or fake track
292 
293  // Check correctness of reconstructed track
294  if (match->GetTrueOverAllHitsRatio() < fQuota) return;
295 
296  CbmTrack* track = static_cast<CbmTrack*>(fMuchTracks->At(trackId));
297  Int_t nofHits = track->GetNofHits();
298  if (nofHits < 1) return; // No hits
299  if (nofHits < fMuchMinNofHits) return; // cut on number of hits in track
300 
301  const CbmLitMCTrack& mcTrack = fMCTrackCreator->GetTrack(mcEventId, mcId);
302 
303  const FairTrackParam* firstParam = track->GetParamFirst();
304  const FairTrackParam* lastParam = track->GetParamLast();
305 
306  FillTrackParamHistogramm("htp_Much_FirstParam", firstParam);
307  FillTrackParamHistogramm("htp_Much_LastParam", lastParam);
308 
309  // Fill histograms for first track parameters
310  const CbmHit* firstHit =
311  static_cast<const CbmHit*>(fMuchPixelHits->At(track->GetHitIndex(0)));
312  // Int_t firstStation = firstHit->GetPlaneId();
313  Int_t firstStation =
315  + 10 * CbmMuchGeoScheme::GetLayerIndex(firstHit->GetAddress())
317  if (mcTrack.GetNofPointsAtStation(ECbmModuleId::kMuch, firstStation) > 0) {
318  const CbmLitMCPoint& firstPoint =
319  mcTrack.GetPointAtStation(ECbmModuleId::kMuch, firstStation, 0);
320  FillResidualsAndPulls(firstParam,
321  &firstPoint,
322  "htf_Much_FirstParam_",
323  nofHits,
325  }
326 
327  // Fill histograms for last track parameters
328  const CbmHit* lastHit = static_cast<const CbmHit*>(
329  fMuchPixelHits->At(track->GetHitIndex(nofHits - 1)));
330  // Int_t lastStation = lastHit->GetPlaneId();
331  Int_t lastStation =
335  if (mcTrack.GetNofPointsAtStation(ECbmModuleId::kMuch, lastStation) > 0) {
336  const CbmLitMCPoint& lastPoint =
337  mcTrack.GetPointAtStation(ECbmModuleId::kMuch, lastStation, 0);
338  FillResidualsAndPulls(lastParam,
339  &lastPoint,
340  "htf_Much_LastParam_",
341  nofHits,
343  }
344 }
345 
346 void CbmLitFitQa::FillResidualsAndPulls(const FairTrackParam* par,
347  const CbmLitMCPoint* mcPoint,
348  const string& histName,
349  Float_t wrongPar,
350  ECbmModuleId detId) {
351  // Residuals
352  Double_t resX = 0., resY = 0., resTx = 0., resTy = 0., resQp = 0.;
353  if (detId == ECbmModuleId::kMvd) { // Use MC average for MVD
354  resX = par->GetX() - mcPoint->GetX();
355  resY = par->GetY() - mcPoint->GetY();
356  resTx = par->GetTx() - mcPoint->GetTx();
357  resTy = par->GetTy() - mcPoint->GetTy();
358  resQp = par->GetQp() - mcPoint->GetQp();
359  } else if (detId == ECbmModuleId::kSts) { // Use MC average for STS
360  resX = par->GetX() - mcPoint->GetX();
361  resY = par->GetY() - mcPoint->GetY();
362  resTx = par->GetTx() - mcPoint->GetTx();
363  resTy = par->GetTy() - mcPoint->GetTy();
364  resQp = par->GetQp() - mcPoint->GetQp();
365  } else if (detId == ECbmModuleId::kTrd) { // Use MC out for TRD
366  resX = par->GetX() - mcPoint->GetXOut();
367  resY = par->GetY() - mcPoint->GetYOut();
368  resTx = par->GetTx() - mcPoint->GetTxOut();
369  resTy = par->GetTy() - mcPoint->GetTyOut();
370  resQp = par->GetQp() - mcPoint->GetQpOut();
371  } else if (detId == ECbmModuleId::kMuch) { // Use MC average for MUCH
372  resX = par->GetX() - mcPoint->GetX();
373  resY = par->GetY() - mcPoint->GetY();
374  resTx = par->GetTx() - mcPoint->GetTx();
375  resTy = par->GetTy() - mcPoint->GetTy();
376  resQp = par->GetQp() - mcPoint->GetQp();
377  }
378  Double_t mcP = mcPoint->GetP();
379 
380  fHM->H2(histName + "Res_X")->Fill(mcP, resX);
381  fHM->H2(histName + "Res_Y")->Fill(mcP, resY);
382  fHM->H2(histName + "Res_Tx")->Fill(mcP, resTx);
383  fHM->H2(histName + "Res_Ty")->Fill(mcP, resTy);
384  fHM->H2(histName + "Res_Qp")->Fill(mcP, resQp);
385 
386  // Pulls
387  Double_t C[15];
388  par->CovMatrix(C);
389 
390  Double_t sigmaX = (C[0] > 0.) ? std::sqrt(C[0]) : -1.;
391  Double_t sigmaY = (C[5] > 0.) ? std::sqrt(C[5]) : -1.;
392  Double_t sigmaTx = (C[9] > 0.) ? std::sqrt(C[9]) : -1.;
393  Double_t sigmaTy = (C[12] > 0.) ? std::sqrt(C[12]) : -1.;
394  Double_t sigmaQp = (C[14] > 0.) ? std::sqrt(C[14]) : -1.;
395  if (sigmaX < 0)
396  fHM->H2(histName + "WrongCov_X")->Fill(mcP, wrongPar);
397  else
398  fHM->H2(histName + "Pull_X")->Fill(mcP, resX / sigmaX);
399  if (sigmaY < 0)
400  fHM->H2(histName + "WrongCov_Y")->Fill(mcP, wrongPar);
401  else
402  fHM->H2(histName + "Pull_Y")->Fill(mcP, resY / sigmaY);
403  if (sigmaTx < 0)
404  fHM->H2(histName + "WrongCov_Tx")->Fill(mcP, wrongPar);
405  else
406  fHM->H2(histName + "Pull_Tx")->Fill(mcP, resTx / sigmaTx);
407  if (sigmaTy < 0)
408  fHM->H2(histName + "WrongCov_Ty")->Fill(mcP, wrongPar);
409  else
410  fHM->H2(histName + "Pull_Ty")->Fill(mcP, resTy / sigmaTy);
411  if (sigmaQp < 0)
412  fHM->H2(histName + "WrongCov_Qp")->Fill(mcP, wrongPar);
413  else
414  fHM->H2(histName + "Pull_Qp")->Fill(mcP, resQp / sigmaQp);
415 }
416 
417 void CbmLitFitQa::FillTrackParamHistogramm(const string& histName,
418  const FairTrackParam* par) {
419  fHM->H1(histName + "_X")->Fill(par->GetX());
420  fHM->H1(histName + "_Y")->Fill(par->GetY());
421  fHM->H1(histName + "_Z")->Fill(par->GetZ());
422  fHM->H1(histName + "_Tx")->Fill(par->GetTx());
423  fHM->H1(histName + "_Ty")->Fill(par->GetTy());
424  fHM->H1(histName + "_Qp")->Fill(par->GetQp());
425  Double_t p = (par->GetQp() != 0) ? std::fabs(1. / par->GetQp()) : 0.;
426  fHM->H1(histName + "_p")->Fill(p);
427  TVector3 mom;
428  par->Momentum(mom);
429  Double_t pt = std::sqrt(mom.X() * mom.X() + mom.Y() * mom.Y());
430  fHM->H1(histName + "_pt")->Fill(pt);
431 }
432 
434  if (fEvents) {
435  Int_t nofEvents = fEvents->GetEntriesFast();
436 
437  for (int i = 0; i < nofEvents; ++i)
438  ProcessTrackParamsAtVertex(static_cast<CbmEvent*>(fEvents->At(i)));
439  } else
441 }
442 
444  Int_t nofTracks = event ? event->GetNofData(ECbmDataType::kStsTrack)
445  : fStsTracks->GetEntriesFast();
446 
447  for (Int_t i = 0; i < nofTracks; ++i) {
448  Int_t iTrack = event ? event->GetIndex(ECbmDataType::kStsTrack, i) : i;
449  CbmStsTrack* track = static_cast<CbmStsTrack*>(fStsTracks->At(iTrack));
450 
451  CbmTrackMatchNew* match =
452  static_cast<CbmTrackMatchNew*>(fStsTrackMatches->At(iTrack));
453  if (match->GetNofLinks() < 1) continue;
454  Int_t mcId = match->GetMatchedLink().GetIndex();
455  Int_t mcEventId = match->GetMatchedLink().GetEntry();
456  if (mcId < 0) continue; // Ghost or fake track
457 
458  const CbmMCTrack* mcTrack =
459  static_cast<const CbmMCTrack*>(fMCTracks->Get(0, mcEventId, mcId));
460  if (mcTrack->GetMotherId() != -1) continue; // only for primaries
461  //if (mcTrack->GetP() < 1.) continue; // only for momentum large 1 GeV/c
462 
463  // Check correctness of reconstructed track
464  if (match->GetTrueOverAllHitsRatio() < fQuota) continue;
465 
466  CbmVertex* primVertex = event ? event->GetVertex() : fPrimVertex;
467  fKFFitter.DoFit(track, 211);
468  Double_t chiPrimary = fKFFitter.GetChiToVertex(track, primVertex);
469  fHM->H1("htf_ChiPrimary")->Fill(chiPrimary);
470 
471  FairTrackParam vtxTrack;
472  fKFFitter.FitToVertex(track, primVertex, &vtxTrack);
473 
474  TVector3 momMC, momRec;
475  mcTrack->GetMomentum(momMC);
476  vtxTrack.Momentum(momRec);
477 
478  Double_t dpp = 100. * (momMC.Mag() - momRec.Mag()) / momMC.Mag();
479  fHM->H1("htf_MomRes_Mom")->Fill(momMC.Mag(), dpp);
480  }
481 }
482 
484  Int_t nofTracks = fGlobalTracks->GetEntriesFast();
485 
486  for (int i = 0; i < nofTracks; ++i) {
487  CbmGlobalTrack* globalTrack =
488  static_cast<CbmGlobalTrack*>(fGlobalTracks->At(i));
489  Int_t stsId = globalTrack->GetStsTrackIndex();
490  CbmStsTrack* stsTrack = static_cast<CbmStsTrack*>(fStsTracks->At(stsId));
491  CbmTrackMatchNew* match =
492  static_cast<CbmTrackMatchNew*>(fStsTrackMatches->At(stsId));
493 
494  if (match->GetNofLinks() < 1) continue;
495 
496  Int_t mcId = match->GetMatchedLink().GetIndex();
497  Int_t mcEventId = match->GetMatchedLink().GetEntry();
498 
499  if (mcId < 0) continue; // Ghost or fake track
500 
501  const CbmMCTrack* mcTrack =
502  static_cast<const CbmMCTrack*>(fMCTracks->Get(0, mcEventId, mcId));
503 
504  if (mcTrack->GetMotherId() != -1) continue; // only for primaries
505 
506  // Check correctness of reconstructed track
507  if (match->GetTrueOverAllHitsRatio() < fQuota) continue;
508 
509  const CbmTrackParam* vtxParam = globalTrack->GetParamVertex();
510  TVector3 momMC;
511  mcTrack->GetMomentum(momMC);
512 
513  fHM->H1("htp_PrimaryVertexResidualPx")->Fill(vtxParam->GetPx() - momMC.x());
514  fHM->H1("htp_PrimaryVertexResidualPy")->Fill(vtxParam->GetPy() - momMC.y());
515  fHM->H1("htp_PrimaryVertexResidualPz")->Fill(vtxParam->GetPz() - momMC.z());
516 
517  fHM->H1("htp_PrimaryVertexPullPx")
518  ->Fill((vtxParam->GetPx() - momMC.x()) / vtxParam->GetDpx());
519  fHM->H1("htp_PrimaryVertexPullPy")
520  ->Fill((vtxParam->GetPy() - momMC.y()) / vtxParam->GetDpy());
521  fHM->H1("htp_PrimaryVertexPullPz")
522  ->Fill((vtxParam->GetPz() - momMC.z()) / vtxParam->GetDpz());
523 
524  fHM->H1("htp_PrimaryVertexResidualTx")
525  ->Fill(vtxParam->GetTx() - mcTrack->GetPx() / mcTrack->GetPz());
526  fHM->H1("htp_PrimaryVertexResidualTy")
527  ->Fill(vtxParam->GetTy() - mcTrack->GetPy() / mcTrack->GetPz());
528  fHM->H1("htp_PrimaryVertexResidualQp")
529  ->Fill(vtxParam->GetQp() - 1 / mcTrack->GetP());
530 
531  Double_t cov[15];
532  vtxParam->CovMatrix(cov);
533  fHM->H1("htp_PrimaryVertexPullTx")
534  ->Fill((vtxParam->GetTx() - mcTrack->GetPx() / mcTrack->GetPz())
535  / TMath::Sqrt(cov[9]));
536  fHM->H1("htp_PrimaryVertexPullTy")
537  ->Fill((vtxParam->GetTy() - mcTrack->GetPy() / mcTrack->GetPz())
538  / TMath::Sqrt(cov[12]));
539  fHM->H1("htp_PrimaryVertexPullQp")
540  ->Fill((vtxParam->GetQp() - 1 / mcTrack->GetP()) / TMath::Sqrt(cov[14]));
541  }
542 }
543 
547 
550 
553 
554  // Momentum resolution vs momwntum
555  fHM->Add("htf_MomRes_Mom",
556  new TH2F("htf_MomRes_Mom",
557  "htf_MomRes_Mom;P [GeV/c];dP/P [%]",
558  fPRangeBins,
559  fPRangeMin,
560  fPRangeMax,
561  100,
562  -3.,
563  3.));
564  fHM->Add(
565  "htf_ChiPrimary",
566  new TH1F("htf_ChiPrimary", "htf_ChiPrimary;#chi_{primary}", 100, 0, 10));
567  fHM->Add("htp_PrimaryVertexResidualPx",
568  new TH1F("htp_PrimaryVertexResidualPx",
569  "htp_PrimaryVertexResidualPx;[cm]",
570  200,
571  -1.,
572  1.));
573  fHM->Add("htp_PrimaryVertexResidualPy",
574  new TH1F("htp_PrimaryVertexResidualPy",
575  "htp_PrimaryVertexResidualPy;[cm]",
576  200,
577  -1.,
578  1.));
579  fHM->Add("htp_PrimaryVertexResidualPz",
580  new TH1F("htp_PrimaryVertexResidualPz",
581  "htp_PrimaryVertexResidualPz;[cm]",
582  200,
583  -1.,
584  1.));
585  fHM->Add(
586  "htp_PrimaryVertexPullPx",
587  new TH1F(
588  "htp_PrimaryVertexPullPx", "htp_PrimaryVertexPullPx;[cm]", 200, -5., 5.));
589  fHM->Add(
590  "htp_PrimaryVertexPullPy",
591  new TH1F(
592  "htp_PrimaryVertexPullPy", "htp_PrimaryVertexPullPy;[cm]", 200, -5., 5.));
593  fHM->Add(
594  "htp_PrimaryVertexPullPz",
595  new TH1F(
596  "htp_PrimaryVertexPullPz", "htp_PrimaryVertexPullPz;[cm]", 200, -5., 5.));
597 
598  fHM->Add("htp_PrimaryVertexResidualTx",
599  new TH1F("htp_PrimaryVertexResidualTx",
600  "htp_PrimaryVertexResidualTx;[cm]",
601  200,
602  -1.,
603  1.));
604  fHM->Add("htp_PrimaryVertexResidualTy",
605  new TH1F("htp_PrimaryVertexResidualTy",
606  "htp_PrimaryVertexResidualTy;[cm]",
607  200,
608  -1.,
609  1.));
610  fHM->Add("htp_PrimaryVertexResidualQp",
611  new TH1F("htp_PrimaryVertexResidualQ[",
612  "htp_PrimaryVertexResidualQp;[cm]",
613  200,
614  -1.,
615  1.));
616  fHM->Add(
617  "htp_PrimaryVertexPullTx",
618  new TH1F(
619  "htp_PrimaryVertexPullTx", "htp_PrimaryVertexPullTx;[cm]", 200, -5., 5.));
620  fHM->Add(
621  "htp_PrimaryVertexPullTy",
622  new TH1F(
623  "htp_PrimaryVertexPullTy", "htp_PrimaryVertexPullTy;[cm]", 200, -5., 5.));
624  fHM->Add(
625  "htp_PrimaryVertexPullQp",
626  new TH1F(
627  "htp_PrimaryVertexPullQp", "htp_PrimaryVertexPullQp;[cm]", 200, -5., 5.));
628 }
629 
631  const string& detName) {
632  if (!fDet.GetDet(detId)) return;
633 
634  // Parameter names of the state vector (x, y, tx, ty, q/p)
635  string parameterNames[] = {"X", "Y", "Tx", "Ty", "Qp"};
636 
637  // Axis titles for residual, pull and wrong covariance histograms
638  string xTitles[] = {"Residual X [cm]",
639  "Residual Y [cm]",
640  "Residual Tx",
641  "Residual Ty",
642  "Residual q/p [(GeV/c)^{-1}]",
643  "Pull X",
644  "Pull Y",
645  "Pull Tx",
646  "Pull Ty",
647  "Pull q/p",
648  "Number of hits",
649  "Number of hits",
650  "Number of hits",
651  "Number of hits",
652  "Number of hits"};
653 
654  // Category names: Residual, Pull, Wrong Covariance
655  string catNames[] = {"Res", "Pull", "WrongCov"};
656 
657  vector<Int_t> bins(15, 200);
658  vector<pair<Float_t, Float_t>> bounds;
659  if (fIsFixedBounds) {
660  vector<pair<Float_t, Float_t>> tmp =
661  boost::assign::list_of(make_pair(-1., 1.)) // X residual
662  (make_pair(-1., 1.)) // Y residual
663  (make_pair(-.01, .01)) // Tx residual
664  (make_pair(-.01, .01)) // Ty residual
665  (make_pair(-.1, .1)) // Qp residual
666  (make_pair(-5., 5.)) // X pull
667  (make_pair(-5., 5.)) // Y pull
668  (make_pair(-5., 5.)) // Tx pull
669  (make_pair(-5., 5.)) // Ty pull
670  (make_pair(-7., 7.)) // Qp pull
671  (make_pair(0., 200.)) // X wrong covariance
672  (make_pair(0., 200.)) // Y wrong covariance
673  (make_pair(0., 200.)) // Tx wrong covariance
674  (make_pair(0., 200.)) // Ty wrong covariance
675  (make_pair(0., 200.)); // Qp wrong covariance
676  bounds = tmp;
677  } else {
678  bounds.assign(15, make_pair(0., 0.));
679  }
680 
681  Double_t momentumMin = 0.;
682  Double_t momentumMax = 10.;
683  Int_t nofMomentumBins = 20;
684  string momentumAxisTitle = "P [GeV/c]";
685 
686  // [0] - for the first track parameter, [1] - for the last track parameter
687  for (Int_t i = 0; i < 2; i++) {
688  string trackParamName = (i == 0) ? "FirstParam" : "LastParam";
689  // [0] - "Res", [1] - "Pull", [2] - "WrongCov"
690  for (Int_t iCat = 0; iCat < 3; iCat++) {
691  for (Int_t iPar = 0; iPar < 5; iPar++) {
692  string histName = "htf_" + detName + "_" + trackParamName + "_"
693  + catNames[iCat] + "_" + parameterNames[iPar];
694  Int_t histId = iCat * 5 + iPar;
695  fHM->Add(histName,
696  new TH2F(histName.c_str(),
697  string(histName + ";" + momentumAxisTitle + +";"
698  + xTitles[histId] + ";Counter")
699  .c_str(),
700  nofMomentumBins,
701  momentumMin,
702  momentumMax,
703  bins[histId],
704  bounds[histId].first,
705  bounds[histId].second));
706  }
707  }
708  }
709 }
710 
712  const string& detName) {
713  if (!fDet.GetDet(detId)) return;
714 
715  // Parameter names of the state vector (x, y, tx, ty, q/p)
716  string parameterNames[] = {"X", "Y", "Z", "Tx", "Ty", "Qp", "p", "pt"};
717 
718  // Axis titles for state vector Components
719  string xTitles[] = {"X [cm]",
720  "Y [cm]",
721  "Z [cm]",
722  "Tx",
723  "Ty",
724  "q/p [(GeV/c)^{-1}]",
725  "Momentum [GeV/c]",
726  "P_{t} [GeV/c]"};
727 
728  vector<Int_t> bins(8, 200);
729  vector<pair<Float_t, Float_t>> bounds;
730  if (fIsFixedBounds) {
731  vector<pair<Float_t, Float_t>> tmp =
732  boost::assign::list_of(make_pair(-500., 500.)) // X
733  (make_pair(-500., 500.)) // Y
734  (make_pair(-500., 500.)) // Z
735  (make_pair(-1., 1.)) // Tx
736  (make_pair(-1., 1.)) // Ty
737  (make_pair(-2., 2.)) // Qp
738  (make_pair(0., 25.)) // Momentum
739  (make_pair(0., 5.)); // Pt
740  bounds = tmp;
741  } else {
742  bounds.assign(8, make_pair(0., 0.));
743  }
744 
745  // [0] - for the first track parameter, [1] - for the last track parameter
746  for (Int_t i = 0; i < 2; i++) {
747  string trackParamName = (i == 0) ? "FirstParam" : "LastParam";
748  for (Int_t iPar = 0; iPar < 8; iPar++) {
749  string histName =
750  "htp_" + detName + "_" + trackParamName + "_" + parameterNames[iPar];
751  fHM->Add(
752  histName,
753  new TH1F(histName.c_str(),
754  string(histName + ";" + xTitles[iPar] + ";Counter").c_str(),
755  bins[iPar],
756  bounds[iPar].first,
757  bounds[iPar].second));
758  }
759  }
760 }
CbmLitFitQa::fPRangeMax
Double_t fPRangeMax
Definition: CbmLitFitQa.h:124
CbmLitFitQa::fQuota
Double_t fQuota
Definition: CbmLitFitQa.h:145
CbmLitFitQa::CreateResidualAndPullHistograms
void CreateResidualAndPullHistograms(ECbmModuleId detId, const string &detName)
Definition: CbmLitFitQa.cxx:630
CbmMCTrack::GetMotherId
Int_t GetMotherId() const
Definition: CbmMCTrack.h:71
CbmMatch::GetMatchedLink
const CbmLink & GetMatchedLink() const
Definition: CbmMatch.h:37
CbmMCDataManager::GetObject
CbmMCDataObject * GetObject(const char *name)
Definition: CbmMCDataManager.cxx:137
CbmMCTrack::GetMomentum
void GetMomentum(TVector3 &momentum) const
Definition: CbmMCTrack.h:172
CbmMCDataManager.h
CbmLitMCTrack::GetNofPointsAtStation
UInt_t GetNofPointsAtStation(ECbmModuleId detId, Int_t stationId) const
Return number of MC points for specified detector ID and station ID.
Definition: CbmLitMCTrack.h:229
CbmLitFitQa::fMuchStripHits
TClonesArray * fMuchStripHits
Definition: CbmLitFitQa.h:141
CbmTrack::GetParamLast
const FairTrackParam * GetParamLast() const
Definition: CbmTrack.h:62
CbmTrackMatchNew::GetTrueOverAllHitsRatio
Double_t GetTrueOverAllHitsRatio() const
Definition: CbmTrackMatchNew.h:35
CbmHistManager::WriteToFile
void WriteToFile()
Write all histograms to current opened file.
Definition: core/base/CbmHistManager.cxx:103
CbmTrdAddress.h
Helper class to convert unique channel ID back and forth.
CbmGlobalTrack::GetMuchTrackIndex
Int_t GetMuchTrackIndex() const
Definition: CbmGlobalTrack.h:40
CbmLitFitQa::~CbmLitFitQa
virtual ~CbmLitFitQa()
Destructor.
Definition: CbmLitFitQa.cxx:66
CbmMatch::GetNofLinks
Int_t GetNofLinks() const
Definition: CbmMatch.h:38
CbmTrack::GetNofHits
virtual Int_t GetNofHits() const
Definition: CbmTrack.h:53
CbmTrackParam::GetPx
Double_t GetPx() const
Definition: CbmTrackParam.h:39
CbmLitFitQa::ProcessGlobalTracks
void ProcessGlobalTracks()
Definition: CbmLitFitQa.cxx:144
CbmLitFitQa::fTrdMinNofHits
Int_t fTrdMinNofHits
Definition: CbmLitFitQa.h:118
CbmStsSetup.h
sqrt
friend F32vec4 sqrt(const F32vec4 &a)
Definition: L1/vectors/P4_F32vec4.h:41
CbmLitMCPoint::GetP
Double_t GetP() const
Definition: CbmLitMCPoint.h:57
CbmLitFitQa::fMuchTrackMatches
TClonesArray * fMuchTrackMatches
Definition: CbmLitFitQa.h:139
CbmLitFitQa::fHM
CbmHistManager * fHM
Definition: CbmLitFitQa.h:127
CbmHistManager::ShrinkEmptyBinsH1ByPattern
void ShrinkEmptyBinsH1ByPattern(const std::string &pattern)
Shrink empty bins in H1.
Definition: core/base/CbmHistManager.cxx:174
CbmStsKFTrackFitter::GetChiToVertex
Double_t GetChiToVertex(CbmStsTrack *track, CbmVertex *vtx=0)
Definition: CbmStsKFTrackFitter.cxx:164
CbmLitFitQa::fStsTrackMatches
TClonesArray * fStsTrackMatches
Definition: CbmLitFitQa.h:132
CbmLitFitQa::fKFFitter
CbmStsKFTrackFitter fKFFitter
Definition: CbmLitFitQa.h:149
CbmLitMCPoint::GetY
Double_t GetY() const
Definition: CbmLitMCPoint.h:50
CbmLitFitQa::fGlobalTracks
TClonesArray * fGlobalTracks
Definition: CbmLitFitQa.h:130
i
int i
Definition: L1/vectors/P4_F32vec4.h:25
CbmLitMCTrackCreator::Instance
static CbmLitMCTrackCreator * Instance()
Singleton instance.
Definition: CbmLitMCTrackCreator.cxx:63
CbmLitFitQa::fPRangeBins
Int_t fPRangeBins
Definition: CbmLitFitQa.h:125
CbmMvdHit::GetStationNr
virtual Int_t GetStationNr() const
Definition: CbmMvdHit.h:61
CbmMCDataManager::InitBranch
CbmMCDataArray * InitBranch(const char *name)
Definition: CbmMCDataManager.cxx:106
CbmLitMCPoint::GetTx
Double_t GetTx() const
Definition: CbmLitMCPoint.h:55
CbmStsSetup::Instance
static CbmStsSetup * Instance()
Definition: CbmStsSetup.cxx:293
CbmLitMCTrackCreator::GetTrack
const CbmLitMCTrack & GetTrack(int mcEventId, int mcId) const
Return CbmLitMCTrack by its index.
Definition: CbmLitMCTrackCreator.h:74
CbmHistManager::H2
TH2 * H2(const std::string &name) const
Return pointer to TH2 histogram.
Definition: CbmHistManager.h:190
CbmLitMCPoint::GetX
Double_t GetX() const
Definition: CbmLitMCPoint.h:49
CbmMCDataArray.h
CbmLitFitQa::fPRangeMin
Double_t fPRangeMin
Definition: CbmLitFitQa.h:123
CbmGlobalTrack.h
ECbmModuleId::kMvd
@ kMvd
Micro-Vertex Detector.
ECbmModuleId
ECbmModuleId
Definition: CbmDefs.h:33
CbmTrdAddress::GetLayerId
static UInt_t GetLayerId(UInt_t address)
Return layer ID from address.
Definition: CbmTrdAddress.h:69
CbmLitFitQa::FillTrackParamHistogramm
void FillTrackParamHistogramm(const string &histName, const FairTrackParam *par)
Definition: CbmLitFitQa.cxx:417
CbmMCTrack::GetPx
Double_t GetPx() const
Definition: CbmMCTrack.h:72
CbmMCTrack::GetPy
Double_t GetPy() const
Definition: CbmMCTrack.h:73
CbmLitMCPoint::GetTxOut
Double_t GetTxOut() const
Definition: CbmLitMCPoint.h:87
CbmLitFitQa::fOutputDir
string fOutputDir
Definition: CbmLitFitQa.h:121
CbmLitFitQa::fTrdTracks
TClonesArray * fTrdTracks
Definition: CbmLitFitQa.h:135
CbmHistManager.h
Histogram manager.
CbmStsKFTrackFitter::DoFit
Int_t DoFit(CbmStsTrack *track, Int_t pidHypo=211)
Definition: CbmStsKFTrackFitter.cxx:79
CbmLitMCPoint::GetTyOut
Double_t GetTyOut() const
Definition: CbmLitMCPoint.h:88
CbmMvdHit
Definition: CbmMvdHit.h:29
CbmGlobalTrack::GetParamVertex
const CbmTrackParam * GetParamVertex() const
Definition: CbmGlobalTrack.h:45
CbmLitFitQa::FillResidualsAndPulls
void FillResidualsAndPulls(const FairTrackParam *par, const CbmLitMCPoint *mcPoint, const string &histName, Float_t wrongPar, ECbmModuleId detId)
Definition: CbmLitFitQa.cxx:346
CbmMuchGeoScheme::GetLayerSideIndex
static Int_t GetLayerSideIndex(Int_t address)
Definition: CbmMuchGeoScheme.h:74
CbmStsTrack::GetNofMvdHits
Int_t GetNofMvdHits() const
Definition: CbmStsTrack.h:84
CbmLitFitQa::ProcessStsTrack
void ProcessStsTrack(Int_t trackId)
Definition: CbmLitFitQa.cxx:155
CbmLitFitQa::fTrdHits
TClonesArray * fTrdHits
Definition: CbmLitFitQa.h:137
CbmEvent.h
CbmStsHit
data class for a reconstructed 3-d hit in the STS
Definition: CbmStsHit.h:31
CbmLitFitQa::ProcessTrackMomentumAtVertex
void ProcessTrackMomentumAtVertex()
Definition: CbmLitFitQa.cxx:483
ECbmDataType::kStsTrack
@ kStsTrack
CbmLitFitQa::fMCTracks
CbmMCDataArray * fMCTracks
Definition: CbmLitFitQa.h:142
CbmTrack
Definition: CbmTrack.h:32
CbmStsTrack.h
Data class for STS tracks.
CbmGlobalTrack::GetStsTrackIndex
Int_t GetStsTrackIndex() const
Definition: CbmGlobalTrack.h:38
CbmMuchGeoScheme::GetLayerIndex
static Int_t GetLayerIndex(Int_t address)
Definition: CbmMuchGeoScheme.h:71
CbmHistManager
Histogram manager.
Definition: CbmHistManager.h:41
CbmMuchGeoScheme::GetStationIndex
static Int_t GetStationIndex(Int_t address)
Definition: CbmMuchGeoScheme.h:68
CbmLitFitQa::fMvdMinNofHits
Int_t fMvdMinNofHits
Definition: CbmLitFitQa.h:116
CbmLitMCPoint
Monte-Carlo point.
Definition: CbmLitMCPoint.h:21
CbmLitFitQa::fStsMinNofHits
Int_t fStsMinNofHits
Definition: CbmLitFitQa.h:117
CbmLitFitQa::ProcessTrdTrack
void ProcessTrdTrack(Int_t trackId)
Definition: CbmLitFitQa.cxx:231
CbmTrack.h
CbmLitFitQa::ProcessMuchTrack
void ProcessMuchTrack(Int_t trackId)
Definition: CbmLitFitQa.cxx:284
CbmLitFitQa::fDet
CbmLitDetectorSetup fDet
Definition: CbmLitFitQa.h:153
CbmStsSetup::GetStationNumber
Int_t GetStationNumber(Int_t address)
Definition: CbmStsSetup.cxx:187
CbmMCDataArray::Get
TObject * Get(const CbmLink *lnk)
Definition: CbmMCDataArray.h:47
CbmLitFitQa::CreateHistograms
void CreateHistograms()
Definition: CbmLitFitQa.cxx:544
CbmTrack::GetHitIndex
Int_t GetHitIndex(Int_t iHit) const
Definition: CbmTrack.h:54
CbmGlobalTrack::GetTrdTrackIndex
Int_t GetTrdTrackIndex() const
Definition: CbmGlobalTrack.h:39
CbmHit::GetAddress
Int_t GetAddress() const
Definition: CbmHit.h:73
CbmTrackMatchNew.h
CbmVertex
Definition: CbmVertex.h:26
CbmLitDetectorSetup::GetDet
bool GetDet(ECbmModuleId detId) const
Return detector presence in setup.
Definition: CbmLitDetectorSetup.cxx:27
CbmSimulationReport::Create
void Create(CbmHistManager *histManager, const std::string &outputDir)
Main function which creates report data.
CbmLitFitQa::Exec
virtual void Exec(Option_t *opt)
Inherited from FairTask.
Definition: CbmLitFitQa.cxx:80
CbmTrackParam::GetPy
Double_t GetPy() const
Definition: CbmTrackParam.h:40
CbmLitFitQa
Track fit QA for track reconstruction.
Definition: CbmLitFitQa.h:38
CbmLitFitQa::fMuchMinNofHits
Int_t fMuchMinNofHits
Definition: CbmLitFitQa.h:119
CbmLitMCPoint::GetQp
Double_t GetQp() const
Definition: CbmLitMCPoint.h:60
CbmLitFitQa::fTrdTrackMatches
TClonesArray * fTrdTrackMatches
Definition: CbmLitFitQa.h:136
CbmLitFitQa::fMvdHits
TClonesArray * fMvdHits
Definition: CbmLitFitQa.h:134
ClassImp
ClassImp(CbmConverterManager) InitStatus CbmConverterManager
Definition: CbmConverterManager.cxx:12
CbmHistManager::H1
TH1 * H1(const std::string &name) const
Return pointer to TH1 histogram.
Definition: CbmHistManager.h:170
CbmLitFitQa::CreateTrackParamHistograms
void CreateTrackParamHistograms(ECbmModuleId detId, const string &detName)
Definition: CbmLitFitQa.cxx:711
CbmTrack::GetParamFirst
const FairTrackParam * GetParamFirst() const
Definition: CbmTrack.h:61
CbmTrackParam::GetDpz
Double_t GetDpz() const
Definition: CbmTrackParam.h:44
CbmLitFitQa::CbmLitFitQa
CbmLitFitQa()
Constructor.
Definition: CbmLitFitQa.cxx:35
CbmLitFitQa::Init
virtual InitStatus Init()
Inherited from FairTask.
Definition: CbmLitFitQa.cxx:70
CbmMCDataManager
Task class creating and managing CbmMCDataArray objects.
Definition: CbmMCDataManager.h:27
CbmStsKFTrackFitter::Init
void Init()
Definition: CbmStsKFTrackFitter.cxx:29
ECbmModuleId::kTrd
@ kTrd
Transition Radiation Detector.
CbmGlobalTrack
Definition: CbmGlobalTrack.h:26
CbmLitFitQa::fPrimVertex
CbmVertex * fPrimVertex
Definition: CbmLitFitQa.h:147
CbmLitMCPoint::GetYOut
Double_t GetYOut() const
Definition: CbmLitMCPoint.h:82
CbmLitFitQa::ProcessTrackParamsAtVertex
void ProcessTrackParamsAtVertex()
Definition: CbmLitFitQa.cxx:433
CbmLitMCTrack
Monte-Carlo track.
Definition: CbmLitMCTrack.h:30
CbmLitMCPoint::GetXOut
Double_t GetXOut() const
Definition: CbmLitMCPoint.h:81
CbmLitFitQa::fEvents
TClonesArray * fEvents
Definition: CbmLitFitQa.h:143
CbmMCTrack.h
CbmLitFitQa::fMCTrackCreator
CbmLitMCTrackCreator * fMCTrackCreator
Definition: CbmLitFitQa.h:151
CbmMCTrack
Definition: CbmMCTrack.h:34
CbmLitFitQa::fMuchPixelHits
TClonesArray * fMuchPixelHits
Definition: CbmLitFitQa.h:140
CbmLitFitQa::fMuchTracks
TClonesArray * fMuchTracks
Definition: CbmLitFitQa.h:138
CbmHit
Definition: CbmHit.h:38
CbmLitFitQa::ReadDataBranches
void ReadDataBranches()
Reads data branches.
Definition: CbmLitFitQa.cxx:105
CbmMvdHit.h
fabs
friend F32vec4 fabs(const F32vec4 &a)
Definition: L1/vectors/P4_F32vec4.h:60
CbmTrackParam
Definition: CbmTrackParam.h:22
CbmLitMCPoint::GetQpOut
Double_t GetQpOut() const
Definition: CbmLitMCPoint.h:92
CbmLitFitQa::fStsHits
TClonesArray * fStsHits
Definition: CbmLitFitQa.h:133
CbmTrackParam::GetPz
Double_t GetPz() const
Definition: CbmTrackParam.h:41
CbmStsKFTrackFitter::FitToVertex
Double_t FitToVertex(CbmStsTrack *track, CbmVertex *vtx, FairTrackParam *v_track)
Definition: CbmStsKFTrackFitter.cxx:200
CbmTrackParam::GetDpy
Double_t GetDpy() const
Definition: CbmTrackParam.h:43
CbmEvent
Class characterising one event by a collection of links (indices) to data objects,...
Definition: CbmEvent.h:30
CbmStsTrack::GetMvdHitIndex
Int_t GetMvdHitIndex(Int_t iHit) const
Definition: CbmStsTrack.h:70
CbmLitMCTrack::GetPointAtStation
const CbmLitMCPoint & GetPointAtStation(ECbmModuleId detId, Int_t stationId, Int_t index) const
Return MC point for specified detector id and point index.
Definition: CbmLitMCTrack.h:217
CbmStsAddress.h
CbmLitFitQa::fStsTracks
TClonesArray * fStsTracks
Definition: CbmLitFitQa.h:131
CbmLitFitQa::fIsFixedBounds
Bool_t fIsFixedBounds
Definition: CbmLitFitQa.h:114
CbmLitFitQa.h
Track fit QA for track reconstruction.
CbmSimulationReport
Base class for simulation reports.
Definition: CbmSimulationReport.h:28
ECbmModuleId::kMuch
@ kMuch
Muon detection system.
CbmTrackMatchNew
Definition: CbmTrackMatchNew.h:19
CbmMuchGeoScheme.h
CbmStsTrack
Definition: CbmStsTrack.h:37
CbmTrackParam::GetDpx
Double_t GetDpx() const
Definition: CbmTrackParam.h:42
ECbmModuleId::kSts
@ kSts
Silicon Tracking System.
CbmLitMCTrackCreator::Create
void Create(Int_t eventNum)
Creates array of CbmLitMCTracks for current event.
Definition: CbmLitMCTrackCreator.cxx:68
CbmLitMCPoint::GetTy
Double_t GetTy() const
Definition: CbmLitMCPoint.h:56
CbmMCTrack::GetP
Double_t GetP() const
Definition: CbmMCTrack.h:100
CbmMCTrack::GetPz
Double_t GetPz() const
Definition: CbmMCTrack.h:74
CbmLitFitQaReport.h
Create report for fit QA.
CbmStsHit.h
Data class for a reconstructed hit in the STS.
CbmLitFitQaReport
Create report for fit QA.
Definition: CbmLitFitQaReport.h:21
CbmLitDetectorSetup::DetermineSetup
void DetermineSetup()
Determines detector presence using TGeoManager.
Definition: CbmLitDetectorSetup.cxx:79
CbmHistManager::Add
void Add(const std::string &name, TNamed *object)
Add new named object to manager.
Definition: CbmHistManager.h:58
CbmStsTrack::GetNofStsHits
Int_t GetNofStsHits() const
Definition: CbmStsTrack.h:90
CbmLitFitQa::Finish
virtual void Finish()
Inherited from FairTask.
Definition: CbmLitFitQa.cxx:90