Belle II Software  release-06-00-14
VXDDedxPIDModule.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 <reconstruction/modules/VXDDedxPID/VXDDedxPIDModule.h>
10 #include <reconstruction/modules/VXDDedxPID/HelixHelper.h>
11 
12 #include <framework/gearbox/Const.h>
13 
14 #include <vxd/geometry/GeoCache.h>
15 
16 #include <genfit/MaterialEffects.h>
17 
18 #include <TH2F.h>
19 
20 #include <memory>
21 #include <cassert>
22 #include <cmath>
23 #include <algorithm>
24 
25 using namespace Belle2;
26 using namespace Dedx;
27 
28 REG_MODULE(VXDDedxPID)
29 
31 {
32  setPropertyFlags(c_ParallelProcessingCertified);
33 
34  //Set module properties
35  setDescription("Extract dE/dx and corresponding log-likelihood from fitted tracks and hits in the SVD and PXD.");
36 
37  //Parameter definitions
38  addParam("useIndividualHits", m_useIndividualHits,
39  "Include PDF value for each hit in likelihood. If false, the truncated mean of dedx values for the detectors will be used.", false);
40  // no need to define highest and lowest truncated value as we always remove higest two dE/dx value from the 8 dE/dx value
41  addParam("onlyPrimaryParticles", m_onlyPrimaryParticles, "Only save data for primary particles (as determined by MC truth)", false);
42  addParam("usePXD", m_usePXD, "Use PXDClusters for dE/dx calculation", false);
43  addParam("useSVD", m_useSVD, "Use SVDClusters for dE/dx calculation", true);
44  addParam("trackDistanceThreshold", m_trackDistanceThreshhold,
45  "Use a faster helix parametrisation, with corrections as soon as the approximation is more than ... cm off.", double(4.0));
46  addParam("enableDebugOutput", m_enableDebugOutput, "Option to write out debugging information to DedxTracks (DataStore objects).",
47  false);
48  addParam("ignoreMissingParticles", m_ignoreMissingParticles, "Ignore particles for which no PDFs are found", false);
49 
50  m_eventID = -1;
51  m_trackID = 0;
52 }
53 
55 
57 {
58  //load dedx:momentum PDFs
59  if (!m_DBDedxPDFs) B2FATAL("No VXD Dedx PDFS available");
60  int nBinsXPXD, nBinsYPXD;
61  double xMinPXD, xMaxPXD, yMinPXD, yMaxPXD;
62  nBinsXPXD = nBinsYPXD = -1;
63  xMinPXD = xMaxPXD = yMinPXD = yMaxPXD = 0.0;
64 
65  int nBinsXSVD, nBinsYSVD;
66  double xMinSVD, xMaxSVD, yMinSVD, yMaxSVD;
67  nBinsXSVD = nBinsYSVD = -1;
68  xMinSVD = xMaxSVD = yMinSVD = yMaxSVD = 0.0;
69 
70  for (unsigned int iPart = 0; iPart < 6; iPart++) {
71  const int pdgCode = Const::chargedStableSet.at(iPart).getPDGCode();
72  const TH2F* svd_pdf = m_DBDedxPDFs->getSVDPDF(iPart, !m_useIndividualHits);
73  const TH2F* pxd_pdf = m_DBDedxPDFs->getPXDPDF(iPart, !m_useIndividualHits);
74 
75  if (pxd_pdf->GetEntries() == 0 || svd_pdf->GetEntries() == 0) {
76  if (m_ignoreMissingParticles)
77  continue;
78  B2FATAL("Couldn't find PDF for PDG " << pdgCode);
79  }
80 
81  //check that PXD PDFs have the same dimensions and same binning
82  const double epsFactor = 1e-5;
83  if (nBinsXPXD == -1 and nBinsYPXD == -1) {
84  nBinsXPXD = pxd_pdf->GetNbinsX();
85  nBinsYPXD = pxd_pdf->GetNbinsY();
86  xMinPXD = pxd_pdf->GetXaxis()->GetXmin();
87  xMaxPXD = pxd_pdf->GetXaxis()->GetXmax();
88  yMinPXD = pxd_pdf->GetYaxis()->GetXmin();
89  yMaxPXD = pxd_pdf->GetYaxis()->GetXmax();
90  } else if (nBinsXPXD != pxd_pdf->GetNbinsX()
91  or nBinsYPXD != pxd_pdf->GetNbinsY()
92  or fabs(xMinPXD - pxd_pdf->GetXaxis()->GetXmin()) > epsFactor * xMaxPXD
93  or fabs(xMaxPXD - pxd_pdf->GetXaxis()->GetXmax()) > epsFactor * xMaxPXD
94  or fabs(yMinPXD - pxd_pdf->GetYaxis()->GetXmin()) > epsFactor * yMaxPXD
95  or fabs(yMaxPXD - pxd_pdf->GetYaxis()->GetXmax()) > epsFactor * yMaxPXD) {
96  B2FATAL("PDF for PDG " << pdgCode << ", PXD has binning/dimensions differing from previous PDF.");
97  }
98 
99  //check that SVD PDFs have the same dimensions and same binning
100  if (nBinsXSVD == -1 and nBinsYSVD == -1) {
101  nBinsXSVD = svd_pdf->GetNbinsX();
102  nBinsYSVD = svd_pdf->GetNbinsY();
103  xMinSVD = svd_pdf->GetXaxis()->GetXmin();
104  xMaxSVD = svd_pdf->GetXaxis()->GetXmax();
105  yMinSVD = svd_pdf->GetYaxis()->GetXmin();
106  yMaxSVD = svd_pdf->GetYaxis()->GetXmax();
107  } else if (nBinsXSVD != svd_pdf->GetNbinsX()
108  or nBinsYSVD != svd_pdf->GetNbinsY()
109  or fabs(xMinSVD - svd_pdf->GetXaxis()->GetXmin()) > epsFactor * xMaxSVD
110  or fabs(xMaxSVD - svd_pdf->GetXaxis()->GetXmax()) > epsFactor * xMaxSVD
111  or fabs(yMinSVD - svd_pdf->GetYaxis()->GetXmin()) > epsFactor * yMaxSVD
112  or fabs(yMaxSVD - svd_pdf->GetYaxis()->GetXmax()) > epsFactor * yMaxSVD) {
113  B2FATAL("PDF for PDG " << pdgCode << ", PXD has binning/dimensions differing from previous PDF.");
114  }
115  }
116 }
117 
119 {
120 
121  // required inputs
122  m_tracks.isRequired();
123  m_recoTracks.isRequired();
124 
125  //optional inputs
126  m_mcparticles.isOptional();
127  m_tracks.optionalRelationTo(m_mcparticles);
128 
129  if (m_useSVD)
130  m_svdClusters.isRequired();
131  else
132  m_svdClusters.isOptional();
133  if (m_usePXD)
134  m_pxdClusters.isRequired();
135  else
136  m_pxdClusters.isOptional();
137 
138  // register optional outputs
139  if (m_enableDebugOutput) {
140  m_dedxTracks.registerInDataStore();
141  m_tracks.registerRelationTo(m_dedxTracks);
142  }
143 
144  // register outputs
145  m_dedxLikelihoods.registerInDataStore();
146  m_tracks.registerRelationTo(m_dedxLikelihoods);
147 
148  m_DBDedxPDFs.addCallback([this]() {checkPDFs();});
149  checkPDFs();
150 
151 
152  // create instances here to not confuse profiling
154 
155  if (!genfit::MaterialEffects::getInstance()->isInitialized()) {
156  B2FATAL("Need to have SetupGenfitExtrapolationModule in path before this one.");
157  }
158 }
159 
161 {
162  // go through Tracks
163  // get fitresult and gftrack and do extrapolations, save corresponding dE/dx and likelihood values
164  // get genfit::TrackCand through genfit::Track::getCand()
165  // get hit indices through genfit::TrackCand::getHit(...)
166  // create one VXDDedxTrack per fitresult/gftrack
167  //create one DedkLikelihood per Track (plus rel)
168  m_eventID++;
169 
170  // inputs
171  const int numMCParticles = m_mcparticles.getEntries();
172 
173  // **************************************************
174  //
175  // LOOP OVER TRACKS
176  //
177  // **************************************************
178 
179  for (const auto& track : m_tracks) {
180  m_trackID++;
181 
182  std::shared_ptr<VXDDedxTrack> dedxTrack = std::make_shared<VXDDedxTrack>();
183  dedxTrack->m_eventID = m_eventID;
184  dedxTrack->m_trackID = m_trackID;
185 
186  // load the pion fit hypothesis or the hypothesis which is the closest in mass to a pion
187  // the tracking will not always successfully fit with a pion hypothesis
188  const TrackFitResult* fitResult = track.getTrackFitResultWithClosestMass(Const::pion);
189  if (!fitResult) {
190  B2WARNING("No related fit for track ...");
191  continue;
192  }
193 
194  if ((m_enableDebugOutput or m_onlyPrimaryParticles) and numMCParticles != 0) {
195  // find MCParticle corresponding to this track
196  const MCParticle* mcpart = track.getRelatedTo<MCParticle>();
197 
198  if (mcpart) {
199  if (m_onlyPrimaryParticles && !mcpart->hasStatus(MCParticle::c_PrimaryParticle)) {
200  continue; //not a primary particle, ignore
201  }
202 
203  //add some MC truths to VXDDedxTrack object
204  dedxTrack->m_pdg = mcpart->getPDG();
205  const MCParticle* mother = mcpart->getMother();
206  dedxTrack->m_motherPDG = mother ? mother->getPDG() : 0;
207 
208  const TVector3 trueMomentum = mcpart->getMomentum();
209  dedxTrack->m_pTrue = trueMomentum.Mag();
210  }
211  }
212 
213  // get momentum (at origin) from fit result
214  const TVector3& trackPos = fitResult->getPosition();
215  const TVector3& trackMom = fitResult->getMomentum();
216  dedxTrack->m_p = trackMom.Mag();
217  dedxTrack->m_cosTheta = trackMom.CosTheta();
218  dedxTrack->m_charge = fitResult->getChargeSign();
219 
220  // dE/dx values will be calculated using associated RecoTrack
221  const RecoTrack* recoTrack = track.getRelatedTo<RecoTrack>();
222  if (!recoTrack) {
223  B2WARNING("No related track for this fit...");
224  continue;
225  }
226 
227  // Check to see if the track is pruned. We use the cardinal representation for this, as we do not expect that one
228  // representation is pruned whereas the other is not.
229  if (recoTrack->getTrackFitStatus()->isTrackPruned()) {
230  B2ERROR("GFTrack is pruned, please run VXDDedxPID only on unpruned tracks! Skipping this track.");
231  continue;
232  }
233 
234  //used for PXD/SVD hits
235  const HelixHelper helixAtOrigin(trackPos, trackMom, dedxTrack->m_charge);
236 
237  if (m_usePXD) {
238  const std::vector<PXDCluster*>& pxdClusters = recoTrack->getPXDHitList();
239  saveSiHits(dedxTrack.get(), helixAtOrigin, pxdClusters);
240  }
241 
242  if (m_useSVD) {
243  const std::vector<SVDCluster*>& svdClusters = recoTrack->getSVDHitList();
244  saveSiHits(dedxTrack.get(), helixAtOrigin, svdClusters);
245  }
246 
247  if (dedxTrack->dedx.empty()) {
248  B2DEBUG(50, "Found track with no hits, ignoring.");
249  continue;
250  }
251 
252  // calculate likelihoods for truncated mean
253  if (!m_useIndividualHits) {
254  for (int detector = 0; detector <= c_SVD; detector++) {
255  if (!detectorEnabled(static_cast<Detector>(detector)))
256  continue; //unwanted detector
257 
258  if (detector == 0) savePXDLogLikelihood(dedxTrack->m_vxdLogl, dedxTrack->m_p, dedxTrack->m_dedxAvgTruncated[detector]);
259  else if (detector == 1) saveSVDLogLikelihood(dedxTrack->m_vxdLogl, dedxTrack->m_p, dedxTrack->m_dedxAvgTruncated[detector]);
260  }
261  }
262 
263  if (m_enableDebugOutput) {
264  // add a few last things to the VXDDedxTrack
265  const int numDedx = dedxTrack->dedx.size();
266  dedxTrack->m_nHits = numDedx;
267  // no need to define lowedgetruncated and highedgetruncated as we always remove the highest 2 dE/dx values from 8 dE/dx value
268 
269  dedxTrack->m_nHitsUsed = numDedx - 2;
270 
271  // now book the information for this track
272  VXDDedxTrack* newVXDDedxTrack = m_dedxTracks.appendNew(*dedxTrack);
273  track.addRelationTo(newVXDDedxTrack);
274  }
275 
276  // save VXDDedxLikelihood
277  VXDDedxLikelihood* likelihoodObj = m_dedxLikelihoods.appendNew(dedxTrack->m_vxdLogl);
278  track.addRelationTo(likelihoodObj);
279 
280  } // end of loop over tracks
281 }
282 
284 {
285 
286  B2DEBUG(50, "VXDDedxPIDModule exiting after processing " << m_trackID <<
287  " tracks in " << m_eventID + 1 << " events.");
288 }
289 // calculateMeans need some change as we always remove highest 2 dE/dx values
290 void VXDDedxPIDModule::calculateMeans(double* mean, double* truncatedMean, double* truncatedMeanErr,
291  const std::vector<double>& dedx) const
292 {
293  // Calculate the truncated average by skipping only highest two value
294  std::vector<double> sortedDedx = dedx;
295  std::sort(sortedDedx.begin(), sortedDedx.end());
296 
297  double truncatedMeanTmp = 0.0;
298  double meanTmp = 0.0;
299  double sumOfSquares = 0.0;
300  const int numDedx = sortedDedx.size();
301 
302 
303  for (int i = 0; i < numDedx; i++) {
304  meanTmp += sortedDedx[i];
305  }
306  if (numDedx != 0) {
307  meanTmp /= numDedx;
308  }
309 
310  for (int i = 0; i < numDedx - 2; i++) {
311  truncatedMeanTmp += sortedDedx[i];
312  sumOfSquares += sortedDedx[i] * sortedDedx[i];
313  }
314  if (numDedx - 2 != 0) {
315  truncatedMeanTmp /= numDedx - 2;
316  }
317 
318  *mean = meanTmp;
319  *truncatedMean = truncatedMeanTmp;
320 
321  if (numDedx - 2 > 1) {
322  *truncatedMeanErr = sqrt(sumOfSquares / double(numDedx - 2) - truncatedMeanTmp * truncatedMeanTmp) / double(
323  (numDedx - 2) - 1);
324  } else {
325  *truncatedMeanErr = 0;
326  }
327 }
328 
330 {
332  const VXD::SensorInfoBase& sensor = geo.get(hit->getSensorID());
333 
334  const TVector3 localPos(hit->getU(), hit->getV(), 0.0); //z-component is height over the center of the detector plane
335  const TVector3& globalPos = sensor.pointToGlobal(localPos);
336  const TVector3& localMomentum = helix->momentum(helix->pathLengthToPoint(globalPos));
337 
338  const TVector3& sensorNormal = sensor.vectorToGlobal(TVector3(0.0, 0.0, 1.0));
339  const double angle = sensorNormal.Angle(localMomentum); //includes theta and phi components
340 
341  //I'm assuming there's only one hit per sensor, there are _very_ rare exceptions to that (most likely curlers)
342  return TMath::Min(sensor.getWidth(), sensor.getThickness() / fabs(cos(angle)));
343 }
344 
345 
347 {
349  const VXD::SensorInfoBase& sensor = geo.get(hit->getSensorID());
350 
351  TVector3 a, b;
352  if (hit->isUCluster()) {
353  const float u = hit->getPosition();
354  a = sensor.pointToGlobal(TVector3(sensor.getBackwardWidth() / sensor.getWidth(0) * u, -0.5 * sensor.getLength(), 0.0));
355  b = sensor.pointToGlobal(TVector3(sensor.getForwardWidth() / sensor.getWidth(0) * u, +0.5 * sensor.getLength(), 0.0));
356  } else {
357  const float v = hit->getPosition();
358  a = sensor.pointToGlobal(TVector3(-0.5 * sensor.getWidth(v), v, 0.0));
359  b = sensor.pointToGlobal(TVector3(+0.5 * sensor.getWidth(v), v, 0.0));
360  }
361  const double pathLength = helix->pathLengthToLine(a, b);
362  const TVector3& localMomentum = helix->momentum(pathLength);
363 
364  const TVector3& sensorNormal = sensor.vectorToGlobal(TVector3(0.0, 0.0, 1.0));
365  const double angle = sensorNormal.Angle(localMomentum); //includes theta and phi components
366 
367  return TMath::Min(sensor.getWidth(), sensor.getThickness() / fabs(cos(angle)));
368 }
369 
370 
371 template <class HitClass> void VXDDedxPIDModule::saveSiHits(VXDDedxTrack* track, const HelixHelper& helix,
372  const std::vector<HitClass*>& hits) const
373 {
374  const int numHits = hits.size();
375  if (numHits == 0)
376  return;
377 
379 
380  //figure out which detector to assign hits to
381  const int currentDetector = geo.get(hits.front()->getSensorID()).getType();
382  assert(currentDetector == VXD::SensorInfoBase::PXD or currentDetector == VXD::SensorInfoBase::SVD);
383  assert(currentDetector <= 1); //used as array index
384 
385  std::vector<double> siliconDedx; //used for averages
386  siliconDedx.reserve(numHits);
387 
388  VxdID prevSensor;
389  for (int i = 0; i < numHits; i++) {
390  const HitClass* hit = hits[i];
391  if (!hit) {
392  B2ERROR("Added hit is a null pointer!");
393  continue;
394  }
395  const VxdID& currentSensor = hit->getSensorID();
396  int layer = -1;
397  if (m_enableDebugOutput) {
398  layer = -currentSensor.getLayerNumber();
399  assert(layer >= -6 && layer < 0);
400  }
401 
402  //active medium traversed, in cm (can traverse one sensor at most)
403  //assumption: Si detectors are close enough to the origin that helix is still accurate
404  const double totalDistance = getTraversedLength(hit, &helix);
405  const float charge = hit->getCharge();
406  const float dedx = charge / totalDistance;
407  if (dedx <= 0) {
408  B2WARNING("dE/dx is " << dedx << " in layer " << layer);
409  } else if (i == 0 or prevSensor != currentSensor) { //only save once per sensor (u and v hits share charge!)
410  prevSensor = currentSensor;
411  //store data
412  siliconDedx.push_back(dedx);
413  track->m_dedxAvg[currentDetector] += dedx;
414  track->addDedx(layer, totalDistance, dedx);
415  if (m_useIndividualHits) {
416  if (currentDetector == 0) savePXDLogLikelihood(track->m_vxdLogl, track->m_p, dedx);
417  else if (currentDetector == 1) saveSVDLogLikelihood(track->m_vxdLogl, track->m_p, dedx);
418  }
419  }
420 
421  if (m_enableDebugOutput) {
422  track->addHit(currentSensor, layer, charge, totalDistance, dedx);
423  }
424  }
425 
426  //save averages averages
427  if (!m_useIndividualHits or m_enableDebugOutput) {
428  calculateMeans(&(track->m_dedxAvg[currentDetector]),
429  &(track->m_dedxAvgTruncated[currentDetector]),
430  &(track->m_dedxAvgTruncatedErr[currentDetector]),
431  siliconDedx);
432  }
433 }
434 
435 void VXDDedxPIDModule::savePXDLogLikelihood(double(&logl)[Const::ChargedStable::c_SetSize], double p, float dedx) const
436 {
437  //all pdfs have the same dimensions
438  const TH2F* pdf = m_DBDedxPDFs->getPXDPDF(0, !m_useIndividualHits);
439  const Int_t binX = pdf->GetXaxis()->FindFixBin(p);
440  const Int_t binY = pdf->GetYaxis()->FindFixBin(dedx);
441 
442  for (unsigned int iPart = 0; iPart < Const::ChargedStable::c_SetSize; iPart++) {
443  pdf = m_DBDedxPDFs->getPXDPDF(iPart, !m_useIndividualHits);
444  if (pdf->GetEntries() == 0) //might be NULL if m_ignoreMissingParticles is set
445  continue;
446  double probability = 0.0;
447 
448  //check if this is still in the histogram, take overflow bin otherwise
449  if (binX < 1 or binX > pdf->GetNbinsX()
450  or binY < 1 or binY > pdf->GetNbinsY()) {
451  probability = pdf->GetBinContent(binX, binY);
452  } else {
453  //in normal histogram range. Of course ROOT has a bug that Interpolate()
454  //is not declared as const but it does not modify the internal state so
455  //fine, const_cast it is.
456  probability = const_cast<TH2F*>(pdf)->Interpolate(p, dedx);
457  }
458 
459  if (probability != probability)
460  B2ERROR("probability NAN for a track with p=" << p << " and dedx=" << dedx);
461 
462  //my pdfs aren't perfect...
463  if (probability == 0.0)
464  probability = m_useIndividualHits ? (1e-5) : (1e-3); //likelihoods for truncated mean are much higher
465 
466  logl[iPart] += log(probability);
467  }
468 }
469 
470 void VXDDedxPIDModule::saveSVDLogLikelihood(double(&logl)[Const::ChargedStable::c_SetSize], double p, float dedx) const
471 {
472  //all pdfs have the same dimensions
473  const TH2F* pdf = m_DBDedxPDFs->getSVDPDF(0, !m_useIndividualHits);
474  const Int_t binX = pdf->GetXaxis()->FindFixBin(p);
475  const Int_t binY = pdf->GetYaxis()->FindFixBin(dedx);
476 
477  for (unsigned int iPart = 0; iPart < Const::ChargedStable::c_SetSize; iPart++) {
478  pdf = m_DBDedxPDFs->getSVDPDF(iPart, !m_useIndividualHits);
479  if (pdf->GetEntries() == 0) //might be NULL if m_ignoreMissingParticles is set
480  continue;
481  double probability = 0.0;
482 
483  //check if this is still in the histogram, take overflow bin otherwise
484  if (binX < 1 or binX > pdf->GetNbinsX()
485  or binY < 1 or binY > pdf->GetNbinsY()) {
486  probability = pdf->GetBinContent(binX, binY);
487  } else {
488  //in normal histogram range. Of course ROOT has a bug that Interpolate()
489  //is not declared as const but it does not modify the internal state so
490  //fine, const_cast it is.
491  probability = const_cast<TH2F*>(pdf)->Interpolate(p, dedx);
492  }
493 
494  if (probability != probability)
495  B2ERROR("probability NAN for a track with p=" << p << " and dedx=" << dedx);
496 
497  //my pdfs aren't perfect...
498  if (probability == 0.0)
499  probability = m_useIndividualHits ? (1e-5) : (1e-3); //likelihoods for truncated mean are much higher
500 
501  logl[iPart] += log(probability);
502  }
503 }
static const unsigned int c_SetSize
Number of elements (for use in array bounds etc.)
Definition: Const.h:496
const ParticleType & at(unsigned int index) const
Return particle at given index, or end() if out of range.
Definition: Const.h:430
int getPDGCode() const
PDG code.
Definition: Const.h:354
static const ParticleSet chargedStableSet
set of charged stable particles
Definition: Const.h:499
static const ChargedStable pion
charged pion particle
Definition: Const.h:542
Helper class representing a helical track.
Definition: HelixHelper.h:28
TVector3 momentum(double s=0) const
momentum of the particle, at the helix point corresponding to a flown path length s (from poca).
Definition: HelixHelper.h:109
double pathLengthToPoint(const TVector3 &p) const
returns the path length (along the helix) to the helix point closest to p.
Definition: HelixHelper.h:68
double pathLengthToLine(const TVector3 &a, const TVector3 &b) const
returns the path length (along the helix) to the helix point closest to the line going through points...
Definition: HelixHelper.h:87
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
bool hasStatus(unsigned short int bitmask) const
Return if specific status bit is set.
Definition: MCParticle.h:129
TVector3 getMomentum() const
Return momentum.
Definition: MCParticle.h:198
int m_pdg
PDG-Code of the particle.
Definition: MCParticle.h:525
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
Base class for Modules.
Definition: Module.h:72
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:76
const genfit::FitStatus * getTrackFitStatus(const genfit::AbsTrackRep *representation=nullptr) const
Return the track fit status for the given representation or for the cardinal one. You are not allowed...
Definition: RecoTrack.h:536
std::vector< Belle2::RecoTrack::UsedSVDHit * > getSVDHitList() const
Return an unsorted list of svd hits.
Definition: RecoTrack.h:448
std::vector< Belle2::RecoTrack::UsedPXDHit * > getPXDHitList() const
Return an unsorted list of pxd hits.
Definition: RecoTrack.h:451
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:28
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.
TVector3 getPosition() const
Getter for vector of position at closest approach of track in r/phi projection.
Container for likelihoods obtained by the VXD dE/dx PID (VXDDedxPIDModule).
Extract dE/dx from fitted tracks.
void saveSVDLogLikelihood(double(&logl)[Const::ChargedStable::c_SetSize], double p, float dedx) const
for all particles in the SVD, save log-likelihood values into 'logl'.
virtual void initialize() override
Initialize the module.
virtual void event() override
This method is called for each event.
void checkPDFs()
Check the pdfs for consistency everytime they change in the database.
virtual void terminate() override
End of the event processing.
void saveSiHits(VXDDedxTrack *track, const HelixHelper &helix, const std::vector< HitClass * > &hits) const
save energy loss and hit information from SVD/PXDHits to track
void calculateMeans(double *mean, double *truncatedMean, double *truncatedMeanErr, const std::vector< double > &dedx) const
Save arithmetic and truncated mean for the 'dedx' values.
virtual ~VXDDedxPIDModule()
Destructor.
static double getTraversedLength(const PXDCluster *hit, const HelixHelper *helix)
returns traversed length through active medium of given PXDCluster.
void savePXDLogLikelihood(double(&logl)[Const::ChargedStable::c_SetSize], double p, float dedx) const
for all particles in the PXD, save log-likelihood values into 'logl'.
Debug output for VXDDedxPID module.
Definition: VXDDedxTrack.h:27
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:213
static const SensorInfoBase & get(Belle2::VxdID id)
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.h:139
Base class to provide Sensor Information for PXD and SVD.
SensorType getType() const
Return the Type of the Sensor.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
bool isTrackPruned() const
Has the track been pruned after the fit?
Definition: FitStatus.h:116
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
MCParticle * getMother() const
Returns a pointer to the mother particle.
Definition: MCParticle.h:582
Abstract base class for different kinds of events.