Belle II Software development
BKLMTrackingModule.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9/* Own header. */
10#include <klm/bklm/modules/bklmTracking/BKLMTrackingModule.h>
11
12/* KLM headers. */
13#include <klm/bklm/geometry/GeometryPar.h>
14#include <klm/bklm/modules/bklmTracking/BKLMTrackFinder.h>
15
16/* Basf2 headers. */
17#include <framework/dataobjects/EventMetaData.h>
18#include <framework/datastore/StoreObjPtr.h>
19#include <framework/logging/Logger.h>
20
21using namespace Belle2;
22using namespace Belle2::bklm;
23using namespace CLHEP;
24
25REG_MODULE(BKLMTracking);
26
28 m_effiYX(nullptr),
29 m_effiYZ(nullptr),
30 m_passYX(nullptr),
31 m_totalYX(nullptr),
32 m_passYZ(nullptr),
33 m_totalYZ(nullptr),
34 m_runTotalEvents(0),
35 m_runTotalEventsWithTracks(0)
36{
37 for (int i = 0; i < 8; ++i) {
38 m_total[0][i] = nullptr;
39 m_total[1][i] = nullptr;
40 m_pass[0][i] = nullptr;
41 m_pass[1][i] = nullptr;
42 m_effiVsLayer[0][i] = nullptr;
43 m_effiVsLayer[1][i] = nullptr;
44 }
45 setDescription("Perform standard-alone straight line tracking for BKLM");
46 addParam("MatchToRecoTrack", m_MatchToRecoTrack, "[bool], whether match BKLMTrack to RecoTrack; (default is false)", false);
47 addParam("MaxAngleRequired", m_maxAngleRequired,
48 "[degree], match BKLMTrack to RecoTrack; angle between them is required to be smaller than (default 10)", double(10.0));
49 addParam("MaxDistance", m_maxDistance,
50 "[cm], During efficiency calculation, distance between track and 2dhit must be smaller than (default 10)", double(10.0));
51 addParam("MaxSigma", m_maxSigma,
52 "[sigma], During efficiency calculation, uncertainty of 2dhit must be smaller than (default 5); ", double(5));
53 addParam("MinHitList", m_minHitList,
54 ", During track finding, a good track after initial seed hits must be larger than is (default 2); ", unsigned(2));
55 addParam("MaxHitList", m_maxHitList,
56 ", During track finding, a good track after initial seed hits must be smaller than is (default 60); ", unsigned(60));
57 addParam("fitGlobalBKLMTrack", m_globalFit,
58 "[bool], do the BKLMTrack fitting in global system (multi-sectors track) or local system (sector by sector) (default is false, local sys.)",
59 false);
60 addParam("StudyEffiMode", m_studyEffi, "[bool], run in efficiency study mode (default is false)", false);
61 addParam("outputName", m_outPath, "[string], output file name containing efficiencies plots ", std::string("bklmEffi.root"));
62}
63
65{
66
67}
68
70{
71
72 hits2D.isRequired();
73 m_storeTracks.registerInDataStore();
74 m_storeTracks.registerRelationTo(hits2D);
75 m_storeTracks.registerRelationTo(recoTracks);
76 recoHitInformation.registerRelationTo(hits2D);
77 hits2D.registerRelationTo(recoTracks);
78
79 if (m_studyEffi)
80 B2INFO("BKLMTrackingModule:: this module is running in efficiency study mode!");
81
82 m_file = new TFile(m_outPath.c_str(), "recreate");
83 TString hname;
84 std::string labelFB[2] = {"BB", "BF"};
85 int Nbin = 16;
86 float gmin = -350;
87 float gmax = 350;
88 int gNbin = 150;
89
90 m_totalYX = new TH2F("totalYX", " denominator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
91 m_passYX = new TH2F("passYX", " numerator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
92 m_totalYZ = new TH2F("totalYZ", " denominator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
93 m_passYZ = new TH2F("passYZ", " numerator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
94 m_effiYX = new TH2F("effiYX", " effi. Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
95 m_effiYZ = new TH2F("effiYZ", " effi. Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
96 m_effiYX->GetXaxis()->SetTitle("x (cm)");
97 m_effiYX->GetYaxis()->SetTitle("y (cm)");
98 m_effiYZ->GetXaxis()->SetTitle("z (cm)");
99 m_effiYZ->GetYaxis()->SetTitle("y (cm)");
100 for (int iF = 0; iF < 2; iF++) {
101 for (int iS = 0; iS < 8; iS++) {
102 hname.Form("effi_%s%i", labelFB[iF].c_str(), iS);
103 m_effiVsLayer[iF][iS] = new TEfficiency(hname, hname, Nbin, 0, 16);
104 hname.Form("total_%s%i", labelFB[iF].c_str(), iS);
105 m_total[iF][iS] = new TH1F(hname, hname, Nbin, 0, 16);
106 hname.Form("pass_%s%i", labelFB[iF].c_str(), iS);
107 m_pass[iF][iS] = new TH1F(hname, hname, Nbin, 0, 16);
108 }
109 }
110
111}
112
114{
115 StoreObjPtr<EventMetaData> eventMetaData("EventMetaData", DataStore::c_Event);
116 m_runNumber.push_back((int)eventMetaData->getRun());
119}
120
122{
123 m_storeTracks.clear();
124 bool thereIsATrack = false;
125
126 if (!m_studyEffi) {
127 runTracking(0, -1, -1, -1);
128 if (m_storeTracks.getEntries() > 0)
129 thereIsATrack = true;
130 } else {
131 for (int iSection = 0; iSection < 2; iSection++) {
132 for (int iSector = 0; iSector < 8; iSector++) {
133 for (int iLayer = 0; iLayer < 15; iLayer++) {
134 runTracking(1, iSection, iSector, iLayer);
135 if (m_storeTracks.getEntries() > 0)
136 thereIsATrack = true;
137 generateEffi(iSection, iSector, iLayer);
138 //clear tracks so prepare for the next layer efficiency study
139 m_storeTracks.clear();
140 }
141 }
142 }
143 }
144
146 if (thereIsATrack)
148}
149
150void BKLMTrackingModule::runTracking(int mode, int iSection, int iSector, int iLayer)
151{
152 m_storeTracks.clear();
153 //std::list<BKLMTrack*> tracks;
154 //tracks.clear();
155
156 BKLMTrackFitter* m_fitter = new BKLMTrackFitter();
157 BKLMTrackFinder* m_finder = new BKLMTrackFinder();
158 m_finder->setGlobalFit(m_globalFit);
159 if (mode == 1)
160 m_finder->setGlobalFit(false);
161 m_finder->registerFitter(m_fitter);
162
163 if (hits2D.getEntries() < 1)
164 return;
165 if (mode == 1) { //efficiency study
166 for (int j = 0; j < hits2D.getEntries(); j++) {
167 if (hits2D[j]->getSubdetector() != KLMElementNumbers::c_BKLM)
168 continue;
169 hits2D[j]->isOnStaTrack(false);
170 }
171 }
172
173 for (int hi = 0; hi < hits2D.getEntries() - 1; ++hi) {
174 if (hits2D[hi]->getSubdetector() != KLMElementNumbers::c_BKLM)
175 continue;
176
177 if (mode == 1 && isLayerUnderStudy(iSection, iSector, iLayer, hits2D[hi]))
178 continue;
179 if (mode == 1 && !isSectorUnderStudy(iSection, iSector, hits2D[hi]))
180 continue;
181 if (hits2D[hi]->isOnStaTrack())
182 continue;
183 if (hits2D[hi]->isOutOfTime())
184 continue;
185 for (int hj = hi + 1; hj < hits2D.getEntries(); ++hj) {
186
187 if (hits2D[hj]->isOnStaTrack())
188 continue;
189 if (hits2D[hj]->isOutOfTime())
190 continue;
191 if (!m_globalFit && !sameSector(hits2D[hi], hits2D[hj]))
192 continue;
193 if (sameSector(hits2D[hi], hits2D[hj]) &&
194 std::abs(hits2D[hi]->getLayer() - hits2D[hj]->getLayer()) < 3)
195 continue;
196
197 std::list<KLMHit2d*> sectorHitList;
198 //sectorHitList.push_back(hits2D[hi]);
199 //sectorHitList.push_back(hits2D[hj]);
200
201 std::list<KLMHit2d*> seed;
202 seed.push_back(hits2D[hi]);
203 seed.push_back(hits2D[hj]);
204
205 for (int ho = 0; ho < hits2D.getEntries(); ++ho) {
206
207 // Exclude seed hits.
208 if (ho == hi || ho == hj)
209 continue;
210 if (mode == 1 && isLayerUnderStudy(iSection, iSector, iLayer, hits2D[hj]))
211 continue;
212 if (mode == 1 && !isSectorUnderStudy(iSection, iSector, hits2D[hj]))
213 continue;
214 if (hits2D[ho]->isOnStaTrack())
215 continue;
216 if (!m_globalFit && !sameSector(hits2D[ho], hits2D[hi]))
217 continue;
218 // if (hits2D[ho]->getLayer() == hits2D[hi]->getLayer() || hits2D[ho]->getLayer() == hits2D[hj]->getLayer())
219 // continue;
220 if (hits2D[ho]->isOutOfTime())
221 continue;
222 sectorHitList.push_back(hits2D[ho]);
223 }
224
225 /* Require at least four hits (minimum for good track, already two as seed, so here we require 2) but
226 * no more than 60 (most likely noise, 60 would be four good tracks).
227 */
228 if (sectorHitList.size() < m_minHitList || sectorHitList.size() > m_maxHitList)
229 continue;
230
231 std::list<KLMHit2d*> m_hits;
232 if (m_finder->filter(seed, sectorHitList, m_hits)) {
233 BKLMTrack* m_track = m_storeTracks.appendNew();
234 m_track->setTrackParam(m_fitter->getTrackParam());
235 m_track->setTrackParamErr(m_fitter->getTrackParamErr());
236 m_track->setLocalTrackParam(m_fitter->getTrackParamSector());
237 m_track->setLocalTrackParamErr(m_fitter->getTrackParamSectorErr());
238 m_track->setTrackChi2(m_fitter->getChi2());
239 m_track->setNumHitOnTrack(m_fitter->getNumHit());
240 m_track->setIsValid(m_fitter->isValid());
241 m_track->setIsGood(m_fitter->isGood());
242 m_hits.sort(sortByLayer);
243 for (KLMHit2d* hit2d : m_hits) {
244 hit2d->isOnStaTrack(true);
245 m_track->addRelationTo(hit2d);
246 }
247 //tracks.push_back(m_track);
248 //m_track->getTrackParam().Print();
249 //m_track->getTrackParamErr().Print();
250 //match BKLMTrack to RecoTrack
251 if (mode == 0) {
252 RecoTrack* closestTrack = nullptr;
253 if (m_MatchToRecoTrack) {
254 if (findClosestRecoTrack(m_track, closestTrack)) {
255 m_track->addRelationTo(closestTrack);
256 for (KLMHit2d* hit2d : m_hits) {
257 unsigned int sortingParameter = closestTrack->getNumberOfTotalHits();
258 closestTrack->addBKLMHit(hit2d, sortingParameter, RecoHitInformation::OriginTrackFinder::c_LocalTrackFinder);
259 }
260 }
261 }//end match
262 }
263 }
264 }
265 }
266
267 delete m_fitter;
268 delete m_finder;
269
270}
271
273{
276}
277
279{
280 for (long unsigned int i = 0; i < m_runNumber.size(); i++) {
281 float ratio = (float)m_totalEventsWithTracks.at(i) / (float)m_totalEvents.at(i);
282 B2INFO("BKLMTrackingModule:: run " << m_runNumber.at(i) << " --> " << ratio * 100 << "% of events has 1+ BKLMTracks");
283 }
284
285 m_file->cd();
286 for (int iF = 0; iF < 2; iF++) {
287 for (int iS = 0; iS < 8; iS++) {
288 m_effiVsLayer[iF][iS]->Write();
289 m_total[iF][iS]->Write();
290 m_pass[iF][iS]->Write();
291 }
292 }
293
294 for (int i = 0; i < m_totalYX->GetNbinsX(); i++) {
295 for (int j = 0; j < m_totalYX->GetNbinsY(); j++) {
296 float num = m_passYX->GetBinContent(i + 1, j + 1);
297 float denom = m_totalYX->GetBinContent(i + 1, j + 1);
298 if (num > 0) {
299 m_effiYX->SetBinContent(i + 1, j + 1, num / denom);
300 m_effiYX->SetBinError(i + 1, j + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
301 } else {
302 m_effiYX->SetBinContent(i + 1, j + 1, 0);
303 m_effiYX->SetBinError(i + 1, j + 1, 0);
304 }
305
306 num = m_passYZ->GetBinContent(i + 1, j + 1);
307 denom = m_totalYZ->GetBinContent(i + 1, j + 1);
308 if (num > 0) {
309 m_effiYZ->SetBinContent(i + 1, j + 1, num / denom);
310 m_effiYZ->SetBinError(i + 1, j + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
311 } else {
312 m_effiYZ->SetBinContent(i + 1, j + 1, 0);
313 m_effiYZ->SetBinError(i + 1, j + 1, 0);
314 }
315 }
316 }
317
318 m_totalYX->Write();
319 m_passYX->Write();
320 m_totalYZ->Write();
321 m_passYZ->Write();
322 m_effiYX->Write();
323 m_effiYZ->Write();
324 m_file->Close();
325
326}
327
329{
330 if (hit1->getSection() == hit2->getSection() && hit1->getSector() == hit2->getSector())
331 return true;
332 else return false;
333}
334
335
337{
338
339 //StoreArray<RecoTrack> recoTracks;
340 RelationVector<KLMHit2d> bklmHits = bklmTrk->getRelationsTo<KLMHit2d> ();
341
342 if (bklmHits.size() < 1) {
343 B2INFO("BKLMTrackingModule::something is wrong! there is BKLMTrack but no bklmHits");
344 return false;
345 }
346 if (recoTracks.getEntries() < 1) {
347 B2INFO("BKLMTrackingModule::there is no recoTrack");
348 return false;
349 }
350 double oldDistanceSq = INFINITY;
351 double oldAngle = INFINITY;
352 closestTrack = nullptr;
353 //bklmHits are already sorted by layer
354 //possible two hits in one layer?
355 //genfit requires TVector3 rather than XYZVector
356 TVector3 firstBKLMHitPosition(bklmHits[0]->getPosition().X(),
357 bklmHits[0]->getPosition().Y(),
358 bklmHits[0]->getPosition().Z());
359
360 // To get direction (angle) below, we have two points on the bklmTrk:
361 // (x1, TrackParam[0]+TrackParam[1]*x1, TrackParam[2]+TrackParam[3]*x1)
362 // (x2, TrackParam[0]+TrackParam[1]*x2, TrackParam[2]+TrackParam[3]*x2)
363 // the difference vector is
364 // (x2-x1, TrackParam[1]*(x2-x1), TrackParam[3]*(x2-x1))
365 // which is proportional to
366 // (1, TrackParam[1], TrackParam[3]).
367 TVector3 bklmTrkVec(1.0, bklmTrk->getTrackParam()[1], bklmTrk->getTrackParam()[3]);
368
369 TMatrixDSym cov(6);
370 TVector3 pos; // initializes to (0,0,0)
371 TVector3 mom; // initializes to (0,0,0)
372
373 for (RecoTrack& track : recoTracks) {
374 try {
375 genfit::MeasuredStateOnPlane state = track.getMeasuredStateOnPlaneFromLastHit();
377 state.getPosMomCov(pos, mom, cov);
378 if (mom.Y() * pos.Y() < 0) {
379 state = track.getMeasuredStateOnPlaneFromFirstHit();
380 }
381 const TVector3& distanceVec = firstBKLMHitPosition - pos;
382 state.extrapolateToPoint(firstBKLMHitPosition);
383 double newDistanceSq = distanceVec.Mag2();
384 double angle = bklmTrkVec.Angle(mom);
385 // choose closest distance or minimum open angle ?
386 // overwrite old distance
387 if (newDistanceSq < oldDistanceSq) {
388 oldDistanceSq = newDistanceSq;
389 closestTrack = &track;
390 oldAngle = angle;
391 }
392
393 /* if(angle<oldAngle)
394 {
395 oldAngle=angle;
396 closestTrack = &track;
397 }
398 */
399 } catch (genfit::Exception& e) {
400 }// try
401 }
402
403 // can not find matched RecoTrack
404 // problem here is the errors of the track parameters are not considered!
405 // best way is the position or vector direction are required within 5/10 sigma ?
406 if (oldAngle > m_maxAngleRequired)
407 return false;
408 // found matched RecoTrack
409 else return true;
410}
411
412void BKLMTrackingModule::generateEffi(int iSection, int iSector, int iLayer)
413{
414
415 std::set<int> m_pointUsed;
416 m_pointUsed.clear();
417 if (m_storeTracks.getEntries() < 1)
418 return;
419 B2DEBUG(10, "BKLMTracking:generateEffi: " << iSection << " " << iSector << " " << iLayer);
420
421
422 for (int it = 0; it < m_storeTracks.getEntries(); it++) {
423 //if(m_storeTracks[it]->getTrackChi2()>10) continue;
424 //if(m_storeTracks[it]->getNumHitOnTrack()<6) continue;
425 int cnt1 = 0;
426 int cnt2 = 0;
427
428 RelationVector<KLMHit2d> relatedHit2D = m_storeTracks[it]->getRelationsTo<KLMHit2d>();
429 for (const KLMHit2d& hit2D : relatedHit2D) {
430 if (hit2D.getLayer() > iLayer + 1)
431 cnt1++;
432 if (hit2D.getLayer() < iLayer + 1)
433 cnt2++;
434 if (hit2D.getLayer() == iLayer + 1) {
435 B2DEBUG(10, "generateEffi: Hit info. Secti/sector/Lay = " << hit2D.getSection()
436 << "/" << hit2D.getSector() - 1 << "/" << hit2D.getLayer() - 1);
437 B2DEBUG(11, "generateEffi: Hit info. x/y/z = " << hit2D.getPositionX()
438 << "/" << hit2D.getPositionY() << "/" << hit2D.getPositionZ());
439 }
440 }
441
442 if (iLayer != 0 && cnt2 < 1)
443 return;
444 if (iLayer != 14 && cnt1 < 1)
445 return;
447 const bklm::Module* module = m_GeoPar->findModule(iSection, iSector + 1, iLayer + 1);
448 int minPhiStrip = module->getPhiStripMin();
449 int maxPhiStrip = module->getPhiStripMax();
450 int minZStrip = module->getZStripMin();
451 int maxZStrip = module->getZStripMax();
452
453 CLHEP::Hep3Vector local = module->getLocalPosition(minPhiStrip, minZStrip);
454 CLHEP::Hep3Vector local2 = module->getLocalPosition(maxPhiStrip, maxZStrip);
455 float minLocalY, maxLocalY;
456 float minLocalZ, maxLocalZ;
457 if (local[1] > local2[1]) {
458 maxLocalY = local[1];
459 minLocalY = local2[1];
460 } else {
461 maxLocalY = local2[1];
462 minLocalY = local[1];
463 }
464 if (local[2] > local2[2]) {
465 maxLocalZ = local[2];
466 minLocalZ = local2[2];
467 } else {
468 maxLocalZ = local2[2];
469 minLocalZ = local[2];
470 }
471
472 TVectorD trkPar = m_storeTracks[it]->getLocalTrackParam();
473
474 //first layer is the reference layer
475 //if (iSection == 1 && (iSector + 1 ) == 5)
476 // cout<<" local X "<<m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1, iLayer + 1) - m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1, 1) << endl;
477 float reflocalX = fabs(m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1,
478 iLayer + 1) - m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1, 1));
479 //if (iSection == 1 && (iSector + 1 ) == 5)
480 // cout<<" local X "<<m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1, iLayer + 1) - m_GeoPar->getActiveMiddleRadius(iSection, iSector + 1, 1) << endl;
481
482 float reflocalY = trkPar[0] + trkPar[1] * reflocalX;
483 float reflocalZ = trkPar[2] + trkPar[3] * reflocalX;
484
485 //reference module is the first layer
486 //module = m_GeoPar->findModule(iSection, iSector + 1, 1);
487 reflocalX = 0.0;
488 Hep3Vector reflocal(reflocalX, reflocalY, reflocalZ);
489 //Hep3Vector global(localX, localY, localZ);
490 Hep3Vector global(0, 0, 0);
491 module = m_GeoPar->findModule(iSection, iSector + 1, iLayer + 1);
492 global = module->localToGlobal(reflocal);
493 //float localX = module->globalToLocal(global)[0];
494 float localY = module->globalToLocal(global)[1];
495 float localZ = module->globalToLocal(global)[2];
496
497 B2DEBUG(10, "BKLM:generateEffi: RefLocal " << reflocalX << " " << reflocalY << " " << reflocalZ);
498 B2DEBUG(10, "BKLM:generateEffi: Global " << global[0] << " " << global[1] << " " << global[2]);
499 B2DEBUG(10, "BKLM:generateEffi: Local " << 0 << " " << localY << " " << localZ);
500
501
502
503 //geometry cut
504 if (localY > minLocalY && localY < maxLocalY && localZ > minLocalZ && localZ < maxLocalZ) {
505
506 bool m_iffound = false;
507 m_total[iSection][iSector]->Fill(iLayer + 1);
508 m_totalYX->Fill(global[0], global[1]);
509 m_totalYZ->Fill(global[2], global[1]);
510
511 for (int he = 0; he < hits2D.getEntries(); ++he) {
512 if (!isLayerUnderStudy(iSection, iSector, iLayer, hits2D[he]))
513 continue;
514 if (hits2D[he]->isOutOfTime())
515 continue;
516 //if already used, skip
517 if (m_pointUsed.find(he) != m_pointUsed.end())
518 continue;
519
520 double error, sigma;
521 float distance = distanceToHit(m_storeTracks[it], hits2D[he], error, sigma);
522 B2DEBUG(10, "BKLM Dist = " << distance << ", error = " << error);
523 if (distance < m_maxDistance && sigma < m_maxSigma) {
524 m_iffound = true;
525 B2DEBUG(10, "BKLMTracking:generateEffi: Hit found!");;
526 }
527 if (m_iffound) {
528 m_pointUsed.insert(he);
529 //global[0] = hits2D[he]->getPosition()[0];
530 //global[1] = hits2D[he]->getPosition()[1];
531 //global[2] = hits2D[he]->getPosition()[2];
532 m_pass[iSection][iSector]->Fill(iLayer + 1);
533 m_passYX->Fill(global[0], global[1]);
534 m_passYZ->Fill(global[2], global[1]);
535 break;
536 }
537 }
538
539 m_effiVsLayer[iSection][iSector]->Fill(m_iffound, iLayer + 1);
540 //cout<<" global "<<global[0]<<", "<< global[1]<<" "<<global[2]<<endl;
541 //m_effiYX->Fill(m_iffound, global[1], global[0]);
542 //m_effiYZ->Fill(m_iffound, global[1], global[2]);
543 //m_effiYX->SetPassedHistogram(*m_passYX);
544 //m_effiYX->SetTotalHistogram(*m_totalYX);
545 //m_effiYZ->SetPassedHistogram(*m_passYZ);
546 //m_effiYZ->SetTotalHistogram(*m_totalYZ);
547 }
548 }//end of loop tracks
549
550}
551
553{
554
555 return hit1->getLayer() < hit2->getLayer();
556
557}
558
559bool BKLMTrackingModule::isLayerUnderStudy(int section, int iSector, int iLayer, KLMHit2d* hit)
560{
561 if (hit->getSection() == section && hit->getSector() == iSector + 1 && hit->getLayer() == iLayer + 1)
562 return true;
563 else return false;
564}
565
566bool BKLMTrackingModule::isSectorUnderStudy(int section, int iSector, KLMHit2d* hit)
567{
568 if (hit->getSection() == section && hit->getSector() == iSector + 1)
569 return true;
570 else return false;
571}
572
574 double& error,
575 double& sigma)
576{
577
578 double x, y, z, dy, dz;
579
580 error = DBL_MAX;
581 sigma = DBL_MAX;
582
583 TVectorD m_SectorPar = track->getLocalTrackParam();
584
586 const Belle2::bklm::Module* refMod = m_GeoPar->findModule(hit->getSection(), hit->getSector(), 1);
587 const Belle2::bklm::Module* corMod = m_GeoPar->findModule(hit->getSection(), hit->getSector(), hit->getLayer());
588
589 CLHEP::Hep3Vector globalPos(hit->getPositionX(), hit->getPositionY(),
590 hit->getPositionZ());
591 CLHEP::Hep3Vector local = refMod->globalToLocal(globalPos);
592
593 x = local[0] ;
594
595 y = m_SectorPar[ 0 ] + x * m_SectorPar[ 1 ];
596 z = m_SectorPar[ 2 ] + x * m_SectorPar[ 3 ];
597
598 dy = y - local[1];
599 dz = z - local[2];
600
601 double distance = sqrt(dy * dy + dz * dz);
602
603 double hit_localPhiErr = corMod->getPhiStripWidth() / sqrt(12);
604 double hit_localZErr = corMod->getZStripWidth() / sqrt(12);
605
606 //error from tracking is ignored here
607 error = sqrt(pow(hit_localPhiErr, 2) +
608 pow(hit_localZErr, 2));
609
610 if (error != 0.0) {
611 sigma = distance / error;
612 } else {
613 sigma = DBL_MAX;
614 }
615
616 return (distance);
617
618}
track finding procedure
bool filter(const std::list< KLMHit2d * > &seed, std::list< KLMHit2d * > &hits, std::list< KLMHit2d * > &track)
find associated hits and do fit.
void registerFitter(BKLMTrackFitter *fitter)
Register a fitter if not constructed with one.
void setGlobalFit(bool localOrGlobal)
set the fitting mode, local system or global system
track fitting procedure
CLHEP::HepVector getTrackParamSector()
Get track parameters in the sector locan system, y = p0 + p1 * x, z = p2 + p3 *x, where the first lay...
float getChi2()
Chi square of the fit.
bool isGood()
Is fit good.
CLHEP::HepSymMatrix getTrackParamSectorErr()
Get invariance matrix of track parameters in the sector local system, where the first layer of the se...
int getNumHit()
number of the hits on this track
CLHEP::HepSymMatrix getTrackParamErr()
Get invariance matrix of track parameters in the global system.
CLHEP::HepVector getTrackParam()
Get track parameters in the global system. y = p0 + p1 * x; y = p2 + p3 * z, if in local sector fit m...
bool isValid()
Is fit valid.
Store one BKLM Track as a ROOT object.
Definition: BKLMTrack.h:35
void setIsValid(const bool valid)
set the fit valid status
Definition: BKLMTrack.h:120
void setTrackChi2(const float chi2)
Set the fitted chi2 of the track.
Definition: BKLMTrack.h:108
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:140
void setTrackParamErr(const CLHEP::HepSymMatrix &trkParErr)
Set invariance matrix of track parameters in the global system.
Definition: BKLMTrack.cc:130
void setNumHitOnTrack(const int NumHit)
Set the number of 2d hits on the track.
Definition: BKLMTrack.h:114
TVectorD getTrackParam()
Get track parameters in the global system. y = p0 + p1 * x; z = p2 + p3 * x.
Definition: BKLMTrack.cc:71
void setIsGood(const bool good)
set the fit good status
Definition: BKLMTrack.h:126
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:150
void setTrackParam(const CLHEP::HepVector &trkPar)
Set track parameters in the global system. y = p0 + p1 * x; z = p2 + p3 * x.
Definition: BKLMTrack.cc:120
bool m_MatchToRecoTrack
whether match BKLMTrack to RecoTrack
TEfficiency * m_effiVsLayer[2][8]
Efficiency of each layer.
std::vector< int > m_runNumber
run number
bool findClosestRecoTrack(BKLMTrack *bklmTrk, RecoTrack *&closestTrack)
find the closest RecoTrack, match BKLMTrack to RecoTrack, if the matched RecoTrack is found,...
TH2F * m_passYZ
passed event at global position Y vs Z
TH2F * m_effiYX
Efficiency at global position Y vs X.
bool m_studyEffi
option for efficiency study mode, in this mode, the layer under study should not be used in tracking
double m_maxSigma
maximum sigma for hit acceptance during efficiency calculation
double m_maxAngleRequired
angle required between RecoTrack and BKLMTrack, if openangle is larger than m_maxAngleRequired,...
double distanceToHit(BKLMTrack *track, KLMHit2d *hit, double &error, double &sigma)
calculate distance from track to hit
void initialize() override
Initialize at start of job.
unsigned int m_minHitList
minimum number of hits in sector for track finder to run (-2 from initial seed)
std::string m_outPath
output file name containing efficiencies plots
void event() override
Unpack one event and create digits.
void endRun() override
end run stuff
StoreArray< RecoTrack > recoTracks
RecoTrack StoreArray.
void runTracking(int mode, int section, int sector, int layer)
run the track finding and fitting
void terminate() override
Terminate at the end of job.
bklm::GeometryPar * m_GeoPar
bklm GeometryPar
TH1F * m_pass[2][8]
Numerator of each layer.
TH2F * m_totalYX
total event at global position Y vs X
StoreArray< RecoHitInformation > recoHitInformation
RecoHitInformation StoreArray.
int m_runTotalEventsWithTracks
total number of processed events in the run with at lease one BKLMTrack
void beginRun() override
begin run stuff
double m_maxDistance
maximum distance required between track and KLMHit2d to be accepted for efficiency calculation
StoreArray< KLMHit2d > hits2D
KLMHit2d StoreArray.
TFile * m_file
TFile that store efficiency plots.
std::vector< int > m_totalEvents
total number of processed events
bool isLayerUnderStudy(int section, int iSector, int iLayer, KLMHit2d *hit)
judge whether the current layer is understudy
TH1F * m_total[2][8]
Denominator of each layer.
StoreArray< BKLMTrack > m_storeTracks
BKLMTrack StoreArray.
bool m_globalFit
do the BKLMTrack fitting in global system (multi-sectors track) or local system (sector by sector)
unsigned int m_maxHitList
max number of hits in sector for track finder to run
std::vector< int > m_totalEventsWithTracks
total number of processed events with at least one BKLMTrack
TH2F * m_effiYZ
Efficiency at global position Y vs Z.
bool isSectorUnderStudy(int section, int iSector, KLMHit2d *hit)
judge whether the hits come from the sctor understudy
static bool sortByLayer(KLMHit2d *hit1, KLMHit2d *hit2)
my defined sort function using layer number
TH2F * m_totalYZ
total event at global position Y vs Z
bool sameSector(KLMHit2d *hit1, KLMHit2d *hit2)
Judge if two hits come from the same sector.
int m_runTotalEvents
total number of processed events in the run
TH2F * m_passYX
passed event at global position Y vs X
void generateEffi(int section, int sector, int layer)
calculate efficiency
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Definition: DataStore.h:59
KLM 2d hit.
Definition: KLMHit2d.h:33
int getLayer() const
Get layer number.
Definition: KLMHit2d.h:132
int getSection() const
Get section number.
Definition: KLMHit2d.h:96
int getSector() const
Get sector number.
Definition: KLMHit2d.h:114
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
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:286
unsigned int getNumberOfTotalHits() const
Return the number of cdc + svd + pxd + bklm + eklm hits.
Definition: RecoTrack.h:436
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
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).
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
const Module * findModule(int section, int sector, int layer) const
Get the pointer to the definition of a module.
Definition: GeometryPar.cc:721
double getActiveMiddleRadius(int section, int sector, int layer) const
Get the radial midpoint of the detector module's active volume of specified layer.
Definition: GeometryPar.cc:607
static GeometryPar * instance(void)
Static method to get a reference to the singleton GeometryPar instance.
Definition: GeometryPar.cc:27
Define the geometry of a BKLM module Each sector [octant] contains Modules.
Definition: Module.h:76
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:339
double getPhiStripWidth() const
Get phi-strip width.
Definition: Module.h:137
double getZStripWidth() const
Get z-strip width.
Definition: Module.h:155
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
ExpRunEvt getPosition(const std::vector< Evt > &events, double tEdge)
Get the exp-run-evt number from the event time [hours].
Definition: Splitter.h:341
Abstract base class for different kinds of events.