Belle II Software development
KLMEventT0EstimatorModule.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 geometry & dataobjects */
12#include <klm/bklm/geometry/GeometryPar.h>
13#include <klm/eklm/geometry/GeometryData.h>
14#include <klm/eklm/geometry/TransformData.h>
15#include <klm/dataobjects/KLMHit2d.h>
16#include <klm/dataobjects/KLMDigit.h>
17#include <klm/dataobjects/bklm/BKLMHit1d.h>
18#include <klm/dataobjects/KLMElementNumbers.h>
19#include <klm/dataobjects/KLMMuidLikelihood.h>
20
21/* KLM conditions */
22#include <klm/dbobjects/KLMChannelStatus.h>
23#include <klm/dbobjects/KLMTimeConstants.h>
24#include <klm/dbobjects/KLMTimeCableDelay.h>
25#include <klm/dbobjects/KLMEventT0HitResolution.h>
26
27/* Tracking / analysis */
28#include <mdst/dataobjects/Track.h>
29#include <tracking/dataobjects/ExtHit.h>
30#include <analysis/dataobjects/Particle.h>
31#include <analysis/dataobjects/ParticleList.h>
32
33/* Framework */
34#include <framework/core/HistoModule.h>
35#include <framework/datastore/StoreObjPtr.h>
36#include <framework/datastore/StoreArray.h>
37#include <framework/database/DBObjPtr.h>
38#include <framework/datastore/RelationVector.h>
39#include <framework/dataobjects/EventT0.h>
40#include <framework/gearbox/Unit.h>
41
42/* ROOT headers. */
43#include <TDirectory.h>
44#include <TH1D.h>
45#include <TH1I.h>
46
47/* C++ headers. */
48#include <cmath>
49#include <vector>
50#include <map>
51#include <string>
52#include <utility>
53
54namespace Belle2 {
59
82 class KLMEventT0EstimatorModule : public HistoModule {
83
84 public:
85 KLMEventT0EstimatorModule();
86 ~KLMEventT0EstimatorModule() override;
87
89 void defineHisto() override;
90
92 void initialize() override;
93
95 void beginRun() override;
96
98 void event() override;
99
101 void endRun() override;
102
104 void terminate() override;
105
106 private:
107 /* Local helpers. */
108
110 using ExtMap = std::multimap<unsigned int, Belle2::ExtHit>;
111
113 using ExtPair = std::pair<Belle2::ExtHit*, Belle2::ExtHit*>;
114
119
128 bool passesADCCut(double charge, int subdetector, int layer, bool inRPC) const;
129
131 void collectExtrapolatedHits(const Track* track, ExtMap& scintMap, ExtMap& rpcMap);
132
134 ExtPair matchExt(unsigned int key, ExtMap& v_ExtHits);
135
145 double getHitSigma(int subdetector, int layer, bool inRPC, int plane = 0) const;
146
152 double& sumW, double& sumWT);
153
158 double& sumW, double& sumWT);
159
167 void accumulateBKLMRPC(RelationVector<KLMHit2d>& klmHit2ds, const ExtMap& rpcMap,
168 double& sumW, double& sumWT);
169
170 /* ---------- Parameters (set via addParam in constructor) ---------- */
171
173 std::string m_MuonListName;
174
177
179 bool m_ignoreBackward{false};
180
182 std::string m_histDirName;
183
185 std::string m_histSubdirUncorr{"uncorrected"};
186
187 /* Geometry and conditions. */
188
191
194
197
200
203
206
207 /* Inputs. */
208
211
214
215 /* Working buffers. */
216
219
222
224 double m_seedT0{0.0};
225
226 /* Monitoring histograms. */
227
229 TH1D* m_hT0Trk_BKLM_Scint{nullptr};
230
232 TH1D* m_hT0Trk_BKLM_RPC{nullptr};
233
235 TH1D* m_hT0Trk_EKLM_Scint{nullptr};
236
239
242
245
247 TH1D* m_hT0Evt_TrkAvg_All{nullptr};
248
251
254
257
260
262 TH1I* m_hFinalSource{nullptr};
263 };
264
266} // namespace Belle2
Class for accessing objects in the database.
Definition DBObjPtr.h:21
EKLM geometry data.
Transformation data.
HistoModule()
Constructor.
Definition HistoModule.h:32
static const KLMElementNumbers & Instance()
Instantiation.
Belle2::EKLM::TransformData * m_transformE
EKLM strip transformation data.
TH1D * m_hT0Evt_TrkAvg_BKLM_RPC_SEM
Per-event T0 track-average for BKLM RPC (SEM) [ns].
TH1D * m_hT0Evt_TrkAvg_EKLM_Scint_SEM
Per-event T0 track-average for EKLM scintillator (SEM) [ns].
TH1D * m_hT0Trk_BKLM_RPC
Per-track T0 for BKLM RPC [ns].
double m_ADCCut_BKLM_Scint_Min
Minimum ADC cut for BKLM scintillator.
TH1D * m_hT0Trk_EKLM_Scint
Per-track T0 for EKLM scintillator [ns].
TH1I * m_hFinalSource
Final EventT0 source selection (7 bins).
std::string m_histDirName
Parent directory inside the ROOT file (HistoManager) for this module.
void accumulateEKLM(const RelationVector< KLMHit2d > &, const ExtMap &, double &sumW, double &sumWT)
Accumulate EKLM scintillator per-digit T0 estimates (weighted).
ExtPair matchExt(unsigned int key, ExtMap &v_ExtHits)
Find earliest (entry) and latest (exit) ExtHits matching a key (channel or module).
const Belle2::EKLM::GeometryData * m_geoParE
EKLM geometry data.
TH1D * m_hT0Trk_BKLM_Scint
Per-track T0 for BKLM scintillator [ns].
bool m_useCDCTemporaryT0
Use CDC temporary EventT0 as a diagnostic seed (not applied to averaging).
void initialize() override
Register inputs/params; get geometry; call REG_HISTOGRAM.
void event() override
Per-event algorithm: collect hits, compute residuals, fill outputs.
std::string m_histSubdirUncorr
Subdirectory name for uncorrected timing histograms.
DBObjPtr< KLMEventT0HitResolution > m_eventT0HitResolution
Per-hit time resolution for EventT0 estimation.
std::multimap< unsigned int, Belle2::ExtHit > ExtMap
Multimap of ExtHit objects keyed by channel or module number.
void endRun() override
Called when the current run ends.
ExtMap m_extRPC
Extrapolated hits keyed by module number (RPC).
void terminate() override
Called at the end of processing.
TH1D * m_hT0Evt_TrkAvg_BKLM_Scint_SEM
Per-event T0 track-average for BKLM scintillator (SEM) [ns].
void collectExtrapolatedHits(const Track *track, ExtMap &scintMap, ExtMap &rpcMap)
Build maps of extrapolated hits for a track (scint: channel key; RPC: module key).
double m_ADCCut_BKLM_Scint_Max
Maximum ADC cut for BKLM scintillator.
DBObjPtr< KLMChannelStatus > m_channelStatus
Channel status (Normal/Dead/etc.).
StoreObjPtr< ParticleList > m_MuonList
Selected muon particle list.
void beginRun() override
Per-run resets if desired (histos remain booked).
StoreArray< Track > m_tracks
Reconstructed tracks.
double m_ADCCut_EKLM_Scint_Max
Maximum ADC cut for EKLM scintillator.
bool passesADCCut(double charge, int subdetector, int layer, bool inRPC) const
Check if a digit passes the ADC charge cut.
double m_ADCCut_EKLM_Scint_Min
Minimum ADC cut for EKLM scintillator.
Belle2::bklm::GeometryPar * m_geoParB
BKLM geometry.
bool m_ignoreBackward
Ignore backward-propagated ExtHits when forming entry/exit pairs.
TH1D * m_hT0Evt_TrkAvg_BKLM_Scint
Per-event T0 track-average for BKLM scintillator (mean) [ns].
std::pair< Belle2::ExtHit *, Belle2::ExtHit * > ExtPair
Pair of entry and exit ExtHit pointers.
double m_seedT0
Optional seed from CDC (for logging only).
std::string m_MuonListName
Input ParticleList (e.g.
double getHitSigma(int subdetector, int layer, bool inRPC, int plane=0) const
Get per-hit sigma for a digit based on detector category.
TH1D * m_hT0Evt_TrkAvg_All
Per-event T0 track-average combined (mean) [ns].
const KLMElementNumbers * m_elementNum
Element numbering helpers.
void accumulateBKLMRPC(RelationVector< KLMHit2d > &klmHit2ds, const ExtMap &rpcMap, double &sumW, double &sumWT)
Accumulate BKLM RPC per-digit T0 estimates (weighted, both readout directions).
TH1D * m_hT0Evt_TrkAvg_All_SEM
Per-event T0 track-average combined (SEM) [ns].
TH1D * m_hT0Evt_TrkAvg_BKLM_RPC
Per-event T0 track-average for BKLM RPC (mean) [ns].
ExtMap m_extScint
Extrapolated hits keyed by channel number (scintillator).
TH1D * m_hT0Evt_TrkAvg_EKLM_Scint
Per-event T0 track-average for EKLM scintillator (mean) [ns].
void accumulateBKLMScint(RelationVector< KLMHit2d > &, const ExtMap &, double &sumW, double &sumWT)
Accumulate BKLM scintillator per-digit T0 estimates (weighted).
void defineHisto() override
Definition of histograms (called once by HistoManager).
Class for type safe access to objects that are referred to in relations.
Accessor to arrays stored in the data store.
Definition StoreArray.h:113
Type-safe access to single objects in the data store.
Definition StoreObjPtr.h:96
Class that bundles various TrackFitResults.
Definition Track.h:25
Provides BKLM geometry parameters for simulation, reconstruction etc (from Gearbox or DataBase)
Definition GeometryPar.h:37
Abstract base class for different kinds of events.