Belle II Software  release-08-01-10
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 using namespace Belle2;
24 
25 //-----------------------------------------------------------------
26 // Register the Module
27 //-----------------------------------------------------------------
28 REG_MODULE(ECLMatchingPerformanceExpert);
29 
30 ECLMatchingPerformanceExpertModule::ECLMatchingPerformanceExpertModule() :
31  Module(), m_trackProperties()
32 {
33  setDescription("Module to test the matching efficiency between tracks and ECLClusters. Additionally, information about the tracks are written into a ROOT file.");
34  setPropertyFlags(c_ParallelProcessingCertified);
35 
36  addParam("outputFileName", m_outputFileName, "Name of output root file.",
37  std::string("ECLMatchingPerformanceOutput.root"));
38  addParam("minimalCalDigitEnergy", m_minCalDigitEnergy, "minimal energy deposition in crystal to distinguish true signal from noise",
39  0.002);
40  addParam("innerDistanceEnergy", m_innerDistanceEnergy, "determine distance to closest crystal with at least this deposited energy",
41  0.01);
42 }
43 
44 void ECLMatchingPerformanceExpertModule::initialize()
45 {
46  // Required modules
48  m_tracks.isRequired();
49  m_trackFitResults.isRequired();
50  m_eclClusters.isRequired();
51  m_extHits.isRequired();
52  m_eclCalDigits.isRequired();
53 
57 
58  m_outputFile = new TFile(m_outputFileName.c_str(), "RECREATE");
59  TDirectory* oldDir = gDirectory;
60  m_outputFile->cd();
61  m_dataTree = new TTree("data", "data");
62  oldDir->cd();
63 
64  setupTree();
65 }
66 
68 {
69  m_iEvent = m_EventMetaData->getEvent();
70  m_iRun = m_EventMetaData->getRun();
71  m_iExperiment = m_EventMetaData->getExperiment();
72  m_trackMultiplicity = m_tracks.getEntries();
73 
74  double distance;
75  ROOT::Math::XYZVector pos_enter, pos_exit, pos_center;
77 
78  for (const Track& track : m_tracks) {
80  const RecoTrack* recoTrack = track.getRelated<RecoTrack>();
81  if (recoTrack) {
82  const TrackFitResult* fitResult = track.getTrackFitResultWithClosestMass(Const::muon);
83  B2ASSERT("Related Belle2 Track has no related track fit result!", fitResult);
84 
85  // write some data to the root tree
86  ROOT::Math::XYZVector mom = fitResult->getMomentum();
87  m_trackProperties.cosTheta = cos(mom.Theta());
88  m_trackProperties.phi = mom.Phi();
89  m_trackProperties.ptot = mom.R();
90  m_trackProperties.pt = mom.Rho();
91  m_trackProperties.px = mom.X();
92  m_trackProperties.py = mom.Y();
93  m_trackProperties.pz = mom.Z();
94  m_trackProperties.x = fitResult->getPosition().X();
95  m_trackProperties.y = fitResult->getPosition().Y();
96  m_trackProperties.z = fitResult->getPosition().Z();
97 
98  m_pValue = fitResult->getPValue();
99  m_charge = (int)fitResult->getChargeSign();
100  m_d0 = fitResult->getD0();
101  m_z0 = fitResult->getZ0();
102  m_ndf = recoTrack->getTrackFitStatus()->getNdf();
103 
104  // Count hits
108 
109  for (auto& eclCluster : track.getRelationsTo<ECLCluster>()) {
110  if (!eclCluster.hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) continue;
111  if (!(eclCluster.isTrack())) continue;
113  m_hypothesisOfMatchedECLCluster = eclCluster.getHypotheses();
114  break;
115  }
116 
117  bool found_first_enter = false;
118  std::set<int> uniqueECLCalDigits;
119  bool inserted;
120  int thetaID, phiID;
121  for (const auto& extHit : track.getRelationsTo<ExtHit>()) {
122  int copyid = extHit.getCopyID();
123  if (copyid == -1) continue;
124  if (extHit.getDetectorID() != Const::EDetector::ECL) continue;
125  ECLCluster* eclClusterNear = extHit.getRelatedFrom<ECLCluster>();
126  if (eclClusterNear) {
127  distance = (eclClusterNear->getClusterPosition() - extHit.getPosition()).R();
128  if (m_distance < 0 || distance < m_distance) {
129  m_distance = distance;
130  }
131  }
132  const int cell = copyid + 1;
133  if (extHit.getStatus() == EXT_ENTER) {
134  m_enter++;
135  inserted = (uniqueECLCalDigits.insert(cell)).second;
136  if (inserted) {
137  for (const auto& eclCalDigit : m_eclCalDigits) {
138  if (eclCalDigit.getCellId() == cell) {
139  if (eclCalDigit.getEnergy() > m_minCalDigitEnergy) {
140  m_deposited_energy += eclCalDigit.getEnergy();
141  }
142  break;
143  }
144  }
145  }
146  if (!found_first_enter) {
147  pos_enter = extHit.getPosition();
148  m_enteringcellid = cell;
149  m_enteringcelltheta = (geometry->GetCrystalPos(cell - 1)).Theta();
150  found_first_enter = true;
151  }
152  //Find ECLCalDigit in cell ID of ExtHit or one of its neighbours
153  if (m_matchedTo1x1Neighbours == 0) {
155  }
156  if (m_matchedTo1x1Neighbours == 1) {
159  }
163  geometry->Mapping(cell - 1);
164  thetaID = geometry->GetThetaID();
165  short int crystalsPerRing = m_eclNeighbours1x1->getCrystalsPerRing(thetaID);
166  phiID = geometry->GetPhiID();
167  const short int phiInc = (((phiID + 1) > (crystalsPerRing - 1)) ? phiID + 1 - crystalsPerRing : phiID + 1);
168  const double fractionalPhiInc = static_cast < double >(phiInc) / crystalsPerRing;
169  const short int phiDec = ((phiID - 1 < 0) ? crystalsPerRing + phiID - 1 : phiID - 1);
170  const double fractionalPhiDec = static_cast < double >(phiDec) / crystalsPerRing;
171  if (m_matchedToDecreasedPhi == 0) {
172  findECLCalDigitMatch(geometry->GetCellID(thetaID, phiDec) + 1, m_matchedToDecreasedPhi);
173  }
174  if (m_matchedToIncreasedPhi == 0) {
175  findECLCalDigitMatch(geometry->GetCellID(thetaID, phiInc) + 1, m_matchedToIncreasedPhi);
176  }
177  if (thetaID < 68) {
178  if (m_matchedToIncreasedTheta == 0) {
179  findECLCalDigitMatch(geometry->GetCellID(thetaID + 1, phiID) + 1, m_matchedToIncreasedTheta);
180  }
182  crystalsPerRing = m_eclNeighbours1x1->getCrystalsPerRing(thetaID + 1);
184  const short int thisPhiInc = std::lround(fractionalPhiInc * crystalsPerRing);
185  findECLCalDigitMatch(geometry->GetCellID(thetaID + 1, thisPhiInc) + 1, m_matchedToIncreasedPhiIncreasedTheta);
186  }
188  short int thisPhiDec = std::lround(fractionalPhiDec * crystalsPerRing);
189  if (thisPhiDec == crystalsPerRing) thisPhiDec = 0;
190  findECLCalDigitMatch(geometry->GetCellID(thetaID + 1, thisPhiDec) + 1, m_matchedToDecreasedPhiIncreasedTheta);
191  }
192  }
193  }
194  if (thetaID > 0) {
195  if (m_matchedToDecreasedTheta == 0) {
196  findECLCalDigitMatch(geometry->GetCellID(thetaID - 1, phiID) + 1, m_matchedToDecreasedTheta);
197  }
199  crystalsPerRing = m_eclNeighbours1x1->getCrystalsPerRing(thetaID - 1);
201  short int thisPhiDec = std::lround(fractionalPhiDec * crystalsPerRing);
202  if (thisPhiDec == crystalsPerRing) thisPhiDec = 0;
203  findECLCalDigitMatch(geometry->GetCellID(thetaID - 1, thisPhiDec) + 1, m_matchedToDecreasedPhiDecreasedTheta);
204  }
206  const short int thisPhiInc = std::lround(fractionalPhiInc * crystalsPerRing);
207  findECLCalDigitMatch(geometry->GetCellID(thetaID - 1, thisPhiInc) + 1, m_matchedToIncreasedPhiDecreasedTheta);
208  }
209  }
210  }
211  }
212  if (m_matchedTo3x3Neighbours == 0) {
214  }
215  if (m_matchedTo3x3Neighbours == 1) {
217  }
218  if (m_matchedTo5x5Neighbours == 0) {
220  }
221  } else if (extHit.getStatus() == EXT_EXIT) {
222  m_exit++;
223  pos_exit = extHit.getPosition();
224  }
225  }
226  m_trackLength = (pos_enter - pos_exit).R();
227  for (const auto& eclCalDigit : m_eclCalDigits) {
228  if (eclCalDigit.getEnergy() < m_innerDistanceEnergy) continue;
229  int cellid = eclCalDigit.getCellId();
230  ROOT::Math::XYZVector cvec = geometry->GetCrystalPos(cellid - 1);
231  ROOT::Math::XYZVector crystalPosition(cvec.X(), cvec.Y(), cvec.Z());
232  distance = (crystalPosition - 0.5 * (pos_enter + pos_exit)).R();
233  if (m_innerdistance < 0 || distance < m_innerdistance) {
234  m_innerdistance = distance;
235  }
236  }
237  }
238  m_dataTree->Fill(); // write data to tree
239  }
240 }
241 
242 
244 {
245  writeData();
246  delete m_eclNeighbours1x1;
247  delete m_eclNeighbours3x3;
248  delete m_eclNeighbours5x5;
249 }
250 
252 {
253  if (m_dataTree == nullptr) {
254  B2FATAL("Data tree was not created.");
255  }
257  addVariableToTree("runNo", m_iRun);
258  addVariableToTree("evtNo", m_iEvent);
260 
262 
264 
268 
272 
274 
276 
277  addVariableToTree("pValue", m_pValue);
278 
279  addVariableToTree("charge", m_charge);
280 
281  addVariableToTree("d0", m_d0);
282  addVariableToTree("z0", m_z0);
283 
284  addVariableToTree("ndf", m_ndf);
285 
289 
292 
293  addVariableToTree("MinDistance", m_distance);
294 
295  addVariableToTree("TrackToOneByOneNeighboursMatch", m_matchedTo1x1Neighbours);
296  addVariableToTree("TrackToThreeByThreeNeighboursMatch", m_matchedTo3x3Neighbours);
297  addVariableToTree("TrackToFiveByFiveNeighboursMatch", m_matchedTo5x5Neighbours);
298 
299  addVariableToTree("TrackToDecreasedPhiNeighbourMatch", m_matchedToDecreasedPhi);
300  addVariableToTree("TrackToIncreasedPhiNeighbourMatch", m_matchedToIncreasedPhi);
301  addVariableToTree("TrackToDecreasedThetaNeighbourMatch", m_matchedToDecreasedTheta);
302  addVariableToTree("TrackToIncreasedThetaNeighbourMatch", m_matchedToIncreasedTheta);
303 
304  addVariableToTree("TrackToDecDecNeighbourMatch", m_matchedToDecreasedPhiDecreasedTheta);
305  addVariableToTree("TrackToIncDecNeighbourMatch", m_matchedToIncreasedPhiDecreasedTheta);
306  addVariableToTree("TrackToDecIncNeighbourMatch", m_matchedToDecreasedPhiIncreasedTheta);
307  addVariableToTree("TrackToIncIncNeighbourMatch", m_matchedToIncreasedPhiIncreasedTheta);
308 
309  addVariableToTree("nECLEnter", m_enter);
310  addVariableToTree("nECLExit", m_exit);
311 
313  addVariableToTree("CrystalTheta", m_enteringcelltheta);
314 
315  addVariableToTree("DepositedEnergy", m_deposited_energy);
316  addVariableToTree("TrackLengthInECL", m_trackLength);
317 
318  addVariableToTree("DistanceToXMeV", m_innerdistance);
319 }
320 
322 {
323  if (m_dataTree != nullptr) {
324  TDirectory* oldDir = gDirectory;
325  if (m_outputFile)
326  m_outputFile->cd();
327  m_dataTree->Write();
328  oldDir->cd();
329  delete m_dataTree;
330  }
331  if (m_outputFile != nullptr) {
332  m_outputFile->Close();
333  delete m_outputFile;
334  }
335 }
336 
338 {
339  m_trackProperties = -999;
340 
341  m_pValue = -999;
342 
343  m_charge = 0;
344 
345  m_d0 = -999;
346 
347  m_ndf = -999;
348 
350 
352 
353  m_distance = -999;
354 
358 
363 
368 
369  m_enter = 0;
370  m_exit = 0;
371 
372  m_deposited_energy = 0;
373  m_trackLength = 0;
374 
375  m_innerdistance = -999;
376 
377  m_enteringcellid = -1;
378  m_enteringcelltheta = -999;
379 }
380 
381 void ECLMatchingPerformanceExpertModule::addVariableToTree(const std::string& varName, double& varReference)
382 {
383  std::stringstream leaf;
384  leaf << varName << "/D";
385  m_dataTree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
386 }
387 
388 void ECLMatchingPerformanceExpertModule::addVariableToTree(const std::string& varName, int& varReference)
389 {
390  std::stringstream leaf;
391  leaf << varName << "/I";
392  m_dataTree->Branch(varName.c_str(), &varReference, leaf.str().c_str());
393 }
394 
396  int& matchedToNeighbours, const int& cell)
397 {
398  const auto& vec_of_neighbouring_cells = eclneighbours->getNeighbours(cell);
399  for (const auto& neighbouringcell : vec_of_neighbouring_cells) {
400  const auto idigit = find_if(m_eclCalDigits.begin(), m_eclCalDigits.end(),
401  [&](const ECLCalDigit & d) { return (d.getCellId() == neighbouringcell && d.getEnergy() > m_minCalDigitEnergy); }
402  );
403  //Found ECLCalDigit close to ExtHit
404  if (idigit != m_eclCalDigits.end()) {
405  matchedToNeighbours = 1;
406  break;
407  }
408  }
409 }
410 
412 {
413  const auto idigit = find_if(m_eclCalDigits.begin(), m_eclCalDigits.end(),
414  [&](const ECLCalDigit & d) { return (d.getCellId() == cell && d.getEnergy() > m_minCalDigitEnergy); }
415  );
416  //Found ECLCalDigit close to ExtHit
417  if (idigit != m_eclCalDigits.end()) {
418  matched = 1;
419  }
420 }
double R
typedef autogenerated by FFTW
static const ChargedStable muon
muon particle
Definition: Const.h:651
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.
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:32
Base class for Modules.
Definition: Module.h:72
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
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
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
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
double getNdf() const
Get the degrees of freedom of the fit.
Definition: FitStatus.h:122
REG_MODULE(arichBtest)
Register the Module.
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