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