Belle II Software  release-06-01-15
ECLTrackClusterMatchingModule.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 <ecl/modules/eclTrackClusterMatching/ECLTrackClusterMatchingModule.h>
10 
11 #include <ecl/dbobjects/ECLTrackClusterMatchingParameterizations.h>
12 #include <ecl/dbobjects/ECLTrackClusterMatchingThresholds.h>
13 #include <ecl/utility/utilityFunctions.h>
14 
15 #include <framework/gearbox/Const.h>
16 #include <framework/logging/Logger.h>
17 #include <mdst/dataobjects/HitPatternCDC.h>
18 
19 using namespace std;
20 using namespace Belle2;
21 
22 REG_MODULE(ECLTrackClusterMatching)
23 
25  m_matchingParameterizations("ECLTrackClusterMatchingParameterizations"),
26  m_matchingThresholds("ECLTrackClusterMatchingThresholds")
27 {
28  setDescription("Creates and saves a Relation between Tracks and ECLCluster in the DataStore. It uses the existing Relation between Tracks and ExtHit as well as the Relation between ECLCluster and ExtHit.");
29  setPropertyFlags(c_ParallelProcessingCertified);
30  addParam("useAngularDistanceMatching", m_angularDistanceMatching,
31  "if true use track cluster matching based on angular distance, if false use matching based on entered crystals", bool(true));
32  addParam("useOptimizedMatchingConsistency", m_useOptimizedMatchingConsistency,
33  "set false if you want to set the matching criterion on your own", bool(true));
34  addParam("matchingConsistency", m_matchingConsistency,
35  "the 2D consistency of Delta theta and Delta phi has to exceed this value for a track to be matched to an ECL cluster", 1e-6);
36  addParam("matchingPTThreshold", m_matchingPTThreshold,
37  "tracks with pt greater than this value will exclusively be matched based on angular distance", 0.3);
38  addParam("brlEdgeTheta", m_brlEdgeTheta,
39  "distance of polar angle from gaps where crystal-entering based matching is applied (in rad)", 0.1);
40  addParam("minimalCDCHits", m_minimalCDCHits,
41  "bad VXD-standalone tracks cause (too) low photon efficiency in end caps, temporarily fixed by requiring minimal number of CDC hits",
42  -1);
43  addParam("skipZeroChargeTracks", m_skipZeroChargeTracks,
44  "switch to exclude tracks with zero charge from track-cluster matching", bool(false));
45 }
46 
47 ECLTrackClusterMatchingModule::~ECLTrackClusterMatchingModule()
48 {
49 }
50 
51 void ECLTrackClusterMatchingModule::initialize()
52 {
53  // Check dependencies
54  m_tracks.isRequired();
55  m_eclClusters.isRequired();
56  m_eclShowers.isRequired();
57  m_eclCalDigits.isRequired();
58  m_extHits.isRequired();
59  m_trackFitResults.isRequired();
60 
61  if (m_angularDistanceMatching) {
62  m_tracks.registerRelationTo(m_eclShowers);
63  m_tracks.registerRelationTo(m_eclClusters);
64  m_tracks.registerRelationTo(m_eclShowers, DataStore::c_Event, DataStore::c_WriteOut, "AngularDistance");
65  m_tracks.registerRelationTo(m_eclClusters, DataStore::c_Event, DataStore::c_WriteOut, "AngularDistance");
66 
67  // function to update parameterization functions if the payload changes
68  auto updateParameterizationFunctions = [this]() {
69  const auto& map = m_matchingParameterizations->getRMSParameterizations();
70 
71  f_phiRMSFWDCROSS = map.at("PhiFWDCROSS");
72  f_phiRMSFWDDL = map.at("PhiFWDDL");
73  f_phiRMSFWDNEAR = map.at("PhiFWDNEAR");
74  f_phiRMSBRLCROSS = map.at("PhiBRLCROSS");
75  f_phiRMSBRLDL = map.at("PhiBRLDL");
76  f_phiRMSBRLNEAR = map.at("PhiBRLNEAR");
77  f_phiRMSBWDCROSS = map.at("PhiBWDCROSS");
78  f_phiRMSBWDDL = map.at("PhiBWDDL");
79  f_phiRMSBWDNEAR = map.at("PhiBWDNEAR");
80  f_thetaRMSFWDCROSS = map.at("ThetaFWDCROSS");
81  f_thetaRMSFWDDL = map.at("ThetaFWDDL");
82  f_thetaRMSFWDNEAR = map.at("ThetaFWDNEAR");
83  f_thetaRMSBRLCROSS = map.at("ThetaBRLCROSS");
84  f_thetaRMSBRLDL = map.at("ThetaBRLDL");
85  f_thetaRMSBRLNEAR = map.at("ThetaBRLNEAR");
86  f_thetaRMSBWDCROSS = map.at("ThetaBWDCROSS");
87  f_thetaRMSBWDDL = map.at("ThetaBWDDL");
88  f_thetaRMSBWDNEAR = map.at("ThetaBWDNEAR");
89  };
90 
91  // function to update matching threshold functions if the payload changes
92  auto updateMatchingThresholds = [this]() {
93  m_matchingThresholdValuesFWD = m_matchingThresholds->getFWDMatchingThresholdValues();
94  m_matchingThresholdValuesBRL = m_matchingThresholds->getBRLMatchingThresholdValues();
95  m_matchingThresholdValuesBWD = m_matchingThresholds->getBWDMatchingThresholdValues();
96  };
97 
98  // Update once right away
99  updateParameterizationFunctions();
100  updateMatchingThresholds();
101  // And register to be called every time the payloads change
102  m_matchingParameterizations.addCallback(updateParameterizationFunctions);
103  m_matchingThresholds.addCallback(updateMatchingThresholds);
104  } else {
105  m_tracks.registerRelationTo(m_eclShowers, DataStore::c_Event, DataStore::c_WriteOut, "EnterCrystal");
106  m_tracks.registerRelationTo(m_eclClusters, DataStore::c_Event, DataStore::c_WriteOut, "EnterCrystal");
107  }
108 }
109 
110 void ECLTrackClusterMatchingModule::event()
111 {
112  for (auto& eclCluster : m_eclClusters) {
113  eclCluster.setIsTrack(false);
114  }
115  for (const Track& track : m_tracks) {
116  const TrackFitResult* fitResult = track.getTrackFitResultWithClosestMass(Const::pion);
117  // TEMPORARY FIX: require minimal number of CDC hits, otherwise exclude tracks from track-cluster matching procedure
118  if (!(fitResult->getHitPatternCDC().getNHits() > m_minimalCDCHits)) continue;
119  if (m_skipZeroChargeTracks && fitResult->getChargeSign() == 0) continue;
120  double theta = TMath::ACos(fitResult->getMomentum().CosTheta());
121  double pt = fitResult->getTransverseMomentum();
122  if (!m_angularDistanceMatching || pt < m_matchingPTThreshold || trackTowardsGap(theta)) {
123 
124  // Unique shower ids related to this track
125  set<int> uniqueShowerIds;
126 
127  // Need to make sure that we match one shower at most
128  set<int> uniquehypothesisIds;
129  vector<int> hypothesisIds;
130  vector<double> energies;
131  vector<int> arrayIndexes;
132 
133  // Find extrapolated track hits in the ECL, considering
134  // only hit points where the track enters the crystal
135  // note that more than one crystal belonging to more than one shower
136  // can be found
137  for (const auto& extHit : track.getRelationsTo<ExtHit>()) {
138  if (!isECLEnterHit(extHit)) continue;
139  const int cell = extHit.getCopyID() + 1;
140 
141  // Find ECLCalDigit with same cell ID as ExtHit
142  const auto idigit = find_if(m_eclCalDigits.begin(), m_eclCalDigits.end(),
143  [&](const ECLCalDigit & d) { return d.getCellId() == cell; }
144  );
145  // Couldn't find ECLCalDigit with same cell ID as the ExtHit
146  if (idigit == m_eclCalDigits.end()) continue;
147 
148  // Save all unique shower IDs of the showers related to idigit
149  for (auto& shower : idigit->getRelationsFrom<ECLShower>()) {
150  bool inserted = (uniqueShowerIds.insert(shower.getUniqueId())).second;
151 
152  // If this track <-> shower relation hasn't been set yet, set it for the shower and the ECLCLuster
153  if (!inserted) continue;
154 
155  hypothesisIds.push_back(shower.getHypothesisId());
156  energies.push_back(shower.getEnergy());
157  arrayIndexes.push_back(shower.getArrayIndex());
158  uniquehypothesisIds.insert(shower.getHypothesisId());
159 
160  B2DEBUG(29, shower.getArrayIndex() << " " << shower.getHypothesisId() << " " << shower.getEnergy() << " " <<
161  shower.getConnectedRegionId());
162 
163  } // end loop on shower related to idigit
164  } // end loop on ExtHit
165 
166  // only set the relation for the highest energetic shower per hypothesis
167  for (auto hypothesisId : uniquehypothesisIds) {
168  double highestEnergy = 0.0;
169  int arrayindex = -1;
170 
171  for (unsigned ix = 0; ix < energies.size(); ix++) {
172  if (hypothesisIds[ix] == hypothesisId and energies[ix] > highestEnergy) {
173  highestEnergy = energies[ix];
174  arrayindex = arrayIndexes[ix];
175  }
176  }
177 
178  // if we find a shower, take that one by directly accessing the store array
179  if (arrayindex > -1) {
180  auto shower = m_eclShowers[arrayindex];
181  shower->setIsTrack(true);
182  if (m_angularDistanceMatching) {
183  track.addRelationTo(shower);
184  track.addRelationTo(shower, 1.0, "AngularDistance");
185  } else {
186  track.addRelationTo(shower, 1.0, "EnterCrystal");
187  }
188  B2DEBUG(29, shower->getArrayIndex() << " " << shower->getIsTrack());
189 
190  // there is a 1:1 relation, just set the relation for the corresponding cluster as well
191  ECLCluster* cluster = shower->getRelatedFrom<ECLCluster>();
192  if (cluster != nullptr) {
193  cluster->setIsTrack(true);
194  if (m_angularDistanceMatching) {
195  track.addRelationTo(cluster);
196  track.addRelationTo(cluster, 1.0, "AngularDistance");
197  } else {
198  track.addRelationTo(cluster, 1.0, "EnterCrystal");
199  }
200  }
201  }
202  } // end loop on hypothesis IDs
203  }
204  if (m_angularDistanceMatching) {
205  // tracks should never be matched to more than one cluster
206  if (track.getRelationsTo<ECLCluster>("", "AngularDistance").size() > 0) continue;
207  // never match tracks pointing towards gaps or adjacent part of barrel using angular distance
208  if (trackTowardsGap(theta)) continue;
209  // for low-pt tracks matching based on the angular distance is only applied if track points towards the FWD
210  ECL::DetectorRegion trackDetectorRegion = ECL::getDetectorRegion(theta);
211  if (pt < m_matchingPTThreshold && trackDetectorRegion != ECL::DetectorRegion::FWD) continue;
212 
213  map<int, pair<double, int>> hypothesisIdBestQualityCROSSArrayIndexMap;
214  map<int, pair<double, int>> hypothesisIdBestQualityDLArrayIndexMap;
215  map<int, pair<double, int>> hypothesisIdBestQualityNEARArrayIndexMap;
216  set<int> uniqueHypothesisIds;
217 
218  // Find extrapolated track hits in the ECL, considering only hit points
219  // that either are on the sphere, closest to, or on radial direction of an
220  // ECLCluster.
221  for (const auto& extHit : track.getRelationsTo<ExtHit>()) {
222  if (!isECLHit(extHit)) continue;
223  ECLCluster* eclCluster = extHit.getRelatedFrom<ECLCluster>();
224  if (!eclCluster) continue;
225  ECLShower* eclShower = eclCluster->getRelatedTo<ECLShower>();
226  if (eclShower != nullptr) {
227  // accept only shower from region matching track direction, exception for gaps
228  int eclDetectorRegion = eclShower->getDetectorRegion();
229  if (abs(eclDetectorRegion - trackDetectorRegion) == 1) continue;
230  // never match low-pt tracks with showers in the barrel
231  if (pt < m_matchingPTThreshold && eclDetectorRegion == ECL::DetectorRegion::BRL) continue;
232  double phiHit = extHit.getPosition().Phi();
233  double phiShower = eclShower->getPhi();
234  double deltaPhi = phiHit - phiShower;
235  if (deltaPhi > M_PI) {
236  deltaPhi = deltaPhi - 2 * M_PI;
237  } else if (deltaPhi < -M_PI) {
238  deltaPhi = deltaPhi + 2 * M_PI;
239  }
240  double thetaHit = extHit.getPosition().Theta();
241  double thetaShower = eclShower->getTheta();
242  double deltaTheta = thetaHit - thetaShower;
243  ExtHitStatus extHitStatus = extHit.getStatus();
244  double quality = showerQuality(deltaPhi, deltaTheta, pt, eclDetectorRegion, extHitStatus);
245  int hypothesisId = eclShower->getHypothesisId();
246  bool inserted = (uniqueHypothesisIds.insert(hypothesisId)).second;
247  if (inserted) {
248  hypothesisIdBestQualityCROSSArrayIndexMap.insert(make_pair(hypothesisId, make_pair(0, -1)));
249  hypothesisIdBestQualityDLArrayIndexMap.insert(make_pair(hypothesisId, make_pair(0, -1)));
250  hypothesisIdBestQualityNEARArrayIndexMap.insert(make_pair(hypothesisId, make_pair(0, -1)));
251  }
252  if (extHitStatus == EXT_ECLCROSS) {
253  if (quality > hypothesisIdBestQualityCROSSArrayIndexMap.at(hypothesisId).first) {
254  hypothesisIdBestQualityCROSSArrayIndexMap[hypothesisId] = make_pair(quality, eclShower->getArrayIndex());
255  }
256  } else if (extHitStatus == EXT_ECLDL) {
257  if (quality > hypothesisIdBestQualityDLArrayIndexMap.at(hypothesisId).first) {
258  hypothesisIdBestQualityDLArrayIndexMap[hypothesisId] = make_pair(quality, eclShower->getArrayIndex());
259  }
260  } else {
261  if (quality > hypothesisIdBestQualityNEARArrayIndexMap.at(hypothesisId).first) {
262  hypothesisIdBestQualityNEARArrayIndexMap[hypothesisId] = make_pair(quality, eclShower->getArrayIndex());
263  }
264  }
265  }
266  } // end loop on ExtHits related to Track
267 
268  vector<map<int, pair<double, int>>> hypothesisIdBestQualityArrayIndexMaps;
269  hypothesisIdBestQualityArrayIndexMaps.push_back(hypothesisIdBestQualityCROSSArrayIndexMap);
270  hypothesisIdBestQualityArrayIndexMaps.push_back(hypothesisIdBestQualityDLArrayIndexMap);
271  hypothesisIdBestQualityArrayIndexMaps.push_back(hypothesisIdBestQualityNEARArrayIndexMap);
272  if (m_useOptimizedMatchingConsistency) optimizedPTMatchingConsistency(theta, pt);
273  for (const auto& uniqueHypothesisId : uniqueHypothesisIds) {
274  for (const auto& hypothesisIdBestQualityArrayIndexMap : hypothesisIdBestQualityArrayIndexMaps) {
275  if (hypothesisIdBestQualityArrayIndexMap.at(uniqueHypothesisId).first > m_matchingConsistency
276  && hypothesisIdBestQualityArrayIndexMap.at(uniqueHypothesisId).second > -1) {
277  auto shower = m_eclShowers[hypothesisIdBestQualityArrayIndexMap.at(uniqueHypothesisId).second];
278  shower->setIsTrack(true);
279  track.addRelationTo(shower);
280  track.addRelationTo(shower, 1.0, "AngularDistance");
281  ECLCluster* cluster = shower->getRelatedFrom<ECLCluster>();
282  if (cluster != nullptr) {
283  cluster->setIsTrack(true);
284  track.addRelationTo(cluster);
285  track.addRelationTo(cluster, 1.0, "AngularDistance");
286  }
287  break;
288  }
289  }
290  }
291  }
292  } // end loop on Tracks
293 } // end event loop
294 
295 void ECLTrackClusterMatchingModule::terminate()
296 {
297 }
298 
299 bool ECLTrackClusterMatchingModule::isECLEnterHit(const ExtHit& extHit) const
300 {
301  if ((extHit.getDetectorID() != Const::EDetector::ECL)) return false;
302  if ((extHit.getStatus() != EXT_ENTER)) return false;
303  if (extHit.getCopyID() == -1) return false;
304  else return true;
305 }
306 
307 bool ECLTrackClusterMatchingModule::isECLHit(const ExtHit& extHit) const
308 {
309  if ((extHit.getDetectorID() != Const::EDetector::ECL)) return false;
310  ExtHitStatus extHitStatus = extHit.getStatus();
311  if (extHitStatus == EXT_ECLCROSS || extHitStatus == EXT_ECLDL || extHitStatus == EXT_ECLNEAR) return true;
312  else return false;
313 }
314 
315 double ECLTrackClusterMatchingModule::showerQuality(double deltaPhi, double deltaTheta, double pt,
316  int eclDetectorRegion, int hitStatus) const
317 {
318  double phi_consistency = phiConsistency(deltaPhi, pt, eclDetectorRegion, hitStatus);
319  double theta_consistency = thetaConsistency(deltaTheta, pt, eclDetectorRegion, hitStatus);
320  return phi_consistency * theta_consistency * (1 - log(phi_consistency * theta_consistency));
321 }
322 
323 double ECLTrackClusterMatchingModule::phiConsistency(double deltaPhi, double pt, int eclDetectorRegion, int hitStatus) const
324 {
325  double phi_RMS;
326  if (eclDetectorRegion == ECL::DetectorRegion::FWD || eclDetectorRegion == ECL::DetectorRegion::FWDGap) {
327  if (hitStatus == EXT_ECLCROSS) {
328  phi_RMS = f_phiRMSFWDCROSS.Eval(pt);
329  } else if (hitStatus == EXT_ECLDL) {
330  phi_RMS = f_phiRMSFWDDL.Eval(pt);
331  } else {
332  phi_RMS = f_phiRMSFWDNEAR.Eval(pt);
333  }
334  } else if (eclDetectorRegion == ECL::DetectorRegion::BRL) {
335  if (hitStatus == EXT_ECLCROSS) {
336  phi_RMS = f_phiRMSBRLCROSS.Eval(pt);
337  } else if (hitStatus == EXT_ECLDL) {
338  phi_RMS = f_phiRMSBRLDL.Eval(pt);
339  } else {
340  phi_RMS = f_phiRMSBRLNEAR.Eval(pt);
341  }
342  } else if (eclDetectorRegion == ECL::DetectorRegion::BWD || eclDetectorRegion == ECL::DetectorRegion::BWDGap) {
343  if (hitStatus == EXT_ECLCROSS) {
344  phi_RMS = f_phiRMSBWDCROSS.Eval(pt);
345  } else if (hitStatus == EXT_ECLDL) {
346  phi_RMS = f_phiRMSBWDDL.Eval(pt);
347  } else {
348  phi_RMS = f_phiRMSBWDNEAR.Eval(pt);
349  }
350  } else { /* ECL cluster below acceptance */
351  return 0;
352  }
353  return erfc(abs(deltaPhi) / phi_RMS);
354 }
355 
356 double ECLTrackClusterMatchingModule::thetaConsistency(double deltaTheta, double pt, int eclDetectorRegion, int hitStatus) const
357 {
358  double theta_RMS;
359  if (eclDetectorRegion == ECL::DetectorRegion::FWD || eclDetectorRegion == ECL::DetectorRegion::FWDGap) {
360  if (hitStatus == EXT_ECLCROSS) {
361  theta_RMS = f_thetaRMSFWDCROSS.Eval(pt);
362  } else if (hitStatus == EXT_ECLDL) {
363  theta_RMS = f_thetaRMSFWDDL.Eval(pt);
364  } else {
365  theta_RMS = f_thetaRMSFWDNEAR.Eval(pt);
366  }
367  } else if (eclDetectorRegion == ECL::DetectorRegion::BRL) {
368  if (hitStatus == EXT_ECLCROSS) {
369  theta_RMS = f_thetaRMSBRLCROSS.Eval(pt);
370  } else if (hitStatus == EXT_ECLDL) {
371  theta_RMS = f_thetaRMSBRLDL.Eval(pt);
372  } else {
373  theta_RMS = f_thetaRMSBRLNEAR.Eval(pt);
374  }
375  } else if (eclDetectorRegion == ECL::DetectorRegion::BWD || eclDetectorRegion == ECL::DetectorRegion::BWDGap) {
376  if (hitStatus == EXT_ECLCROSS) {
377  theta_RMS = f_thetaRMSBWDCROSS.Eval(pt);
378  } else if (hitStatus == EXT_ECLDL) {
379  theta_RMS = f_thetaRMSBWDDL.Eval(pt);
380  } else {
381  theta_RMS = f_thetaRMSBWDNEAR.Eval(pt);
382  }
383  } else { /* ECL cluster below acceptance */
384  return 0;
385  }
386  return erfc(abs(deltaTheta) / theta_RMS);
387 }
388 
389 bool ECLTrackClusterMatchingModule::trackTowardsGap(double theta) const
390 {
391  if (ECL::getDetectorRegion(theta) == ECL::DetectorRegion::BRL) {
392  if (ECL::getDetectorRegion(theta - m_brlEdgeTheta) != ECL::DetectorRegion::BRL) return true;
393  else if (ECL::getDetectorRegion(theta + m_brlEdgeTheta) != ECL::DetectorRegion::BRL) return true;
394  else return false;
395  } else return false;
396 }
397 
398 void ECLTrackClusterMatchingModule::optimizedPTMatchingConsistency(double theta, double pt)
399 {
400  if (ECL::getDetectorRegion(theta) == ECL::DetectorRegion::FWD || ECL::getDetectorRegion(theta) == ECL::DetectorRegion::FWDGap) {
401  for (const auto& matchingThresholdPair : m_matchingThresholdValuesFWD) {
402  if (pt < matchingThresholdPair.first) {
403  m_matchingConsistency = matchingThresholdPair.second;
404  break;
405  }
406  }
407  } else if (ECL::getDetectorRegion(theta) == ECL::DetectorRegion::BRL) {
408  for (const auto& matchingThresholdPair : m_matchingThresholdValuesBRL) {
409  if (theta < matchingThresholdPair.first && pt < matchingThresholdPair.second.first) {
410  m_matchingConsistency = matchingThresholdPair.second.second;
411  break;
412  }
413  }
414  } else if (ECL::getDetectorRegion(theta) == ECL::DetectorRegion::BWD
415  || ECL::getDetectorRegion(theta) == ECL::DetectorRegion::BWDGap) {
416  for (const auto& matchingThresholdPair : m_matchingThresholdValuesBWD) {
417  if (pt < matchingThresholdPair.first) {
418  m_matchingConsistency = matchingThresholdPair.second;
419  break;
420  }
421  }
422  }
423 }
Class to store calibrated ECLDigits: ECLCalDigits.
Definition: ECLCalDigit.h:23
ECL cluster data.
Definition: ECLCluster.h:27
Class to store ECL Showers.
Definition: ECLShower.h:25
int getHypothesisId() const
Get Hypothesis Id.
Definition: ECLShower.h:259
double getPhi() const
Get Phi.
Definition: ECLShower.h:284
int getDetectorRegion() const
Return detector region: 0: below acceptance, 1: FWD, 2: BRL, 3: BWD, 11: FWDGAP, 13: BWDGAP.
Definition: ECLShower.h:447
double getTheta() const
Get Theta.
Definition: ECLShower.h:279
The module creates and saves a Relation between Tracks and ECLCluster in the DataStore.
Class to hold the parameterizations of the RMS for the difference in polar and azimuthal angle betwee...
Class to hold the matching thresholds for the track-ECLCluster matching.
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:30
ExtHitStatus getStatus() const
Get state of extrapolation at this hit.
Definition: ExtHit.h:128
int getCopyID() const
Get detector-element ID of sensitive element within detector.
Definition: ExtHit.h:124
Const::EDetector getDetectorID() const
Get detector ID of this extrapolation hit.
Definition: ExtHit.h:120
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
Base class for Modules.
Definition: Module.h:72
int getArrayIndex() const
Returns this object's array index (in StoreArray), or -1 if not found.
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
Values of the result of a track fit with a given particle hypothesis.
short getChargeSign() const
Return track charge (1 or -1).
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
double getTransverseMomentum() const
Getter for the absolute value of the transverse momentum at the perigee.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
Class that bundles various TrackFitResults.
Definition: Track.h:25
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
ExtHitStatus
Define state of extrapolation for each recorded hit.
Definition: ExtHit.h:25
Abstract base class for different kinds of events.