Belle II Software development
KLMTimeCollectorModule.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/KLMTimeCollector/KLMTimeCollectorModule.h>
11
12/* KLM headers. */
13#include <klm/dataobjects/bklm/BKLMHit1d.h>
14#include <klm/dataobjects/KLMDigit.h>
15#include <klm/dataobjects/KLMElementNumbers.h>
16#include <klm/dataobjects/KLMMuidLikelihood.h>
17
18/* Basf2 headers. */
19#include <analysis/dataobjects/Particle.h>
20#include <analysis/dataobjects/ParticleList.h>
21#include <framework/dataobjects/EventT0.h>
22#include <framework/gearbox/Unit.h>
23#include <framework/logging/Logger.h>
24#include <mdst/dataobjects/Track.h>
25#include <tracking/dataobjects/ExtHit.h>
26
27/* ROOT headers. */
28#include <TTree.h>
29
30/* C++ headers. */
31#include <utility>
32
33using namespace Belle2;
34using namespace Belle2::bklm;
35using namespace Belle2::EKLM;
36
37REG_MODULE(KLMTimeCollector);
38
41 m_geoParB(nullptr),
42 m_geoParE(nullptr),
43 m_TransformData(nullptr),
44 m_HeventT0_0(nullptr),
45 m_HeventT0_1(nullptr),
46 m_HnumTrack(nullptr),
47 m_HnKLMHit2dOfTrack(nullptr),
48 m_HpositionDiff{nullptr},
49 m_HpositionXDiff(nullptr),
50 m_HpositionYDiff(nullptr),
51 m_HpositionZDiff(nullptr),
52 m_HflyTimeB(nullptr),
53 m_HflyTimeE(nullptr),
54 m_HnumDigit_scint_end(nullptr),
55 m_HnumDigit_scint(nullptr),
56 m_HnumDigit_rpc(nullptr),
57 m_outTree(nullptr)
58{
59 setDescription("Module for KLM time calibration (data collection).");
61
62 addParam("inputParticleList", m_inputListName, "input particle list.", std::string("mu+:cali"));
63 addParam("useEventT0", m_useEvtT0, "Use event T0 or not.", true);
64 addParam("IgnoreBackwardPropagation", m_IgnoreBackwardPropagation,
65 "Whether to ignore ExtHits with backward propagation.", false);
66
68}
69
71{
72}
73
75{
76 /* Initialize geometry. */
80
81 /* Require input dataobjects. */
82 if (m_useEvtT0)
83 m_eventT0.isRequired("EventT0");
84 m_tracks.isRequired();
85
86 /* Initialization of output tree. */
87 m_outTree = new TTree("time_calibration_data", "");
88 m_outTree->Branch("t0", &m_Event.t0, "t0/D");
89 m_outTree->Branch("flyTime", &m_Event.flyTime, "flyTime/D");
90 m_outTree->Branch("recTime", &m_Event.recTime, "recTime/D");
91 m_outTree->Branch("dist", &m_Event.dist, "dist/D");
92 m_outTree->Branch("diffDistX", &m_Event.diffDistX, "diffDistX/D");
93 m_outTree->Branch("diffDistY", &m_Event.diffDistY, "diffDistY/D");
94 m_outTree->Branch("diffDistZ", &m_Event.diffDistZ, "diffDistZ/D");
95 m_outTree->Branch("eDep", &m_Event.eDep, "eDep/D");
96 m_outTree->Branch("nPE", &m_Event.nPE, "nPE/D");
97 m_outTree->Branch("channelId", &m_Event.channelId, "channelId/I");
98 m_outTree->Branch("inRPC", &m_Event.inRPC, "inRPC/O");
99 m_outTree->Branch("isFlipped", &m_Event.isFlipped, "isFlipped/O");
100
101 registerObject<TTree>("time_calibration_data", m_outTree);
102
103 /* Initialization of histograms. */
104 m_HeventT0_0 = new TH1D("m_HeventT0_0", "collision time before track number request;t0[ns]", 200, -100, 100);
105 m_HeventT0_1 = new TH1D("m_HeventT0_1", "collision time after track number request;t0[ns]", 200, -100, 100);
106 m_HnumTrack = new TH1I("m_HnnumTrack", "Number of Track;nTrack", 30, 0, 30);
107
108 m_HnKLMHit2dOfTrack = new TH1I("m_HnKLMHit2dOfTrack", "Number of KLMHit2d belong to recTrack;num of KLMHit2d", 20, 0, 20);
109
110 m_HpositionDiff = new TH1D("m_HpositionDiff", "Dist between extHit and KLMHit2d;dist", 160, 0, 8);
111 m_HpositionXDiff = new TH1D("m_HpositionXDiff", "DistX between extHit and KLMHit2d;distX", 100, 0, 5);
112 m_HpositionYDiff = new TH1D("m_HpositionYDiff", "DistY between extHit and KLMHit2d;distY", 100, 0, 5);
113 m_HpositionZDiff = new TH1D("m_HpositionZDiff", "DistZ between extHit and KLMHit2d;distZ", 100, 0, 5);
114
115 m_HflyTimeB = new TH2D("m_HflyTimeB", "flyTime;flyTime(wet)/ns;Layer", 40, 0, 20, 20, 0, 20);
116 m_HflyTimeE = new TH2D("m_HflyTimeE", "flyTime;flyTime(wet)/ns;Layer", 40, 0, 20, 20, 0, 20);
117
118 m_HnumDigit_rpc = new TH1I("m_HnumDigit_rpc", "Number of digit per bklmHit1d (RPC);number of digit", 15, 0, 15);
119 m_HnumDigit_scint = new TH1I("m_HnumDigit_scint", "Number of digit per bklmHit1d (scint);number of digit", 15, 0, 15);
120 m_HnumDigit_scint_end = new TH1I("m_HnumDigit_scint_end", "Number of eklmDigit per eklmHit2d (scint);number of eklmDigit", 15, 0,
121 15);
122
123 registerObject<TH1D>("m_HevtT0_0", m_HeventT0_0);
124 registerObject<TH1D>("m_HevtT0_1", m_HeventT0_1);
125 registerObject<TH1I>("m_HnTrack", m_HnumTrack);
126 registerObject<TH1I>("m_HnKLMHit2d", m_HnKLMHit2dOfTrack);
127
128 registerObject<TH1D>("m_HposiDiff", m_HpositionDiff);
129 registerObject<TH1D>("m_HposiXDiff", m_HpositionXDiff);
130 registerObject<TH1D>("m_HposiYDiff", m_HpositionYDiff);
131 registerObject<TH1D>("m_HposiZDiff", m_HpositionZDiff);
132
133 registerObject<TH2D>("m_HfTimeB", m_HflyTimeB);
134 registerObject<TH2D>("m_HfTimeE", m_HflyTimeE);
135
136 registerObject<TH1I>("m_HnDigit_rpc", m_HnumDigit_rpc);
137 registerObject<TH1I>("m_HnDigit_scint", m_HnumDigit_scint);
138 registerObject<TH1I>("m_HnDigit_scint_end", m_HnumDigit_scint_end);
139}
140
142{
143 StoreObjPtr<EventMetaData> eventMetaData("EventMetaData", DataStore::c_Event);
144
145 m_Event.t0 = 0.0;
146 /* Require event T0 determined from CDC. */
147 if (m_useEvtT0) {
148 if (!m_eventT0.isValid())
149 return;
150 if (!m_eventT0->hasTemporaryEventT0(Const::EDetector::CDC))
151 return;
152 const auto bestCDCEvtT0C = m_eventT0->getBestCDCTemporaryEventT0();
153 m_Event.t0 = bestCDCEvtT0C->eventT0;
154 }
155
156 /* Read event metadata. */
157 int runId = eventMetaData->getRun();
158 int evtId = eventMetaData->getEvent();
159
160 getObjectPtr<TH1D>("m_HevtT0_0")->Fill(m_Event.t0);
161
163 unsigned n_track = inputList->getListSize();
164 /* Return if there are no tracks. */
165 if (n_track < 1) {
166 B2DEBUG(20, "No necessary tracks in" << LogVar("run", runId) << LogVar("event", evtId));
167 return;
168 }
169
170 B2DEBUG(20, "debug info for" << LogVar("run", runId) << LogVar("event", evtId) << LogVar("number of rec tracks", n_track));
171
172 getObjectPtr<TH1D>("m_HevtT0_1")->Fill(m_Event.t0);
173 getObjectPtr<TH1I>("m_HnTrack")->Fill(n_track);
174
175 /* Here begins the ext track sequence */
176 B2DEBUG(20, "Collect :: Track loop begins.");
177
178 /* Main loop */
179 for (unsigned iT = 0; iT < n_track; ++iT) {
180 // Good track selection
181 const Particle* particle = inputList->getParticle(iT);
182 const Track* track = particle->getTrack();
183
184 // Find data objects related to track
185 RelationVector<ExtHit> extHits = track->getRelationsTo<ExtHit>();
186 RelationVector<KLMHit2d> klmHit2ds = track->getRelationsTo<KLMHit2d>();
187
188 getObjectPtr<TH1I>("m_HnKLMHit2d")->Fill(int(klmHit2ds.size()));
189
190 B2DEBUG(20, "Track" << LogVar("exthits", extHits.size())
191 << LogVar("KLMHit2d", klmHit2ds.size()));
192 if (klmHit2ds.size() < 2)
193 continue;
194
195 // Loop for extroplate hits
196 m_vExtHits.clear();
197 m_vExtHits_RPC.clear();
198 B2DEBUG(20, "Collect :: extHits loop begins.");
199 for (const ExtHit& eHit : extHits) {
200 if (eHit.getStatus() != EXT_EXIT)
201 continue;
202
204 if (eHit.isBackwardPropagated())
205 continue;
206 }
207
208 bool bklmCover = (eHit.getDetectorID() == Const::EDetector::BKLM);
209 bool eklmCover = (eHit.getDetectorID() == Const::EDetector::EKLM);
210 if (!bklmCover && !eklmCover)
211 continue;
212
213 int copyId = eHit.getCopyID();
214 int tFor, tSec, tLay, tPla, tStr;
215 int tSub = -1;
216
217 if (eklmCover) {
219 EKLMElementNumbers::Instance().stripNumberToElementNumbers(copyId, &tFor, &tLay, &tSec, &tPla, &tStr);
220 }
221 if (bklmCover) {
223 BKLMElementNumbers::channelNumberToElementNumbers(copyId, &tFor, &tSec, &tLay, &tPla, &tStr);
224 }
225 if (tSub < 0)
226 continue;
227 B2DEBUG(20, "Collect :: Assign elementNum based on copyId for extHits." << LogVar("Sub from elementNumber",
228 tSub) << LogVar("bklmCover", bklmCover) << LogVar("eklmCover", eklmCover));
229
230 bool crossed = false; // should be only once ?
231 KLMMuidLikelihood* muidLikelihood = track->getRelatedTo<KLMMuidLikelihood>();
232 if (bklmCover) {
233 crossed = muidLikelihood->isExtrapolatedBarrelLayerCrossed(tLay - 1);
234 if (crossed) {
235 if (tLay > 2) {
236 unsigned int tModule = m_elementNum->moduleNumber(tSub, tFor, tSec, tLay);
237 m_vExtHits_RPC.insert(std::pair<unsigned int, ExtHit>(tModule, eHit));
238 } else {
239 unsigned int tChannel = m_elementNum->channelNumber(tSub, tFor, tSec, tLay, tPla, tStr);
240 if (m_channelStatus->getChannelStatus(tChannel) != KLMChannelStatus::c_Normal)
241 continue;
242 m_vExtHits.insert(std::pair<unsigned int, ExtHit>(tChannel, eHit));
243 }
244 }
245 }
246 if (eklmCover) {
247 crossed = muidLikelihood->isExtrapolatedEndcapLayerCrossed(tLay - 1);
248 if (crossed) {
249 unsigned int tChannel = m_elementNum->channelNumber(tSub, tFor, tSec, tLay, tPla, tStr);
250 if (m_channelStatus->getChannelStatus(tChannel) != KLMChannelStatus::c_Normal)
251 continue;
252 m_vExtHits.insert(std::pair<unsigned int, ExtHit>(tChannel, eHit));
253 }
254 }
255 }
256 if (m_vExtHits.size() > 0) {
257 collectScintEnd(klmHit2ds);
258 collectScint(klmHit2ds);
259 }
260 if (m_vExtHits_RPC.size() > 0) {
261 collectRPC(klmHit2ds);
262 }
263 }
264}
265
267{
268 const HepGeom::Transform3D* tr;
269 HepGeom::Point3D<double> hitGlobal_extHit, hitLocal_extHit;
270 double l;
271
272 for (const KLMHit2d& hit2d : klmHit2ds) {
273 if (hit2d.getSubdetector() != KLMElementNumbers::c_EKLM)
274 continue;
275 RelationVector<KLMDigit> digits = hit2d.getRelationsTo<KLMDigit>();
276 unsigned nDigit = digits.size();
277 getObjectPtr<TH1I>("m_HnDigit_scint_end")->Fill(nDigit);
278 if (nDigit != 2)
279 B2FATAL("Wrong number of related KLMDigits.");
280
281 if (!digits[0]->isGood() || !digits[1]->isGood())
282 continue;
283
284 for (KLMDigit& digitHit : digits) {
285 m_Event.channelId = digitHit.getUniqueChannelID();
286 m_Event.inRPC = 0;
288 continue;
289
290 std::pair<ExtHit*, ExtHit*> pair_extHits = matchExt(m_Event.channelId, m_vExtHits);
291 if (pair_extHits.first == nullptr || pair_extHits.second == nullptr)
292 continue;
293 ExtHit& entryHit = *(pair_extHits.first);
294 ExtHit& exitHit = *(pair_extHits.second);
295
296 m_Event.flyTime = 0.5 * (entryHit.getTOF() + exitHit.getTOF());
297
298 ROOT::Math::XYZVector positionGlobal_extHit = 0.5 * (entryHit.getPosition() + exitHit.getPosition());
299 ROOT::Math::XYZVector positionGlobal_digit = hit2d.getPosition();
300 ROOT::Math::XYZVector positionGlobal_diff = positionGlobal_extHit - positionGlobal_digit;
301
302 storeDistDiff(positionGlobal_diff);
303
304 l = m_geoParE->getStripLength(digitHit.getStrip()) / CLHEP::mm * Unit::mm;
305 hitGlobal_extHit.setX(positionGlobal_extHit.X() / Unit::mm * CLHEP::mm);
306 hitGlobal_extHit.setY(positionGlobal_extHit.Y() / Unit::mm * CLHEP::mm);
307 hitGlobal_extHit.setZ(positionGlobal_extHit.Z() / Unit::mm * CLHEP::mm);
308 tr = m_TransformData->getStripGlobalToLocal(&digitHit);
309 hitLocal_extHit = (*tr) * hitGlobal_extHit;
310
311 m_Event.dist = 0.5 * l - hitLocal_extHit.x() / CLHEP::mm * Unit::mm;
312 m_Event.recTime = digitHit.getTime();
313 m_Event.eDep = digitHit.getEnergyDeposit();
314 m_Event.nPE = digitHit.getNPhotoelectrons();
315
316 getObjectPtr<TH2D>("m_HfTimeE")->Fill(m_Event.flyTime, digitHit.getLayer());
317 getObjectPtr<TTree>("time_calibration_data")->Fill();
318 }
319 }
320}
321
323{
324 double stripWidtm_HZ, stripWidtm_HPhi;
325 for (KLMHit2d& hit2d : klmHit2ds) {
326 if (hit2d.getSubdetector() != KLMElementNumbers::c_BKLM)
327 continue;
328 if (hit2d.inRPC())
329 continue;
330 if (hit2d.isOutOfTime())
331 continue;
332
333 RelationVector<BKLMHit1d> bklmHit1ds = hit2d.getRelationsTo<BKLMHit1d>();
334 if (bklmHit1ds.size() != 2)
335 continue;
336
337 ROOT::Math::XYZVector positionGlobal_hit2d = hit2d.getPosition();
338 const bklm::Module* corMod = m_geoParB->findModule(hit2d.getSection(), hit2d.getSector(), hit2d.getLayer());
339 stripWidtm_HZ = corMod->getZStripWidth();
340 stripWidtm_HPhi = corMod->getPhiStripWidth();
341
342 for (const BKLMHit1d& hit1d : bklmHit1ds) {
343 RelationVector<KLMDigit> digits = hit1d.getRelationsTo<KLMDigit>();
344 getObjectPtr<TH1I>("m_HnDigit_scint")->Fill(digits.size());
345 if (digits.size() > 3)
346 continue;
347
348 for (const KLMDigit& digitHit : digits) {
349 if (digitHit.inRPC() || !digitHit.isGood())
350 continue;
351
352 unsigned int channelId_digit = digitHit.getUniqueChannelID();
353 if (m_channelStatus->getChannelStatus(channelId_digit) != KLMChannelStatus::c_Normal)
354 continue;
355
356 std::pair<ExtHit*, ExtHit*> pair_extHits = matchExt(channelId_digit, m_vExtHits);
357 if (pair_extHits.first == nullptr || pair_extHits.second == nullptr)
358 continue;
359 ExtHit& entryHit = *(pair_extHits.first);
360 ExtHit& exitHit = *(pair_extHits.second);
361
362 m_Event.inRPC = digitHit.inRPC();
363 m_Event.flyTime = 0.5 * (entryHit.getTOF() + exitHit.getTOF());
364
365 ROOT::Math::XYZVector positionGlobal_extHit = 0.5 * (entryHit.getPosition() + exitHit.getPosition());
366 ROOT::Math::XYZVector positionGlobal_diff = positionGlobal_extHit - positionGlobal_hit2d;
367 B2DEBUG(20, LogVar("Distance between digit and hit2d", positionGlobal_diff.R()));
368
369 storeDistDiff(positionGlobal_diff);
370 const CLHEP::Hep3Vector positionLocal_extHit = corMod->globalToLocal(CLHEP::Hep3Vector(
371 positionGlobal_extHit.X(), positionGlobal_extHit.Y(), positionGlobal_extHit.Z()), true);
372 const CLHEP::Hep3Vector positionLocal_hit2d = corMod->globalToLocal(CLHEP::Hep3Vector(
373 positionGlobal_hit2d.X(), positionGlobal_hit2d.Y(), positionGlobal_hit2d.Z()), true);
374 const CLHEP::Hep3Vector diffLocal = positionLocal_extHit - positionLocal_hit2d;
375 if (fabs(diffLocal.z()) > stripWidtm_HZ || fabs(diffLocal.y()) > stripWidtm_HPhi)
376 continue;
377
378 double propaLength = corMod->getPropagationDistance(positionLocal_extHit, digitHit.getStrip(), digitHit.isPhiReadout());
379
380 m_Event.isFlipped = corMod->isFlipped();
381 m_Event.recTime = digitHit.getTime();
382 m_Event.dist = propaLength;
383 m_Event.eDep = digitHit.getEnergyDeposit();
384 m_Event.nPE = digitHit.getNPhotoelectrons();
385 m_Event.channelId = channelId_digit;
386
387 getObjectPtr<TH2D>("m_HfTimeB")->Fill(m_Event.flyTime, digitHit.getLayer());
388 getObjectPtr<TTree>("time_calibration_data")->Fill();
389 }
390 }
391 }
392}
393
395{
396 for (KLMHit2d& hit2d : klmHit2ds) {
397 if (hit2d.getSubdetector() != KLMElementNumbers::c_BKLM)
398 continue;
399 if (!hit2d.inRPC())
400 continue;
401 if (hit2d.isOutOfTime())
402 continue;
403 RelationVector<BKLMHit1d> bklmHit1ds = hit2d.getRelationsTo<BKLMHit1d>();
404 if (bklmHit1ds.size() != 2)
405 continue;
406
407 ROOT::Math::XYZVector positionGlobal_hit2d = hit2d.getPosition();
408 const bklm::Module* corMod = m_geoParB->findModule(hit2d.getSection(), hit2d.getSector(), hit2d.getLayer());
409 double stripWidtm_HZ = corMod->getZStripWidth();
410 double stripWidtm_HPhi = corMod->getPhiStripWidth();
411
412 for (const BKLMHit1d& hit1d : bklmHit1ds) {
413 RelationVector<KLMDigit> digits = hit1d.getRelationsTo<KLMDigit>();
414 getObjectPtr<TH1I>("m_HnDigit_rpc")->Fill(digits.size());
415 if (digits.size() > 5)
416 continue;
417
418 for (const KLMDigit& digitHit : digits) {
419 unsigned int channelId_digit = digitHit.getUniqueChannelID();
420 m_Event.inRPC = digitHit.inRPC();
421 if (m_channelStatus->getChannelStatus(channelId_digit) != KLMChannelStatus::c_Normal)
422 continue;
423 unsigned int tModule = m_elementNum->moduleNumber(digitHit.getSubdetector(), digitHit.getSection(), digitHit.getSector(),
424 digitHit.getLayer());
425
426 std::pair<ExtHit*, ExtHit*> pair_extHits = matchExt(tModule, m_vExtHits_RPC);
427 if (pair_extHits.first == nullptr || pair_extHits.second == nullptr)
428 continue;
429 ExtHit& entryHit = *(pair_extHits.first);
430 ExtHit& exitHit = *(pair_extHits.second);
431
432 m_Event.flyTime = 0.5 * (entryHit.getTOF() + exitHit.getTOF());
433 ROOT::Math::XYZVector positionGlobal_extHit = 0.5 * (entryHit.getPosition() + exitHit.getPosition());
434 ROOT::Math::XYZVector positionGlobal_diff = positionGlobal_extHit - positionGlobal_hit2d;
435 B2DEBUG(20, LogVar("Distance between digit and hit2d", positionGlobal_diff.R()));
436
437 storeDistDiff(positionGlobal_diff);
438 const CLHEP::Hep3Vector positionLocal_extHit = corMod->globalToLocal(CLHEP::Hep3Vector(
439 positionGlobal_extHit.X(), positionGlobal_extHit.Y(), positionGlobal_extHit.Z()), true);
440 const CLHEP::Hep3Vector positionLocal_hit2d = corMod->globalToLocal(CLHEP::Hep3Vector(
441 positionGlobal_hit2d.X(), positionGlobal_hit2d.Y(), positionGlobal_hit2d.Z()), true);
442 const CLHEP::Hep3Vector diffLocal = positionLocal_extHit - positionLocal_hit2d;
443 if (fabs(diffLocal.z()) > stripWidtm_HZ || fabs(diffLocal.y()) > stripWidtm_HPhi)
444 continue;
445
446 const CLHEP::Hep3Vector propaLengthV3 = corMod->getPropagationDistance(positionLocal_extHit);
447 double propaLength = propaLengthV3[2 - int(hit1d.isPhiReadout())];
448
449 m_Event.isFlipped = corMod->isFlipped();
450 m_Event.recTime = digitHit.getTime();
451 m_Event.dist = propaLength;
452 m_Event.eDep = digitHit.getEnergyDeposit();
453 m_Event.nPE = digitHit.getNPhotoelectrons();
454 m_Event.channelId = channelId_digit;
455
456 getObjectPtr<TH2D>("m_HfTimeB")->Fill(m_Event.flyTime, digitHit.getLayer());
457 getObjectPtr<TTree>("time_calibration_data")->Fill();
458 }
459 }
460 }
461}
462
463std::pair<ExtHit*, ExtHit*> KLMTimeCollectorModule::matchExt(
464 KLMChannelNumber channelID, std::multimap<unsigned int, ExtHit>& v_ExtHits)
465{
466 ExtHit* entryHit = nullptr;
467 ExtHit* exitHit = nullptr;
468 std::multimap<unsigned int, ExtHit>::iterator it, itlow, itup;
469 itlow = v_ExtHits.lower_bound(channelID);
470 itup = v_ExtHits.upper_bound(channelID);
471 for (it = itlow; it != itup; ++it) {
472 if (entryHit == nullptr) {
473 entryHit = &(it->second);
474 } else if ((it->second).getTOF() < entryHit->getTOF()) {
475 entryHit = &(it->second);
476 }
477 if (exitHit == nullptr) {
478 exitHit = &(it->second);
479 } else if ((it->second).getTOF() > exitHit->getTOF()) {
480 exitHit = &(it->second);
481 }
482 }
483 std::pair<ExtHit*, ExtHit*> p_extHits = std::make_pair(entryHit, exitHit);
484 return p_extHits;
485}
486
487void KLMTimeCollectorModule::storeDistDiff(ROOT::Math::XYZVector& pDiff)
488{
489 double diffM = pDiff.R();
490 double diffX = pDiff.X();
491 double diffY = pDiff.Y();
492 double diffZ = pDiff.Z();
493 getObjectPtr<TH1D>("m_HposiDiff")->Fill(diffM);
494 getObjectPtr<TH1D>("m_HposiXDiff")->Fill(diffX);
495 getObjectPtr<TH1D>("m_HposiYDiff")->Fill(diffY);
496 getObjectPtr<TH1D>("m_HposiZDiff")->Fill(diffZ);
497 m_Event.diffDistX = diffX;
498 m_Event.diffDistY = diffY;
499 m_Event.diffDistZ = diffZ;
500}
501
503{
504}
static void channelNumberToElementNumbers(KLMChannelNumber channel, int *section, int *sector, int *layer, int *plane, int *strip)
Get element numbers by channel number.
Store one reconstructed BKLM 1D hit as a ROOT object.
Definition: BKLMHit1d.h:30
Calibration collector module base class.
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Definition: DataStore.h:59
static const EKLMElementNumbers & Instance()
Instantiation.
void stripNumberToElementNumbers(int stripGlobal, int *section, int *layer, int *sector, int *plane, int *strip) const
Get element numbers by strip global number.
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Definition: GeometryData.cc:33
double getStripLength(int strip) const
Get strip length.
Definition: GeometryData.h:71
Transformation data.
Definition: TransformData.h:35
@ c_None
Displacement is not used.
Definition: TransformData.h:43
const HepGeom::Transform3D * getStripGlobalToLocal(KLMDigit *hit) const
Get strip global to local transformation by hit.
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:31
ROOT::Math::XYZVector getPosition() const
Get position of this extrapolation hit.
Definition: ExtHit.h:144
double getTOF() const
Get time of flight from the point of closest approach near the origin to this hit.
Definition: ExtHit.h:136
@ c_Normal
Normally operating channel.
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:29
unsigned int getUniqueChannelID() const override
Get unique channel identifier.
Definition: KLMDigit.cc:83
KLMChannelNumber channelNumber(int subdetector, int section, int sector, int layer, int plane, int strip) const
Get channel number.
static const KLMElementNumbers & Instance()
Instantiation.
KLMModuleNumber moduleNumber(int subdetector, int section, int sector, int layer) const
Get module number.
KLM 2d hit.
Definition: KLMHit2d.h:33
Class to store the likelihoods from KLM with additional information related to the extrapolation.
bool isExtrapolatedEndcapLayerCrossed(int layer) const
Check whether the given EKLM layer is crossed during extrapolation.
bool isExtrapolatedBarrelLayerCrossed(int layer) const
Check whether the given BKLM layer is crossed during extrapolation.
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.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
Class to store reconstructed particles.
Definition: Particle.h:75
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
Class that bundles various TrackFitResults.
Definition: Track.h:25
static const double mm
[millimeters]
Definition: Unit.h:70
const Module * findModule(int section, int sector, int layer) const
Get the pointer to the definition of a module.
Definition: GeometryPar.cc:721
static GeometryPar * instance(void)
Static method to get a reference to the singleton GeometryPar instance.
Definition: GeometryPar.cc:27
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 getPropagationDistance(const CLHEP::Hep3Vector &) const
Convert local coordinates to signal-propagation distance (cm).
Definition: Module.cc:290
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
Class to store variables with their name which were sent to the logging service.
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:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
uint16_t KLMChannelNumber
Channel number.
Abstract base class for different kinds of events.
double flyTime
Particle flying time.
double diffDistX
Global position difference between klmHit2d and ExtHit (X).
double t0
EventT0 for the digit.
int channelId
Unique channel id Barral and endcap merged.
bool inRPC
BKLM RPC flag, used for testing and not necessary.
double eDep
Collect energy eV.
double diffDistY
Global position difference between klmHit2d and ExtHit (Y).
double diffDistZ
Global position difference between klmHit2d and ExtHit (Z).
double dist
Propagation distance from hit to FEE.
double nPE
Number of photon electron.
double recTime
Recosntruction time respected to the trigger time.
bool isFlipped
If phi and z plane flipped, used for testing and not necessary.