Belle II Software development
KLMTrackingModule.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/modules/KLMTracking/KLMTrackingModule.h>
11
12/* KLM headers. */
13#include <klm/bklm/geometry/GeometryPar.h>
14#include <klm/modules/KLMTracking/KLMTrackFinder.h>
15#include <klm/eklm/geometry/TransformDataGlobalAligned.h> //TODO: Is this the right module
16
17/* Basf2 headers. */
18#include <framework/dataobjects/EventMetaData.h>
19#include <framework/datastore/StoreObjPtr.h>
20#include <framework/datastore/StoreArray.h>
21#include <framework/logging/Logger.h>
22#include <tracking/dataobjects/RecoHitInformation.h>
23
24/* C++ standard libraries*/
25#include <set>
26#include <iostream>
27
28using namespace Belle2;
29using namespace Belle2::KLM;
30using namespace CLHEP;
31
32REG_MODULE(KLMTracking);
33
35 m_effiYX(nullptr),
36 m_effiYZ(nullptr),
37 m_passYX(nullptr),
38 m_totalYX(nullptr),
39 m_passYZ(nullptr),
40 m_totalYZ(nullptr),
41 m_runTotalEvents(0),
42 m_runTotalEventsWithTracks(0)
43{
44 for (int i = 0; i < 8; ++i) {
45 m_total[0][i] = nullptr;
46 m_total[1][i] = nullptr;
47 m_pass[0][i] = nullptr;
48 m_pass[1][i] = nullptr;
49 m_effiVsLayer[0][i] = nullptr;
50 m_effiVsLayer[1][i] = nullptr;
51 }
52 setDescription("Perform standard-alone straight line tracking for KLM. ");
53 addParam("MatchToRecoTrack", m_MatchToRecoTrack, "[bool], whether match KLMTrack to RecoTrack; (default is false)", false);
54 addParam("MaxAngleRequired", m_maxAngleRequired,
55 "[degree], match KLMTrack to RecoTrack; angle between them is required to be smaller than (default 10)", double(10.0));
56 addParam("MaxDistance", m_maxDistance,
57 "[cm], During efficiency calculation, distance between track and 2dhit must be smaller than (default 10)", double(10.0));
58 addParam("MaxSigma", m_maxSigma,
59 "[sigma], During efficiency calculation, uncertainty of 2dhit must be smaller than (default 5); ", double(5));
60 addParam("MinHitList", m_minHitList,
61 ", During track finding, a good track after initial seed hits must be larger than is (default 2); ", unsigned(2));
62 addParam("MaxHitList", m_maxHitList,
63 ", During track finding, a good track after initial seed hits must be smaller than is (default 60); ", unsigned(60));
64 addParam("MinNLayer", m_minNLayer,
65 ", Only look at tracks with more than n number of layers; ", int(4));
66 addParam("StudyEffiMode", m_studyEffi, "[bool], run in efficieny study mode (default is false)", false);
67 addParam("outputName", m_outPath, "[string], output file name containing efficiencies plots ",
68 std::string("standaloneKLMEffi.root"));
69}
70
72{
73
74}
75
77{
78
79 hits2D.isRequired();
80 m_storeTracks.registerInDataStore();
81 m_storeTracks.registerRelationTo(hits2D);
82 m_storeTracks.registerRelationTo(recoTracks);
83 recoHitInformation.registerRelationTo(hits2D);
84 hits2D.registerRelationTo(recoTracks);
85
86 if (m_studyEffi)
87 B2INFO("KLMTrackingModule::initialize this module is running in efficiency study mode!");
88
89 m_file = new TFile(m_outPath.c_str(), "recreate");
90 TString hname;
91 std::string labelFB[2] = {"BB", "BF"};
92 int Nbin = 16;
93 float gmin = -350;
94 float gmax = 350;
95 int gNbin = 150;
96
97 //TODO: Extend to include EKLM...
98 m_totalYX = new TH2F("totalYX", " denominator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
99 m_passYX = new TH2F("passYX", " numerator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
100 m_totalYZ = new TH2F("totalYZ", " denominator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
101 m_passYZ = new TH2F("passYZ", " numerator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
102 m_effiYX = new TH2F("effiYX", " effi. Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
103 m_effiYZ = new TH2F("effiYZ", " effi. Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
104 m_effiYX->GetXaxis()->SetTitle("x (cm)");
105 m_effiYX->GetYaxis()->SetTitle("y (cm)");
106 m_effiYZ->GetXaxis()->SetTitle("z (cm)");
107 m_effiYZ->GetYaxis()->SetTitle("y (cm)");
108 for (int iF = 0; iF < 2; iF++) {
109 for (int iS = 0; iS < 8; iS++) {
110 hname.Form("effi_%s%i", labelFB[iF].c_str(), iS);
111 m_effiVsLayer[iF][iS] = new TEfficiency(hname, hname, Nbin, 0, 16);
112 hname.Form("total_%s%i", labelFB[iF].c_str(), iS);
113 m_total[iF][iS] = new TH1F(hname, hname, Nbin, 0, 16);
114 hname.Form("pass_%s%i", labelFB[iF].c_str(), iS);
115 m_pass[iF][iS] = new TH1F(hname, hname, Nbin, 0, 16);
116 }
117 }
118
119 //EKLM Plots TODO: Not tested yet
120 /*
121 m_totalYXE = new TH2F("totalYX", " denominator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
122 m_passYXE = new TH2F("passYX", " numerator Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
123 m_totalYZE = new TH2F("totalYZ", " denominator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
124 m_passYZE = new TH2F("passYZ", " numerator Y vs. Z", gNbin, gmin, gmax, gNbin, gmin, gmax);
125 m_effiYXE = new TH2F("effiYX", " effi. Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
126 m_effiYZE = new TH2F("effiYZ", " effi. Y vs. X", gNbin, gmin, gmax, gNbin, gmin, gmax);
127 m_effiYXE->GetXaxis()->SetTitle("x (cm)");
128 m_effiYXE->GetYaxis()->SetTitle("y (cm)");
129 m_effiYZE->GetXaxis()->SetTitle("z (cm)");
130 m_effiYZE->GetYaxis()->SetTitle("y (cm)");
131 for (int iF = 0; iF < 2; iF++) {
132 for (int iS = 0; iS < 8; iS++) {
133 hname.Form("effi_%s%i", labelFB[iF].c_str(), iS);
134 m_effiVsLayer[iF][iS] = new TEfficiency(hname, hname, Nbin, 0, 16);
135 hname.Form("total_%s%i", labelFB[iF].c_str(), iS);
136 m_total[iF][iS] = new TH1F(hname, hname, Nbin, 0, 16);
137 hname.Form("pass_%s%i", labelFB[iF].c_str(), iS);
138 m_pass[iF][iS] = new TH1F(hname, hname, Nbin, 0, 16);
139 }
140 } //end of EKLM layer info
141 */
142
143}
144
146{
147 StoreObjPtr<EventMetaData> eventMetaData("EventMetaData", DataStore::c_Event);
148 m_runNumber.push_back((int)eventMetaData->getRun());
151}
152
154{
155 m_storeTracks.clear();
156 bool thereIsATrack = false;
157
158 if (!m_studyEffi) {
161 if (m_storeTracks.getEntries() > 0)
162 thereIsATrack = true;
163 } else {
164 for (int iSection = 0; iSection < 2; iSection++) {
165 for (int iSector = 0; iSector < 8; iSector++) {
166 for (int iLayer = 0; iLayer < 15; iLayer++) {
167 runTracking(1, KLMElementNumbers::c_BKLM, iSection, iSector, iLayer);
168 if (m_storeTracks.getEntries() > 0)
169 thereIsATrack = true;
170 generateEffi(KLMElementNumbers::c_BKLM, iSection, iSector, iLayer);
171 //clear tracks so prepare for the next layer efficieny study
172 m_storeTracks.clear();
173 }
174 }
175 }
176 }
177
179 if (thereIsATrack)
181}
182
183void KLMTrackingModule::runTracking(int mode, int iSubdetector, int iSection, int iSector, int iLayer)
184{
185 //m_storeTracks.clear(); //done in event stage
186
187 KLMTrackFitter* m_fitter = new KLMTrackFitter();
188 KLMTrackFinder* m_finder = new KLMTrackFinder();
189 m_finder->registerFitter(m_fitter);
190
191 if (hits2D.getEntries() < 1)
192 return;
193 if (mode == 1) { //efficieny study
194 for (int j = 0; j < hits2D.getEntries(); j++) {
195 if (hits2D[j]->getSubdetector() != iSubdetector) //TODO: shoud we kee?
196 continue;
197 hits2D[j]->isOnStaTrack(false);
198 }
199 }
200
201 for (int hi = 0; hi < hits2D.getEntries() - 1; ++hi) {
202 if (hits2D[hi]->getSubdetector() != iSubdetector) //TODO: Should we keep?
203 continue;
204
205 if (mode == 1 && isLayerUnderStudy(iSection, iSector, iLayer, hits2D[hi]))
206 continue;
207 if (mode == 1 && !isSectorUnderStudy(iSection, iSector, hits2D[hi]))
208 continue;
209 if (hits2D[hi]->isOnStaTrack())
210 continue;
211 if (hits2D[hi]->isOutOfTime())
212 continue;
213 for (int hj = hi + 1; hj < hits2D.getEntries(); ++hj) {
214
215 if (hits2D[hj]->isOnStaTrack())
216 continue;
217 if (hits2D[hj]->isOutOfTime())
218 continue;
219 // at least for track seed, hits should remain in the same subdetector
220 if (hits2D[hi]->getSubdetector() != hits2D[hj]->getSubdetector())
221 continue;
222 if (sameSector(hits2D[hi], hits2D[hj]) &&
223 std::abs(hits2D[hi]->getLayer() - hits2D[hj]->getLayer()) < 3)
224 continue;
225
226 std::list<KLMHit2d*> sectorHitList;
227
228
229 std::list<KLMHit2d*> seed;
230 seed.push_back(hits2D[hi]);
231 seed.push_back(hits2D[hj]);
232
233 for (int ho = 0; ho < hits2D.getEntries(); ++ho) {
234
235 // Exclude seed hits.
236 if (ho == hi || ho == hj)
237 continue;
238 if (mode == 1 && (hits2D[ho]->getSubdetector() != iSubdetector))
239 continue;
240 if (mode == 1 && isLayerUnderStudy(iSection, iSector, iLayer, hits2D[hj]))
241 continue;
242 if (mode == 1 && !isSectorUnderStudy(iSection, iSector, hits2D[hj]))
243 continue;
244 if (hits2D[ho]->isOnStaTrack())
245 continue;
246 //TODO: consider removing the commented lines below
247 if (mode == 1 && !sameSector(hits2D[ho], hits2D[hi]))
248 continue;
249 if (hits2D[ho]->isOutOfTime())
250 continue;
251 sectorHitList.push_back(hits2D[ho]);
252 }
253
254 /* Require at least four hits (minimum for good track, already two as seed, so here we require 2) but
255 * no more than 60 (most likely noise, 60 would be four good tracks).
256 * TODO: Should be tuned since we have EKLM hits now. 60 was from BKLMTracking
257 */
258 if (sectorHitList.size() < m_minHitList || sectorHitList.size() > m_maxHitList)
259 continue;
260
261 std::list<KLMHit2d*> m_hits;
262
263 if (m_finder->filter(seed, sectorHitList, m_hits, iSubdetector)) {
264 KLMTrack* m_track = m_storeTracks.appendNew();
265 m_track->setTrackParam(m_fitter->getTrackParam());
266 m_track->setTrackParamErr(m_fitter->getTrackParamErr());
267 m_track->setTrackChi2(m_fitter->getChi2());
268 m_track->setNumHitOnTrack(m_fitter->getNumHit());
269 m_track->setIsValid(m_fitter->isValid());
270 m_track->setIsGood(m_fitter->isGood());
271 std::list<KLMHit2d*>::iterator j;
272 m_hits.sort(sortByLayer);
273 int nBKLM = 0; int nEKLM = 0;
274 for (j = m_hits.begin(); j != m_hits.end(); ++j) {
275 (*j)->isOnStaTrack(true);
276 m_track->addRelationTo((*j));
277 if ((*j)->getSubdetector() == KLMElementNumbers::c_BKLM)
278 nBKLM += 1;
279 else if ((*j)->getSubdetector() == KLMElementNumbers::c_EKLM)
280 nEKLM += 1;
281 } //end of klmhit2d loop
282 B2DEBUG(31, "KLMTracking::runTracking totalHit " << m_hits.size() << ", nBKLM " << nBKLM << ", nEKLM " << nEKLM);
283 m_track->setInSubdetector(nBKLM, nEKLM);
284
285 //match KLMTrack to RecoTrack
286 if (mode == 0) {
287 RecoTrack* closestTrack = nullptr;
288 B2DEBUG(30, "KLMTracking::runTracking started RecoTrack matching");
289 if (m_MatchToRecoTrack) {
290 if (findClosestRecoTrack(m_track, closestTrack)) {
291 B2DEBUG(30, "KLMTracking::runTracking was able to find ClosestRecoTrack");
292 m_track->addRelationTo(closestTrack);
293 for (j = m_hits.begin(); j != m_hits.end(); ++j) {
294 unsigned int sortingParameter = closestTrack->getNumberOfTotalHits();
295 if ((*j)->getSubdetector() == KLMElementNumbers::c_BKLM)
296 closestTrack->addBKLMHit((*j), sortingParameter, RecoHitInformation::OriginTrackFinder::c_LocalTrackFinder);
297 else if ((*j)->getSubdetector() == KLMElementNumbers::c_EKLM) {
298 for (const EKLMAlignmentHit& alignmentHit : (*j)->getRelationsFrom<EKLMAlignmentHit>()) {
299 closestTrack->addEKLMHit(&(alignmentHit), sortingParameter,
300 RecoHitInformation::OriginTrackFinder::c_LocalTrackFinder);
301 }
302 }
303 } //end of hit loop
304 } //end of closest recotrack
305 }
306 } //end of KLMTrack to RecoTrack
307 } //end of finder loop
308 }
309 }
310
311 delete m_fitter;
312 delete m_finder;
313
314}
315
317{
320}
321
323{
324 for (long unsigned int i = 0; i < m_runNumber.size(); i++) {
325 float ratio = (float)m_totalEventsWithTracks.at(i) / (float)m_totalEvents.at(i);
326 B2INFO("KLMTrackingModule::terminate run " << m_runNumber.at(i) << " --> " << ratio * 100 << "% of events has 1+ KLMTracks");
327 }
328
329 m_file->cd();
330 for (int iF = 0; iF < 2; iF++) {
331 for (int iS = 0; iS < 8; iS++) {
332 m_effiVsLayer[iF][iS]->Write();
333 m_total[iF][iS]->Write();
334 m_pass[iF][iS]->Write();
335 }
336 }
337
338 for (int i = 0; i < m_totalYX->GetNbinsX(); i++) {
339 for (int j = 0; j < m_totalYX->GetNbinsY(); j++) {
340 float num = m_passYX->GetBinContent(i + 1, j + 1);
341 float denom = m_totalYX->GetBinContent(i + 1, j + 1);
342 if (num > 0) {
343 m_effiYX->SetBinContent(i + 1, j + 1, num / denom);
344 m_effiYX->SetBinError(i + 1, j + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
345 } else {
346 m_effiYX->SetBinContent(i + 1, j + 1, 0);
347 m_effiYX->SetBinError(i + 1, j + 1, 0);
348 }
349
350 num = m_passYZ->GetBinContent(i + 1, j + 1);
351 denom = m_totalYZ->GetBinContent(i + 1, j + 1);
352 if (num > 0) {
353 m_effiYZ->SetBinContent(i + 1, j + 1, num / denom);
354 m_effiYZ->SetBinError(i + 1, j + 1, sqrt(num * (denom - num) / (denom * denom * denom)));
355 } else {
356 m_effiYZ->SetBinContent(i + 1, j + 1, 0);
357 m_effiYZ->SetBinError(i + 1, j + 1, 0);
358 }
359 }
360 }
361
362
363 m_totalYX->SetOption("colz");
364 m_passYX->SetOption("colz");
365 m_totalYZ->SetOption("colz");
366 m_passYZ->SetOption("colz");
367 m_effiYX->SetOption("colz");
368 m_effiYZ->SetOption("colz");
369
370 m_totalYX->Write();
371 m_passYX->Write();
372 m_totalYZ->Write();
373 m_passYZ->Write();
374 m_effiYX->Write();
375 m_effiYZ->Write();
376 m_file->Close();
377
378}
379
381{
382 if (hit1->getSection() == hit2->getSection() && hit1->getSector() == hit2->getSector())
383 return true;
384 else return false;
385}
386
387
389{
390
391 //StoreArray<RecoTrack> recoTracks;
393
394 if (klmHits.size() < 1) {
395 B2INFO("KLMTrackingModule::findClosestRecoTrack, something is wrong! there is a KLMTrack but no klmHits");
396 return false;
397 }
398 if (recoTracks.getEntries() < 1) {
399 B2DEBUG(20, "KLMTrackingModule::findClosestRecoTrack, there is no recoTrack");
400 return false;
401 }
402 double oldDistanceSq = INFINITY;
403 double oldAngle = INFINITY;
404 closestTrack = nullptr;
405 //klmHits are already sorted by layer
406 //possible two hits in one layer?
407 //genfit requires TVector3 rather than XYZVector
408
409 TVector3 firstKLMHitPosition(klmHits[0]->getPosition().X(),
410 klmHits[0]->getPosition().Y(),
411 klmHits[0]->getPosition().Z());
412
413 // To get direction (angle) below, we have two points on the klmTrk:
414 // (x1, TrackParam[0]+TrackParam[1]*x1, TrackParam[2]+TrackParam[3]*x1)
415 // (x2, TrackParam[0]+TrackParam[1]*x2, TrackParam[2]+TrackParam[3]*x2)
416 // the difference vector is
417 // (x2-x1, TrackParam[1]*(x2-x1), TrackParam[3]*(x2-x1))
418 // which is proportional to
419 // (1, TrackParam[1], TrackParam[3]).
420 TVector3 klmTrkVec(1.0, klmTrk->getTrackParam()[1], klmTrk->getTrackParam()[3]);
421
422 TMatrixDSym cov(6);
423 TVector3 pos; // initializes to (0,0,0)
424 TVector3 mom; // initializes to (0,0,0)
425
426 for (RecoTrack& track : recoTracks) {
427 if (track.wasFitSuccessful()) {
428 try {
429 genfit::MeasuredStateOnPlane state = track.getMeasuredStateOnPlaneFromLastHit();
430 B2DEBUG(30, "KLMTracking::findClosestRecoTrack, finished MSOP from last hit");
432 state.getPosMomCov(pos, mom, cov);
433 if (mom.Y() * pos.Y() < 0) {
434 state = track.getMeasuredStateOnPlaneFromFirstHit();
435 }
436 const TVector3& distanceVec = firstKLMHitPosition - pos;
437 state.extrapolateToPoint(firstKLMHitPosition);
438 double newDistanceSq = distanceVec.Mag2();
439 double angle = klmTrkVec.Angle(mom);
440 // choose closest distance or minimum open angle ?
441 // overwrite old distance
442 if (newDistanceSq < oldDistanceSq) {
443 oldDistanceSq = newDistanceSq;
444 closestTrack = &track;
445 oldAngle = angle;
446 }
447
448 /* if(angle<oldAngle)
449 {
450 oldAngle=angle;
451 closestTrack = &track;
452 }
453 */
454 B2DEBUG(30, "KLMTracking::findClosestRecoTrack, step one done");
455 } catch (genfit::Exception& e) {
456 }// try
457 }
458 }
459
460 // can not find matched RecoTrack
461 // problem here is the errors of the track parameters are not considered!
462 // best way is the positon or vector direction are required within 5/10 sigma ?
463 if (oldAngle > m_maxAngleRequired)
464 return false;
465 // found matched RecoTrack
466 else {
467 B2DEBUG(28, "KLMTrackingModule::findClosestRecoTrack RecoTrack found! ");
468 return true;
469 }
470
471}
472
473//TODO: GENERALIZE ME [HARDEST PART!? :'( ]
474void KLMTrackingModule::generateEffi(int iSubdetector, int iSection, int iSector, int iLayer)
475{
476 //TODO: let's comment out during testing. remove this later
477
478 std::set<int> m_pointUsed;
479 std::set<int> layerList;
480 m_pointUsed.clear();
481 if (m_storeTracks.getEntries() < 1)
482 return;
483 B2DEBUG(10, "KLMTrackingModule:generateEffi: " << iSection << " " << iSector << " " << iLayer);
484
485 for (int it = 0; it < m_storeTracks.getEntries(); it++) {
486 //if(m_storeTracks[it]->getTrackChi2()>10) continue;
487 //if(m_storeTracks[it]->getNumHitOnTrack()<6) continue;
488 int cnt1 = 0;
489 int cnt2 = 0;
490
491 layerList.clear();
492
493
494 RelationVector<KLMHit2d> relatedHit2D = m_storeTracks[it]->getRelationsTo<KLMHit2d>();
495 for (const KLMHit2d& hit2D : relatedHit2D) {
496 if (hit2D.getSubdetector() != iSubdetector)
497 continue;
498 if (hit2D.getLayer() > iLayer + 1)
499 {cnt1++; layerList.insert(hit2D.getLayer());}
500 if (hit2D.getLayer() < iLayer + 1)
501 {cnt2++; layerList.insert(hit2D.getLayer());}
502 if (hit2D.getLayer() == iLayer + 1) {
503 B2DEBUG(10, "generateEffi: Hit info. Secti/sector/Lay = " << hit2D.getSection()
504 << "/" << hit2D.getSector() - 1 << "/" << hit2D.getLayer() - 1);
505 B2DEBUG(11, "generateEffi: Hit info. x/y/z = " << hit2D.getPositionX()
506 << "/" << hit2D.getPositionY() << "/" << hit2D.getPositionZ());
507 }
508 }
509
510 if ((int)layerList.size() < m_minNLayer)
511 continue;
512
513 if (iLayer != 0 && cnt2 < 1)
514 return;
515 if (iLayer != 14 && cnt1 < 1)
516 return;
517 //TODO: Extend to includ EKLM?
518 //m_GeoPar = GeometryPar::instance(); w/ geometry cuts.
519
520 if (iSubdetector == KLMElementNumbers::c_BKLM) {
521
522 const bklm::GeometryPar* bklmGeo = m_GeoPar->BarrelInstance();
523 const bklm::Module* module = bklmGeo->findModule(iSection, iSector + 1, iLayer + 1);
524 const bklm::Module* refmodule = bklmGeo->findModule(iSection, iSector + 1, 1);
525 int minPhiStrip = module->getPhiStripMin();
526 int maxPhiStrip = module->getPhiStripMax();
527 int minZStrip = module->getZStripMin();
528 int maxZStrip = module->getZStripMax();
529
530 CLHEP::Hep3Vector local = module->getLocalPosition(minPhiStrip, minZStrip);
531 CLHEP::Hep3Vector local2 = module->getLocalPosition(maxPhiStrip, maxZStrip);
532 float minLocalY, maxLocalY;
533 float minLocalZ, maxLocalZ;
534 if (local[1] > local2[1]) {
535 maxLocalY = local[1];
536 minLocalY = local2[1];
537 } else {
538 maxLocalY = local2[1];
539 minLocalY = local[1];
540 }
541 if (local[2] > local2[2]) {
542 maxLocalZ = local[2];
543 minLocalZ = local2[2];
544 } else {
545 maxLocalZ = local2[2];
546 minLocalZ = local[2];
547 }
548
549 //in global coordinates
550 TVectorD trkPar = m_storeTracks[it]->getTrackParam();
551
552 //line in local coordinates
553 Hep3Vector point1(0, trkPar[0], trkPar[2]);
554 Hep3Vector point2(1, trkPar[0] + trkPar[1], trkPar[2] + trkPar[3]);
555
556 Hep3Vector refPoint1(0., 0., 0.); Hep3Vector refPoint2(0., 0., 0.);
557 refPoint1 = refmodule->globalToLocal(point1);
558 refPoint2 = refmodule->globalToLocal(point2);
559
560 Hep3Vector refSlope(refPoint2[0] - refPoint1[0], refPoint2[1] - refPoint1[1], refPoint2[2] - refPoint1[2]);
561
562
563
564 //defined in coordinates relative to layer 1 of this sector.
565 float reflocalX = fabs(bklmGeo->getActiveMiddleRadius(iSection, iSector + 1,
566 iLayer + 1) - bklmGeo->getActiveMiddleRadius(iSection, iSector + 1, 1));
567 if (refmodule->isFlipped())
568 reflocalX = -reflocalX;
569 float X_coord = (reflocalX - refPoint1[0]) / refSlope[0];
570 float reflocalY = refPoint1[1] + refSlope[1] * X_coord;
571 float reflocalZ = refPoint1[2] + refSlope[2] * X_coord;
572
573 //Hep3Vector global(globalX, globalY, globalZ);
574 Hep3Vector reflocal(reflocalX, reflocalY, reflocalZ);
575 Hep3Vector global(0., 0., 0.);
576 global = refmodule->localToGlobal(reflocal);
577
578
579 float localX = module->globalToLocal(global)[0];
580 float localY = module->globalToLocal(global)[1];
581 float localZ = module->globalToLocal(global)[2];
582
583 B2DEBUG(10, "KLMTrackingModule:generateEffi: RefLocal " << reflocalX << " " << reflocalY << " " << reflocalZ);
584 B2DEBUG(10, "KLMTrackingModule:generateEffi: Global " << global[0] << " " << global[1] << " " << global[2]);
585 B2DEBUG(10, "KLMTrackingModule:generateEffi: Local " << localX << " " << localY << " " << localZ);
586
587
588
589 //geometry cut
590 if (localY > minLocalY && localY < maxLocalY && localZ > minLocalZ && localZ < maxLocalZ) {
591
592 bool m_iffound = false;
593 m_total[iSection][iSector]->Fill(iLayer + 1);
594 m_totalYX->Fill(global[0], global[1]);
595 m_totalYZ->Fill(global[2], global[1]);
596
597 for (int he = 0; he < hits2D.getEntries(); ++he) {
598 if (!isLayerUnderStudy(iSection, iSector, iLayer, hits2D[he])) {
599 B2DEBUG(11, "not isLayerUnderStudy");
600 continue;
601 }
602 if (hits2D[he]->isOutOfTime()) {
603 B2DEBUG(11, "hit isOutOfTime");
604 continue;
605 }
606 //if alreday used, skip
607 if (m_pointUsed.find(he) != m_pointUsed.end()) {
608 B2DEBUG(11, "passed unused");
609 continue;
610 }
611 B2DEBUG(11, "KLMTrackingModule:generateEffi: Reached Distance Check");
612 double error, sigma;
613 float distance = distanceToHit(m_storeTracks[it], hits2D[he], error, sigma);
614 float hitX = hits2D[he]->getPositionX();
615 float hitY = hits2D[he]->getPositionY();
616 float hitZ = hits2D[he]->getPositionZ();
617 float deltaX = hitX - global[0]; float deltaY = hitY - global[1]; float deltaZ = hitZ - global[2];
618 float dist = sqrt(deltaX * deltaX + deltaY * deltaY + deltaZ * deltaZ);
619 B2DEBUG(10, "dist w/ hit = " << dist << ", dist func = " << distance << ", error = " << error);
620 if (distance < m_maxDistance && sigma < m_maxSigma) {
621 m_iffound = true;
622 B2DEBUG(10, "KLMTrackingModule:generateEffi: Hit found!");
623 }
624 if (m_iffound) {
625 m_pointUsed.insert(he);
626 m_pass[iSection][iSector]->Fill(iLayer + 1);
627 m_passYX->Fill(global[0], global[1]);
628 m_passYZ->Fill(global[2], global[1]);
629 break;
630 }
631 }
632
633 m_effiVsLayer[iSection][iSector]->Fill(m_iffound, iLayer + 1);
634 //efficiencies will be defined at terminate stage
635 } //end of BKLM geometry cut
636
637 } //end of BKLM section
638
639
640 }//end of loop tracks
641}
642
643
645{
646
647 return hit1->getLayer() < hit2->getLayer();
648
649}
650
651bool KLMTrackingModule::isLayerUnderStudy(int section, int iSector, int iLayer, KLMHit2d* hit)
652{
653 if (hit->getSection() == section && hit->getSector() == iSector + 1 && hit->getLayer() == iLayer + 1)
654 return true;
655 else return false;
656}
657
658bool KLMTrackingModule::isSectorUnderStudy(int section, int iSector, KLMHit2d* hit)
659{
660 if (hit->getSection() == section && hit->getSector() == iSector + 1)
661 return true;
662 else return false;
663}
664
666 double& error,
667 double& sigma)
668{
669
670 double x, y, z, dx, dy, dz, distance;
671
672 error = DBL_MAX;
673 sigma = DBL_MAX;
674
675 TVectorD m_GlobalPar = track->getTrackParam();
676
677
678 if (hit->getSubdetector() == KLMElementNumbers::c_BKLM) {
679
680 const bklm::GeometryPar* bklmGeo = m_GeoPar->BarrelInstance();
681 const Belle2::bklm::Module* corMod = bklmGeo->findModule(hit->getSection(), hit->getSector(), hit->getLayer());
682
683 //x = hit->getPositionX();
684 //y = m_GlobalPar[ 0 ] + x * m_GlobalPar[ 1 ];
685 //z = m_GlobalPar[ 2 ] + x * m_GlobalPar[ 3 ];
686
687 //since there are z-planes, let's exploit this fact.
688 z = hit->getPositionZ();
689 x = (z - m_GlobalPar[ 2 ]) / m_GlobalPar[ 3 ];
690 y = m_GlobalPar[ 0 ] + x * m_GlobalPar[ 1 ];
691
692 dx = x - hit->getPositionX() ;
693 dy = y - hit->getPositionY();
694 dz = z - hit->getPositionZ();
695
696 double x2 = hit->getPositionX();
697 double y2 = m_GlobalPar[ 0 ] + x2 * m_GlobalPar[ 1 ];
698 double z2 = m_GlobalPar[ 2 ] + x2 * m_GlobalPar[ 3 ];
699
700 double dx2 = x2 - hit->getPositionX();
701 double dy2 = y2 - hit->getPositionY();
702 double dz2 = z2 - hit->getPositionZ();
703
704 double dist2 = sqrt(dx2 * dx2 + dy2 * dy2 + dz2 * dz2);
705
706
707 distance = sqrt(dx * dx + dy * dy + dz * dz);
708 double hit_localPhiErr = corMod->getPhiStripWidth() / sqrt(12);
709 double hit_localZErr = corMod->getZStripWidth() / sqrt(12);
710
711 //error from tracking is ignored here
712 error = sqrt(pow(hit_localPhiErr, 2) +
713 pow(hit_localZErr, 2));
714 B2DEBUG(11, "Dist = " << distance << ", error = " << error);
715 B2DEBUG(11, "Dist2 = " << dist2 << ", error = " << error);
716 } //end of BKLM section
717
718 else if (hit->getSubdetector() == KLMElementNumbers::c_EKLM) {
719
720 const EKLM::GeometryData* eklmGeo = m_GeoPar->EndcapInstance();
721
722
723 // use z coordinate as main point of interest
724 // should be close enough to distance of closest appraoch
725 z = hit->getPositionZ();
726 x = (z - m_GlobalPar[ 2 ]) / m_GlobalPar[ 3 ];
727 y = m_GlobalPar[ 0 ] + x * m_GlobalPar[ 1 ];
728
729 dx = x - hit->getPositionX();
730 dy = y - hit->getPositionY();
731 dz = 0.;
732
733 distance = sqrt(dx * dx + dy * dy + dz * dz);
734
735
736 //here get the resolustion of a hit, repeated several times, ugly. should we store this in KLMHit2d object ?
737 double hit_xErr = (eklmGeo->getStripGeometry()->getWidth()) * (Unit::cm / CLHEP::cm) *
738 (hit->getXStripMax() - hit->getXStripMin()) / sqrt(12);
739 double hit_yErr = (eklmGeo->getStripGeometry()->getWidth()) * (Unit::cm / CLHEP::cm) *
740 (hit->getYStripMax() - hit->getYStripMin()) / sqrt(12);
741
742
743 //error from tracking is ignored here
744 error = sqrt(pow(hit_xErr, 2) +
745 pow(hit_yErr, 2));
746 } //end of EKLM section
747 else {
748 B2WARNING("KLMTracking::distanceToHit Received KLMHit2d that's not from E/B-KLM. Setting distance to -1");
749 distance = -1.;
750 }
751
752
753 if (error != 0.0) {
754 sigma = distance / error;
755 } else {
756 sigma = DBL_MAX;
757 }
758
759 return distance;
760
761}
762
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Definition: DataStore.h:59
This dataobject is used only for EKLM alignment.
double getWidth() const
Get width.
const StripGeometry * getStripGeometry() const
Get strip geometry data.
EKLM geometry data.
Definition: GeometryData.h:38
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
track finding procedure
bool filter(const std::list< KLMHit2d * > &seed, std::list< KLMHit2d * > &hits, std::list< KLMHit2d * > &track, int iSubdetector)
find associated hits and do fit.
void registerFitter(KLMTrackFitter *fitter)
Register a fitter if not constructed with one.
track fitting procedure
float getChi2()
Chi square of the fit.
bool isGood()
Is fit good.
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; z = p2 + p3 * x.
bool isValid()
Is fit valid.
Store one KLM Track as a ROOT object.
Definition: KLMTrack.h:38
void setIsValid(const bool valid)
set the fit valid status
Definition: KLMTrack.h:112
void setTrackChi2(const float chi2)
Set the fitted chi2 of the track.
Definition: KLMTrack.h:93
void setTrackParamErr(const CLHEP::HepSymMatrix &trkParErr)
Set invariance matrix of track parameters in the global system.
Definition: KLMTrack.cc:102
void setNumHitOnTrack(const int NumHit)
Set the number of 2d hits on the track.
Definition: KLMTrack.h:99
TVectorD getTrackParam()
Get track parameters in the global system. y = p0 + p1 * x; z = p2 + p3 * x.
Definition: KLMTrack.cc:66
void setIsGood(const bool good)
set the fit good status
Definition: KLMTrack.h:118
void setInSubdetector(int nBKLM, int nEKLM)
setting whether track passes through E/B-KLM
Definition: KLMTrack.h:105
void setTrackParam(const CLHEP::HepVector &trkPar)
Set track parameters in the global system. y = p0 + p1 * x; z = p2 + p3 * x.
Definition: KLMTrack.cc:92
void generateEffi(int iSubdetector, int section, int sector, int layer)
calculate efficiency
bool m_MatchToRecoTrack
whether match KLMTrack to RecoTrack
TEfficiency * m_effiVsLayer[2][8]
Efficieny of each layer.
std::vector< int > m_runNumber
run number
TH2F * m_passYZ
passed event at global position Y vs Z
TH2F * m_effiYX
Efficieny at global position Y vs X.
bool m_studyEffi
option for efficieny 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
void runTracking(int mode, int iSubdetector, int section, int sector, int layer)
run the track finding and fitting
double m_maxAngleRequired
angle required between RecoTrack and KLMTrack, if openangle is larger than m_maxAngleRequired,...
int m_minNLayer
minimum number of layers for track finder to run
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.
Belle2::KLM::KLMGeometryPar * m_GeoPar
KLMGeometryPar to call on B/E-KLM.
bool findClosestRecoTrack(KLMTrack *klmTrk, RecoTrack *&closestTrack)
find the closest RecoTrack, match KLMTrack to RecoTrack, if the matched RecoTrack is found,...
void endRun() override
end run stuff
StoreArray< RecoTrack > recoTracks
RecoTrack StoreArray.
void terminate() override
Terminate at the end of job.
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< KLMTrack > m_storeTracks
KLMTrack StoreArray.
StoreArray< KLMHit2d > hits2D
KLMHit2d StoreArray.
TFile * m_file
TFile that store efficieny 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.
double distanceToHit(KLMTrack *track, KLMHit2d *hit, double &error, double &sigma)
calculate distance from track to hit
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
Efficieny at global position Y vs Z.
bool isSectorUnderStudy(int section, int iSector, KLMHit2d *hit)
judge whether the hits come from the sctor under study
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
static const bklm::GeometryPar * BarrelInstance()
Return a pointer to the bklm::GeometryPar instance.
static const EKLM::GeometryData * EndcapInstance()
Return a pointer to the EKLM::GeometryData instance.
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
bool addEKLMHit(const UsedEKLMHit *eklmHit, const unsigned int sortingParameter, OriginTrackFinder foundByTrackFinder=OriginTrackFinder::c_undefinedTrackFinder)
Adds an eklm hit with the given information to the reco track.
Definition: RecoTrack.h:300
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< FROM > getRelationsFrom(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from another store array to this object.
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
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47
Provides BKLM geometry parameters for simulation, reconstruction etc (from Gearbox or DataBase)
Definition: GeometryPar.h:37
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
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
const CLHEP::Hep3Vector localToGlobal(const CLHEP::Hep3Vector &v, bool reco=false) const
Transform space-point within this module from local to global coordinates.
Definition: Module.cc:326
double getZStripWidth() const
Get z-strip width.
Definition: Module.h:155
bool isFlipped() const
Determine if this module is flipped by 180 degrees about z axis within its air gap.
Definition: Module.h:113
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.