Belle II Software  release-05-02-19
BKLMTrackingModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Yinghui GUAN *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 /* Own header. */
12 #include <klm/bklm/modules/bklmTracking/BKLMTrackingModule.h>
13 
14 /* KLM headers. */
15 #include <klm/bklm/geometry/GeometryPar.h>
16 #include <klm/bklm/modules/bklmTracking/BKLMTrackFinder.h>
17 
18 /* Belle 2 headers. */
19 #include <framework/dataobjects/EventMetaData.h>
20 #include <framework/datastore/StoreObjPtr.h>
21 #include <framework/logging/Logger.h>
22 
23 using namespace std;
24 using namespace Belle2;
25 using namespace Belle2::bklm;
26 using namespace CLHEP;
27 
28 REG_MODULE(BKLMTracking)
29 
31  m_effiYX(nullptr),
32  m_effiYZ(nullptr),
33  m_passYX(nullptr),
34  m_totalYX(nullptr),
35  m_passYZ(nullptr),
36  m_totalYZ(nullptr),
37  m_runTotalEvents(0),
38  m_runTotalEventsWithTracks(0)
39 {
40  for (int i = 0; i < 8; ++i) {
41  m_total[0][i] = nullptr;
42  m_total[1][i] = nullptr;
43  m_pass[0][i] = nullptr;
44  m_pass[1][i] = nullptr;
45  m_effiVsLayer[0][i] = nullptr;
46  m_effiVsLayer[1][i] = nullptr;
47  }
48  setDescription("Perform standard-alone straight line tracking for BKLM");
49  addParam("MatchToRecoTrack", m_MatchToRecoTrack, "[bool], whether match BKLMTrack to RecoTrack; (default is false)", false);
50  addParam("MaxAngleRequired", m_maxAngleRequired,
51  "[degree], match BKLMTrack to RecoTrack; angle between them is required to be smaller than (default 10)", double(10.0));
52  addParam("fitGlobalBKLMTrack", m_globalFit,
53  "[bool], do the BKLMTrack fitting in global system (multi-sectors track) or local system (sector by sector) (default is false, local sys.)",
54  false);
55  addParam("StudyEffiMode", m_studyEffi, "[bool], run in efficieny study mode (default is false)", false);
56  addParam("outputName", m_outPath , "[string], output file name containing efficiencies plots ", string("bklmEffi.root"));
57 }
58 
59 BKLMTrackingModule::~BKLMTrackingModule()
60 {
61 
62 }
63 
64 void BKLMTrackingModule::initialize()
65 {
66 
67  hits2D.isRequired();
68  m_storeTracks.registerInDataStore();
69  m_storeTracks.registerRelationTo(hits2D);
70  m_storeTracks.registerRelationTo(recoTracks);
71  recoHitInformation.registerRelationTo(hits2D);
72  hits2D.registerRelationTo(recoTracks);
73 
74  if (m_studyEffi)
75  B2INFO("BKLMTrackingModule:: this module is running in efficiency study mode!");
76 
77  m_file = new TFile(m_outPath.c_str(), "recreate");
78  TString hname;
79  std::string labelFB[2] = {"BB", "BF"};
80  int Nbin = 16;
81  float gmin = -350;
82  float gmax = 350;
83  int gNbin = 150;
84 
85  m_totalYX = new TH2F("totalYX", " denominator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
86  m_passYX = new TH2F("passYX", " numerator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
87  m_totalYZ = new TH2F("totalYZ", " denominator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
88  m_passYZ = new TH2F("passYZ", " numerator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
89  m_effiYX = new TH2F("effiYX", " effi. Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
90  m_effiYZ = new TH2F("effiYZ", " effi. Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
91  m_effiYX->GetXaxis()->SetTitle("x (cm)");
92  m_effiYX->GetYaxis()->SetTitle("y (cm)");
93  m_effiYZ->GetXaxis()->SetTitle("z (cm)");
94  m_effiYZ->GetYaxis()->SetTitle("y (cm)");
95  for (int iF = 0; iF < 2; iF++) {
96  for (int iS = 0; iS < 8; iS++) {
97  hname.Form("effi_%s%i", labelFB[iF].c_str(), iS);
98  m_effiVsLayer[iF][iS] = new TEfficiency(hname, hname, Nbin, 0, 16);
99  hname.Form("total_%s%i", labelFB[iF].c_str(), iS);
100  m_total[iF][iS] = new TH1F(hname, hname, Nbin, 0, 16);
101  hname.Form("pass_%s%i", labelFB[iF].c_str(), iS);
102  m_pass[iF][iS] = new TH1F(hname, hname, Nbin, 0, 16);
103  }
104  }
105 
106 }
107 
108 void BKLMTrackingModule::beginRun()
109 {
110  StoreObjPtr<EventMetaData> eventMetaData("EventMetaData", DataStore::c_Event);
111  m_runNumber.push_back((int)eventMetaData->getRun());
112  m_runTotalEvents = 0;
113  m_runTotalEventsWithTracks = 0;
114 }
115 
116 void BKLMTrackingModule::event()
117 {
118  m_storeTracks.clear();
119  bool thereIsATrack = false;
120 
121  if (!m_studyEffi) {
122  runTracking(0, -1, -1, -1);
123  if (m_storeTracks.getEntries() > 0)
124  thereIsATrack = true;
125  } else {
126  for (int iSection = 0; iSection < 2; iSection++) {
127  for (int iSector = 0; iSector < 8; iSector++) {
128  for (int iLayer = 0; iLayer < 15; iLayer++) {
129  runTracking(1, iSection, iSector , iLayer);
130  if (m_storeTracks.getEntries() > 0)
131  thereIsATrack = true;
132  generateEffi(iSection, iSector, iLayer);
133  //clear tracks so prepare for the next layer efficieny study
134  m_storeTracks.clear();
135  }
136  }
137  }
138  }
139 
140  m_runTotalEvents++;
141  if (thereIsATrack)
142  m_runTotalEventsWithTracks++;
143 }
144 
145 void BKLMTrackingModule::runTracking(int mode, int iSection, int iSector, int iLayer)
146 {
147  m_storeTracks.clear();
148  //std::list<BKLMTrack*> tracks;
149  //tracks.clear();
150 
151  BKLMTrackFitter* m_fitter = new BKLMTrackFitter();
152  BKLMTrackFinder* m_finder = new BKLMTrackFinder();
153  m_finder->setGlobalFit(m_globalFit);
154  if (mode == 1) m_finder->setGlobalFit(false);
155  m_finder->registerFitter(m_fitter);
156 
157  if (hits2D.getEntries() < 1) return;
158  if (mode == 1) { //efficieny study
159  for (int j = 0; j < hits2D.getEntries(); j++) {
160  hits2D[j]->isOnStaTrack(false);
161  }
162  }
163 
164  for (int hi = 0; hi < hits2D.getEntries() - 1; ++hi) {
165 
166  if (mode == 1 && isLayerUnderStudy(iSection, iSector, iLayer, hits2D[hi])) continue;
167  if (mode == 1 && !isSectorUnderStudy(iSection, iSector, hits2D[hi])) continue;
168  if (hits2D[hi]->isOnStaTrack()) continue;
169  if (hits2D[hi]->isOutOfTime()) continue;
170  for (int hj = hi + 1; hj < hits2D.getEntries(); ++hj) {
171 
172  if (hits2D[hj]->isOnStaTrack()) { continue; }
173  if (hits2D[hj]->isOutOfTime()) continue;
174  if (!m_globalFit && !sameSector(hits2D[hi], hits2D[hj])) { continue; }
175  if (sameSector(hits2D[hi], hits2D[hj]) && abs(hits2D[hi]->getLayer() - hits2D[hj]->getLayer()) < 3) { continue;}
176 
177  std::list<BKLMHit2d*> sectorHitList;
178  //sectorHitList.push_back(hits2D[hi]);
179  //sectorHitList.push_back(hits2D[hj]);
180 
181  std::list<BKLMHit2d*> seed;
182  seed.push_back(hits2D[hi]);
183  seed.push_back(hits2D[hj]);
184 
185  for (int ho = 0; ho < hits2D.getEntries(); ++ho) {
186 
187  if (ho == hi || ho == hj) continue; //exclude seed hits
188  if (mode == 1 && isLayerUnderStudy(iSection, iSector, iLayer, hits2D[hj])) continue;
189  if (mode == 1 && !isSectorUnderStudy(iSection, iSector, hits2D[hj])) continue;
190  if (hits2D[ho]->isOnStaTrack()) continue;
191  if (!m_globalFit && !sameSector(hits2D[ho], hits2D[hi])) continue;
192  // if (hits2D[ho]->getLayer() == hits2D[hi]->getLayer() || hits2D[ho]->getLayer() == hits2D[hj]->getLayer()) continue;
193  if (hits2D[ho]->isOutOfTime()) continue;
194  sectorHitList.push_back(hits2D[ho]);
195  }
196 
197  /* Require at least four hits (minimum for good track, already two as seed, so here we require 2) but
198  * no more than 60 (most likely noise, 60 would be four good tracks).
199  */
200  if (sectorHitList.size() < 2 || sectorHitList.size() > 60) continue;
201 
202  std::list<BKLMHit2d*> m_hits;
203  if (m_finder->filter(seed, sectorHitList, m_hits)) {
204  BKLMTrack* m_track = m_storeTracks.appendNew();
205  m_track->setTrackParam(m_fitter->getTrackParam());
206  m_track->setTrackParamErr(m_fitter->getTrackParamErr());
207  m_track->setLocalTrackParam(m_fitter->getTrackParamSector());
208  m_track->setLocalTrackParamErr(m_fitter->getTrackParamSectorErr());
209  m_track->setTrackChi2(m_fitter->getChi2());
210  m_track->setNumHitOnTrack(m_fitter->getNumHit());
211  m_track->setIsValid(m_fitter->isValid());
212  m_track->setIsGood(m_fitter->isGood());
213  std::list<BKLMHit2d*>::iterator j;
214  m_hits.sort(sortByLayer);
215  for (j = m_hits.begin(); j != m_hits.end(); ++j) {
216  (*j)->isOnStaTrack(true);
217  m_track->addRelationTo((*j));
218  }
219  //tracks.push_back(m_track);
220  //m_track->getTrackParam().Print();
221  //m_track->getTrackParamErr().Print();
222  //match BKLMTrack to RecoTrack
223  if (mode == 0) {
224  RecoTrack* closestTrack = nullptr;
225  if (m_MatchToRecoTrack) {
226  if (findClosestRecoTrack(m_track, closestTrack)) {
227  m_track->addRelationTo(closestTrack);
228  for (j = m_hits.begin(); j != m_hits.end(); ++j) {
229  unsigned int sortingParameter = closestTrack->getNumberOfTotalHits();
230  closestTrack->addBKLMHit((*j), sortingParameter, RecoHitInformation::OriginTrackFinder::c_LocalTrackFinder);
231  }
232  }
233  }//end match
234  }
235  }
236  }
237  }
238 
239  delete m_fitter;
240  delete m_finder;
241 
242 }
243 
244 void BKLMTrackingModule::endRun()
245 {
246  m_totalEvents.push_back(m_runTotalEvents);
247  m_totalEventsWithTracks.push_back(m_runTotalEventsWithTracks);
248 }
249 
250 void BKLMTrackingModule::terminate()
251 {
252  for (long unsigned int i = 0; i < m_runNumber.size(); i++) {
253  float ratio = (float)m_totalEventsWithTracks.at(i) / (float)m_totalEvents.at(i);
254  B2INFO("BKLMTrackingModule:: run " << m_runNumber.at(i) << " --> " << ratio * 100 << "% of events has 1+ BKLMTracks");
255  }
256 
257  m_file->cd();
258  for (int iF = 0; iF < 2; iF++) {
259  for (int iS = 0; iS < 8; iS++) {
260  m_effiVsLayer[iF][iS]->Write();
261  m_total[iF][iS]->Write();
262  m_pass[iF][iS]->Write();
263  }
264  }
265 
266  for (int i = 0; i < m_totalYX->GetNbinsX(); i++) {
267  for (int j = 0; j < m_totalYX->GetNbinsY(); j++) {
268  float num = m_passYX->GetBinContent(i + 1, j + 1);
269  float denom = m_totalYX->GetBinContent(i + 1, j + 1);
270  if (num > 0) {
271  m_effiYX->SetBinContent(i + 1, j + 1, num / denom);
272  m_effiYX->SetBinError(i + 1, j + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
273  } else {
274  m_effiYX->SetBinContent(i + 1, j + 1, 0);
275  m_effiYX->SetBinError(i + 1, j + 1, 0);
276  }
277 
278  num = m_passYZ->GetBinContent(i + 1, j + 1);
279  denom = m_totalYZ->GetBinContent(i + 1, j + 1);
280  if (num > 0) {
281  m_effiYZ->SetBinContent(i + 1, j + 1, num / denom);
282  m_effiYZ->SetBinError(i + 1, j + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
283  } else {
284  m_effiYZ->SetBinContent(i + 1, j + 1, 0);
285  m_effiYZ->SetBinError(i + 1, j + 1, 0);
286  }
287  }
288  }
289 
290  m_totalYX->Write();
291  m_passYX->Write();
292  m_totalYZ->Write();
293  m_passYZ->Write();
294  m_effiYX->Write();
295  m_effiYZ->Write();
296  m_file->Close();
297 
298 }
299 
300 bool BKLMTrackingModule::sameSector(BKLMHit2d* hit1, BKLMHit2d* hit2)
301 {
302  if (hit1->getSection() == hit2->getSection() && hit1->getSector() == hit2->getSector()) return true;
303  else return false;
304 }
305 
306 
307 bool BKLMTrackingModule::findClosestRecoTrack(BKLMTrack* bklmTrk, RecoTrack*& closestTrack)
308 {
309 
310  //StoreArray<RecoTrack> recoTracks;
311  RelationVector<BKLMHit2d> bklmHits = bklmTrk->getRelationsTo<BKLMHit2d> ();
312 
313  if (bklmHits.size() < 1) { B2INFO("BKLMTrackingModule::something is wrong! there is BKLMTrack but no bklmHits"); return false;}
314  if (recoTracks.getEntries() < 1) { B2INFO("BKLMTrackingModule::there is no recoTrack"); return false;}
315  double oldDistance = INFINITY;
316  double oldAngle = INFINITY;
317  closestTrack = nullptr;
318  //TVector3 poca = TVector3(0, 0, 0);
319  TVector3 firstBKLMHitPosition(0, 0, 0);
320  //bklmHits are already sorted by layer
321  //possible two hits in one layer?
322  firstBKLMHitPosition = bklmHits[0]->getGlobalPosition();
323 
324  TMatrixDSym cov(6);
325  TVector3 pos(0, 0, 0);
326  TVector3 mom(0, 0, 0);
327 
328  for (RecoTrack& track : recoTracks) {
329  try {
330  genfit::MeasuredStateOnPlane state = track.getMeasuredStateOnPlaneFromLastHit();
332  state.getPosMomCov(pos, mom, cov);
333  if (mom.Y() * pos.Y() < 0)
334  { state = track.getMeasuredStateOnPlaneFromFirstHit(); }
335  //pos.Print(); mom.Print();
336  const TVector3& distanceVec = firstBKLMHitPosition - pos;
337  state.extrapolateToPoint(firstBKLMHitPosition);
338  double newDistance = distanceVec.Mag2();
339  // two points on the track, (x1,TrkParam[0]+TrkParam[1]*x1, TrkParam[2]+TrkParam[3]*x1),
340  // and (x2,TrkParam[0]+TrkParam[1]*x2, TrkParam[2]+TrkParam[3]*x2),
341  // then we got the vector (x2-x1,....), that is same with (1,TrkParam[1], TrkParam[3]).
342  TVector3 trkVec(1, bklmTrk->getTrackParam()[1], bklmTrk->getTrackParam()[3]);
343  double angle = trkVec.Angle(mom);
344  // choose closest distance or minimum open angle ?
345  // overwrite old distance
346  if (newDistance < oldDistance) {
347  oldDistance = newDistance;
348  closestTrack = &track;
349  //poca = pos;
350  oldAngle = angle;
351  }
352 
353  /* if(angle<oldAngle)
354  {
355  oldAngle=angle;
356  closestTrack = &track;
357  }
358  */
359  } catch (genfit::Exception& e) {
360  }// try
361  }
362 
363  // can not find matched RecoTrack
364  // problem here is the errors of the track parameters are not considered!
365  // best way is the positon or vector direction are required within 5/10 sigma ?
366  if (oldAngle > m_maxAngleRequired) return false;
367  // found matched RecoTrack
368  else return true;
369 }
370 
371 void BKLMTrackingModule::generateEffi(int iSection, int iSector, int iLayer)
372 {
373 
374  set<int> m_pointUsed;
375  m_pointUsed.clear();
376  if (m_storeTracks.getEntries() < 1) return;
377 
378  for (int it = 0; it < m_storeTracks.getEntries(); it++) {
379  //if(m_storeTracks[it]->getTrackChi2()>10) continue;
380  //if(m_storeTracks[it]->getNumHitOnTrack()<6) continue;
381  int cnt1 = 0;
382  int cnt2 = 0;
383 
384  RelationVector<BKLMHit2d> relatedHit2D = m_storeTracks[it]->getRelationsTo<BKLMHit2d>();
385  for (const BKLMHit2d& hit2D : relatedHit2D) {
386  if (hit2D.getLayer() > iLayer + 1) cnt1 ++;
387  if (hit2D.getLayer() < iLayer + 1) cnt2 ++;
388  }
389 
390  if (iLayer != 0 && cnt2 < 1) return;
391  if (iLayer != 14 && cnt1 < 1) return;
392  m_GeoPar = GeometryPar::instance();
393  const bklm::Module* module = m_GeoPar->findModule(iSection, iSector + 1, iLayer + 1);
394  int minPhiStrip = module->getPhiStripMin();
395  int maxPhiStrip = module->getPhiStripMax();
396  int minZStrip = module->getZStripMin();
397  int maxZStrip = module->getZStripMax();
398 
399  CLHEP::Hep3Vector local = module->getLocalPosition(minPhiStrip, minZStrip);
400  CLHEP::Hep3Vector local2 = module->getLocalPosition(maxPhiStrip, maxZStrip);
401  float minLocalY, maxLocalY;
402  float minLocalZ, maxLocalZ;
403  if (local[1] > local2[1]) { maxLocalY = local[1]; minLocalY = local2[1]; } else { maxLocalY = local2[1]; minLocalY = local[1];}
404  if (local[2] > local2[2]) { maxLocalZ = local[2]; minLocalZ = local2[2]; } else { maxLocalZ = local2[2]; minLocalZ = local[2];}
405 
406  TVectorD trkPar = m_storeTracks[it]->getLocalTrackParam();
407 
408  //first layer is the reference layer
409  //if (iSection == 1 && (iSector + 1 ) ==5) cout<<" local X "<<m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1, iLayer + 1) - m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1, 1)<<endl;
410  float reflocalX = fabs(m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1,
411  iLayer + 1) - m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1, 1));
412  //if (iSection == 1 && (iSector + 1 ) ==5) cout<<" local X "<<m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1, iLayer + 1) - m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1, 1)<<endl;
413 
414  float reflocalY = trkPar[0] + trkPar[1] * reflocalX;
415  float reflocalZ = trkPar[2] + trkPar[3] * reflocalX;
416 
417  //reference module is the first layer
418  //module = m_GeoPar->findModule(iSection, iSector + 1, 1);
419  reflocalX = 0.0;
420  Hep3Vector reflocal(reflocalX, reflocalY, reflocalZ);
421  //Hep3Vector global(localX, localY, localZ);
422  Hep3Vector global(0, 0, 0);
423  module = m_GeoPar->findModule(iSection, iSector + 1, iLayer + 1);
424  global = module->localToGlobal(reflocal);
425  //float localX = module->globalToLocal(global)[0];
426  float localY = module->globalToLocal(global)[1];
427  float localZ = module->globalToLocal(global)[2];
428 
429 
430  //geometry cut
431  if (localY > minLocalY && localY < maxLocalY && localZ > minLocalZ && localZ < maxLocalZ) {
432 
433  bool m_iffound = false;
434  m_total[iSection][iSector]->Fill(iLayer + 1);
435  m_totalYX->Fill(global[0], global[1]);
436  m_totalYZ->Fill(global[2], global[1]);
437 
438  for (int he = 0; he < hits2D.getEntries(); ++he) {
439  if (!isLayerUnderStudy(iSection, iSector, iLayer, hits2D[he])) continue;
440  if (hits2D[he]->isOutOfTime()) continue;
441  //if alreday used, skip
442  if (m_pointUsed.find(he) != m_pointUsed.end()) continue;
443 
444  double error, sigma;
445  float distance = distanceToHit(m_storeTracks[it], hits2D[he], error, sigma);
446 
447  if (distance < 10 && sigma < 5) { m_iffound = true; }
448  if (m_iffound) {
449  m_pointUsed.insert(he);
450  //global[0] = hits2D[he]->getGlobalPosition()[0];
451  //global[1] = hits2D[he]->getGlobalPosition()[1];
452  //global[2] = hits2D[he]->getGlobalPosition()[2];
453  m_pass[iSection][iSector]->Fill(iLayer + 1);
454  m_passYX->Fill(global[0], global[1]);
455  m_passYZ->Fill(global[2], global[1]);
456  break;
457  }
458  }
459 
460  m_effiVsLayer[iSection][iSector]->Fill(m_iffound, iLayer + 1);
461  //cout<<" global "<<global[0]<<", "<< global[1]<<" "<<global[2]<<endl;
462  //m_effiYX->Fill(m_iffound, global[1], global[0]);
463  //m_effiYZ->Fill(m_iffound, global[1], global[2]);
464  //m_effiYX->SetPassedHistogram(*m_passYX);
465  //m_effiYX->SetTotalHistogram(*m_totalYX);
466  //m_effiYZ->SetPassedHistogram(*m_passYZ);
467  //m_effiYZ->SetTotalHistogram(*m_totalYZ);
468  }
469  }//end of loop tracks
470 
471 }
472 
473 bool BKLMTrackingModule::sortByLayer(BKLMHit2d* hit1, BKLMHit2d* hit2)
474 {
475 
476  return hit1->getLayer() < hit2->getLayer();
477 
478 }
479 
480 bool BKLMTrackingModule::isLayerUnderStudy(int section, int iSector, int iLayer, BKLMHit2d* hit)
481 {
482  if (hit->getSection() == section && hit->getSector() == iSector + 1 && hit->getLayer() == iLayer + 1) return true;
483  else return false;
484 }
485 
486 bool BKLMTrackingModule::isSectorUnderStudy(int section, int iSector, BKLMHit2d* hit)
487 {
488  if (hit->getSection() == section && hit->getSector() == iSector + 1) return true;
489  else return false;
490 }
491 
492 double BKLMTrackingModule::distanceToHit(BKLMTrack* track, BKLMHit2d* hit,
493  double& error,
494  double& sigma)
495 {
496 
497  double x, y, z, dy, dz;
498 
499  error = DBL_MAX;
500  sigma = DBL_MAX;
501 
502  TVectorD m_SectorPar = track->getLocalTrackParam();
503 
504  m_GeoPar = GeometryPar::instance();
505  const Belle2::bklm::Module* refMod = m_GeoPar->findModule(hit->getSection(), hit->getSector(), 1);
506  const Belle2::bklm::Module* corMod = m_GeoPar->findModule(hit->getSection(), hit->getSector(), hit->getLayer());
507 
508  CLHEP::Hep3Vector globalPos(hit->getGlobalPosition()[0], hit->getGlobalPosition()[1], hit->getGlobalPosition()[2]);
509  CLHEP::Hep3Vector local = refMod->globalToLocal(globalPos);
510 
511  x = local[0] ;
512 
513  y = m_SectorPar[ 0 ] + x * m_SectorPar[ 1 ];
514  z = m_SectorPar[ 2 ] + x * m_SectorPar[ 3 ];
515 
516  dy = y - local[1];
517  dz = z - local[2];
518 
519  double distance = sqrt(dy * dy + dz * dz);
520 
521  double hit_localPhiErr = corMod->getPhiStripWidth() / sqrt(12);
522  double hit_localZErr = corMod->getZStripWidth() / sqrt(12);
523 
524  //error from tracking is ignored here
525  error = sqrt(pow(hit_localPhiErr, 2) +
526  pow(hit_localZErr, 2));
527 
528  if (error != 0.0) {
529  sigma = distance / error;
530  } else {
531  sigma = DBL_MAX;
532  }
533 
534  return (distance);
535 
536 }
Belle2::RelationVector::size
size_t size() const
Get number of relations.
Definition: RelationVector.h:98
genfit::Exception
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Definition: Exception.h:48
Belle2::BKLMTrackingModule
This module perform straight line track finding and fitting for BKLM.
Definition: BKLMTrackingModule.h:47
Belle2::BKLMTrackFitter
track fitting procedure
Definition: BKLMTrackFitter.h:39
Belle2::BKLMTrack::setIsValid
void setIsValid(const bool valid)
set the fit valid status
Definition: BKLMTrack.h:122
genfit::MeasuredStateOnPlane
#StateOnPlane with additional covariance matrix.
Definition: MeasuredStateOnPlane.h:39
Belle2::BKLMTrack::setIsGood
void setIsGood(const bool good)
set the fit good status
Definition: BKLMTrack.h:128
Belle2::BKLMTrackFinder::filter
bool filter(const std::list< BKLMHit2d * > &seed, std::list< BKLMHit2d * > &hits, std::list< BKLMHit2d * > &track)
find associated hits and do fit.
Definition: BKLMTrackFinder.cc:48
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::BKLMTrackFitter::getTrackParamSectorErr
CLHEP::HepSymMatrix getTrackParamSectorErr()
Get invariance matrix of track parameters in the sector local system, where the first layer of the se...
Definition: BKLMTrackFitter.h:93
Belle2::BKLMHit2d::getSection
int getSection() const
Get section number.
Definition: BKLMHit2d.h:74
Belle2::BKLMTrack::setTrackParam
void setTrackParam(const CLHEP::HepVector &trkPar)
Set track parameters in the global system. y = p0 + p1 * x; z = p2 + p3 * x.
Definition: BKLMTrack.cc:122
Belle2::RelationsInterface::addRelationTo
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
Definition: RelationsObject.h:144
Belle2::bklm::Module::getPhiStripWidth
double getPhiStripWidth() const
Get phi-strip width.
Definition: Module.h:147
Belle2::BKLMTrackFinder::registerFitter
void registerFitter(BKLMTrackFitter *fitter)
Register a fitter if not constructed with one.
Definition: BKLMTrackFinder.cc:41
Belle2::BKLMTrackFitter::isValid
bool isValid()
Is fit valid.
Definition: BKLMTrackFitter.h:99
Belle2::bklm::Module
Define the geometry of a BKLM module Each sector [octant] contains Modules.
Definition: Module.h:86
Belle2::RelationsInterface::getRelationsTo
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
Definition: RelationsObject.h:199
Belle2::BKLMTrackFitter::getTrackParamSector
CLHEP::HepVector getTrackParamSector()
Get track parameters in the sector locan system, y = p0 + p1 * x, z = p2 + p3 *x, where the first lay...
Definition: BKLMTrackFitter.h:87
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::BKLMTrack::setLocalTrackParamErr
void setLocalTrackParamErr(const CLHEP::HepSymMatrix &trkParErr)
Set invariance matrix of track parameters in the sector local system, where the first layer of the se...
Definition: BKLMTrack.cc:152
Belle2::BKLMTrackFitter::getChi2
float getChi2()
Chi square of the fit.
Definition: BKLMTrackFitter.h:111
Belle2::RecoTrack
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:78
Belle2::BKLMTrack::getTrackParam
TVectorD getTrackParam()
Get track parameters in the global system. y = p0 + p1 * x; z = p2 + p3 * x.
Definition: BKLMTrack.cc:73
Belle2::RelationVector
Class for type safe access to objects that are referred to in relations.
Definition: DataStore.h:38
Belle2::BKLMHit2d::getSector
int getSector() const
Get sector number.
Definition: BKLMHit2d.h:81
Belle2::bklm::Module::globalToLocal
const CLHEP::Hep3Vector globalToLocal(const CLHEP::Hep3Vector &v, bool reco=false) const
Transform space-point within this module from global to local coordinates.
Definition: Module.cc:320
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
Belle2::RecoTrack::addBKLMHit
bool addBKLMHit(const UsedBKLMHit *bklmHit, const unsigned int sortingParameter, OriginTrackFinder foundByTrackFinder=OriginTrackFinder::c_undefinedTrackFinder)
Adds a bklm hit with the given information to the reco track.
Definition: RecoTrack.h:283
Belle2::BKLMTrack::setNumHitOnTrack
void setNumHitOnTrack(const int NumHit)
Set the number of 2d hits on the track.
Definition: BKLMTrack.h:116
Belle2::BKLMTrackFitter::getTrackParam
CLHEP::HepVector getTrackParam()
Get track parameters in the global system. y = p0 + p1 * x; y = p2 + p3 * z, if in local sector fit m...
Definition: BKLMTrackFitter.h:75
Belle2::BKLMTrackFitter::getTrackParamErr
CLHEP::HepSymMatrix getTrackParamErr()
Get invariance matrix of track parameters in the global system.
Definition: BKLMTrackFitter.h:81
Belle2::BKLMTrack
Store one BKLM Track as a ROOT object.
Definition: BKLMTrack.h:37
Belle2::BKLMTrackFitter::getNumHit
int getNumHit()
number of the hits on this track
Definition: BKLMTrackFitter.h:117
Belle2::RecoTrack::getNumberOfTotalHits
unsigned int getNumberOfTotalHits() const
Return the number of cdc + svd + pxd + bklm + eklm hits.
Definition: RecoTrack.h:432
Belle2::BKLMTrack::setTrackChi2
void setTrackChi2(const float chi2)
Set the fitted chi2 of the track.
Definition: BKLMTrack.h:110
Belle2::bklm::Module::getZStripWidth
double getZStripWidth() const
Get z-strip width.
Definition: Module.h:165
Belle2::BKLMTrackFitter::isGood
bool isGood()
Is fit good.
Definition: BKLMTrackFitter.h:105
Belle2::BKLMTrackFinder
track finding procedure
Definition: BKLMTrackFinder.h:35
Belle2::BKLMTrack::setTrackParamErr
void setTrackParamErr(const CLHEP::HepSymMatrix &trkParErr)
Set invariance matrix of track parameters in the global system.
Definition: BKLMTrack.cc:132
Belle2::BKLMHit2d
Store one BKLM strip hit as a ROOT object.
Definition: BKLMHit2d.h:42
Belle2::BKLMHit2d::getLayer
int getLayer() const
Get layer number.
Definition: BKLMHit2d.h:88
Belle2::BKLMTrackFinder::setGlobalFit
void setGlobalFit(bool localOrGlobal)
set the fitting mode, local system or global system
Definition: BKLMTrackFinder.h:57
Belle2::BKLMTrack::setLocalTrackParam
void setLocalTrackParam(const CLHEP::HepVector &trkPar)
Set track parameters in the sector local system, where the first layer of the sector is used as refer...
Definition: BKLMTrack.cc:142