Belle II Software development
ECLMatchingPerformanceExpertModule.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/* Own header. */
10#include <ecl/modules/eclMatchingPerformance/ECLMatchingPerformanceExpertModule.h>
11
12/* ECL headers. */
13#include <ecl/geometry/ECLGeometryPar.h>
14
15/* Basf2 headers. */
16#include <framework/datastore/RelationVector.h>
17
18/* ROOT headers. */
19#include <Math/Vector3D.h>
20#include <TFile.h>
21#include <TTree.h>
22
23/* C++ headers. */
24#include <algorithm>
25
26using namespace Belle2;
27
28//-----------------------------------------------------------------
29// Register the Module
30//-----------------------------------------------------------------
31REG_MODULE(ECLMatchingPerformanceExpert);
32
33ECLMatchingPerformanceExpertModule::ECLMatchingPerformanceExpertModule() :
34 Module(), m_trackProperties()
35{
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);
38
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",
42 0.002);
43 addParam("innerDistanceEnergy", m_innerDistanceEnergy, "determine distance to closest crystal with at least this deposited energy",
44 0.01);
45}
46
48{
49 // Required modules
51 m_tracks.isRequired();
52 m_trackFitResults.isRequired();
53 m_eclClusters.isRequired();
54 m_extHits.isRequired();
55 m_eclCalDigits.isRequired();
56
60
61 m_outputFile = new TFile(m_outputFileName.c_str(), "RECREATE");
62 TDirectory* oldDir = gDirectory;
63 m_outputFile->cd();
64 m_dataTree = new TTree("data", "data");
65 oldDir->cd();
66
67 setupTree();
68}
69
71{
72 m_iEvent = m_EventMetaData->getEvent();
73 m_iRun = m_EventMetaData->getRun();
74 m_iExperiment = m_EventMetaData->getExperiment();
75 m_trackMultiplicity = m_tracks.getEntries();
76
77 double distance;
78 ROOT::Math::XYZVector pos_enter, pos_exit, pos_center;
80
81 for (const Track& track : m_tracks) {
83 const RecoTrack* recoTrack = track.getRelated<RecoTrack>();
84 if (recoTrack) {
85 const TrackFitResult* fitResult = track.getTrackFitResultWithClosestMass(Const::muon);
86 B2ASSERT("Related Belle2 Track has no related track fit result!", fitResult);
87
88 // write some data to the root tree
89 ROOT::Math::XYZVector mom = fitResult->getMomentum();
90 m_trackProperties.cosTheta = cos(mom.Theta());
91 m_trackProperties.phi = mom.Phi();
92 m_trackProperties.ptot = mom.R();
93 m_trackProperties.pt = mom.Rho();
94 m_trackProperties.px = mom.X();
95 m_trackProperties.py = mom.Y();
96 m_trackProperties.pz = mom.Z();
97 m_trackProperties.x = fitResult->getPosition().X();
98 m_trackProperties.y = fitResult->getPosition().Y();
99 m_trackProperties.z = fitResult->getPosition().Z();
100
101 m_pValue = fitResult->getPValue();
102 m_charge = (int)fitResult->getChargeSign();
103 m_d0 = fitResult->getD0();
104 m_z0 = fitResult->getZ0();
105 m_ndf = recoTrack->getTrackFitStatus()->getNdf();
106
107 // Count hits
111
112 for (auto& eclCluster : track.getRelationsTo<ECLCluster>()) {
113 if (!eclCluster.hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) continue;
114 if (!(eclCluster.isTrack())) continue;
116 m_hypothesisOfMatchedECLCluster = eclCluster.getHypotheses();
117 break;
118 }
119
120 bool found_first_enter = false;
121 std::set<int> uniqueECLCalDigits;
122 bool inserted;
123 int thetaID, phiID;
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;
128 ECLCluster* eclClusterNear = extHit.getRelatedFrom<ECLCluster>();
129 if (eclClusterNear) {
130 distance = (eclClusterNear->getClusterPosition() - extHit.getPosition()).R();
131 if (m_distance < 0 || distance < m_distance) {
132 m_distance = distance;
133 }
134 }
135 const int cell = copyid + 1;
136 if (extHit.getStatus() == EXT_ENTER) {
137 m_enter++;
138 inserted = (uniqueECLCalDigits.insert(cell)).second;
139 if (inserted) {
140 for (const auto& eclCalDigit : m_eclCalDigits) {
141 if (eclCalDigit.getCellId() == cell) {
142 if (eclCalDigit.getEnergy() > m_minCalDigitEnergy) {
143 m_deposited_energy += eclCalDigit.getEnergy();
144 }
145 break;
146 }
147 }
148 }
149 if (!found_first_enter) {
150 pos_enter = extHit.getPosition();
151 m_enteringcellid = cell;
152 m_enteringcelltheta = (geometry->GetCrystalPos(cell - 1)).Theta();
153 found_first_enter = true;
154 }
155 //Find ECLCalDigit in cell ID of ExtHit or one of its neighbours
156 if (m_matchedTo1x1Neighbours == 0) {
158 }
159 if (m_matchedTo1x1Neighbours == 1) {
162 }
166 geometry->Mapping(cell - 1);
167 thetaID = geometry->GetThetaID();
168 short int crystalsPerRing = m_eclNeighbours1x1->getCrystalsPerRing(thetaID);
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;
174 if (m_matchedToDecreasedPhi == 0) {
175 findECLCalDigitMatch(geometry->GetCellID(thetaID, phiDec) + 1, m_matchedToDecreasedPhi);
176 }
177 if (m_matchedToIncreasedPhi == 0) {
178 findECLCalDigitMatch(geometry->GetCellID(thetaID, phiInc) + 1, m_matchedToIncreasedPhi);
179 }
180 if (thetaID < 68) {
181 if (m_matchedToIncreasedTheta == 0) {
182 findECLCalDigitMatch(geometry->GetCellID(thetaID + 1, phiID) + 1, m_matchedToIncreasedTheta);
183 }
185 crystalsPerRing = m_eclNeighbours1x1->getCrystalsPerRing(thetaID + 1);
187 const short int thisPhiInc = std::lround(fractionalPhiInc * crystalsPerRing);
188 findECLCalDigitMatch(geometry->GetCellID(thetaID + 1, thisPhiInc) + 1, m_matchedToIncreasedPhiIncreasedTheta);
189 }
191 short int thisPhiDec = std::lround(fractionalPhiDec * crystalsPerRing);
192 if (thisPhiDec == crystalsPerRing) thisPhiDec = 0;
193 findECLCalDigitMatch(geometry->GetCellID(thetaID + 1, thisPhiDec) + 1, m_matchedToDecreasedPhiIncreasedTheta);
194 }
195 }
196 }
197 if (thetaID > 0) {
198 if (m_matchedToDecreasedTheta == 0) {
199 findECLCalDigitMatch(geometry->GetCellID(thetaID - 1, phiID) + 1, m_matchedToDecreasedTheta);
200 }
202 crystalsPerRing = m_eclNeighbours1x1->getCrystalsPerRing(thetaID - 1);
204 short int thisPhiDec = std::lround(fractionalPhiDec * crystalsPerRing);
205 if (thisPhiDec == crystalsPerRing) thisPhiDec = 0;
206 findECLCalDigitMatch(geometry->GetCellID(thetaID - 1, thisPhiDec) + 1, m_matchedToDecreasedPhiDecreasedTheta);
207 }
209 const short int thisPhiInc = std::lround(fractionalPhiInc * crystalsPerRing);
210 findECLCalDigitMatch(geometry->GetCellID(thetaID - 1, thisPhiInc) + 1, m_matchedToIncreasedPhiDecreasedTheta);
211 }
212 }
213 }
214 }
215 if (m_matchedTo3x3Neighbours == 0) {
217 }
218 if (m_matchedTo3x3Neighbours == 1) {
220 }
221 if (m_matchedTo5x5Neighbours == 0) {
223 }
224 } else if (extHit.getStatus() == EXT_EXIT) {
225 m_exit++;
226 pos_exit = extHit.getPosition();
227 }
228 }
229 m_trackLength = (pos_enter - pos_exit).R();
230 for (const auto& eclCalDigit : m_eclCalDigits) {
231 if (eclCalDigit.getEnergy() < m_innerDistanceEnergy) continue;
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();
236 if (m_innerdistance < 0 || distance < m_innerdistance) {
237 m_innerdistance = distance;
238 }
239 }
240 }
241 m_dataTree->Fill(); // write data to tree
242 }
243}
244
245
247{
248 writeData();
249 delete m_eclNeighbours1x1;
250 delete m_eclNeighbours3x3;
251 delete m_eclNeighbours5x5;
252}
253
255{
256 if (m_dataTree == nullptr) {
257 B2FATAL("Data tree was not created.");
258 }
260 addVariableToTree("runNo", m_iRun);
261 addVariableToTree("evtNo", m_iEvent);
263
265
267
271
275
277
279
280 addVariableToTree("pValue", m_pValue);
281
282 addVariableToTree("charge", m_charge);
283
284 addVariableToTree("d0", m_d0);
285 addVariableToTree("z0", m_z0);
286
287 addVariableToTree("ndf", m_ndf);
288
292
295
296 addVariableToTree("MinDistance", m_distance);
297
298 addVariableToTree("TrackToOneByOneNeighboursMatch", m_matchedTo1x1Neighbours);
299 addVariableToTree("TrackToThreeByThreeNeighboursMatch", m_matchedTo3x3Neighbours);
300 addVariableToTree("TrackToFiveByFiveNeighboursMatch", m_matchedTo5x5Neighbours);
301
302 addVariableToTree("TrackToDecreasedPhiNeighbourMatch", m_matchedToDecreasedPhi);
303 addVariableToTree("TrackToIncreasedPhiNeighbourMatch", m_matchedToIncreasedPhi);
304 addVariableToTree("TrackToDecreasedThetaNeighbourMatch", m_matchedToDecreasedTheta);
305 addVariableToTree("TrackToIncreasedThetaNeighbourMatch", m_matchedToIncreasedTheta);
306
307 addVariableToTree("TrackToDecDecNeighbourMatch", m_matchedToDecreasedPhiDecreasedTheta);
308 addVariableToTree("TrackToIncDecNeighbourMatch", m_matchedToIncreasedPhiDecreasedTheta);
309 addVariableToTree("TrackToDecIncNeighbourMatch", m_matchedToDecreasedPhiIncreasedTheta);
310 addVariableToTree("TrackToIncIncNeighbourMatch", m_matchedToIncreasedPhiIncreasedTheta);
311
312 addVariableToTree("nECLEnter", m_enter);
313 addVariableToTree("nECLExit", m_exit);
314
316 addVariableToTree("CrystalTheta", m_enteringcelltheta);
317
318 addVariableToTree("DepositedEnergy", m_deposited_energy);
319 addVariableToTree("TrackLengthInECL", m_trackLength);
320
321 addVariableToTree("DistanceToXMeV", m_innerdistance);
322}
323
325{
326 if (m_dataTree != nullptr) {
327 TDirectory* oldDir = gDirectory;
328 if (m_outputFile)
329 m_outputFile->cd();
330 m_dataTree->Write();
331 oldDir->cd();
332 delete m_dataTree;
333 }
334 if (m_outputFile != nullptr) {
335 m_outputFile->Close();
336 delete m_outputFile;
337 }
338}
339
341{
342 m_trackProperties = -999;
343
344 m_pValue = -999;
345
346 m_charge = 0;
347
348 m_d0 = -999;
349
350 m_ndf = -999;
351
353
355
356 m_distance = -999;
357
361
366
371
372 m_enter = 0;
373 m_exit = 0;
374
376 m_trackLength = 0;
377
378 m_innerdistance = -999;
379
380 m_enteringcellid = -1;
381 m_enteringcelltheta = -999;
382}
383
384void ECLMatchingPerformanceExpertModule::addVariableToTree(const std::string& varName, double& varReference)
385{
386 std::stringstream leaf;
387 leaf << varName << "/D";
388 m_dataTree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
389}
390
391void ECLMatchingPerformanceExpertModule::addVariableToTree(const std::string& varName, int& varReference)
392{
393 std::stringstream leaf;
394 leaf << varName << "/I";
395 m_dataTree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
396}
397
399 int& matchedToNeighbours, const int& cell)
400{
401 const auto& vec_of_neighbouring_cells = eclneighbours->getNeighbours(cell);
402 for (const auto& neighbouringcell : vec_of_neighbouring_cells) {
403 const auto idigit = std::find_if(m_eclCalDigits.begin(), m_eclCalDigits.end(),
404 [&](const ECLCalDigit & d) { return (d.getCellId() == neighbouringcell && d.getEnergy() > m_minCalDigitEnergy); }
405 );
406 //Found ECLCalDigit close to ExtHit
407 if (idigit != m_eclCalDigits.end()) {
408 matchedToNeighbours = 1;
409 break;
410 }
411 }
412}
413
415{
416 const auto idigit = std::find_if(m_eclCalDigits.begin(), m_eclCalDigits.end(),
417 [&](const ECLCalDigit & d) { return (d.getCellId() == cell && d.getEnergy() > m_minCalDigitEnergy); }
418 );
419 //Found ECLCalDigit close to ExtHit
420 if (idigit != m_eclCalDigits.end()) {
421 matched = 1;
422 }
423}
double R
typedef autogenerated by FFTW
static const ChargedStable muon
muon particle
Definition: Const.h:660
Class to store calibrated ECLDigits: ECLCalDigits.
Definition: ECLCalDigit.h:23
ECL cluster data.
Definition: ECLCluster.h:27
ROOT::Math::XYZVector getClusterPosition() const
Return ROOT::Math::XYZVector on cluster position (x,y,z)
Definition: ECLCluster.cc:44
@ c_nPhotons
CR is split into n photons (N1)
double m_minCalDigitEnergy
noise threshold of deposited energy in crystal
ECL::ECLNeighbours * m_eclNeighbours3x3
Neighbour map of 9 crystals.
double m_d0
signed distance of the track to the IP in the r-phi plane
int m_matchedToDecreasedPhiIncreasedTheta
boolean for match between track and neighbouring cell with lower phi value
ECL::ECLNeighbours * m_eclNeighbours1x1
Neighbour map of 1 crystal.
void initialize() override
Register the needed StoreArrays and open the output TFile.
int m_matchedToIncreasedPhiDecreasedTheta
boolean for match between track and neighbouring cell with higher phi value
int m_matchedToECLCluster
boolean for match between track and ECL cluster
void event() override
Fill the tree with the event data.
void findECLCalDigitMatchInNeighbouringCell(ECL::ECLNeighbours *eclneighbours, int &matchedToNeighbours, const int &cell)
find a match between crystals in which energy was deposited and the cell or its neighbors that a trac...
int m_matchedToIncreasedTheta
boolean for match between track and neighbouring cell with higher phi value
int m_matchedToIncreasedPhi
boolean for match between track and neighbouring cell with higher phi value
double m_enteringcelltheta
theta of first crystal that is entered by track
StoreArray< TrackFitResult > m_trackFitResults
Required input array of TrackFitResults.
void terminate() override
Write the tree into the opened root file.
double m_distance
minimal distance between track and ECLCluster
ECL::ECLNeighbours * m_eclNeighbours5x5
Neighbour map of 25 crystals.
int m_matchedToDecreasedPhiDecreasedTheta
boolean for match between track and neighbouring cell with lower phi value
ParticleProperties m_trackProperties
properties of a reconstructed track
void setVariablesToDefaultValue()
Sets all variables to the default value, here -999.
double m_innerdistance
minimal distance between track at center of ECL and ECLCalDigit with at least 10 MeV
void findECLCalDigitMatch(const int &cell, int &matched)
determine whether energy has been deposited in crystal with ID cell and write result to matched
int m_ndf
number of degrees of freedom of the track (should be #CDC hits - 5 (helix parameters))
StoreArray< Track > m_tracks
Required input array of Tracks.
int m_matchedTo5x5Neighbours
boolean for match between track and one of 25 ECLCalDigit neighbouring cells
StoreArray< ECLCluster > m_eclClusters
Required input array of ECLClusters.
void addVariableToTree(const std::string &varName, double &varReference)
add a variable with double format
int m_matchedTo3x3Neighbours
boolean for match between track and one of 9 ECLCalDigit neighbouring cells
double m_innerDistanceEnergy
minimal deposited energy for search of closest crystal
int m_hypothesisOfMatchedECLCluster
hypothesis of matched ECL cluster
StoreArray< RecoTrack > m_recoTracks
Required input array of RecoTracks.
StoreArray< ExtHit > m_extHits
Required input array of ExtHits.
int m_matchedToDecreasedPhi
boolean for match between track and neighbouring cell with lower phi value
int m_matchedToDecreasedTheta
boolean for match between track and neighbouring cell with lower phi value
int m_matchedToIncreasedPhiIncreasedTheta
boolean for match between track and neighbouring cell with higher phi value
double m_z0
distance of the track to the IP along the beam axis
void writeData()
write root tree to output file and close the file
int m_enteringcellid
cell id of first crystal that is entered by track
int m_matchedTo1x1Neighbours
boolean for match between track and ECLCalDigit cell
StoreObjPtr< EventMetaData > m_EventMetaData
Event metadata.
StoreArray< ECLCalDigit > m_eclCalDigits
Required input array of ECLCalDigits.
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.
Definition: ECLNeighbours.h:25
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
Definition: ECLNeighbours.h:39
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:31
Base class for Modules.
Definition: Module.h:72
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
unsigned int getNumberOfSVDHits() const
Return the number of svd hits.
Definition: RecoTrack.h:424
unsigned int getNumberOfCDCHits() const
Return the number of cdc hits.
Definition: RecoTrack.h:427
unsigned int getNumberOfPXDHits() const
Return the number of pxd hits.
Definition: RecoTrack.h:421
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:621
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.
Definition: Track.h:25
#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
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