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