Belle II Software development
KLMTimeCollectorModule.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/calibration/KLMTimeAlgorithm.h>
13#include <klm/bklm/geometry/GeometryPar.h>
14#include <klm/dataobjects/KLMHit2d.h>
15#include <klm/dbobjects/KLMChannelStatus.h>
16#include <klm/eklm/geometry/GeometryData.h>
17#include <klm/eklm/geometry/TransformData.h>
18
19/* Basf2 headers */
20#include <calibration/CalibrationCollectorModule.h>
21#include <framework/dataobjects/EventT0.h>
22#include <framework/datastore/StoreArray.h>
23#include <framework/datastore/StoreObjPtr.h>
24#include <mdst/dataobjects/Track.h>
25#include <tracking/dataobjects/ExtHit.h>
26
27/* ROOT headers */
28#include <TH1D.h>
29#include <TH1I.h>
30#include <TH2D.h>
31#include <TTree.h>
32
33namespace Belle2 {
43
44 public:
45
50
55
59 void prepare() override;
60
64 void collect() override;
65
69 void finish() override;
70
71 private:
72
77
82
87
91 std::pair<ExtHit*, ExtHit*> matchExt(
92 KLMChannelNumber channelID, std::multimap<unsigned int, ExtHit>&);
93
95 void storeDistDiff(ROOT::Math::XYZVector&);
96
99
102
104 std::string m_inputListName;
105
107 std::multimap<unsigned int, ExtHit> m_vExtHits_RPC;
108
110 std::multimap<unsigned int, ExtHit> m_vExtHits;
111
114
117
120
123
126
129
132
135
136 /* Monitor histograms. */
137
140
143
146
149
152
155
158
161
164
167
168 /* Sum number of digits collected. */
169
172
175
178
180 TTree* m_outTree;
181
182 };
184}
185
Calibration collector module base class.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
EKLM geometry data.
Definition: GeometryData.h:38
Transformation data.
Definition: TransformData.h:35
KLM element numbers.
Collect hit information for KLM time calibration with CAF.
TH1D * m_HeventT0_1
Event T0 distribution after track selection.
StoreObjPtr< EventT0 > m_eventT0
Event T0 array.
const bklm::GeometryPar * m_geoParB
BKLM geometry parameters.
std::pair< ExtHit *, ExtHit * > matchExt(KLMChannelNumber channelID, std::multimap< unsigned int, ExtHit > &)
Match KLM hit and extHit.
std::multimap< unsigned int, ExtHit > m_vExtHits
Map for handle the extHit related to scint.
EKLM::TransformData * m_TransformData
Transformation data.
TTree * m_outTree
Data collection tree.
TH1I * m_HnKLMHit2dOfTrack
Number of KLM hits related to track.
TH2D * m_HflyTimeB
Particle flying time versus detector layers (for BKLM).
TH1I * m_HnumDigit_scint
BKLM scitillator part.
TH1D * m_HpositionZDiff
Difference between global and local position (Z).
DBObjPtr< KLMChannelStatus > m_channelStatus
Channel status.
void collect() override
Event action, collect information for calibration.
virtual ~KLMTimeCollectorModule()
Destructor.
KLMTimeAlgorithm::Event m_Event
Time calibration data event.
void collectRPC(RelationVector< KLMHit2d > &)
Collect hits information for RPC of BKLM.
StoreArray< Track > m_tracks
Global tracks.
bool m_IgnoreBackwardPropagation
Whether to ignore ExtHits with backward propagation.
TH1D * m_HpositionXDiff
Difference between global and local position (X).
std::multimap< unsigned int, ExtHit > m_vExtHits_RPC
Map for handle the extHit related to RPC.
void collectScint(RelationVector< KLMHit2d > &)
Collect hits information for scintillator of BKLM.
void prepare() override
Initializes the module.
const EKLM::GeometryData * m_geoParE
EKLM geometry parameters.
void finish() override
Termination action.
TH1D * m_HeventT0_0
Event T0 distribution before track selection.
void storeDistDiff(ROOT::Math::XYZVector &)
Save position difference between matched kLMHit and ExtHit.
const KLMElementNumbers * m_elementNum
KLM element numbers.
TH1D * m_HpositionYDiff
Difference between global and local position (Y).
bool m_useEvtT0
Use event T0 or not.
void collectScintEnd(const RelationVector< KLMHit2d > &)
Collect hits information for scintillator of EKLM.
TH1D * m_HpositionDiff
Difference between global and local position.
TH2D * m_HflyTimeE
Particle flying time versus detector layers (for EKLM).
std::string m_inputListName
Input partilce list name.
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
Provides BKLM geometry parameters for simulation, reconstruction etc (from Gearbox or DataBase)
Definition: GeometryPar.h:37
uint16_t KLMChannelNumber
Channel number.
Abstract base class for different kinds of events.