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
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("t0_uc", &m_Event.t0_uc, "t0_uc/D");
90 m_outTree->Branch("flyTime", &m_Event.flyTime, "flyTime/D");
91 m_outTree->Branch("recTime", &m_Event.recTime, "recTime/D");
92 m_outTree->Branch("dist", &m_Event.dist, "dist/D");
93 m_outTree->Branch("diffDistX", &m_Event.diffDistX, "diffDistX/D");
94 m_outTree->Branch("diffDistY", &m_Event.diffDistY, "diffDistY/D");
95 m_outTree->Branch("diffDistZ", &m_Event.diffDistZ, "diffDistZ/D");
96 m_outTree->Branch("eDep", &m_Event.eDep, "eDep/D");
97 m_outTree->Branch("nPE", &m_Event.nPE, "nPE/D");
98 m_outTree->Branch("channelId", &m_Event.channelId, "channelId/I");
99 m_outTree->Branch("inRPC", &m_Event.inRPC, "inRPC/O");
100 m_outTree->Branch("isFlipped", &m_Event.isFlipped, "isFlipped/O");
101 m_outTree->Branch("isGood", &m_Event.isGood, "isGood/O");
102 m_outTree->Branch("getADCcount", &m_Event.getADCcount, "getADCcount/s");
103
104 m_outTree->Branch("Run", &m_Event.Run, "Run/D");
105 m_outTree->Branch("Event", &m_Event.Events, "Event/D");
106 m_outTree->Branch("nTrack", &m_Event.nTrack, "nTrack/D");
107 m_outTree->Branch("Track_Charge", &m_Event.Track_Charge, "Track_Charge/D");
108
109 registerObject<TTree>("time_calibration_data", m_outTree);
110
111 /* Initialization of histograms. */
112 m_HeventT0_0 = new TH1D("m_HeventT0_0", "collision time before track number request;t0[ns]", 200, -100, 100);
113 m_HeventT0_1 = new TH1D("m_HeventT0_1", "collision time after track number request;t0[ns]", 200, -100, 100);
114 m_HnumTrack = new TH1I("m_HnnumTrack", "Number of Track;nTrack", 30, 0, 30);
115
116 m_HnKLMHit2dOfTrack = new TH1I("m_HnKLMHit2dOfTrack", "Number of KLMHit2d belong to recTrack;num of KLMHit2d", 20, 0, 20);
117
118 m_HpositionDiff = new TH1D("m_HpositionDiff", "Dist between extHit and KLMHit2d;dist", 160, 0, 8);
119 m_HpositionXDiff = new TH1D("m_HpositionXDiff", "DistX between extHit and KLMHit2d;distX", 100, 0, 5);
120 m_HpositionYDiff = new TH1D("m_HpositionYDiff", "DistY between extHit and KLMHit2d;distY", 100, 0, 5);
121 m_HpositionZDiff = new TH1D("m_HpositionZDiff", "DistZ between extHit and KLMHit2d;distZ", 100, 0, 5);
122
123 m_HflyTimeB = new TH2D("m_HflyTimeB", "flyTime;flyTime(wet)/ns;Layer", 40, 0, 20, 20, 0, 20);
124 m_HflyTimeE = new TH2D("m_HflyTimeE", "flyTime;flyTime(wet)/ns;Layer", 40, 0, 20, 20, 0, 20);
125
126 m_HnumDigit_rpc = new TH1I("m_HnumDigit_rpc", "Number of digit per bklmHit1d (RPC);number of digit", 15, 0, 15);
127 m_HnumDigit_scint = new TH1I("m_HnumDigit_scint", "Number of digit per bklmHit1d (scint);number of digit", 15, 0, 15);
128 m_HnumDigit_scint_end = new TH1I("m_HnumDigit_scint_end", "Number of eklmDigit per eklmHit2d (scint);number of eklmDigit", 15, 0,
129 15);
130
131 registerObject<TH1D>("m_HevtT0_0", m_HeventT0_0);
132 registerObject<TH1D>("m_HevtT0_1", m_HeventT0_1);
133 registerObject<TH1I>("m_HnTrack", m_HnumTrack);
135
137 registerObject<TH1D>("m_HposiXDiff", m_HpositionXDiff);
138 registerObject<TH1D>("m_HposiYDiff", m_HpositionYDiff);
139 registerObject<TH1D>("m_HposiZDiff", m_HpositionZDiff);
140
141 registerObject<TH2D>("m_HfTimeB", m_HflyTimeB);
142 registerObject<TH2D>("m_HfTimeE", m_HflyTimeE);
143
144 registerObject<TH1I>("m_HnDigit_rpc", m_HnumDigit_rpc);
145 registerObject<TH1I>("m_HnDigit_scint", m_HnumDigit_scint);
146 registerObject<TH1I>("m_HnDigit_scint_end", m_HnumDigit_scint_end);
147}
148
150{
151 StoreObjPtr<EventMetaData> eventMetaData("EventMetaData", DataStore::c_Event);
152
153 m_Event.t0 = 0.0;
154 /* Require event T0 determined from CDC. */
155 if (m_useEvtT0) {
156 if (!m_eventT0.isValid())
157 return;
158 if (!m_eventT0->hasTemporaryEventT0(Const::EDetector::CDC))
159 return;
160 const auto bestCDCEvtT0C = m_eventT0->getBestCDCTemporaryEventT0();
161 m_Event.t0 = bestCDCEvtT0C->eventT0;
162
163 double Uncertainty = m_eventT0->getEventT0Uncertainty();
164 m_Event.t0_uc = Uncertainty;
165 }
166
167 /* Read event metadata. */
168 int runId = eventMetaData->getRun();
169 int evtId = eventMetaData->getEvent();
170 m_Event.Run = runId;
171 m_Event.Events = evtId;
172
173 getObjectPtr<TH1D>("m_HevtT0_0")->Fill(m_Event.t0);
174
176 unsigned n_track = inputList->getListSize();
177 /* Return if there are no tracks. */
178 if (n_track < 1) {
179 B2DEBUG(20, "No necessary tracks in" << LogVar("run", runId) << LogVar("event", evtId));
180 return;
181 }
182
183 B2DEBUG(20, "debug info for" << LogVar("run", runId) << LogVar("event", evtId) << LogVar("number of rec tracks", n_track));
184
185 getObjectPtr<TH1D>("m_HevtT0_1")->Fill(m_Event.t0);
186 getObjectPtr<TH1I>("m_HnTrack")->Fill(n_track);
187
188 /* Here begins the ext track sequence */
189 B2DEBUG(20, "Collect :: Track loop begins.");
190
191 /* Main loop */
192 for (unsigned iT = 0; iT < n_track; ++iT) {
193 // Good track selection
194 const Particle* particle = inputList->getParticle(iT);
195 const Track* track = particle->getTrack();
196
197 m_Event.nTrack = iT;
198 double T_Charge = particle->getCharge();
199 m_Event.Track_Charge = T_Charge;
200
201 // Find data objects related to track
202 RelationVector<ExtHit> extHits = track->getRelationsTo<ExtHit>();
203 RelationVector<KLMHit2d> klmHit2ds = track->getRelationsTo<KLMHit2d>();
204
205 getObjectPtr<TH1I>("m_HnKLMHit2d")->Fill(int(klmHit2ds.size()));
206
207 B2DEBUG(20, "Track" << LogVar("exthits", extHits.size())
208 << LogVar("KLMHit2d", klmHit2ds.size()));
209 if (klmHit2ds.size() < 2)
210 continue;
211
212 // Loop for extroplate hits
213 m_vExtHits.clear();
214 m_vExtHits_RPC.clear();
215 B2DEBUG(20, "Collect :: extHits loop begins.");
216 for (const ExtHit& eHit : extHits) {
217 if (eHit.getStatus() != EXT_EXIT)
218 continue;
219
221 if (eHit.isBackwardPropagated())
222 continue;
223 }
224
225 bool bklmCover = (eHit.getDetectorID() == Const::EDetector::BKLM);
226 bool eklmCover = (eHit.getDetectorID() == Const::EDetector::EKLM);
227 if (!bklmCover && !eklmCover)
228 continue;
229
230 int copyId = eHit.getCopyID();
231 int tFor, tSec, tLay, tPla, tStr;
232 int tSub = -1;
233
234 if (eklmCover) {
236 EKLMElementNumbers::Instance().stripNumberToElementNumbers(copyId, &tFor, &tLay, &tSec, &tPla, &tStr);
237 }
238 if (bklmCover) {
240 BKLMElementNumbers::channelNumberToElementNumbers(copyId, &tFor, &tSec, &tLay, &tPla, &tStr);
241 }
242 if (tSub < 0)
243 continue;
244 B2DEBUG(20, "Collect :: Assign elementNum based on copyId for extHits." << LogVar("Sub from elementNumber",
245 tSub) << LogVar("bklmCover", bklmCover) << LogVar("eklmCover", eklmCover));
246
247 bool crossed = false; // should be only once ?
248 KLMMuidLikelihood* muidLikelihood = track->getRelatedTo<KLMMuidLikelihood>();
249 if (bklmCover) {
250 crossed = muidLikelihood->isExtrapolatedBarrelLayerCrossed(tLay - 1);
251 if (crossed) {
252 if (tLay > 2) {
253 unsigned int tModule = m_elementNum->moduleNumber(tSub, tFor, tSec, tLay);
254 m_vExtHits_RPC.insert(std::pair<unsigned int, ExtHit>(tModule, eHit));
255 } else {
256 unsigned int tChannel = m_elementNum->channelNumber(tSub, tFor, tSec, tLay, tPla, tStr);
257 if (m_channelStatus->getChannelStatus(tChannel) != KLMChannelStatus::c_Normal)
258 continue;
259 m_vExtHits.insert(std::pair<unsigned int, ExtHit>(tChannel, eHit));
260 }
261 }
262 }
263 if (eklmCover) {
264 crossed = muidLikelihood->isExtrapolatedEndcapLayerCrossed(tLay - 1);
265 if (crossed) {
266 unsigned int tChannel = m_elementNum->channelNumber(tSub, tFor, tSec, tLay, tPla, tStr);
267 if (m_channelStatus->getChannelStatus(tChannel) != KLMChannelStatus::c_Normal)
268 continue;
269 m_vExtHits.insert(std::pair<unsigned int, ExtHit>(tChannel, eHit));
270 }
271 }
272 }
273 if (m_vExtHits.size() > 0) {
274 collectScintEnd(klmHit2ds);
275 collectScint(klmHit2ds);
276 }
277 if (m_vExtHits_RPC.size() > 0) {
278 collectRPC(klmHit2ds);
279 }
280 }
281}
282
284{
285 const HepGeom::Transform3D* tr;
286 HepGeom::Point3D<double> hitGlobal_extHit, hitLocal_extHit;
287 double l;
288
289 for (const KLMHit2d& hit2d : klmHit2ds) {
290 if (hit2d.getSubdetector() != KLMElementNumbers::c_EKLM)
291 continue;
292 RelationVector<KLMDigit> digits = hit2d.getRelationsTo<KLMDigit>();
293 unsigned nDigit = digits.size();
294 getObjectPtr<TH1I>("m_HnDigit_scint_end")->Fill(nDigit);
295 if (nDigit != 2)
296 B2FATAL("Wrong number of related KLMDigits.");
297
298 if (!digits[0]->isGood() || !digits[1]->isGood())
299 continue;
300
301 for (KLMDigit& digitHit : digits) {
302 m_Event.channelId = digitHit.getUniqueChannelID();
303 m_Event.inRPC = 0;
304 if (m_channelStatus->getChannelStatus(m_Event.channelId) != KLMChannelStatus::c_Normal)
305 continue;
306
307 std::pair<ExtHit*, ExtHit*> pair_extHits = matchExt(m_Event.channelId, m_vExtHits);
308 if (pair_extHits.first == nullptr || pair_extHits.second == nullptr)
309 continue;
310 ExtHit& entryHit = *(pair_extHits.first);
311 ExtHit& exitHit = *(pair_extHits.second);
312
313 m_Event.flyTime = 0.5 * (entryHit.getTOF() + exitHit.getTOF());
314
315 ROOT::Math::XYZVector positionGlobal_extHit = 0.5 * (entryHit.getPosition() + exitHit.getPosition());
316 ROOT::Math::XYZVector positionGlobal_digit = hit2d.getPosition();
317 ROOT::Math::XYZVector positionGlobal_diff = positionGlobal_extHit - positionGlobal_digit;
318
319 storeDistDiff(positionGlobal_diff);
320
321 l = m_geoParE->getStripLength(digitHit.getStrip()) / CLHEP::mm * Unit::mm;
322 hitGlobal_extHit.setX(positionGlobal_extHit.X() / Unit::mm * CLHEP::mm);
323 hitGlobal_extHit.setY(positionGlobal_extHit.Y() / Unit::mm * CLHEP::mm);
324 hitGlobal_extHit.setZ(positionGlobal_extHit.Z() / Unit::mm * CLHEP::mm);
325 tr = m_TransformData->getStripGlobalToLocal(&digitHit);
326 hitLocal_extHit = (*tr) * hitGlobal_extHit;
327
328 m_Event.dist = 0.5 * l - hitLocal_extHit.x() / CLHEP::mm * Unit::mm;
329 m_Event.recTime = digitHit.getTime();
330 m_Event.eDep = digitHit.getEnergyDeposit();
331 m_Event.nPE = digitHit.getNPhotoelectrons();
332 m_Event.isGood = digitHit.isGood();
333 m_Event.getADCcount = digitHit.getCharge();
334
335 getObjectPtr<TH2D>("m_HfTimeE")->Fill(m_Event.flyTime, digitHit.getLayer());
336 getObjectPtr<TTree>("time_calibration_data")->Fill();
337 }
338 }
339}
340
342{
343 double stripWidtm_HZ, stripWidtm_HPhi;
344 for (KLMHit2d& hit2d : klmHit2ds) {
345 if (hit2d.getSubdetector() != KLMElementNumbers::c_BKLM)
346 continue;
347 if (hit2d.inRPC())
348 continue;
349 if (hit2d.isOutOfTime())
350 continue;
351
352 RelationVector<BKLMHit1d> bklmHit1ds = hit2d.getRelationsTo<BKLMHit1d>();
353 if (bklmHit1ds.size() != 2)
354 continue;
355
356 ROOT::Math::XYZVector positionGlobal_hit2d = hit2d.getPosition();
357 const bklm::Module* corMod = m_geoParB->findModule(hit2d.getSection(), hit2d.getSector(), hit2d.getLayer());
358 stripWidtm_HZ = corMod->getZStripWidth();
359 stripWidtm_HPhi = corMod->getPhiStripWidth();
360
361 for (const BKLMHit1d& hit1d : bklmHit1ds) {
362 RelationVector<KLMDigit> digits = hit1d.getRelationsTo<KLMDigit>();
363 getObjectPtr<TH1I>("m_HnDigit_scint")->Fill(digits.size());
364 if (digits.size() > 3)
365 continue;
366
367 for (const KLMDigit& digitHit : digits) {
368 if (digitHit.inRPC() || !digitHit.isGood())
369 continue;
370
371 unsigned int channelId_digit = digitHit.getUniqueChannelID();
372 if (m_channelStatus->getChannelStatus(channelId_digit) != KLMChannelStatus::c_Normal)
373 continue;
374
375 std::pair<ExtHit*, ExtHit*> pair_extHits = matchExt(channelId_digit, m_vExtHits);
376 if (pair_extHits.first == nullptr || pair_extHits.second == nullptr)
377 continue;
378 ExtHit& entryHit = *(pair_extHits.first);
379 ExtHit& exitHit = *(pair_extHits.second);
380
381 m_Event.inRPC = digitHit.inRPC();
382 m_Event.flyTime = 0.5 * (entryHit.getTOF() + exitHit.getTOF());
383
384 ROOT::Math::XYZVector positionGlobal_extHit = 0.5 * (entryHit.getPosition() + exitHit.getPosition());
385 ROOT::Math::XYZVector positionGlobal_diff = positionGlobal_extHit - positionGlobal_hit2d;
386 B2DEBUG(20, LogVar("Distance between digit and hit2d", positionGlobal_diff.R()));
387
388 storeDistDiff(positionGlobal_diff);
389 const CLHEP::Hep3Vector positionLocal_extHit = corMod->globalToLocal(CLHEP::Hep3Vector(
390 positionGlobal_extHit.X(), positionGlobal_extHit.Y(), positionGlobal_extHit.Z()), true);
391 const CLHEP::Hep3Vector positionLocal_hit2d = corMod->globalToLocal(CLHEP::Hep3Vector(
392 positionGlobal_hit2d.X(), positionGlobal_hit2d.Y(), positionGlobal_hit2d.Z()), true);
393 const CLHEP::Hep3Vector diffLocal = positionLocal_extHit - positionLocal_hit2d;
394 if (fabs(diffLocal.z()) > stripWidtm_HZ || fabs(diffLocal.y()) > stripWidtm_HPhi)
395 continue;
396
397 double propaLength = corMod->getPropagationDistance(positionLocal_extHit, digitHit.getStrip(), digitHit.isPhiReadout());
398
399 m_Event.isFlipped = corMod->isFlipped();
400 m_Event.recTime = digitHit.getTime();
401 m_Event.dist = propaLength;
402 m_Event.eDep = digitHit.getEnergyDeposit();
403 m_Event.nPE = digitHit.getNPhotoelectrons();
404 m_Event.channelId = channelId_digit;
405 m_Event.isGood = digitHit.isGood();
406 m_Event.getADCcount = digitHit.getCharge();
407
408 getObjectPtr<TH2D>("m_HfTimeB")->Fill(m_Event.flyTime, digitHit.getLayer());
409 getObjectPtr<TTree>("time_calibration_data")->Fill();
410 }
411 }
412 }
413}
414
416{
417 for (KLMHit2d& hit2d : klmHit2ds) {
418 if (hit2d.getSubdetector() != KLMElementNumbers::c_BKLM)
419 continue;
420 if (!hit2d.inRPC())
421 continue;
422 if (hit2d.isOutOfTime())
423 continue;
424 RelationVector<BKLMHit1d> bklmHit1ds = hit2d.getRelationsTo<BKLMHit1d>();
425 if (bklmHit1ds.size() != 2)
426 continue;
427
428 ROOT::Math::XYZVector positionGlobal_hit2d = hit2d.getPosition();
429 const bklm::Module* corMod = m_geoParB->findModule(hit2d.getSection(), hit2d.getSector(), hit2d.getLayer());
430 double stripWidtm_HZ = corMod->getZStripWidth();
431 double stripWidtm_HPhi = corMod->getPhiStripWidth();
432
433 for (const BKLMHit1d& hit1d : bklmHit1ds) {
434 RelationVector<KLMDigit> digits = hit1d.getRelationsTo<KLMDigit>();
435 getObjectPtr<TH1I>("m_HnDigit_rpc")->Fill(digits.size());
436 if (digits.size() > 5)
437 continue;
438
439 for (const KLMDigit& digitHit : digits) {
440 unsigned int channelId_digit = digitHit.getUniqueChannelID();
441 m_Event.inRPC = digitHit.inRPC();
442 if (m_channelStatus->getChannelStatus(channelId_digit) != KLMChannelStatus::c_Normal)
443 continue;
444 unsigned int tModule = m_elementNum->moduleNumber(digitHit.getSubdetector(), digitHit.getSection(), digitHit.getSector(),
445 digitHit.getLayer());
446
447 std::pair<ExtHit*, ExtHit*> pair_extHits = matchExt(tModule, m_vExtHits_RPC);
448 if (pair_extHits.first == nullptr || pair_extHits.second == nullptr)
449 continue;
450 ExtHit& entryHit = *(pair_extHits.first);
451 ExtHit& exitHit = *(pair_extHits.second);
452
453 m_Event.flyTime = 0.5 * (entryHit.getTOF() + exitHit.getTOF());
454 ROOT::Math::XYZVector positionGlobal_extHit = 0.5 * (entryHit.getPosition() + exitHit.getPosition());
455 ROOT::Math::XYZVector positionGlobal_diff = positionGlobal_extHit - positionGlobal_hit2d;
456 B2DEBUG(20, LogVar("Distance between digit and hit2d", positionGlobal_diff.R()));
457
458 storeDistDiff(positionGlobal_diff);
459 const CLHEP::Hep3Vector positionLocal_extHit = corMod->globalToLocal(CLHEP::Hep3Vector(
460 positionGlobal_extHit.X(), positionGlobal_extHit.Y(), positionGlobal_extHit.Z()), true);
461 const CLHEP::Hep3Vector positionLocal_hit2d = corMod->globalToLocal(CLHEP::Hep3Vector(
462 positionGlobal_hit2d.X(), positionGlobal_hit2d.Y(), positionGlobal_hit2d.Z()), true);
463 const CLHEP::Hep3Vector diffLocal = positionLocal_extHit - positionLocal_hit2d;
464 if (fabs(diffLocal.z()) > stripWidtm_HZ || fabs(diffLocal.y()) > stripWidtm_HPhi)
465 continue;
466
467 const CLHEP::Hep3Vector propaLengthV3 = corMod->getPropagationDistance(positionLocal_extHit);
468 double propaLength = propaLengthV3[2 - int(hit1d.isPhiReadout())];
469
470 m_Event.isFlipped = corMod->isFlipped();
471 m_Event.recTime = digitHit.getTime();
472 m_Event.dist = propaLength;
473 m_Event.eDep = digitHit.getEnergyDeposit();
474 m_Event.nPE = digitHit.getNPhotoelectrons();
475 m_Event.channelId = channelId_digit;
476 m_Event.isGood = 0;
477 m_Event.getADCcount = 0;
478
479 getObjectPtr<TH2D>("m_HfTimeB")->Fill(m_Event.flyTime, digitHit.getLayer());
480 getObjectPtr<TTree>("time_calibration_data")->Fill();
481 }
482 }
483 }
484}
485
486std::pair<ExtHit*, ExtHit*> KLMTimeCollectorModule::matchExt(
487 KLMChannelNumber channelID, std::multimap<unsigned int, ExtHit>& v_ExtHits)
488{
489 ExtHit* entryHit = nullptr;
490 ExtHit* exitHit = nullptr;
491 std::multimap<unsigned int, ExtHit>::iterator it, itlow, itup;
492 itlow = v_ExtHits.lower_bound(channelID);
493 itup = v_ExtHits.upper_bound(channelID);
494 for (it = itlow; it != itup; ++it) {
495 if (entryHit == nullptr) {
496 entryHit = &(it->second);
497 } else if ((it->second).getTOF() < entryHit->getTOF()) {
498 entryHit = &(it->second);
499 }
500 if (exitHit == nullptr) {
501 exitHit = &(it->second);
502 } else if ((it->second).getTOF() > exitHit->getTOF()) {
503 exitHit = &(it->second);
504 }
505 }
506 std::pair<ExtHit*, ExtHit*> p_extHits = std::make_pair(entryHit, exitHit);
507 return p_extHits;
508}
509
510void KLMTimeCollectorModule::storeDistDiff(ROOT::Math::XYZVector& pDiff)
511{
512 double diffM = pDiff.R();
513 double diffX = pDiff.X();
514 double diffY = pDiff.Y();
515 double diffZ = pDiff.Z();
516 getObjectPtr<TH1D>("m_HposiDiff")->Fill(diffM);
517 getObjectPtr<TH1D>("m_HposiXDiff")->Fill(diffX);
518 getObjectPtr<TH1D>("m_HposiYDiff")->Fill(diffY);
519 getObjectPtr<TH1D>("m_HposiZDiff")->Fill(diffZ);
520 m_Event.diffDistX = diffX;
521 m_Event.diffDistY = diffY;
522 m_Event.diffDistZ = diffZ;
523}
524
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
void registerObject(std::string name, T *obj)
Register object with a name, takes ownership, do not access the pointer beyond prepare()
CalibrationCollectorModule()
Constructor. Sets the default prefix for calibration dataobjects.
T * getObjectPtr(std::string name)
Calls the CalibObjManager to get the requested stored collector data.
@ 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.
Transformation data.
@ c_None
Displacement is not used.
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
static const KLMElementNumbers & Instance()
Instantiation.
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.
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:76
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
static GeometryPar * instance(void)
Static method to get a reference to the singleton GeometryPar instance.
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
uint16_t KLMChannelNumber
Channel number.
Abstract base class for different kinds of events.