Belle II Software development
TrackFinderMCTruthRecoTracksModule.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#include <tracking/modules/trackFinderMCTruth/TrackFinderMCTruthRecoTracksModule.h>
10
11#include <framework/datastore/RelationArray.h>
12#include <framework/datastore/StoreObjPtr.h>
13#include <framework/dataobjects/EventMetaData.h>
14#include <framework/gearbox/Const.h>
15#include <cdc/geometry/CDCGeometryPar.h>
16#include <top/dataobjects/TOPBarHit.h>
17#include <simulation/monopoles/MonopoleConstants.h>
18
19#include <framework/geometry/BFieldManager.h>
20#include <framework/geometry/B2Vector3.h>
21#include <framework/dataobjects/Helix.h>
22
23#include <TRandom.h>
24
25#include <vector>
26#include <string>
27#include <sstream>
28#include <cmath>
29
30
31using namespace Belle2;
32
33//-----------------------------------------------------------------
34// Register the Module
35//-----------------------------------------------------------------
36REG_MODULE(TrackFinderMCTruthRecoTracks);
37
38//-----------------------------------------------------------------
39// Implementation
40//-----------------------------------------------------------------
41
43{
44 //Set module properties
45 setDescription("Uses the MC information to create genfit::TrackCandidates for primary MCParticles and Relations between them. "
46 "Fills the created genfit::TrackCandidates with all information (start values, hit indices) needed for the fitting.");
48
49 // choose which hits to use, all hits assigned to the track candidate will be used in the fit
50 addParam("UsePXDHits",
52 "Set true if PXDHits or PXDClusters should be used",
53 true);
54 addParam("UseSVDHits",
56 "Set true if SVDHits or SVDClusters should be used",
57 true);
58 addParam("UseCDCHits",
60 "Set true if CDCHits should be used",
61 true);
62 addParam("UseOnlyAxialCDCHits",
64 "Set true if only the axial CDCHits should be used",
65 false);
66 addParam("UseNLoops",
68 "Set the number of loops whose hits will be marked as priortiy hits. All other hits "
69 "will be marked as auxiliary and therfore not considered for efficiency computations "
70 "By default, all hits will be priority hits.",
71 INFINITY);
72 addParam("UseOnlyBeforeTOP",
74 "Mark hits as auxiliary after the track left the CDC and touched the TOP detector.",
75 false);
76 addParam("UseReassignedHits",
78 "Include hits reassigned from discarded seconardy daughters in the tracks.",
79 false);
80 addParam("UseSecondCDCHits",
82 "Also includes the CDC 2nd hit information in the MC tracks.",
83 false);
84 addParam("MinPXDHits",
86 "Minimum number of PXD hits needed to allow the created of a track candidate",
87 0);
88 addParam("MinSVDHits",
90 "Minimum number of SVD hits needed to allow the created of a track candidate",
91 0);
92 addParam("MinCDCAxialHits",
94 "Minimum number of CDC hits form an axial wire needed to allow the created of a track candidate",
95 0);
96 addParam("MinCDCStereoHits",
98 "Minimum number of CDC hits form a stereo wire needed to allow the created of a track candidate",
99 0);
100 addParam("AllowFirstCDCSuperLayerOnly",
102 "Allow tracks to pass the stereo hit requirement if they touched only the first (axial) CDC layer",
103 false);
104
105 addParam("MinimalNDF",
107 "Minimum number of total hits needed to allow the creation of a track candidate. "
108 "It is called NDF (number of degrees of freedom) because it counts the dimensionality. "
109 "2D hits are counted as 2",
110 5);
111
112 //choose for which particles a track candidate should be created
113 //this is just an attempt to find out what is the most suitable way to select particles, if you have other/better ideas, communicate it to the tracking group...
114 addParam("WhichParticles",
116 "List of keywords to mark what properties particles must have to get a track candidate. "
117 "If several properties are given all of them must be true: "
118 "\"primary\" particle must come from the generator, "
119 "\"PXD\", \"SVD\", \"CDC\", \"TOP\", \"ARICH\", \"ECL\" or \"KLM\" particle must have hits in the subdetector with that name. "
120 "\"is:X\" where X is a PDG code: particle must have this code. "
121 "\"from:X\" any of the particles's ancestors must have this (X) code",
122 std::vector<std::string>(1, "primary"));
123 addParam("onlyCheckDirectParentPdgCode",
125 "To be used together with WhichParticles to select the ancestor and daughters. "
126 "If true, only check the direct parent to be contained in the list of possible ancestors. "
127 "If false, check all ancestors in in the list of possible ancestors. "
128 "This could be used to e.g. only create MCRecoTracks for slow pions from D* decays instead of creating an MCRecoTrack "
129 "for every pion as long as the D* is one of its ancestors.",
131
132 addParam("EnergyCut",
134 "Track candidates are only created for MCParticles with energy larger than this cut ",
135 double(0.0));
136
137 addParam("Neutrals",
139 "Set true if track candidates should be created also for neutral particles",
140 bool(false));
141
142 addParam("MergeDecayInFlight",
144 "Merge decay in flights that produce a single charged particle to the parent particle",
145 bool(false));
146
147 addParam("SetTimeSeed",
149 "Set true to forward the production time as time seed of the particles to the RecoTrack",
150 true);
151
152 //smearing of MCMomentum
153 addParam("Smearing",
155 "Smearing of MCMomentum/MCVertex prior to storing it in genfit::TrackCandidate (in %). "
156 "A negative value will switch off smearing. This is also the default.",
157 -1.0);
158
159 addParam("SmearingCov",
161 "Covariance matrix used to smear the true pos and mom before passed to track candidate. "
162 "This matrix will also passed to Genfit as the initial covarance matrix. "
163 "If any diagonal value is negative this feature will not be used. "
164 "OFF DIAGONAL ELEMENTS DO NOT HAVE AN EFFECT AT THE MOMENT",
165 std::vector<double>(36, -1.0));
166
167 addParam("SplitAfterDeltaT",
169 "Minimal time delay between two sim hits (in ns) after which MC reco track will be "
170 "split into seperate tracks. If < 0, don't do splitting."
171 "This feature was designed to be used in MC cosmics reconstruction to get two MCRecoTracks"
172 "when track pass through empty SVD region, so that number of MCRecoTracks can be compared with"
173 "number of non merged reco tracks. ",
174 -1.0);
175
176 // names of output containers
177 addParam("RecoTracksStoreArrayName",
179 "Name of store array holding the RecoTracks (output)",
180 std::string(""));
181
182 addParam("TrueHitMustExist",
184 "If set true only cluster hits that have a relation to a TrueHit will be included in the track candidate",
185 false);
186
187 addParam("discardAuxiliaryHits",
189 "If set true hits marked as auxiliary (e.g. hits in higher loops) will not be included in the track candidate.",
191
192 addParam("useCDCSuperLayers",
194 "List of CDC super layers to be used.",
196
197 addParam("useCDCLayers",
199 "List of layers to be used. "
200 "If a layer number is given in 'ignoreLayers', too, this will result on a B2FATAL. "
201 "Either use 'useLayers' and provide a set of layers to be used, "
202 "or use 'ignoreLayers' and provide a list of layers to be ignored, but not both at the same time.",
204
205 addParam("ignoreCDCLayers",
207 "List of layers to be ignored. "
208 "If a layer number is given in 'useLayers', too, this will result on a B2FATAL. "
209 "Either use 'useLayers' and provide a set of layers to be used, "
210 "or use 'ignoreLayers' and provide a list of layers to be ignored, but not both at the same time.",
212
213
214
215}
216
218{
221
222 //output store arrays have to be registered in initialize()
224
226
228
229 // build a bit mask with all properties a MCParticle should have to lead to the creation of a track candidate
231 int aPdgCode = 0;
232 const int nProperties = m_whichParticles.size();
233 for (int i = 0; i not_eq nProperties; ++i) {
234 if (m_whichParticles[i] == "primary") {
236 } else if (m_whichParticles[i] == "PXD") {
238 } else if (m_whichParticles[i] == "SVD") {
240 } else if (m_whichParticles[i] == "CDC") {
242 } else if (m_whichParticles[i] == "TOP") {
244 } else if (m_whichParticles[i] == "ARICH") {
246 } else if (m_whichParticles[i] == "ECL") {
248 } else if (m_whichParticles[i] == "KLM") {
250 } else if (m_whichParticles[i].substr(0, 3) == "is:") {
251 std::string pdgCodeString = m_whichParticles[i].substr(3);
252 std::stringstream(pdgCodeString) >> aPdgCode;
253 B2DEBUG(20, "PDG code added to m_particlePdgCodes " << aPdgCode << " *******");
254 m_particlePdgCodes.push_back(aPdgCode);
255 } else if (m_whichParticles[i].substr(0, 5) == "from:") {
256 std::string pdgCodeString = m_whichParticles[i].substr(5);
257 std::stringstream(pdgCodeString) >> aPdgCode;
258 B2DEBUG(20, "PDG code added to m_fromPdgCodes " << aPdgCode << " *******");
259 m_fromPdgCodes.push_back(aPdgCode);
260 } else {
261 B2FATAL("Invalid values were given to the MCTrackFinder parameter WhichParticles");
262 }
263 }
264
265
266
267 //transfom the smearingCov vector into a TMatrixD
268 //first check if it can be transformed into a 6x6 matrix
269 if (m_smearingCov.size() != 36) {
270 B2FATAL("SmearingCov does not have exactly 36 elements. So 6x6 covariance matrix can be formed from it");
271 }
272 m_initialCov.ResizeTo(6, 6);
273 m_initialCov = TMatrixDSym(6, &m_smearingCov[0]);
274 for (int i = 0; i != 6; ++i) {
275 if (m_initialCov(i, i) < 0.0) {
276 m_initialCov(0, 0) = -1.0; // if first element of matrix is negative this using this matrix will be switched off
277 }
278 }
279
280 if (m_smearing > 0.0 && m_initialCov(0, 0) > 0.0) {
281 B2FATAL("Both relative smearing (Smearing) and using a smearing cov (SmearingCov) is activated but only one of both can be used");
282 }
283 }
284
285 if (not m_param_useCDCSuperLayers.empty()) {
286 for (const uint useSuperLayer : m_param_useCDCSuperLayers) {
287 m_useCDCSuperLayers.at(useSuperLayer) = true;
288 }
289 } else {
290 m_useCDCSuperLayers.fill(true);
291 }
292
293 // Check for common value in the two vectors for using / ignoring layers
294 for (const uint useLayer : m_param_useCDCLayers) {
295 for (const uint ingoreLayer : m_param_ignoreCDCLayers) {
296 if (useLayer == ingoreLayer) {
297 B2FATAL("You chose to use and ignore CDC layer " << useLayer << " at the same time. "
298 "Please decide to either use or to ignore the layer.");
299 }
300 }
301 }
302
303 // fill all layers that should be used
304 if (not m_param_useCDCLayers.empty()) {
305 for (const uint layer : m_param_useCDCLayers) {
306 m_useCDCLayers.at(layer) = true;
307 }
308 } else {
309 m_useCDCLayers.fill(true);
310 }
311 // set layers that should be ignored to false
312 if (not m_param_ignoreCDCLayers.empty()) {
313 for (const uint layer : m_param_ignoreCDCLayers) {
314 m_useCDCLayers.at(layer) = false;
315 }
316 }
317
318}
319
321{
324 m_nRecoTracks = 0;
325}
326
327
329{
330 // Skip in the case there are no MC particles present.
331 if (not m_mcParticlesPresent) {
332 B2DEBUG(20, "Skipping MC Track Finder as there are no MC Particles registered in the DataStore.");
333 return;
334 }
335
336 StoreObjPtr<EventMetaData> eventMetaDataPtr("EventMetaData", DataStore::c_Event);
337 const int eventCounter = eventMetaDataPtr->getEvent();
338 B2DEBUG(20, "******* MCTrackFinderModule processing event number: " << eventCounter << " *******");
339
340 //all the input containers. First: MCParticles
341 const int nMcParticles = m_MCParticles.getEntries();
342 B2DEBUG(20, "MCTrackFinder: total Number of MCParticles: " << nMcParticles);
343
344 //PXD trueHits
345 const int nPXDHits = m_PXDTrueHits.getEntries();
346 B2DEBUG(20, "MCTrackFinder: Number of PXDTrueHits: " << nPXDHits);
347
348 RelationArray mcPartToPXDTrueHits(m_MCParticles, m_PXDTrueHits);
349 const int nMcPartToPXDHits = mcPartToPXDTrueHits.getEntries();
350 B2DEBUG(20, "MCTrackFinder: Number of relations between MCParticles and PXDHits: " << nMcPartToPXDHits);
351
352 //PXD clusters
353 const int nPXDClusters = m_PXDClusters.getEntries();
354 B2DEBUG(20, "MCTrackFinder: Number of PXDClusters: " << nPXDClusters);
355
356 RelationArray pxdClusterToMCParticle(m_PXDClusters, m_MCParticles);
357 const int nPxdClusterToMCPart = pxdClusterToMCParticle.getEntries();
358 B2DEBUG(20, "MCTrackFinder: Number of relations between PXDCluster and MCParticles: " << nPxdClusterToMCPart);
359
360 //SVD truehits
361 const int nSVDHits = m_SVDTrueHits.getEntries();
362 B2DEBUG(20, "MCTrackFinder: Number of SVDDHits: " << nSVDHits);
363
364 RelationArray mcPartToSVDTrueHits(m_MCParticles, m_SVDTrueHits);
365 const int nMcPartToSVDHits = mcPartToSVDTrueHits.getEntries();
366 B2DEBUG(20, "MCTrackFinder: Number of relations between MCParticles and SVDHits: " << nMcPartToSVDHits);
367
368 //SVD clusters
369 const int nSVDClusters = m_SVDClusters.getEntries();
370 B2DEBUG(20, "MCTrackFinder: Number of SVDClusters: " << nSVDClusters);
371
372 RelationArray svdClusterToMCParticle(m_SVDClusters, m_MCParticles);
373 const int nSvdClusterToMCPart = svdClusterToMCParticle.getEntries();
374 B2DEBUG(20, "MCTrackFinder: Number of relations between SVDCluster and MCParticles: " << nSvdClusterToMCPart);
375
376 //CDC
377 const int nCDCHits = m_CDCHits.getEntries();
378 B2DEBUG(20, "MCTrackFinder: Number of CDCHits: " << nCDCHits);
379
380 RelationArray mcPartToCDCHits(m_MCParticles, m_CDCHits);
381 const int nMcPartToCDCHits = mcPartToCDCHits.getEntries();
382 B2DEBUG(20, "MCTrackFinder: Number of relations between MCParticles and CDCHits: " << nMcPartToCDCHits);
383
384 const int nCDCSimHits = m_CDCSimHits.getEntries();
385 B2DEBUG(20, "MCTrackFinder: Number of CDCSimHits: " << nCDCSimHits);
386
387 RelationArray cdcSimHitToHitRel(m_CDCSimHits, m_CDCHits);
388 const int nCdcSimHitToHitRel = cdcSimHitToHitRel.getEntries();
389 B2DEBUG(20, "MCTrackFinder: Number of relations between CDCSimHit and CDCHits: " << nCdcSimHitToHitRel);
390
391 // prepare rejection of CDC/PXD/SVD hits from higher order loops
392 const double Bz = BFieldManager::getField(0, 0, 0).Z() / Unit::T;
393
394 // loop over MCParticles. And check several user selected properties. Make a track candidate only if MCParticle has properties wanted by user options.
395 std::set<int> alreadyConsumedMCParticles;
396 for (int iPart = 0; iPart < nMcParticles; ++iPart) {
397 if (alreadyConsumedMCParticles.count(iPart)) continue;
398 alreadyConsumedMCParticles.insert(iPart);
399
400 MCParticle* aMcParticlePtr = m_MCParticles[iPart];
401 // Ignore particles that didn't propagate significantly, they cannot make tracks.
402 if ((aMcParticlePtr->getDecayVertex() - aMcParticlePtr->getProductionVertex()).R() < 1 * Unit::cm) {
403 B2DEBUG(20, "Particle that did not propagate significantly cannot make track.");
404 continue;
405 }
406
407 //set the property mask for this particle and compare it to the one generated from user input
408 int mcParticleProperties = 0;
409 if (aMcParticlePtr->hasStatus(MCParticle::c_PrimaryParticle)) {
410 mcParticleProperties += 1;
411 }
412 if (aMcParticlePtr->hasSeenInDetector(Const::PXD)) {
413 mcParticleProperties += 2;
414 }
415 if (aMcParticlePtr->hasSeenInDetector(Const::SVD)) {
416 mcParticleProperties += 4;
417 }
418 if (aMcParticlePtr->hasSeenInDetector(Const::CDC)) {
419 mcParticleProperties += 8;
420 }
421 if (aMcParticlePtr->hasSeenInDetector(Const::TOP)) {
422 mcParticleProperties += 16;
423 }
424 if (aMcParticlePtr->hasSeenInDetector(Const::ARICH)) {
425 mcParticleProperties += 32;
426 }
427 if (aMcParticlePtr->hasSeenInDetector(Const::ECL)) {
428 mcParticleProperties += 64;
429 }
430 if (aMcParticlePtr->hasSeenInDetector(Const::KLM)) {
431 mcParticleProperties += 128;
432 }
433 // check all "seen in" properties that the mcparticle should have in one line.
434 if ((mcParticleProperties bitand m_particleProperties) != m_particleProperties) {
435 B2DEBUG(20, "PDG: " << aMcParticlePtr->getPDG() << " | property mask of particle " << mcParticleProperties <<
436 " demanded property mask " << m_particleProperties);
437 continue; //goto next mcParticle, do not make track candidate
438 }
439 //make links only for interesting MCParticles: energy cut
440 if (aMcParticlePtr->getEnergy() < m_energyCut) {
441 B2DEBUG(20, "particle energy too low. MC particle will be skipped");
442 continue; //goto next mcParticle, do not make track candidate
443 }
444
445 //check if particle has the pdg code the user wants to have. If user did not set any pdg code every code is fine for track candidate creation
446
447 const int nPdgCodes = m_particlePdgCodes.size();
448 if (nPdgCodes not_eq 0) {
449 const int currentPdgCode = aMcParticlePtr->getPDG();
450 int nFalsePdgCodes = 0;
451 for (int i = 0; i not_eq nPdgCodes; ++i) {
452 if (m_particlePdgCodes[i] not_eq currentPdgCode) {
453 ++nFalsePdgCodes;
454 }
455 }
456 if (nFalsePdgCodes == nPdgCodes) {
457 B2DEBUG(20, "particle does not have one of the user provided pdg codes and will therefore be skipped");
458 continue; //goto next mcParticle, do not make track candidate
459 }
460 }
461
462
463 //check if particle has an ancestor selected by the user. If user did not set any pdg code every code is fine for track candidate creation
464 const int nFromPdgCodes = m_fromPdgCodes.size();
465 if (nFromPdgCodes not_eq 0) {
466 MCParticle* currentParent = aMcParticlePtr->getMother();
467 int nFalsePdgCodes = 0;
468 int nAncestor = 0;
469 bool foundParent = false;
470 while (currentParent not_eq nullptr) {
471 int currentParentPdgCode = currentParent->getPDG();
472 for (int i = 0; i not_eq nFromPdgCodes; ++i) {
473 if (m_fromPdgCodes[i] not_eq currentParentPdgCode) {
474 ++nFalsePdgCodes;
475 } else {
476 foundParent = true;
477 // if (m_onlyCheckDirectParentPdgCode) {
478 // break;
479 // }
480 }
481 }
482
484 currentParent = nullptr;
485 } else {
486 currentParent = currentParent->getMother();
487 ++nAncestor;
488 }
489 }
490 if (nFalsePdgCodes == (nAncestor * nFromPdgCodes) or not(m_onlyCheckDirectParentPdgCode and foundParent)) {
491 B2DEBUG(20, "particle does not have and ancestor with one of the user provided pdg codes and will therefore be skipped");
492 continue; //goto next mcParticle, do not make track candidate
493 }
494 }
495
496 // Ignore baryons, except for deuteron. The purpose is mainly to
497 // avoid an error message when getCharge() is called below.
498 if (abs(aMcParticlePtr->getPDG()) > 1000000000
499 && abs(aMcParticlePtr->getPDG()) != Const::deuteron.getPDGCode()) {
500 B2DEBUG(20, "Skipped Baryon.");
501 continue; //goto next mcParticle, do not make track candidate
502
503 }
504
505 // ignore neutrals (unless requested) (and unless monopoles)
506 if (!m_neutrals && (aMcParticlePtr->getCharge() == 0 && abs(aMcParticlePtr->getPDG()) != Monopoles::c_monopolePDGCode)) {
507 B2DEBUG(20, "particle does not have the right charge. MC particle will be skipped");
508 continue; //goto next mcParticle, do not make track candidate
509 }
510
511
512 B2DEBUG(20, "Build a track for the MCParticle with index: " << iPart << " (PDG: " << aMcParticlePtr->getPDG() << ")");
513
514 std::vector<const MCParticle*> trackMCParticles;
515 trackMCParticles.push_back(aMcParticlePtr);
516
517 if (m_mergeDecayInFlight and aMcParticlePtr->getCharge() != 0) {
518 std::vector<MCParticle*> daughters = aMcParticlePtr->getDaughters();
519 std::vector<MCParticle*> decayInFlightParticles;
520 for (MCParticle* daughter : daughters) {
521 if (daughter->getSecondaryPhysicsProcess() > 200 and daughter->getCharge() != 0) {
522 decayInFlightParticles.push_back(daughter);
523 }
524 }
525 if (decayInFlightParticles.size() == 1) {
526 trackMCParticles.push_back(decayInFlightParticles.front());
527 alreadyConsumedMCParticles.insert(trackMCParticles.back()->getArrayIndex());
528 }
529 }
530
531 //assign indices of the Hits from all detectors.
532 // entry 0: global event of each hit
533 // entry 1: indices in the respective hit arrays (PXD, SVD, CDC, ...)
534 // entry 2: the category of the hit: either a high priority and should be found or not so important
535 // entry 3: identifier for the detector type
536 typedef std::tuple<double, int, RecoHitInformation::OriginTrackFinder, Const::EDetector> TimeHitIDDetector;
537 std::vector<TimeHitIDDetector> hitsWithTimeAndDetectorInformation;
538
539 int ndf = 0; // count the ndf of one track candidate
540
541 // create a list containing the indices to the PXDHits that belong to one track
542 if (m_usePXDHits) {
543 int hitCounter = 0;
544 for (const MCParticle* trackMCParticle : trackMCParticles) {
545 const RelationVector<PXDCluster>& relatedClusters = trackMCParticle->getRelationsFrom<PXDCluster>();
546
547 for (size_t i = 0; i < relatedClusters.size(); ++i) {
548 bool isReassigned = relatedClusters.weight(i) < 0;
549 if (!m_useReassignedHits && isReassigned) continue; // skip hits from secondary particles
550
551 // currently only priority Hits for PXD
552 auto mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderPriorityHit;
553
554 // Reassigned hits are auxiliary
555 if (isReassigned) {
556 mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit;
557 }
558
559 // Mark higher order hits as auxiliary, if m_useNLoops has been set
560 if (std::isfinite(m_useNLoops) and not isWithinNLoops<PXDCluster, PXDTrueHit>(Bz, relatedClusters.object(i), m_useNLoops)) {
561 mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit;
562 }
563
564 // check if the particle left hits in TOP and add mark hits after that as auxiliary.
565 if (m_useOnlyBeforeTOP and didParticleExitCDC<PXDCluster, PXDTrueHit>(relatedClusters.object(i))) {
566 mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit;
567 }
568
569 // if flag is set discard all auxiliary hits:
570 if (m_discardAuxiliaryHits and mcFinder == RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit) continue;
571
572 const PXDCluster* pxdCluster = relatedClusters.object(i);
573 const RelationVector<PXDTrueHit>& relatedTrueHits = pxdCluster->getRelationsTo<PXDTrueHit>();
574
575 if (relatedTrueHits.size() == 0 and m_enforceTrueHit) {
576 // there is not trueHit! throw away hit because there is no time information for sorting
578 continue;
579 }
580
581 float time = NAN;
582 for (const PXDTrueHit& pxdTrueHit : relatedTrueHits) {
583 // Make sure only a true hit is taken that really comes from the current mcParticle.
584 // This must be carefully checked because several trueHits from different real tracks can be melted into one cluster
585 const RelationVector<MCParticle>& relatedMCParticles = pxdTrueHit.getRelationsFrom<MCParticle>();
586 if (std::find_if(relatedMCParticles.begin(),
587 relatedMCParticles.end(),
588 [trackMCParticle](const MCParticle & mcParticle) {
589 return &mcParticle == trackMCParticle;
590 }) != relatedMCParticles.end()) {
591 time = pxdTrueHit.getGlobalTime();
592 break;
593 }
594 }
595 if (not std::isnan(time)) {
596 hitsWithTimeAndDetectorInformation.emplace_back(time, pxdCluster->getArrayIndex(), mcFinder, Const::PXD);
597 ++hitCounter;
598 ndf += 2;
599 }
600 }
601
602 B2DEBUG(20, " add " << hitCounter << " PXDClusters. " << relatedClusters.size() - hitCounter <<
603 " PXDClusters were not added because they do not have a corresponding PXDTrueHit");
604 } // next trackMCParticle
605
606 if (hitCounter < m_minPXDHits) {
608 continue; //goto next mcParticle, do not make track candidate
609 }
610 } // end if m_usePXDHits
611
612 // create a list containing the indices to the SVDHits that belong to one track
613 if (m_useSVDHits) {
614 int hitCounter = 0;
615 for (const MCParticle* trackMCParticle : trackMCParticles) {
616 const RelationVector<SVDCluster>& relatedClusters = trackMCParticle->getRelationsFrom<SVDCluster>();
617
618 for (size_t i = 0; i < relatedClusters.size(); ++i) {
619 bool isReassigned = relatedClusters.weight(i) < 0;
620 if (!m_useReassignedHits && isReassigned) continue; // skip hits from secondary particles
621
622 auto mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderPriorityHit;
623
624 // Reassigned hits are auxiliary
625 if (isReassigned) {
626 mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit;
627 }
628
629 // Mark higher order hits as auxiliary, if m_useNLoops has been set
630 if (std::isfinite(m_useNLoops) and not isWithinNLoops<SVDCluster, SVDTrueHit>(Bz, relatedClusters.object(i), m_useNLoops)) {
631 mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit;
632 }
633
634 // check if the particle left hits in TOP and add mark hits after that as auxiliary.
635 if (m_useOnlyBeforeTOP and didParticleExitCDC<SVDCluster, SVDTrueHit>(relatedClusters.object(i))) {
636 mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit;
637 }
638
639 // if flag is set discard all auxiliary hits:
640 if (m_discardAuxiliaryHits and mcFinder == RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit) continue;
641
642 const SVDCluster* svdCluster = relatedClusters.object(i);
643 const RelationVector<SVDTrueHit>& relatedTrueHits = svdCluster->getRelationsTo<SVDTrueHit>();
644
645 if (relatedTrueHits.size() == 0 and m_enforceTrueHit) {
646 // there is not trueHit! throw away hit because there is no time information for sorting
648 continue;
649 }
650
651 float time = NAN;
652 for (const SVDTrueHit& svdTrueHit : relatedTrueHits) {
653 // Make sure only a true hit is taken that really comes from the current mcParticle.
654 // This must be carefully checked because several trueHits from different real tracks can be melted into one cluster
655 const RelationVector<MCParticle>& relatedMCParticles = svdTrueHit.getRelationsFrom<MCParticle>();
656 if (std::find_if(relatedMCParticles.begin(),
657 relatedMCParticles.end(),
658 [trackMCParticle](const MCParticle & mcParticle) {
659 return &mcParticle == trackMCParticle;
660 }) != relatedMCParticles.end()) {
661 time = svdTrueHit.getGlobalTime();
662 break;
663 }
664 }
665 if (not std::isnan(time)) {
666 hitsWithTimeAndDetectorInformation.emplace_back(time, svdCluster->getArrayIndex(), mcFinder, Const::SVD);
667 ++hitCounter;
668 ndf += 1;
669 }
670 }
671
672 B2DEBUG(20, " add " << hitCounter << " SVDClusters. " << relatedClusters.size() - hitCounter <<
673 " SVDClusters were not added because they do not have a corresponding SVDTrueHit");
674 }
675
676 if (hitCounter < m_minSVDHits) {
678 continue; //goto next mcParticle, do not make track candidate
679 }
680 } // end if m_useSVDHits
681
682 if (m_useCDCHits) {
683 // create a list containing the indices to the CDCHits that belong to one track
684 int nAxialHits = 0;
685 int nStereoHits = 0;
686 std::array<int, 9> nHitsBySuperLayerId{};
687
688 for (const MCParticle* trackMCParticle : trackMCParticles) {
689 const RelationVector<CDCHit>& relatedHits = trackMCParticle->getRelationsTo<CDCHit>();
690
691 for (size_t i = 0; i < relatedHits.size(); ++i) {
692 bool isReassigned = relatedHits.weight(i) < 0;
693 if (!m_useReassignedHits && isReassigned) continue; // skip hits from secondary particles
694
695 const CDCHit* cdcHit = relatedHits.object(i);
696
697 unsigned short superLayerId = cdcHit->getISuperLayer();
698 if (not m_useCDCSuperLayers[superLayerId]) continue;
699 unsigned short layerID = cdcHit->getICLayer();
700 if (not m_useCDCLayers[layerID]) continue;
701
702 // continue if this is a 2nd CDC hit information and we do not want to use it
703 if (!m_useSecondCDCHits && cdcHit->is2ndHit()) {
704 continue;
705 }
706
707 auto mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderPriorityHit;
708
709 // Reassigned hits are auxiliary
710 if (isReassigned) {
711 mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit;
712 }
713
714 // Mark higher order hits as auxiliary, if m_useNLoops has been set
715 if (std::isfinite(m_useNLoops) and not isWithinNLoops<CDCHit, CDCSimHit>(Bz, cdcHit, m_useNLoops)) {
716 mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit;
717 }
718
719 // check if the particle left hits in TOP and add mark hits after that as auxiliary.
720 if (m_useOnlyBeforeTOP and didParticleExitCDC<CDCHit, CDCSimHit>(cdcHit)) {
721 mcFinder = RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit;
722 }
723
724 // if flag is set discard all auxiliary hits:
725 if (m_discardAuxiliaryHits and mcFinder == RecoHitInformation::OriginTrackFinder::c_MCTrackFinderAuxiliaryHit) continue;
726
727
728 const CDCSimHit* aCDCSimHitPtr = cdcHit->getRelatedFrom<CDCSimHit>();
729 if (not aCDCSimHitPtr) {
730 B2DEBUG(20, " Skipping CDCHit without related CDCSimHit.");
731 continue;
732 }
733 double time = aCDCSimHitPtr->getFlightTime();
734
735 // Here it is hardcoded what superlayer has axial wires and what has stereo wires.
736 // Maybe it would be better if the WireId would know this
737 if (superLayerId % 2 == 0) {
738 hitsWithTimeAndDetectorInformation.emplace_back(time, cdcHit->getArrayIndex(), mcFinder, Const::CDC);
739 ndf += 1;
740 ++nAxialHits;
741 ++nHitsBySuperLayerId[superLayerId];
742 } else {
743 if (not m_useOnlyAxialCDCHits) {
744 hitsWithTimeAndDetectorInformation.emplace_back(time, cdcHit->getArrayIndex(), mcFinder, Const::CDC);
745 ndf += 1;
746 ++nStereoHits;
747 ++nHitsBySuperLayerId[superLayerId];
748 }
749 }
750 }
751 } // next trackMCParticle
752 B2DEBUG(20, " added " << nAxialHits << " axial and " << nStereoHits << " stereo CDCHits");
753 if (nAxialHits < m_minCDCAxialHits) {
754 // Not enough axial hits. Next MCParticle.
756 continue;
757 }
758
759 if (not m_useOnlyAxialCDCHits and (nStereoHits < m_minCDCStereoHits)) {
761 nHitsBySuperLayerId[0] == nAxialHits and
762 (nAxialHits + nStereoHits >= m_minCDCAxialHits + m_minCDCStereoHits)) {
763 // Special rule for low momentum tracks that only touch the first axial superlayer
764 // If there still enough hits combined from the axial and stereo hit limit -> keep the track.
765 } else {
766 // Not enough stereo hits. Next MCParticle.
768 continue;
769 }
770 }
771 } // end if m_useCDCHits
772
773 if (m_initialCov(0, 0) > 0.0) { //using a user set initial cov and corresponding smearing of inital state adds information
774 ndf += 5;
775 }
776 if (ndf < m_minimalNdf) {
778 continue; //goto next mcParticle, do not make track candidate
779 }
780
781 std::sort(hitsWithTimeAndDetectorInformation.begin(), hitsWithTimeAndDetectorInformation.end(),
782 [](const TimeHitIDDetector & rhs, const TimeHitIDDetector & lhs) {
783 return std::get<0>(rhs) < std::get<0>(lhs);
784 });
785
786 // Now create vectors of vectors, which will either contain hits the full vector WithTimeanddetectorInformation
787 // or slices of it, when cutting SplitAfterDeltaT is positive, each slice corresponding to a (sub-) reco track
788 // created by the MC particle.
789 std::vector< std::vector<TimeHitIDDetector> > hitsWithTimeAndDetectorInformationVectors;
790
791 if (m_splitAfterDeltaT < 0.0) { // no splitting, vector will only contain a single hitInformation vector
792 hitsWithTimeAndDetectorInformationVectors.push_back(hitsWithTimeAndDetectorInformation);
793 } else { // split on delta t
794
795 std::vector<TimeHitIDDetector>::size_type splitFromIdx = 0; // whenever splitting subtrack, start slice from this index
796 for (std::vector<TimeHitIDDetector>::size_type i = 1; i != hitsWithTimeAndDetectorInformation.size(); i++) {
797
798 double delta_t = (std::get<0>(hitsWithTimeAndDetectorInformation[i])
799 - std::get<0>(hitsWithTimeAndDetectorInformation[i - 1]));
800
801 if (delta_t > m_splitAfterDeltaT) {
802 // push slice of `hitsWithTimeAndDetectorInformation' between splitFromidx and previous index
803 hitsWithTimeAndDetectorInformationVectors
804 .emplace_back(hitsWithTimeAndDetectorInformation.begin() + splitFromIdx,
805 hitsWithTimeAndDetectorInformation.begin() + i);
806 splitFromIdx = i;
807 }
808 }
809 // add subtrack after last splitting to list of tracks
810 hitsWithTimeAndDetectorInformationVectors
811 .emplace_back(hitsWithTimeAndDetectorInformation.begin() + splitFromIdx,
812 hitsWithTimeAndDetectorInformation.end());
813 }
814
815 //Now create TrackCandidate
816 int counter = m_RecoTracks.getEntries();
817 B2DEBUG(20, "We came pass all filter of the MCPartile and hit properties. TrackCandidate " << counter <<
818 " will be created from the MCParticle with index: " << iPart << " (PDG: " << aMcParticlePtr->getPDG() << ")");
819
820
821
822 //set track parameters from MCParticle information
823 ROOT::Math::XYZVector positionTrue = aMcParticlePtr->getProductionVertex();
824 ROOT::Math::XYZVector momentumTrue = aMcParticlePtr->getMomentum();
825 double timeTrue = aMcParticlePtr->getProductionTime();
826
827 // if no kind of smearing is activated the initial values (seeds) for track fit will be the simulated truth
828 ROOT::Math::XYZVector momentum = momentumTrue;
829 ROOT::Math::XYZVector position = positionTrue;
830 double time = timeTrue;
831 TVectorD stateSeed(6); //this will
832 TMatrixDSym covSeed(6);
833 covSeed.Zero(); // just to be save
834 covSeed(0, 0) = 1; covSeed(1, 1) = 1; covSeed(2, 2) = 2 * 2;
835 covSeed(3, 3) = 0.1 * 0.1; covSeed(4, 4) = 0.1 * 0.1; covSeed(5, 5) = 0.2 * 0.2;
836 //it may have positive effect on the fit not to start with exactly precise true values (or it may be just interesting to study this)
837 //one can smear the starting momentum values with a gaussian
838 //this calculation is always performed, but with the default value of m_smearing = 0 it has no effect on momentum and position (true values are taken)
839
840 if (m_smearing > 0.0) {
841 double smearing = m_smearing / 100.0; //the module parameter m_smearing goes from 0 to 100, smearing should go from 0 to 1
842
843 double smearedX = gRandom->Gaus(positionTrue.x(), smearing * positionTrue.x());
844 double smearedY = gRandom->Gaus(positionTrue.y(), smearing * positionTrue.y());
845 double smearedZ = gRandom->Gaus(positionTrue.z(), smearing * positionTrue.z());
846 position.SetXYZ(smearedX, smearedY, smearedZ);
847 double smearedPX = gRandom->Gaus(momentumTrue.x(), smearing * momentumTrue.x());
848 double smearedPY = gRandom->Gaus(momentumTrue.y(), smearing * momentumTrue.y());
849 double smearedPZ = gRandom->Gaus(momentumTrue.z(), smearing * momentumTrue.z());
850 momentum.SetXYZ(smearedPX, smearedPY, smearedPZ);
851 }
852
853 //Errors for the position/momentum values can also be passed to genfit::TrackCandidate
854 //Default values in Genfit are (1.,1.,1.,), they seem to be not good!!
855 //The best way to set the 'correct' errors has to be investigated....
856 if (m_initialCov(0, 0) > 0.0) { // alternative seamring with according to a covariance matrix
857 double smearedX = gRandom->Gaus(positionTrue.x(), sqrt(m_initialCov(0, 0)));
858 double smearedY = gRandom->Gaus(positionTrue.y(), sqrt(m_initialCov(1, 1)));
859 double smearedZ = gRandom->Gaus(positionTrue.z(), sqrt(m_initialCov(2, 2)));
860 position.SetXYZ(smearedX, smearedY, smearedZ);
861 double smearedPX = gRandom->Gaus(momentumTrue.x(), sqrt(m_initialCov(3, 3)));
862 double smearedPY = gRandom->Gaus(momentumTrue.y(), sqrt(m_initialCov(4, 4)));
863 double smearedPZ = gRandom->Gaus(momentumTrue.z(), sqrt(m_initialCov(5, 5)));
864 momentum.SetXYZ(smearedPX, smearedPY, smearedPZ);
865 covSeed = m_initialCov;
866 }
867
868 // Finally create RecoTracks for MC particle.
869 // Either only one or multiple tracks, if SplitAfterDeltaT is positive
870 for (const auto& hitInformationVector : hitsWithTimeAndDetectorInformationVectors) {
871
872 // TODO: In former times, the track candidate also stored the PDG code!!!
873 short int charge = static_cast<short int>(aMcParticlePtr->getCharge());
874
875 // extrapolate the position and momentum to the point of the first hit on the helix
876 if (hitInformationVector.size() != 0) {
877 // reset the time to the time of the first hit (assumes time > production time)
878 time = std::get<0>(hitInformationVector.at(0));
879 const double deltaT = time - aMcParticlePtr->getProductionTime();
880 const double energy = sqrt(momentum.Mag2() + aMcParticlePtr->get4Vector().M() * aMcParticlePtr->get4Vector().M());
881 const double beta_xy = momentum.Rho() / energy;
882 // calculate arclength in 2D of the track
883 const double arclength2D = beta_xy * Const::speedOfLight * deltaT;
884
885 // get the position and momentum from the helix at the acrlength corresponding to the first hit
886 Belle2::Helix helix = Belle2::Helix(position, momentum, charge, Bz);
887 momentum = helix.getMomentumAtArcLength2D(arclength2D, Bz);
888 position = helix.getPositionAtArcLength2D(arclength2D);
889 }
890
891
892 RecoTrack* newRecoTrack = m_RecoTracks.appendNew(B2Vector3D(position), B2Vector3D(momentum), charge);
893 if (m_setTimeSeed) {
894 newRecoTrack->setTimeSeed(time);
895 }
897
898 //create relation between the track candidates and the mcParticle (redundant to saving the MCId)
899 newRecoTrack->addRelationTo(aMcParticlePtr);
900 B2DEBUG(20, " --- Create relation between genfit::TrackCand " << counter << " and MCParticle " << iPart);
901
902 int hitCounter = 0;
903 for (const TimeHitIDDetector& hitInformation : hitInformationVector) {
904
905 const Const::EDetector& detectorInformation = std::get<3>(hitInformation);
906 const int hitID = std::get<1>(hitInformation);
907 const auto hitOriginMCFinderType = std::get<2>(hitInformation);
908
909
910 if (detectorInformation == Const::CDC) {
911 const CDCHit* cdcHit = m_CDCHits[hitID];
912 const CDCSimHit* aCDCSimHitPtr = cdcHit->getRelatedFrom<CDCSimHit>();
913
914 //now determine the correct sign to resolve the left right ambiguity in the fitter
915 ROOT::Math::XYZVector simHitPos = aCDCSimHitPtr->getPosTrack();
916 ROOT::Math::XYZVector simMom = aCDCSimHitPtr->getMomentum();
917 ROOT::Math::XYZVector simHitPosOnWire = aCDCSimHitPtr->getPosWire();
918
920 const unsigned short isRightHit = cdcGeometry.getNewLeftRightRaw(simHitPosOnWire, simHitPos, simMom);
921
922 if (isRightHit) {
923 newRecoTrack->addCDCHit(cdcHit, hitCounter, RecoHitInformation::RightLeftInformation::c_right, hitOriginMCFinderType);
924 } else {
925 newRecoTrack->addCDCHit(cdcHit, hitCounter, RecoHitInformation::RightLeftInformation::c_left, hitOriginMCFinderType);
926 }
927 B2DEBUG(20, "CDC hit " << hitID << " has reft/right sign " << isRightHit);
928 } else if (detectorInformation == Const::PXD) {
929 const PXDCluster* pxdCluster = m_PXDClusters[hitID];
930 newRecoTrack->addPXDHit(pxdCluster, hitCounter, hitOriginMCFinderType);
931 } else if (detectorInformation == Const::SVD) {
932 const SVDCluster* svdCluster = m_SVDClusters[hitID];
933 newRecoTrack->addSVDHit(svdCluster, hitCounter, hitOriginMCFinderType);
934 }
935 ++hitCounter;
936 } // end loop over hits
937
938 B2DEBUG(20, "New RecoTrack: #PXDHits: " << newRecoTrack->getPXDHitList().size() <<
939 " #SVDHits: " << newRecoTrack->getSVDHitList().size() <<
940 " #CDCHits: " << newRecoTrack->getCDCHitList().size());
941
942 } // end loop over vector
943 }//end loop over MCParticles
944}
945
946
947
948template< class THit, class TSimHit>
949bool TrackFinderMCTruthRecoTracksModule::isWithinNLoops(double Bz, const THit* aHit, double nLoops)
950{
951 // for SVD there are cases with more than one simhit attached
952 const RelationVector<TSimHit>& relatedSimHits = aHit->template getRelationsWith<TSimHit>();
953
954 // take the first best simhit with mcParticle attached
955 const MCParticle* mcParticle = nullptr;
956 const TSimHit* aSimHit = nullptr;
957 for (const auto& thisSimHit : relatedSimHits) {
958 mcParticle = thisSimHit.template getRelated<MCParticle>();
959 aSimHit = &thisSimHit;
960 if (mcParticle) break;
961 }
962 if (not mcParticle or not aSimHit) {
963 return false;
964 }
965
966 // subtract the production time here in order for this classification to also work
967 // for particles produced at times t' > t0
968 const double tof = aSimHit->getGlobalTime() - mcParticle->getProductionTime();
969 const double speed = mcParticle->get4Vector().Beta() * Const::speedOfLight;
970 const float absMom3D = mcParticle->getMomentum().R();
971
972 const double loopLength = 2 * M_PI * absMom3D / (Bz * 0.00299792458);
973 const double loopTOF = loopLength / speed;
974 if (tof > loopTOF * nLoops) {
975 return false;
976 } else {
977 return true;
978 }
979}
980
981template< class THit, class TSimHit>
983{
984 const TSimHit* simHit = ahit->template getRelated<TSimHit>();
985 if (not simHit) return false;
986
987 const MCParticle* mcParticle = simHit->template getRelated<MCParticle>();
988 if (not mcParticle) return false;
989 if (not mcParticle->hasSeenInDetector(Const::TOP)) return false;
990
992 if (topHits.size() == 0) return false;
993
994 // Get hit with the smallest time.
995 auto lessTime = [](const TOPBarHit & lhs, const TOPBarHit & rhs) {
996 return lhs.getTime() < rhs.getTime();
997 };
998 auto itFirstTopHit = std::min_element(topHits.begin(), topHits.end(), lessTime);
999
1000 return simHit->getGlobalTime() > itFirstTopHit->getTime();
1001};
1002
1003
1004
1006{
1007 if (m_notEnoughtHitsCounter != 0) {
1008 B2WARNING(m_notEnoughtHitsCounter << " tracks had not enough hits to have at least " << m_minimalNdf <<
1009 " number of degrees of freedom (NDF). No Track Candidates were created from them so they will not be passed to the track fitter");
1010 }
1011 if (m_noTrueHitCounter != 0) {
1012 B2WARNING(m_noTrueHitCounter <<
1013 " cluster hits did not have a relation to a true hit and were therefore not included in a track candidate");
1014 }
1015 B2INFO("The MCTrackFinder created a total of " << m_nRecoTracks << " track candidates");
1016}
Class containing the result of the unpacker in raw data and the result of the digitizer in simulation...
Definition: CDCHit.h:40
bool is2ndHit() const
Getter for 2nd hit flag.
Definition: CDCHit.h:205
unsigned short getICLayer() const
Getter for iCLayer (0-55).
Definition: CDCHit.h:178
unsigned short getISuperLayer() const
Getter for iSuperLayer.
Definition: CDCHit.h:184
Example Detector.
Definition: CDCSimHit.h:21
double getFlightTime() const
The method to get flight time.
Definition: CDCSimHit.h:184
B2Vector3D getPosWire() const
The method to get position on wire.
Definition: CDCSimHit.h:199
B2Vector3D getPosTrack() const
The method to get position on the track.
Definition: CDCSimHit.h:217
B2Vector3D getMomentum() const
The method to get momentum.
Definition: CDCSimHit.h:193
The Class for CDC Geometry Parameters.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
unsigned short getNewLeftRightRaw(const B2Vector3D &posOnWire, const B2Vector3D &posOnTrack, const B2Vector3D &momentum) const
Returns new left/right_raw.
Helix parameter class.
Definition: Helix.h:48
int getPDGCode() const
PDG code.
Definition: Const.h:473
EDetector
Enum for identifying the detector components (detector and subdetector).
Definition: Const.h:42
static const double speedOfLight
[cm/ns]
Definition: Const.h:695
static const ChargedStable deuteron
deuteron particle
Definition: Const.h:664
@ c_ErrorIfAlreadyRegistered
If the object/array was already registered, produce an error (aborting initialisation).
Definition: DataStore.h:72
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Definition: DataStore.h:59
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
float getEnergy() const
Return particle energy in GeV.
Definition: MCParticle.h:147
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
ROOT::Math::XYZVector getDecayVertex() const
Return decay vertex.
Definition: MCParticle.h:219
bool hasSeenInDetector(Const::DetectorSet set) const
Return if the seen-in flag for a specific subdetector is set or not.
Definition: MCParticle.h:310
std::vector< Belle2::MCParticle * > getDaughters() const
Get vector of all daughter particles, empty vector if none.
Definition: MCParticle.cc:52
bool hasStatus(unsigned short int bitmask) const
Return if specific status bit is set.
Definition: MCParticle.h:129
ROOT::Math::XYZVector getProductionVertex() const
Return production vertex position.
Definition: MCParticle.h:189
float getCharge() const
Return the particle charge defined in TDatabasePDG.
Definition: MCParticle.cc:36
ROOT::Math::PxPyPzEVector get4Vector() const
Return 4Vector of particle.
Definition: MCParticle.h:207
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
float getProductionTime() const
Return production time in ns.
Definition: MCParticle.h:159
ROOT::Math::XYZVector getMomentum() const
Return momentum.
Definition: MCParticle.h:198
Base class for Modules.
Definition: Module.h:72
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
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
Class PXDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: PXDTrueHit.h:31
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
bool addCDCHit(const UsedCDCHit *cdcHit, const unsigned int sortingParameter, RightLeftInformation rightLeftInformation=RightLeftInformation::c_undefinedRightLeftInformation, OriginTrackFinder foundByTrackFinder=OriginTrackFinder::c_undefinedTrackFinder)
Adds a cdc hit with the given information to the reco track.
Definition: RecoTrack.h:243
std::vector< Belle2::RecoTrack::UsedPXDHit * > getPXDHitList() const
Return an unsorted list of pxd hits.
Definition: RecoTrack.h:449
std::vector< Belle2::RecoTrack::UsedSVDHit * > getSVDHitList() const
Return an unsorted list of svd hits.
Definition: RecoTrack.h:452
bool addPXDHit(const UsedPXDHit *pxdHit, const unsigned int sortingParameter, OriginTrackFinder foundByTrackFinder=OriginTrackFinder::c_undefinedTrackFinder)
Adds a pxd hit with the given information to the reco track.
Definition: RecoTrack.h:258
std::vector< Belle2::RecoTrack::UsedCDCHit * > getCDCHitList() const
Return an unsorted list of cdc hits.
Definition: RecoTrack.h:455
void setTimeSeed(const double timeSeed)
Set the time seed. ATTENTION: This is not the fitted time.
Definition: RecoTrack.h:604
static void registerRequiredRelations(StoreArray< RecoTrack > &recoTracks, std::string const &pxdHitsStoreArrayName="", std::string const &svdHitsStoreArrayName="", std::string const &cdcHitsStoreArrayName="", std::string const &bklmHitsStoreArrayName="", std::string const &eklmHitsStoreArrayName="", std::string const &recoHitInformationStoreArrayName="")
Convenience method which registers all relations required to fully use a RecoTrack.
Definition: RecoTrack.cc:53
bool addSVDHit(const UsedSVDHit *svdHit, const unsigned int sortingParameter, OriginTrackFinder foundByTrackFinder=OriginTrackFinder::c_undefinedTrackFinder)
Adds a svd hit with the given information to the reco track.
Definition: RecoTrack.h:272
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:62
int getEntries() const
Get the number of elements.
Class for type safe access to objects that are referred to in relations.
T * object(int index) const
Get object with index.
size_t size() const
Get number of relations.
iterator end()
Return iterator to last entry +1.
iterator begin()
Return iterator to first entry.
float weight(int index) const
Get weight with index.
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
int getArrayIndex() const
Returns this object's array index (in StoreArray), or -1 if not found.
RelationVector< T > getRelationsWith(const std::string &name="", const std::string &namedRelation="") const
Get the relations between this object and another store array.
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:29
Class SVDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: SVDTrueHit.h:33
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
T * appendNew()
Construct a new T object at the end of the array.
Definition: StoreArray.h:246
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
bool registerRelationTo(const StoreArray< TO > &toArray, DataStore::EDurability durability=DataStore::c_Event, DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut, const std::string &namedRelation="") const
Register a relation to the given StoreArray.
Definition: StoreArray.h:140
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
Class to store track parameters of incoming MC particles relation to MCParticle filled in top/simulat...
Definition: TOPBarHit.h:27
double getTime() const
Returns time of impact.
Definition: TOPBarHit.h:128
bool m_useCDCHits
Boolean to select if CDCHits should be used.
bool didParticleExitCDC(const THit *aHit)
helper function which returns true if the current hit has a creation time after the MCParticle left t...
bool m_enforceTrueHit
If set true only cluster hits that have a relation to a TrueHit will be included in the track candida...
std::vector< int > m_fromPdgCodes
if size() is not 0, only for particles having an ancestor (parent or parent of parent etc) with PDG c...
std::vector< int > m_particlePdgCodes
if size() is not 0, only for particles with PDG codes same as in this vector a track candidate will b...
StoreArray< SVDCluster > m_SVDClusters
SVDTrueHits StoreArray.
double m_energyCut
Create track candidates only for MCParticles with energy above this cut.
std::vector< std::string > m_whichParticles
List of keywords to mark what properties particles must have to get a track candidate .
bool m_discardAuxiliaryHits
if true hits marked as auxiliary will not be included in the RecoTrack
void event() override
This method is the core of the module.
int m_minimalNdf
Minimum number of total hits per track to allow track candidate creation.
StoreArray< PXDTrueHit > m_PXDTrueHits
PXDTrueHits StoreArray.
void endRun() override
This method is called if the current run ends.
int m_particleProperties
Internal encoding of m_whichParticles to avoid string comparisons.
bool m_useReassignedHits
Boolean to select the inclusion of hits form discarded secondary daughters.
std::vector< double > m_smearingCov
Covariance matrix used to smear the true pos and mom before passed to track candidate.
bool m_mcParticlesPresent
This flag is set to false if there are no MC Particles in the data store (probably data run?...
bool m_allowFirstCDCSuperLayerOnly
Boolean to allow tracks to pass the stereo hit requirement if they touched only the first (axial) CDC...
bool m_onlyCheckDirectParentPdgCode
To be used together with WhichParticles to select the ancestor and daughters.
bool m_mergeDecayInFlight
Boolean to merge decay in flight chains that involve a single charged particle.
bool m_setTimeSeed
Boolean to forward the production time as seed time.
float m_useNLoops
Number of loops to include in the MC tracks.
double m_splitAfterDeltaT
Minimal time delay between two sim hits (in ns) after which MC reco track will be split into seperate...
bool m_useSecondCDCHits
Also includes the CDC 2nd hit information in the mc tracks.
StoreArray< CDCSimHit > m_CDCSimHits
CDCSimHits StoreArray.
int m_minCDCStereoHits
Minimum number of CDC hits from stereo wires per track to allow track candidate creation.
StoreArray< PXDCluster > m_PXDClusters
PXDClusters StoreArray.
bool m_usePXDHits
Boolean to select if PXDHits should be used.
void beginRun() override
Called when entering a new run.
std::vector< uint > m_param_useCDCLayers
List of layers to be used.
bool m_useOnlyBeforeTOP
Boolean to select if hits after TOP detector are discarded.
std::vector< uint > m_param_ignoreCDCLayers
List of layers to be ignored in tracking.
std::vector< uint > m_param_useCDCSuperLayers
List of super layers to be used.
int m_minCDCAxialHits
Minimum number of CDC hits from axial wires per track to allow track candidate creation.
int m_noTrueHitCounter
will hold number of cluster hits that do not have a corresponding true hit
int m_minSVDHits
Minimum number of SVD hits per track to allow track candidate creation.
int m_notEnoughtHitsCounter
will hold number of tracks that do not have enough hits to form a track candidate (total NDF less tha...
bool m_useOnlyAxialCDCHits
Boolean to select if only axial CDCHits should be used.
StoreArray< RecoTrack > m_RecoTracks
RecoTracks StoreArray.
std::array< bool, 56 > m_useCDCLayers
Bits for the used layers ATTENTION: hardcoded value for number of layers.
double m_smearing
Smearing of MCMomentum and MCVertex in %.
std::array< bool, 9 > m_useCDCSuperLayers
Bits for the used super layers ATTENTION: hardcoded value for number of super layers.
std::string m_recoTracksStoreArrayName
RecoTracks StoreArray name.
int m_minPXDHits
Minimum number of PXD hits per track to allow track candidate creation.
bool m_neutrals
Boolean to mark if track candidates should also be created for neutral particles.
TMatrixDSym m_initialCov
The std::vector m_smearingCov will be translated into this TMatrixD.
StoreArray< SVDTrueHit > m_SVDTrueHits
SVDTrueHits StoreArray.
StoreArray< MCParticle > m_MCParticles
MCParticles StoreArray.
int m_nRecoTracks
will hold the total number of created track candidates
bool isWithinNLoops(double Bz, const THit *aHit, double nLoops)
helper function which returns true if the current hit is within n loops the template give the hit typ...
bool m_useSVDHits
Boolean to select if SVDHits should be used.
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47
static const double T
[tesla]
Definition: Unit.h:120
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
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
static void getField(const double *pos, double *field)
return the magnetic field at a given position.
Definition: BFieldManager.h:91
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
MCParticle * getMother() const
Returns a pointer to the mother particle.
Definition: MCParticle.h:600
Abstract base class for different kinds of events.