9#include <ecl/modules/eclTrackClusterMatching/ECLTrackClusterMatchingPerformanceModule.h>
11#include <framework/datastore/RelationVector.h>
12#include <mdst/dataobjects/HitPatternCDC.h>
13#include <mdst/dataobjects/HitPatternVXD.h>
14#include <framework/gearbox/Const.h>
16#include <root/TFile.h>
17#include <root/TTree.h>
26ECLTrackClusterMatchingPerformanceModule::ECLTrackClusterMatchingPerformanceModule() :
27 Module(), m_trackProperties()
29 setDescription(
"Module to test the track cluster matching efficiency. Writes information about the tracks and MCParticles in a ROOT file.");
30 setPropertyFlags(c_ParallelProcessingCertified);
32 addParam(
"outputFileName", m_outputFileName,
"Name of output root file.",
33 std::string(
"ECLTrackClusterMatchingPerformance.root"));
34 addParam(
"minClusterEnergy", m_minClusterEnergy,
"Minimal cluster energy in units of particle's true energy for MC match",
36 addParam(
"minWeight", m_minWeight,
"Fraction of cluster energy required to originate from particle to be matched to cluster.",
38 addParam(
"trackClusterRelationName", m_trackClusterRelationName,
"Name of relation array between tracks and ECL clusters",
51 TDirectory* oldDir = gDirectory;
68 bool found_photon =
false, found_mcmatch =
false, found_charged_stable =
false;
70 const auto& relatedMCParticles = eclCluster.getRelationsWith<
MCParticle>();
71 for (
unsigned int index = 0; index < relatedMCParticles.size() && (!found_photon || !found_charged_stable); ++index) {
72 const MCParticle* relatedMCParticle = relatedMCParticles.object(index);
77 const double weight = relatedMCParticles.weight(index);
87 if (eclCluster.getEnergy(hypo) >= 0.5 * relatedMCParticle->
getEnergy()) {
89 found_charged_stable =
true;
98 if (eclCluster.isTrack()) {
131 int pdgCode = mcParticle.getPDG();
132 B2DEBUG(29,
"Primary MCParticle has PDG code " << pdgCode);
146 const Track* b2Track =
nullptr;
147 double maximumWeight = -2;
150 for (
unsigned int index = 0; index < relatedTracks.size(); ++index) {
151 const Track* relatedTrack = relatedTracks.object(index);
152 const double weight = relatedTracks.weight(index);
154 if (weight > maximumWeight) {
155 maximumWeight = weight;
156 b2Track = relatedTrack;
162 B2ASSERT(
"Related Belle2 Track has no related track fit result!", fitResult);
165 ROOT::Math::XYZVector mom = fitResult->
getMomentum();
189 const ECLCluster* eclCluster_PhotonHypothesis_track_related =
nullptr;
190 const ECLCluster* eclCluster_HadronHypothesis_track_related =
nullptr;
192 if (!(eclCluster.isTrack()))
continue;
200 eclCluster_PhotonHypothesis_track_related = &eclCluster;
209 eclCluster_HadronHypothesis_track_related = &eclCluster;
214 const ECLCluster* eclCluster_matchedBestToMCParticle =
nullptr;
216 for (
unsigned int index = 0; index < relatedECLClusters.size(); ++index) {
217 const ECLCluster* relatedECLCluster = relatedECLClusters.object(index);
219 const double weight = relatedECLClusters.weight(index);
220 if (weight > maximumWeight) {
221 eclCluster_matchedBestToMCParticle = relatedECLCluster;
222 maximumWeight = weight;
225 if (eclCluster_matchedBestToMCParticle !=
nullptr) {
231 if ((eclCluster_PhotonHypothesis_track_related == eclCluster_matchedBestToMCParticle
233 || (eclCluster_HadronHypothesis_track_related == eclCluster_matchedBestToMCParticle
264 B2FATAL(
"Data tree or event tree was not created.");
362 TDirectory* oldDir = gDirectory;
457 std::stringstream leaf;
458 leaf << varName <<
"/D";
459 tree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
464 std::stringstream leaf;
465 leaf << varName <<
"/I";
466 tree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
Provides a type-safe way to pass members of the chargedStableSet set.
const ParticleType & find(int pdg) const
Returns particle in set with given PDG code, or invalidParticle if not found.
int getPDGCode() const
PDG code.
static const ParticleSet chargedStableSet
set of charged stable particles
static const ParticleType invalidParticle
Invalid particle, used internally.
static const ParticleType photon
photon particle
bool hasHypothesis(EHypothesisBit bitmask) const
Return if specific hypothesis bit is set.
double getPhi() const
Return Corrected Phi of Shower (radian).
int getDetectorRegion() const
Return detector region: 0: below acceptance, 1: FWD, 2: BRL, 3: BWD, 11: FWDGAP, 13: BWDGAP.
double getTheta() const
Return Corrected Theta of Shower (radian).
EHypothesisBit
The hypothesis bits for this ECLCluster (Connected region (CR) is split using this hypothesis.
@ c_nPhotons
CR is split into n photons (N1)
@ c_neutralHadron
CR is reconstructed as a neutral hadron (N2)
short getLastLayer() const
Returns the index of the last layer with a hit.
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
unsigned short getNPXDHits() const
Get total number of hits in the PXD.
unsigned short getNSVDHits() const
Get total number of hits in the SVD.
A Class to store the Monte Carlo particle information.
float getEnergy() const
Return particle energy in GeV.
@ c_PrimaryParticle
bit 0: Particle is primary particle.
@ c_StableInGenerator
bit 1: Particle is stable, i.e., not decaying in the generator.
bool hasStatus(unsigned short int bitmask) const
Return if specific status bit is set.
int getPDG() const
Return PDG code of particle.
RelationVector< FROM > getRelationsFrom(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from another store array to this object.
RelationVector< T > getRelationsWith(const std::string &name="", const std::string &namedRelation="") const
Get the relations between this object and another store array.
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Values of the result of a track fit with a given particle hypothesis.
short getChargeSign() const
Return track charge (1 or -1).
double getPValue() const
Getter for Chi2 Probability of the track fit.
double getD0() const
Getter for d0.
double getZ0() const
Getter for z0.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
ROOT::Math::XYZVector getPosition() const
Getter for vector of position at closest approach of track in r/phi projection.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
HitPatternVXD getHitPatternVXD() const
Getter for the hit pattern in the VXD;.
Class that bundles various TrackFitResults.
const TrackFitResult * getTrackFitResultWithClosestMass(const Const::ChargedStable &requestedType) const
Return the track fit for a fit hypothesis with the closest mass.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
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
double phi_gen
azimuthal angle of generated momentum vector
int pdg_gen
PDG code of generated particle.
double pt_gen
generated transverse momentum
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
double z
measured z value of position
double y
measured y value of position
double ptot_gen
generated total momentum
double py
measured momentum in y direction
double phi
azimuthal angle of measured momentum vector
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