Belle II Software development
BKLMTrackingModule.h
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#pragma once
10
11/* KLM headers. */
12#include <klm/dataobjects/bklm/BKLMTrack.h>
13#include <klm/dataobjects/KLMHit2d.h>
14#include <klm/bklm/geometry/GeometryPar.h>
15
16/* Basf2 headers. */
17#include <framework/core/Module.h>
18#include <framework/datastore/StoreArray.h>
19#include <tracking/dataobjects/RecoTrack.h>
20
21/* ROOT headers. */
22#include <TEfficiency.h>
23#include <TFile.h>
24#include <TH1F.h>
25#include <TH2F.h>
26
27namespace Belle2 {
34 class BKLMTrackingModule: public Module {
35
36 public:
37
40
43
45 void initialize() override;
46
48 void beginRun() override;
49
51 void event() override;
52
54 void endRun() override;
55
57 void terminate() override;
58
60 bool sameSector(KLMHit2d* hit1, KLMHit2d* hit2);
61
63 bool findClosestRecoTrack(BKLMTrack* bklmTrk, RecoTrack*& closestTrack);
64
65
66 protected:
67
70
73
75 double m_maxAngleRequired = 10;
76
78 double m_maxDistance = 10;
79
81 double m_maxSigma = 5;
82
84 unsigned int m_minHitList = 2;
85
87 unsigned int m_maxHitList = 60;
88
91
93 std::string m_outPath = "bklmEffi.root";
94
95 private:
96
99
101 TFile* m_file = nullptr;
102
104 TH1F* m_total[2][8];
105
107 TH1F* m_pass[2][8];
108
110 TEfficiency* m_effiVsLayer[2][8];
111
113 //TEfficiency* m_effiYX;
114 TH2F* m_effiYX;
115
117 //TEfficiency* m_effiYZ;
118 TH2F* m_effiYZ;
119
121 TH2F* m_passYX;
122
125
127 TH2F* m_passYZ;
128
133
136
139
142
144 void runTracking(int mode, int section, int sector, int layer);
145
147 void generateEffi(int section, int sector, int layer);
148
150 static bool sortByLayer(KLMHit2d* hit1, KLMHit2d* hit2);
151
153 bool isLayerUnderStudy(int section, int iSector, int iLayer, KLMHit2d* hit);
154
156 bool isSectorUnderStudy(int section, int iSector, KLMHit2d* hit);
157
159 double distanceToHit(BKLMTrack* track, KLMHit2d* hit,
160 double& error,
161 double& sigma);
162
164 std::vector<int> m_runNumber;
165
168
170 std::vector<int> m_totalEvents;
171
174
176 std::vector<int> m_totalEventsWithTracks;
177 };
179} // end namespace Belle2
Store one BKLM Track as a ROOT object.
Definition: BKLMTrack.h:35
This module perform straight line track finding and fitting for BKLM.
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
KLM 2d hit.
Definition: KLMHit2d.h:33
Base class for Modules.
Definition: Module.h:72
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Provides BKLM geometry parameters for simulation, reconstruction etc (from Gearbox or DataBase)
Definition: GeometryPar.h:37
Abstract base class for different kinds of events.