Belle II Software  release-06-01-15
StandardTrackingPerformanceModule.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/standardTrackingPerformance/StandardTrackingPerformanceModule.h>
10 
11 #include <framework/datastore/StoreArray.h>
12 #include <framework/datastore/StoreObjPtr.h>
13 #include <framework/dataobjects/EventMetaData.h>
14 #include <framework/datastore/RelationVector.h>
15 #include <framework/gearbox/Const.h>
16 
17 #include <tracking/dataobjects/RecoTrack.h>
18 #include <genfit/TrackPoint.h>
19 #include <genfit/KalmanFitterInfo.h>
20 
21 #include <mdst/dataobjects/MCParticle.h>
22 #include <mdst/dataobjects/Track.h>
23 
24 #include <pxd/reconstruction/PXDRecoHit.h>
25 #include <svd/reconstruction/SVDRecoHit.h>
26 #include <svd/reconstruction/SVDRecoHit2D.h>
27 #include <cdc/dataobjects/CDCRecoHit.h>
28 
29 #include <root/TFile.h>
30 #include <root/TTree.h>
31 
32 using namespace Belle2;
33 
34 //-----------------------------------------------------------------
35 // Register the Module
36 //-----------------------------------------------------------------
37 REG_MODULE(StandardTrackingPerformance)
38 
39 StandardTrackingPerformanceModule::StandardTrackingPerformanceModule() :
40  Module(), m_outputFile(nullptr), m_dataTree(nullptr), m_trackProperties(-999), m_pValue(-999),
41  m_nGeneratedChargedStableMcParticles(-999), m_nReconstructedChargedStableTracks(-999),
42  m_nFittedChargedStabletracks(-999)
43 {
44  setDescription("Module to test the tracking efficiency. Writes information about the tracks and MCParticles in a ROOT file.");
45  setPropertyFlags(c_ParallelProcessingCertified);
46 
47  addParam("outputFileName", m_outputFileName, "Name of output root file.",
48  std::string("StandardTrackingPerformanceOutput.root"));
49  addParam("recoTracksStoreArrayName", m_recoTracksStoreArrayName, "Name of the RecoTracks StoreArray.",
50  std::string(""));
51  addParam("daughterPDGs", m_signalDaughterPDGs, "PDG codes of B daughters.",
52  std::vector<int>(0));
53 }
54 
55 void StandardTrackingPerformanceModule::initialize()
56 {
57  // MCParticles and Tracks needed for this module
58  StoreArray<MCParticle> mcParticles;
59  mcParticles.isRequired();
60  StoreArray<Track> tracks;
61  tracks.isRequired();
62 
64  recoTracks.isRequired();
65 
66  StoreArray<TrackFitResult> trackFitResults;
67  trackFitResults.isRequired();
68 
69  m_outputFile = new TFile(m_outputFileName.c_str(), "RECREATE");
70  TDirectory* oldDir = gDirectory;
71  m_outputFile->cd();
72  m_dataTree = new TTree("data", "data");
73  oldDir->cd();
74 
75  setupTree();
76 }
77 
79 {
80  StoreObjPtr<EventMetaData> eventMetaData("EventMetaData", DataStore::c_Event);
81  unsigned long eventNumber = eventMetaData->getEvent();
82  unsigned long runNumber = eventMetaData->getRun();
83  unsigned long expNumber = eventMetaData->getExperiment();
84 
85  B2DEBUG(99,
86  "Processes experiment " << expNumber << " run " << runNumber << " event " << eventNumber);
87 
89  StoreArray<MCParticle> mcParticles;
90 
94 
95  for (MCParticle& mcParticle : mcParticles) {
96  // check status of mcParticle
97  if (isPrimaryMcParticle(mcParticle) && isChargedStable(mcParticle) && mcParticle.hasStatus(MCParticle::c_StableInGenerator)) {
99 
100  int pdgCode = mcParticle.getPDG();
101  B2DEBUG(99, "Primary MCParticle has PDG code " << pdgCode);
102 
104 
105  m_trackProperties.cosTheta_gen = mcParticle.getMomentum().CosTheta();
106  m_trackProperties.ptot_gen = mcParticle.getMomentum().Mag();
107  m_trackProperties.pt_gen = mcParticle.getMomentum().Pt();
108  m_trackProperties.px_gen = mcParticle.getMomentum().Px();
109  m_trackProperties.py_gen = mcParticle.getMomentum().Py();
110  m_trackProperties.pz_gen = mcParticle.getMomentum().Pz();
111  m_trackProperties.x_gen = mcParticle.getVertex().X();
112  m_trackProperties.y_gen = mcParticle.getVertex().Y();
113  m_trackProperties.z_gen = mcParticle.getVertex().Z();
114 
115 
116  const RecoTrack* recoTrack = nullptr;
117  double maximumWeight = -2;
118  /*// MIMIK OLD MODULE, which took always the last one (which is probably a CDC track).
119  for(const Track& relatedTrack : mcParticle.getRelationsWith<Track>()) {
120  b2Track = &relatedTrack;
121  }*/
122  // find highest rated Track
123  const auto& relatedRecoTracks = mcParticle.getRelationsWith<RecoTrack>(m_recoTracksStoreArrayName);
124  for (unsigned int index = 0; index < relatedRecoTracks.size(); ++index) {
125  const RecoTrack* relatedRecoTrack = relatedRecoTracks.object(index);
126  const double weight = relatedRecoTracks.weight(index);
127 
128  const unsigned int numberOfRelatedTracks = relatedRecoTrack->getRelationsWith<Track>().size();
129  B2ASSERT("B2Track <-> RecoTrack is not a 1:1 relation as expected!", numberOfRelatedTracks <= 1);
130  // use only the fitted reco tracks
131  if (numberOfRelatedTracks == 1) {
132  if (weight > maximumWeight) {
133  maximumWeight = weight;
134  recoTrack = relatedRecoTrack;
135  }
136  }
137  }
138 
139  if (recoTrack) {
140  const Track* b2Track = recoTrack->getRelated<Track>();
141  const TrackFitResult* fitResult = b2Track->getTrackFitResultWithClosestMass(Const::pion);
142  B2ASSERT("Related Belle2 Track has no related track fit result for pion!", fitResult);
143 
145  // write some data to the root tree
146  TVector3 mom = fitResult->getMomentum();
147  m_trackProperties.cosTheta = mom.CosTheta();
148  m_trackProperties.ptot = mom.Mag();
149  m_trackProperties.pt = mom.Pt();
150  m_trackProperties.px = mom.Px();
151  m_trackProperties.py = mom.Py();
152  m_trackProperties.pz = mom.Pz();
153  m_trackProperties.x = fitResult->getPosition().X();
154  m_trackProperties.y = fitResult->getPosition().Y();
155  m_trackProperties.z = fitResult->getPosition().Z();
156 
157  m_pValue = fitResult->getPValue();
158 
159  // Count hits
164  for (genfit::TrackPoint* tp : recoTrack->getHitPointsWithMeasurement()) {
165  for (genfit::AbsMeasurement* m : tp->getRawMeasurements()) {
166  if (dynamic_cast<PXDRecoHit*>(m))
168  else if (dynamic_cast<SVDRecoHit*>(m))
170  else if (dynamic_cast<SVDRecoHit2D*>(m))
172  else if (dynamic_cast<CDCRecoHit*>(m))
174  else
175  B2ERROR("Unknown AbsMeasurement in track.");
176 
177  std::vector<double> weights;
178  genfit::KalmanFitterInfo* kalmanInfo = tp->getKalmanFitterInfo();
179  if (kalmanInfo)
180  weights = kalmanInfo->getWeights();
181 
182  for (size_t i = 0;
183  (i < weights.size()
185  ++i) {
187  }
188  }
189  }
190  }
191 
192  m_dataTree->Fill(); // write data to tree
193  }
194  }
195 }
196 
197 
199 {
200  writeData();
201 }
202 
204 {
205  return mcParticle.hasStatus(MCParticle::c_PrimaryParticle);
206 }
207 
209 {
210  return Const::chargedStableSet.find(abs(mcParticle.getPDG()))
212 }
213 
215 {
216  if (m_dataTree == nullptr) {
217  B2FATAL("Data tree was not created.");
218  }
219 
222 
225 
228 
231 
234 
237 
240 
243 
246 
249 
250  addVariableToTree("pValue", m_pValue);
251 
255 
256  m_dataTree->Branch("nWeights", &m_trackProperties.nWeights, "nWeights/I");
257  m_dataTree->Branch("weights", &m_trackProperties.weights, "weights[nWeights]/F");
258 }
259 
261 {
262  if (m_dataTree != nullptr) {
263  TDirectory* oldDir = gDirectory;
264  if (m_outputFile)
265  m_outputFile->cd();
266  m_dataTree->Write();
267  oldDir->cd();
268  }
269  if (m_outputFile != nullptr) {
270  m_outputFile->Close();
271  }
272 }
273 
275 {
276  std::sort(m_signalDaughterPDGs.begin(), m_signalDaughterPDGs.end());
277 
278  for (const MCParticle& mcParticle : mcParticles) {
279  // continue if mcParticle is not a B meson
280  if (abs(mcParticle.getPDG()) != 511 && abs(mcParticle.getPDG()) != 521)
281  continue;
282 
283  if (isSignalDecay(mcParticle)) {
284  addChargedStable(mcParticle);
285  break;
286  }
287  }
288 }
289 
291 {
292  std::vector<int> daughterPDGs;
293  std::vector<MCParticle*> daughterMcParticles = mcParticle.getDaughters();
294 
295  // remove photons from list
296  daughterMcParticles = removeFinalStateRadiation(daughterMcParticles);
297 
298  for (auto daughterMCParticle : daughterMcParticles) {
299  daughterPDGs.push_back(daughterMCParticle->getPDG());
300  }
301 
302  std::sort(daughterPDGs.begin(), daughterPDGs.end());
303 
304  bool isSignal = (daughterPDGs == m_signalDaughterPDGs);
305 
306  if (isSignal)
307  m_signalMCParticles = daughterMcParticles;
308 
309  return isSignal;
310 }
311 
314 std::vector<MCParticle*> StandardTrackingPerformanceModule::removeFinalStateRadiation(const std::vector<MCParticle*>& in_daughters)
315 {
316  std::vector<MCParticle*> daughtersWOFSR;
317  for (unsigned int iDaughter = 0; iDaughter < in_daughters.size(); iDaughter++)
318  if (abs(in_daughters[iDaughter]->getPDG()) != Const::photon.getPDGCode())
319  daughtersWOFSR.push_back(in_daughters[iDaughter]);
320 
321  return daughtersWOFSR;
322 }
323 
325 {
326 
327  // mcparticle is not a charged stable decays into daughters
328  // loop over daughters and add the charged stable particles recursively to the vector
329  std::vector<MCParticle*> daughters = mcParticle.getDaughters();
330  if (daughters.size() != 0) {
331  for (auto daughterMcParticle : daughters) {
332  addChargedStable(*daughterMcParticle);
333  }
334  }
335 
336  // charged stable particle is added to the interesting particle vector
337  if (isChargedStable(mcParticle) && isPrimaryMcParticle(mcParticle)
339 // B2DEBUG(99,
340 // "Found a charged stable particle. Add it to interesting MCParticles. PDG(" << mcParticle->getPDG() << ").");
341  m_interestingChargedStableMcParcticles.push_back(&mcParticle);
342  return;
343  }
344 
345  return;
346 }
347 
349 {
350  m_trackProperties = -999;
351 
352  m_pValue = -999;
353 }
354 
355 void StandardTrackingPerformanceModule::addVariableToTree(const std::string& varName, double& varReference)
356 {
357  std::stringstream leaf;
358  leaf << varName << "/D";
359  m_dataTree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
360 }
361 
362 void StandardTrackingPerformanceModule::addVariableToTree(const std::string& varName, int& varReference)
363 {
364  std::stringstream leaf;
365  leaf << varName << "/I";
366  m_dataTree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
367 }
This class is used to transfer CDC information to the track fit.
Definition: CDCRecoHit.h:32
const ParticleType & find(int pdg) const
Returns particle in set with given PDG code, or invalidParticle if not found.
Definition: Const.h:452
static const ParticleSet chargedStableSet
set of charged stable particles
Definition: Const.h:499
static const ChargedStable pion
charged pion particle
Definition: Const.h:542
static const ParticleType invalidParticle
Invalid particle, used internally.
Definition: Const.h:561
static const ParticleType photon
photon particle
Definition: Const.h:554
@ 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
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
@ c_StableInGenerator
bit 1: Particle is stable, i.e., not decaying in the generator.
Definition: MCParticle.h:49
std::vector< Belle2::MCParticle * > getDaughters() const
Get vector of all daughter particles, empty vector if none.
Definition: MCParticle.cc:50
bool hasStatus(unsigned short int bitmask) const
Return if specific status bit is set.
Definition: MCParticle.h:129
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
Base class for Modules.
Definition: Module.h:72
PXDRecoHit - an extended form of PXDCluster containing geometry information.
Definition: PXDRecoHit.h:49
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:76
const std::vector< genfit::TrackPoint * > & getHitPointsWithMeasurement() const
Return a list of measurements and track points, which can be used e.g. to extrapolate....
Definition: RecoTrack.h:606
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.
SVDRecoHit - an extended form of SVDHit containing geometry information.
Definition: SVDRecoHit2D.h:46
SVDRecoHit - an extended form of SVDHit containing geometry information.
Definition: SVDRecoHit.h:47
double m_nReconstructedChargedStableTracks
total number of reconstructed track candidates
void event() override
Fill the tree with the event data.
std::vector< const MCParticle * > m_interestingChargedStableMcParcticles
vector with all interesting charged stable MCParticles in the event
std::vector< MCParticle * > m_signalMCParticles
vector with all MCParticles of the searched signal decay
double m_nGeneratedChargedStableMcParticles
total number of genrated charged stable MCParticles
void terminate() override
Write the tree into the opened root file.
void addChargedStable(const MCParticle &mcParticle)
Add all charged stable particles to a vector which originate from.
ParticleProperties m_trackProperties
properties of a reconstructed track
bool isChargedStable(const MCParticle &mcParticle)
Tests if MCPArticle is a charged stable particle.
void setVariablesToDefaultValue()
Sets all variables to the default value, here -999.
bool isPrimaryMcParticle(const MCParticle &mcParticle)
Tests if MCParticle is a primary one.
void addVariableToTree(const std::string &varName, double &varReference)
add a variable with double format
TTree * m_dataTree
root tree with all output data.
std::vector< MCParticle * > removeFinalStateRadiation(const std::vector< MCParticle * > &in_daughters)
Remove all photons from a MCParticle vector.
double m_nFittedChargedStabletracks
total number of fitted tracks
void writeData()
write root tree to output file and close the file
void findSignalMCParticles(const StoreArray< MCParticle > &mcParticles)
Find a MCParticle of a decay chain specified by the user (not implemented yet).
bool isSignalDecay(const MCParticle &mcParticle)
Tests if mcParticle has the searched decay chain.
std::string m_recoTracksStoreArrayName
genfit::Track collection name
std::vector< int > m_signalDaughterPDGs
PDG codes of the B daughters of the interesting decay channel.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
Values of the result of a track fit with a given particle hypothesis.
Class that bundles various TrackFitResults.
Definition: Track.h:25
const TrackFitResult * getTrackFitResultWithClosestMass(const Const::ChargedStable &requestedType) const
Return the track fit for a fit hypothesis with the closest mass.
Definition: Track.cc:68
Contains the measurement and covariance in raw detector coordinates.
Collects information needed and produced by a AbsKalmanFitter implementations and is specific to one ...
std::vector< double > getWeights() const
Get weights of measurements.
Object containing AbsMeasurement and AbsFitterInfo objects.
Definition: TrackPoint.h:46
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.
double pt
measured transverse momentum
int nCDChits
Number of CDC hits in reconstructed track
double pz
measured momentum in z direction
double x_gen
x value of generated position
float weights[maxNweights]
Weights of the hits in sequence
double pt_gen
generated transverse momentum
double mass_gen
generated mass
double y_gen
y value of generated position
int nPXDhits
Number of PXD hits in reconstructed track
int nSVDhits
Number of SVD hits in reconstructed track
double px
measured momentum in x direction
static const int maxNweights
the maximum number of stored weights
double z
measured z value of position
int nWeights
Number of entries in weights array
double y
measured y value of position
double ptot_gen
generated total momentum
double py
measured momentum in y direction
double cosTheta
polar angle of measured momentum vector
double pz_gen
generated momentum in z direction
double py_gen
generated momentum in y direction
double cosTheta_gen
polar angle of generated momentum vector
double px_gen
generated momentum in x direction
double z_gen
z value of generated position
double x
measured x value of position
double ptot
measured total momentum