10#include <ecl/modules/eclMatchingPerformance/ECLMatchingPerformanceExpertModule.h>
13#include <ecl/geometry/ECLGeometryPar.h>
16#include <framework/datastore/RelationVector.h>
19#include <Math/Vector3D.h>
33ECLMatchingPerformanceExpertModule::ECLMatchingPerformanceExpertModule() :
34 Module(), m_trackProperties()
36 setDescription(
"Module to test the matching efficiency between tracks and ECLClusters. Additionally, information about the tracks are written into a ROOT file.");
37 setPropertyFlags(c_ParallelProcessingCertified);
39 addParam(
"outputFileName", m_outputFileName,
"Name of output root file.",
40 std::string(
"ECLMatchingPerformanceOutput.root"));
41 addParam(
"minimalCalDigitEnergy", m_minCalDigitEnergy,
"minimal energy deposition in crystal to distinguish true signal from noise",
43 addParam(
"innerDistanceEnergy", m_innerDistanceEnergy,
"determine distance to closest crystal with at least this deposited energy",
62 TDirectory* oldDir = gDirectory;
78 ROOT::Math::XYZVector pos_enter, pos_exit, pos_center;
86 B2ASSERT(
"Related Belle2 Track has no related track fit result!", fitResult);
89 ROOT::Math::XYZVector mom = fitResult->
getMomentum();
112 for (
auto& eclCluster : track.getRelationsTo<
ECLCluster>()) {
114 if (!(eclCluster.isTrack()))
continue;
120 bool found_first_enter =
false;
121 std::set<int> uniqueECLCalDigits;
124 for (
const auto& extHit : track.getRelationsTo<
ExtHit>()) {
125 int copyid = extHit.getCopyID();
126 if (copyid == -1)
continue;
127 if (extHit.getDetectorID() != Const::EDetector::ECL)
continue;
129 if (eclClusterNear) {
135 const int cell = copyid + 1;
136 if (extHit.getStatus() == EXT_ENTER) {
138 inserted = (uniqueECLCalDigits.insert(cell)).second;
141 if (eclCalDigit.getCellId() == cell) {
149 if (!found_first_enter) {
150 pos_enter = extHit.getPosition();
153 found_first_enter =
true;
166 geometry->Mapping(cell - 1);
167 thetaID = geometry->GetThetaID();
169 phiID = geometry->GetPhiID();
170 const short int phiInc = (((phiID + 1) > (crystalsPerRing - 1)) ? phiID + 1 - crystalsPerRing : phiID + 1);
171 const double fractionalPhiInc =
static_cast < double >(phiInc) / crystalsPerRing;
172 const short int phiDec = ((phiID - 1 < 0) ? crystalsPerRing + phiID - 1 : phiID - 1);
173 const double fractionalPhiDec =
static_cast < double >(phiDec) / crystalsPerRing;
187 const short int thisPhiInc = std::lround(fractionalPhiInc * crystalsPerRing);
191 short int thisPhiDec = std::lround(fractionalPhiDec * crystalsPerRing);
192 if (thisPhiDec == crystalsPerRing) thisPhiDec = 0;
204 short int thisPhiDec = std::lround(fractionalPhiDec * crystalsPerRing);
205 if (thisPhiDec == crystalsPerRing) thisPhiDec = 0;
209 const short int thisPhiInc = std::lround(fractionalPhiInc * crystalsPerRing);
224 }
else if (extHit.getStatus() == EXT_EXIT) {
226 pos_exit = extHit.getPosition();
232 int cellid = eclCalDigit.getCellId();
233 ROOT::Math::XYZVector cvec = geometry->GetCrystalPos(cellid - 1);
234 ROOT::Math::XYZVector crystalPosition(cvec.X(), cvec.Y(), cvec.Z());
235 distance = (crystalPosition - 0.5 * (pos_enter + pos_exit)).R();
257 B2FATAL(
"Data tree was not created.");
327 TDirectory* oldDir = gDirectory;
386 std::stringstream leaf;
387 leaf << varName <<
"/D";
388 m_dataTree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
393 std::stringstream leaf;
394 leaf << varName <<
"/I";
395 m_dataTree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
399 int& matchedToNeighbours,
const int& cell)
401 const auto& vec_of_neighbouring_cells = eclneighbours->
getNeighbours(cell);
402 for (
const auto& neighbouringcell : vec_of_neighbouring_cells) {
404 [&](
const ECLCalDigit & d) { return (d.getCellId() == neighbouringcell && d.getEnergy() > m_minCalDigitEnergy); }
408 matchedToNeighbours = 1;
417 [&](
const ECLCalDigit & d) { return (d.getCellId() == cell && d.getEnergy() > m_minCalDigitEnergy); }
static const ChargedStable muon
muon particle
Class to store calibrated ECLDigits: ECLCalDigits.
ROOT::Math::XYZVector getClusterPosition() const
Return ROOT::Math::XYZVector on cluster position (x,y,z)
@ c_nPhotons
CR is split into n photons (N1)
The Class for ECL Geometry Parameters.
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
Class to get the neighbours for a given cell id.
const std::vector< short int > & getNeighbours(short int cid) const
Return the neighbours for a given cell ID.
short int getCrystalsPerRing(const short int thetaid) const
return number of crystals in a given theta ring
Store one Ext hit as a ROOT object.
This is the Reconstruction Event-Data Model Track.
unsigned int getNumberOfSVDHits() const
Return the number of svd hits.
unsigned int getNumberOfCDCHits() const
Return the number of cdc hits.
unsigned int getNumberOfPXDHits() const
Return the number of pxd hits.
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...
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
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.
Class that bundles various TrackFitResults.
#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
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 py
measured momentum in y direction
double phi
azimuthal angle of measured momentum vector
double cosTheta
polar angle of measured momentum vector
double x
measured x value of position
double ptot
measured total momentum