Belle II Software  release-05-02-19
BKLMTrackingModule.h
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 #pragma once
12 
13 /* KLM headers. */
14 #include <klm/dataobjects/bklm/BKLMHit2d.h>
15 #include <klm/dataobjects/bklm/BKLMTrack.h>
16 #include <klm/bklm/geometry/GeometryPar.h>
17 
18 /* Belle 2 headers. */
19 #include <framework/core/Module.h>
20 #include <framework/datastore/StoreArray.h>
21 #include <tracking/dataobjects/RecoTrack.h>
22 
23 /* ROOT headers. */
24 #include <TEfficiency.h>
25 #include <TFile.h>
26 #include <TH1F.h>
27 #include <TH2F.h>
28 
29 namespace Belle2 {
35 
39  class BKLMTrackingModule: public Module {
40 
41  public:
42 
45 
48 
50  virtual void initialize() override;
51 
53  virtual void beginRun() override;
54 
56  virtual void event() override;
57 
59  virtual void endRun() override;
60 
62  virtual void terminate() override;
63 
65  bool sameSector(BKLMHit2d* hit1, BKLMHit2d* hit2);
66 
68  bool findClosestRecoTrack(BKLMTrack* bklmTrk, RecoTrack*& closestTrack);
69 
70 
71  protected:
72 
74  bool m_studyEffi;
75 
77  bool m_MatchToRecoTrack;
78 
80  double m_maxAngleRequired = 10;
81 
83  bool m_globalFit;
84 
86  std::string m_outPath = "bklmEffi.root";
87 
88  private:
89 
92 
94  TFile* m_file = nullptr;
95 
97  TH1F* m_total[2][8];
98 
100  TH1F* m_pass[2][8];
101 
103  TEfficiency* m_effiVsLayer[2][8];
104 
106  //TEfficiency* m_effiYX;
107  TH2F* m_effiYX;
108 
110  //TEfficiency* m_effiYZ;
111  TH2F* m_effiYZ;
112 
114  TH2F* m_passYX;
115 
117  TH2F* m_totalYX;
118 
120  TH2F* m_passYZ;
121 
123  TH2F* m_totalYZ;
126 
129 
132 
135 
137  void runTracking(int mode, int section, int sector, int layer);
138 
140  void generateEffi(int section, int sector, int layer);
141 
143  static bool sortByLayer(BKLMHit2d* hit1, BKLMHit2d* hit2);
144 
146  bool isLayerUnderStudy(int section, int iSector, int iLayer, BKLMHit2d* hit);
147 
149  bool isSectorUnderStudy(int section, int iSector, BKLMHit2d* hit);
150 
152  double distanceToHit(BKLMTrack* track, BKLMHit2d* hit,
153  double& error,
154  double& sigma);
155 
157  std::vector<int> m_runNumber;
158 
160  int m_runTotalEvents;
161 
163  std::vector<int> m_totalEvents;
164 
167 
169  std::vector<int> m_totalEventsWithTracks;
170  };
172 } // end namespace Belle2
Belle2::BKLMTrackingModule::isLayerUnderStudy
bool isLayerUnderStudy(int section, int iSector, int iLayer, BKLMHit2d *hit)
judge whether the current layer is understudy
Definition: BKLMTrackingModule.cc:480
Belle2::BKLMTrackingModule::m_effiYZ
TH2F * m_effiYZ
Efficieny at global position Y vs Z.
Definition: BKLMTrackingModule.h:119
Belle2::BKLMTrackingModule::m_runTotalEventsWithTracks
int m_runTotalEventsWithTracks
total number of processed events in the run with at lease one BKLMTrack
Definition: BKLMTrackingModule.h:174
Belle2::BKLMTrackingModule::distanceToHit
double distanceToHit(BKLMTrack *track, BKLMHit2d *hit, double &error, double &sigma)
calculate distance from track to hit
Definition: BKLMTrackingModule.cc:492
Belle2::BKLMTrackingModule::m_pass
TH1F * m_pass[2][8]
Numerator of each layer.
Definition: BKLMTrackingModule.h:108
Belle2::BKLMTrackingModule::beginRun
virtual void beginRun() override
begin run stuff
Definition: BKLMTrackingModule.cc:108
Belle2::BKLMTrackingModule::m_runNumber
std::vector< int > m_runNumber
run number
Definition: BKLMTrackingModule.h:165
Belle2::BKLMTrackingModule::BKLMTrackingModule
BKLMTrackingModule()
Constructor.
Definition: BKLMTrackingModule.cc:30
Belle2::BKLMTrackingModule::m_globalFit
bool m_globalFit
do the BKLMTrack fitting in global system (multi-sectors track) or local system (sector by sector)
Definition: BKLMTrackingModule.h:91
Belle2::BKLMTrackingModule::m_GeoPar
bklm::GeometryPar * m_GeoPar
bklm GeometryPar
Definition: BKLMTrackingModule.h:99
Belle2::BKLMTrackingModule::m_effiYX
TH2F * m_effiYX
Efficieny at global position Y vs X.
Definition: BKLMTrackingModule.h:115
Belle2::BKLMTrackingModule::m_totalEventsWithTracks
std::vector< int > m_totalEventsWithTracks
total number of processed events with at least one BKLMTrack
Definition: BKLMTrackingModule.h:177
Belle2::BKLMTrackingModule::isSectorUnderStudy
bool isSectorUnderStudy(int section, int iSector, BKLMHit2d *hit)
judge whether the hits come from the sctor understudy
Definition: BKLMTrackingModule.cc:486
Belle2::BKLMTrackingModule::sameSector
bool sameSector(BKLMHit2d *hit1, BKLMHit2d *hit2)
Judge if two hits come from the same sector.
Definition: BKLMTrackingModule.cc:300
Belle2::BKLMTrackingModule::m_totalYZ
TH2F * m_totalYZ
total event at global position Y vs Z
Definition: BKLMTrackingModule.h:131
Belle2::BKLMTrackingModule::hits2D
StoreArray< BKLMHit2d > hits2D
BKLMHit2d StoreArray.
Definition: BKLMTrackingModule.h:136
Belle2::BKLMTrackingModule::recoTracks
StoreArray< RecoTrack > recoTracks
RecoTrack StoreArray.
Definition: BKLMTrackingModule.h:139
Belle2::BKLMTrackingModule::m_runTotalEvents
int m_runTotalEvents
total number of processed events in the run
Definition: BKLMTrackingModule.h:168
Belle2::bklm::GeometryPar
Provides BKLM geometry parameters for simulation, reconstruction etc (from Gearbox or DataBase)
Definition: GeometryPar.h:48
Belle2::RecoTrack
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:78
Belle2::BKLMTrackingModule::m_total
TH1F * m_total[2][8]
Denominator of each layer.
Definition: BKLMTrackingModule.h:105
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::BKLMTrackingModule::m_totalEvents
std::vector< int > m_totalEvents
total number of processed events
Definition: BKLMTrackingModule.h:171
Belle2::BKLMTrackingModule::m_file
TFile * m_file
TFile that store efficieny plots.
Definition: BKLMTrackingModule.h:102
Belle2::BKLMTrackingModule::m_maxAngleRequired
double m_maxAngleRequired
angle required between RecoTrack and BKLMTrack, if openangle is larger than m_maxAngleRequired,...
Definition: BKLMTrackingModule.h:88
Belle2::BKLMTrackingModule::runTracking
void runTracking(int mode, int section, int sector, int layer)
run the track finding and fitting
Definition: BKLMTrackingModule.cc:145
Belle2::BKLMTrackingModule::sortByLayer
static bool sortByLayer(BKLMHit2d *hit1, BKLMHit2d *hit2)
my defined sort function using layer number
Definition: BKLMTrackingModule.cc:473
Belle2::BKLMTrackingModule::m_studyEffi
bool m_studyEffi
option for efficieny study mode, in this mode, the layer under study should not be used in tracking
Definition: BKLMTrackingModule.h:82
Belle2::BKLMTrack
Store one BKLM Track as a ROOT object.
Definition: BKLMTrack.h:37
Belle2::BKLMTrackingModule::m_storeTracks
StoreArray< BKLMTrack > m_storeTracks
BKLMTrack StoreArray.
Definition: BKLMTrackingModule.h:133
Belle2::BKLMTrackingModule::generateEffi
void generateEffi(int section, int sector, int layer)
calculate efficiency
Definition: BKLMTrackingModule.cc:371
Belle2::BKLMTrackingModule::initialize
virtual void initialize() override
Initialize at start of job.
Definition: BKLMTrackingModule.cc:64
Belle2::BKLMTrackingModule::recoHitInformation
StoreArray< RecoHitInformation > recoHitInformation
RecoHitInformation StoreArray.
Definition: BKLMTrackingModule.h:142
Belle2::BKLMTrackingModule::event
virtual void event() override
Unpack one event and create digits.
Definition: BKLMTrackingModule.cc:116
Belle2::BKLMTrackingModule::endRun
virtual void endRun() override
end run stuff
Definition: BKLMTrackingModule.cc:244
Belle2::BKLMTrackingModule::~BKLMTrackingModule
virtual ~BKLMTrackingModule()
Destructor.
Definition: BKLMTrackingModule.cc:59
Belle2::BKLMTrackingModule::findClosestRecoTrack
bool findClosestRecoTrack(BKLMTrack *bklmTrk, RecoTrack *&closestTrack)
find the closest RecoTrack, match BKLMTrack to RecoTrack, if the matched RecoTrack is found,...
Definition: BKLMTrackingModule.cc:307
Belle2::StoreArray
Accessor to arrays stored in the data store.
Definition: ECLMatchingPerformanceExpertModule.h:33
Belle2::BKLMTrackingModule::m_outPath
std::string m_outPath
output file name containing efficiencies plots
Definition: BKLMTrackingModule.h:94
Belle2::BKLMTrackingModule::m_effiVsLayer
TEfficiency * m_effiVsLayer[2][8]
Efficieny of each layer.
Definition: BKLMTrackingModule.h:111
Belle2::BKLMTrackingModule::m_passYX
TH2F * m_passYX
passed event at global position Y vs X
Definition: BKLMTrackingModule.h:122
Belle2::BKLMTrackingModule::m_passYZ
TH2F * m_passYZ
passed event at global position Y vs Z
Definition: BKLMTrackingModule.h:128
Belle2::BKLMTrackingModule::m_MatchToRecoTrack
bool m_MatchToRecoTrack
whether match BKLMTrack to RecoTrack
Definition: BKLMTrackingModule.h:85
Belle2::BKLMTrackingModule::m_totalYX
TH2F * m_totalYX
total event at global position Y vs X
Definition: BKLMTrackingModule.h:125
Belle2::BKLMHit2d
Store one BKLM strip hit as a ROOT object.
Definition: BKLMHit2d.h:42
Belle2::BKLMTrackingModule::terminate
virtual void terminate() override
Terminate at the end of job.
Definition: BKLMTrackingModule.cc:250