Belle II Software  release-06-00-14
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 #include <ecl/modules/eclMatchingPerformance/ECLMatchingPerformanceExpertModule.h>
10 
11 #include <ecl/geometry/ECLGeometryPar.h>
12 #include <framework/datastore/RelationVector.h>
13 
14 #include <root/TFile.h>
15 #include <root/TTree.h>
16 
17 using namespace Belle2;
18 
19 //-----------------------------------------------------------------
20 // Register the Module
21 //-----------------------------------------------------------------
22 REG_MODULE(ECLMatchingPerformanceExpert)
23 
24 ECLMatchingPerformanceExpertModule::ECLMatchingPerformanceExpertModule() :
25  Module(), m_trackProperties()
26 {
27  setDescription("Module to test the matching efficiency between tracks and ECLClusters. Additionally, information about the tracks are written into a ROOT file.");
28  setPropertyFlags(c_ParallelProcessingCertified);
29 
30  addParam("outputFileName", m_outputFileName, "Name of output root file.",
31  std::string("ECLMatchingPerformanceOutput.root"));
32  addParam("minimalCalDigitEnergy", m_minCalDigitEnergy, "minimal energy deposition in crystal to distinguish true signal from noise",
33  0.002);
34  addParam("innerDistanceEnergy", m_innerDistanceEnergy, "determine distance to closest crystal with at least this deposited energy",
35  0.01);
36 }
37 
38 void ECLMatchingPerformanceExpertModule::initialize()
39 {
40  // Required modules
41  m_recoTracks.isRequired();
42  m_tracks.isRequired();
43  m_trackFitResults.isRequired();
44  m_eclClusters.isRequired();
45  m_extHits.isRequired();
46  m_eclCalDigits.isRequired();
47 
51 
52  m_outputFile = new TFile(m_outputFileName.c_str(), "RECREATE");
53  TDirectory* oldDir = gDirectory;
54  m_outputFile->cd();
55  m_dataTree = new TTree("data", "data");
56  oldDir->cd();
57 
58  setupTree();
59 }
60 
62 {
63  m_iEvent = m_EventMetaData->getEvent();
64  m_iRun = m_EventMetaData->getRun();
65  m_iExperiment = m_EventMetaData->getExperiment();
66  m_trackMultiplicity = m_tracks.getEntries();
67 
68  double distance;
69  TVector3 pos_enter, pos_exit, pos_center;
71 
72  for (const Track& track : m_tracks) {
74  const RecoTrack* recoTrack = track.getRelated<RecoTrack>();
75  if (recoTrack) {
76  const TrackFitResult* fitResult = track.getTrackFitResultWithClosestMass(Const::muon);
77  B2ASSERT("Related Belle2 Track has no related track fit result!", fitResult);
78 
79  // write some data to the root tree
80  TVector3 mom = fitResult->getMomentum();
81  m_trackProperties.cosTheta = mom.CosTheta();
82  m_trackProperties.phi = mom.Phi();
83  m_trackProperties.ptot = mom.Mag();
84  m_trackProperties.pt = mom.Pt();
85  m_trackProperties.px = mom.Px();
86  m_trackProperties.py = mom.Py();
87  m_trackProperties.pz = mom.Pz();
88  m_trackProperties.x = fitResult->getPosition().X();
89  m_trackProperties.y = fitResult->getPosition().Y();
90  m_trackProperties.z = fitResult->getPosition().Z();
91 
92  m_pValue = fitResult->getPValue();
93  m_charge = (int)fitResult->getChargeSign();
94  m_d0 = fitResult->getD0();
95  m_z0 = fitResult->getZ0();
96  m_ndf = recoTrack->getTrackFitStatus()->getNdf();
97 
98  // Count hits
102 
103  for (auto& eclCluster : track.getRelationsTo<ECLCluster>()) {
104  if (!eclCluster.hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) continue;
105  if (!(eclCluster.isTrack())) continue;
107  m_hypothesisOfMatchedECLCluster = eclCluster.getHypotheses();
108  break;
109  }
110 
111  bool found_first_enter = false;
112  std::set<int> uniqueECLCalDigits;
113  bool inserted;
114  int thetaID, phiID;
115  for (const auto& extHit : track.getRelationsTo<ExtHit>()) {
116  int copyid = extHit.getCopyID();
117  if (copyid == -1) continue;
118  if (extHit.getDetectorID() != Const::EDetector::ECL) continue;
119  ECLCluster* eclClusterNear = extHit.getRelatedFrom<ECLCluster>();
120  if (eclClusterNear) {
121  distance = (eclClusterNear->getClusterPosition() - extHit.getPosition()).Mag();
122  if (m_distance < 0 || distance < m_distance) {
123  m_distance = distance;
124  }
125  }
126  const int cell = copyid + 1;
127  if (extHit.getStatus() == EXT_ENTER) {
128  m_enter++;
129  inserted = (uniqueECLCalDigits.insert(cell)).second;
130  if (inserted) {
131  for (const auto& eclCalDigit : m_eclCalDigits) {
132  if (eclCalDigit.getCellId() == cell) {
133  if (eclCalDigit.getEnergy() > m_minCalDigitEnergy) {
134  m_deposited_energy += eclCalDigit.getEnergy();
135  }
136  break;
137  }
138  }
139  }
140  if (!found_first_enter) {
141  pos_enter = extHit.getPosition();
142  m_enteringcellid = cell;
143  m_enteringcelltheta = (geometry->GetCrystalPos(cell - 1)).Theta();
144  found_first_enter = true;
145  }
146  //Find ECLCalDigit in cell ID of ExtHit or one of its neighbours
147  if (m_matchedTo1x1Neighbours == 0) {
149  }
150  if (m_matchedTo1x1Neighbours == 1) {
153  }
157  geometry->Mapping(cell - 1);
158  thetaID = geometry->GetThetaID();
159  short int crystalsPerRing = m_eclNeighbours1x1->getCrystalsPerRing(thetaID);
160  phiID = geometry->GetPhiID();
161  const short int phiInc = (((phiID + 1) > (crystalsPerRing - 1)) ? phiID + 1 - crystalsPerRing : phiID + 1);
162  const double fractionalPhiInc = static_cast < double >(phiInc) / crystalsPerRing;
163  const short int phiDec = ((phiID - 1 < 0) ? crystalsPerRing + phiID - 1 : phiID - 1);
164  const double fractionalPhiDec = static_cast < double >(phiDec) / crystalsPerRing;
165  if (m_matchedToDecreasedPhi == 0) {
166  findECLCalDigitMatch(geometry->GetCellID(thetaID , phiDec) + 1, m_matchedToDecreasedPhi);
167  }
168  if (m_matchedToIncreasedPhi == 0) {
169  findECLCalDigitMatch(geometry->GetCellID(thetaID , phiInc) + 1, m_matchedToIncreasedPhi);
170  }
171  if (thetaID < 68) {
172  if (m_matchedToIncreasedTheta == 0) {
173  findECLCalDigitMatch(geometry->GetCellID(thetaID + 1, phiID) + 1, m_matchedToIncreasedTheta);
174  }
176  crystalsPerRing = m_eclNeighbours1x1->getCrystalsPerRing(thetaID + 1);
178  const short int thisPhiInc = std::lround(fractionalPhiInc * crystalsPerRing);
179  findECLCalDigitMatch(geometry->GetCellID(thetaID + 1, thisPhiInc) + 1, m_matchedToIncreasedPhiIncreasedTheta);
180  }
182  short int thisPhiDec = std::lround(fractionalPhiDec * crystalsPerRing);
183  if (thisPhiDec == crystalsPerRing) thisPhiDec = 0;
184  findECLCalDigitMatch(geometry->GetCellID(thetaID + 1, thisPhiDec) + 1, m_matchedToDecreasedPhiIncreasedTheta);
185  }
186  }
187  }
188  if (thetaID > 0) {
189  if (m_matchedToDecreasedTheta == 0) {
190  findECLCalDigitMatch(geometry->GetCellID(thetaID - 1, phiID) + 1, m_matchedToDecreasedTheta);
191  }
193  crystalsPerRing = m_eclNeighbours1x1->getCrystalsPerRing(thetaID - 1);
195  short int thisPhiDec = std::lround(fractionalPhiDec * crystalsPerRing);
196  if (thisPhiDec == crystalsPerRing) thisPhiDec = 0;
197  findECLCalDigitMatch(geometry->GetCellID(thetaID - 1, thisPhiDec) + 1, m_matchedToDecreasedPhiDecreasedTheta);
198  }
200  const short int thisPhiInc = std::lround(fractionalPhiInc * crystalsPerRing);
201  findECLCalDigitMatch(geometry->GetCellID(thetaID - 1, thisPhiInc) + 1, m_matchedToIncreasedPhiDecreasedTheta);
202  }
203  }
204  }
205  }
206  if (m_matchedTo3x3Neighbours == 0) {
208  }
209  if (m_matchedTo3x3Neighbours == 1) {
211  }
212  if (m_matchedTo5x5Neighbours == 0) {
214  }
215  } else if (extHit.getStatus() == EXT_EXIT) {
216  m_exit++;
217  pos_exit = extHit.getPosition();
218  }
219  }
220  m_trackLength = (pos_enter - pos_exit).Mag();
221  for (const auto& eclCalDigit : m_eclCalDigits) {
222  if (eclCalDigit.getEnergy() < m_innerDistanceEnergy) continue;
223  int cellid = eclCalDigit.getCellId();
224  TVector3 cvec = geometry->GetCrystalPos(cellid - 1);
225  distance = (cvec - 0.5 * (pos_enter + pos_exit)).Mag();
226  if (m_innerdistance < 0 || distance < m_innerdistance) {
227  m_innerdistance = distance;
228  }
229  }
230  }
231  m_dataTree->Fill(); // write data to tree
232  }
233 }
234 
235 
237 {
238  writeData();
239  delete m_eclNeighbours1x1;
240  delete m_eclNeighbours3x3;
241  delete m_eclNeighbours5x5;
242 }
243 
245 {
246  if (m_dataTree == nullptr) {
247  B2FATAL("Data tree was not created.");
248  }
250  addVariableToTree("runNo", m_iRun);
251  addVariableToTree("evtNo", m_iEvent);
253 
255 
257 
261 
265 
267 
269 
270  addVariableToTree("pValue", m_pValue);
271 
272  addVariableToTree("charge", m_charge);
273 
274  addVariableToTree("d0", m_d0);
275  addVariableToTree("z0", m_z0);
276 
277  addVariableToTree("ndf", m_ndf);
278 
282 
285 
286  addVariableToTree("MinDistance", m_distance);
287 
288  addVariableToTree("TrackToOneByOneNeighboursMatch", m_matchedTo1x1Neighbours);
289  addVariableToTree("TrackToThreeByThreeNeighboursMatch", m_matchedTo3x3Neighbours);
290  addVariableToTree("TrackToFiveByFiveNeighboursMatch", m_matchedTo5x5Neighbours);
291 
292  addVariableToTree("TrackToDecreasedPhiNeighbourMatch", m_matchedToDecreasedPhi);
293  addVariableToTree("TrackToIncreasedPhiNeighbourMatch", m_matchedToIncreasedPhi);
294  addVariableToTree("TrackToDecreasedThetaNeighbourMatch", m_matchedToDecreasedTheta);
295  addVariableToTree("TrackToIncreasedThetaNeighbourMatch", m_matchedToIncreasedTheta);
296 
297  addVariableToTree("TrackToDecDecNeighbourMatch", m_matchedToDecreasedPhiDecreasedTheta);
298  addVariableToTree("TrackToIncDecNeighbourMatch", m_matchedToIncreasedPhiDecreasedTheta);
299  addVariableToTree("TrackToDecIncNeighbourMatch", m_matchedToDecreasedPhiIncreasedTheta);
300  addVariableToTree("TrackToIncIncNeighbourMatch", m_matchedToIncreasedPhiIncreasedTheta);
301 
302  addVariableToTree("nECLEnter", m_enter);
303  addVariableToTree("nECLExit", m_exit);
304 
306  addVariableToTree("CrystalTheta", m_enteringcelltheta);
307 
308  addVariableToTree("DepositedEnergy", m_deposited_energy);
309  addVariableToTree("TrackLengthInECL", m_trackLength);
310 
311  addVariableToTree("DistanceToXMeV", m_innerdistance);
312 }
313 
315 {
316  if (m_dataTree != nullptr) {
317  TDirectory* oldDir = gDirectory;
318  if (m_outputFile)
319  m_outputFile->cd();
320  m_dataTree->Write();
321  oldDir->cd();
322  delete m_dataTree;
323  }
324  if (m_outputFile != nullptr) {
325  m_outputFile->Close();
326  delete m_outputFile;
327  }
328 }
329 
331 {
332  m_trackProperties = -999;
333 
334  m_pValue = -999;
335 
336  m_charge = 0;
337 
338  m_d0 = -999;
339 
340  m_ndf = -999;
341 
343 
345 
346  m_distance = -999;
347 
351 
356 
361 
362  m_enter = 0;
363  m_exit = 0;
364 
365  m_deposited_energy = 0;
366  m_trackLength = 0;
367 
368  m_innerdistance = -999;
369 
370  m_enteringcellid = -1;
371  m_enteringcelltheta = -999;
372 }
373 
374 void ECLMatchingPerformanceExpertModule::addVariableToTree(const std::string& varName, double& varReference)
375 {
376  std::stringstream leaf;
377  leaf << varName << "/D";
378  m_dataTree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
379 }
380 
381 void ECLMatchingPerformanceExpertModule::addVariableToTree(const std::string& varName, int& varReference)
382 {
383  std::stringstream leaf;
384  leaf << varName << "/I";
385  m_dataTree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
386 }
387 
389  int& matchedToNeighbours, const int& cell)
390 {
391  const auto& vec_of_neighbouring_cells = eclneighbours->getNeighbours(cell);
392  for (const auto& neighbouringcell : vec_of_neighbouring_cells) {
393  const auto idigit = find_if(m_eclCalDigits.begin(), m_eclCalDigits.end(),
394  [&](const ECLCalDigit & d) { return (d.getCellId() == neighbouringcell && d.getEnergy() > m_minCalDigitEnergy); }
395  );
396  //Found ECLCalDigit close to ExtHit
397  if (idigit != m_eclCalDigits.end()) {
398  matchedToNeighbours = 1;
399  break;
400  }
401  }
402 }
403 
405 {
406  const auto idigit = find_if(m_eclCalDigits.begin(), m_eclCalDigits.end(),
407  [&](const ECLCalDigit & d) { return (d.getCellId() == cell && d.getEnergy() > m_minCalDigitEnergy); }
408  );
409  //Found ECLCalDigit close to ExtHit
410  if (idigit != m_eclCalDigits.end()) {
411  matched = 1;
412  }
413 }
static const ChargedStable muon
muon particle
Definition: Const.h:541
Class to store calibrated ECLDigits: ECLCalDigits.
Definition: ECLCalDigit.h:23
ECL cluster data.
Definition: ECLCluster.h:27
TVector3 getClusterPosition() const
Return TVector3 on cluster position (x,y,z)
Definition: ECLCluster.cc:40
@ 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.
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:23
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:37
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:30
Base class for Modules.
Definition: Module.h:72
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
unsigned int getNumberOfSVDHits() const
Return the number of svd hits.
Definition: RecoTrack.h:419
unsigned int getNumberOfCDCHits() const
Return the number of cdc hits.
Definition: RecoTrack.h:416
unsigned int getNumberOfPXDHits() const
Return the number of pxd hits.
Definition: RecoTrack.h:422
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
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.
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
double getD0() const
Getter for d0.
TVector3 getPosition() const
Getter for vector of position at closest approach of track in r/phi projection.
double getZ0() const
Getter for z0.
Class that bundles various TrackFitResults.
Definition: Track.h:25
double getNdf() const
Get the degrees of freedom of the fit.
Definition: FitStatus.h:122
#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