Belle II Software  release-08-01-10
harvester.py
1 
8 
9 from tracking.harvest.harvesting import HarvestingModule
10 from tracking.harvest import refiners
11 from ROOT import Belle2
12 from ROOT.Math import XYZVector
13 import numpy as np
14 
15 from tracking.ipython_tools.wrapper import QueueHarvester
16 
17 
19 
20  """ The base class with all static methods to use. """
21 
22  @staticmethod
23  def do_for_each_hit_type(cluster, svd_function, pxd_function):
24  """Apply a PXD selection to a PXDCluster or an SVD selection to an SVDCluster"""
25  cluster_type = cluster.__class__.__name__
26  if cluster_type == "Belle2::PXDCluster":
27  return pxd_function(cluster)
28  elif cluster_type == "Belle2::SVDCluster":
29  return svd_function(cluster)
30  else:
31  raise TypeError("Unknown hit type")
32 
33  @staticmethod
35  """Return lists of charges and path lengths for the clusters associated with an MCParticle"""
36  charge_list = []
37  path_length_list = []
38 
39  tools = Belle2.VXDMomentumEstimationTools(clusters[0].__class__.__name__).getInstance()
40 
41  momentum = mc_particle.getMomentum()
42  position = mc_particle.getProductionVertex()
43  charge = mc_particle.getCharge()
44  b_field = Belle2.BFieldManager.getField(position).Z() / Belle2.Unit.T
45  helix = Belle2.Helix(position, momentum, charge, b_field)
46 
47  for cluster in clusters:
48 
49  calibrated_charge = tools.getCalibratedCharge(cluster)
50  path_length = tools.getPathLength(cluster, helix)
51 
52  charge_list.append(calibrated_charge)
53  path_length_list.append(path_length)
54 
55  return charge_list, path_length_list, list(np.divide(charge_list, path_length_list))
56 
57  @staticmethod
58  def generate_truncated(charge_list):
59  """Sort then truncate a list to all but the last 2 entries or the first 4 entries or the first 6 entries"""
60  sorted_list = sorted(charge_list)
61  if len(sorted_list) > 2:
62  return sorted_list[:-2], sorted_list[:4], sorted_list[:6]
63  else:
64  return sorted_list, sorted_list[:4], sorted_list[:6]
65 
66 
68 
69  """ A harvester to check for the positions of the track points in the MCParticleTrajectories"""
70 
71  def __init__(self):
72  """Constructor"""
73  HarvestingModule.__init__(self, foreach="MCParticleTrajectorys", output_file_name="mc_trajectory.root")
74 
75  def peel(self, mc_particle_trajectory):
76  """Aggregate track-point position information for the trajectory"""
77  for track_point in mc_particle_trajectory:
78  yield {"x": track_point.x, "y": track_point.y, "z": track_point.z, "index": self.counter}
79 
80 
82  save_tree = refiners.SaveTreeRefiner()
83 
84 
86 
87  """ A harvester to redo parts of the analysis in the Belle II Paper by Robert """
88 
89  def __init__(self):
90  """Constructor"""
91  HarvestingModule.__init__(self, foreach="MCParticles", output_file_name="mc_particle.root")
92 
93  def pick(self, mc_particle):
94  """Select the MCParticle if it is a primary pion and has some PXD and/or SVD clusters"""
95  pxd_clusters = mc_particle.getRelationsFrom("PXDClusters")
96  svd_clusters = mc_particle.getRelationsFrom("SVDClusters")
97  return (mc_particle.hasStatus(Belle2.MCParticle.c_PrimaryParticle) and
98  abs(mc_particle.getPDG()) == 211 and
99  len(pxd_clusters) + len(svd_clusters) > 0)
100 
101  def generate_cluster_dicts(self, charge_list, path_length_list, normalized_charge_list, name):
102  """Create a dictionary from the lists of charges, normalized charges, and path lengths of the
103  clusters associated with an MCParticle"""
104  result = dict()
105 
106  truncated, first4, first6 = VXDMomentumEnergyEstimator.generate_truncated(normalized_charge_list)
107 
108  result.update({"sum_%s_charges" % name: sum(charge_list)})
109  result.update({"mean_%s_charges" % name: np.mean(charge_list)})
110  result.update({"sum_%s_normalized_charges" % name: sum(normalized_charge_list)})
111  result.update({"mean_%s_normalized_charges" % name: np.mean(normalized_charge_list)})
112 
113  result.update({"sum_%s_normalized_charges_truncated" % name: sum(truncated)})
114  result.update({"mean_%s_normalized_charges_truncated" % name: np.mean(truncated)})
115  result.update({"sum_%s_normalized_charges_first4" % name: sum(first4)})
116  result.update({"mean_%s_normalized_charges_first4" % name: np.mean(first4)})
117  result.update({"sum_%s_normalized_charges_first6" % name: sum(first6)})
118  result.update({"mean_%s_normalized_charges_first6" % name: np.mean(first6)})
119 
120  return result
121 
122  def peel(self, mc_particle):
123  """Aggregate the PXD and SVD cluster information for an MCParticle"""
124  result = dict()
125 
126  pxd_clusters = mc_particle.getRelationsFrom("PXDClusters")
127  svd_clusters = mc_particle.getRelationsFrom("SVDClusters")
128 
129  pxd_results = VXDMomentumEnergyEstimator.calculate_charges_and_path_lengths_for_one_type(pxd_clusters, mc_particle)
130  svd_results = VXDMomentumEnergyEstimator.calculate_charges_and_path_lengths_for_one_type(svd_clusters, mc_particle)
131 
132  pxd_cluster_dicts = self.generate_cluster_dictsgenerate_cluster_dicts(*pxd_results, name="pxd")
133  pxd_charges, pxd_path_length, pxd_normalized_charges = pxd_results
134 
135  svd_charges, svd_path_length, svd_normalized_charges = svd_results
136  svd_cluster_dicts = self.generate_cluster_dictsgenerate_cluster_dicts(*svd_results, name="svd")
137 
138  combined_cluster_dicts = self.generate_cluster_dictsgenerate_cluster_dicts(pxd_charges + svd_charges,
139  pxd_path_length + svd_path_length,
140  pxd_normalized_charges + svd_normalized_charges,
141  name="combined")
142 
143  result.update(pxd_cluster_dicts)
144  result.update(svd_cluster_dicts)
145  result.update(combined_cluster_dicts)
146 
147  return result
148 
149 
151  save_tree = refiners.SaveTreeRefiner()
152 
153 
154 class VXDHarvester(QueueHarvester):
155 
156  """ A base class for the VXD hitwise analysis. Collect dE/dX and the correct p of each hit of the MC particles. """
157 
158  def __init__(self, clusters, detector, output_file_name, use_mc_info=True):
159  """Constructor"""
160  HarvestingModule.__init__(self, foreach="TrackCands", output_file_name=output_file_name)
161 
162 
163  self.svd_toolssvd_tools = Belle2.VXDMomentumEstimationTools("Belle2::SVDCluster").getInstance()
164 
165  self.pxd_toolspxd_tools = Belle2.VXDMomentumEstimationTools("Belle2::PXDCluster").getInstance()
166 
167 
168  self.clustersclusters = clusters
169 
170  self.detectordetector = detector
171 
172 
173  self.use_mc_infouse_mc_info = use_mc_info
174 
175  def is_valid_cluster(self, cluster):
176  """Determine if a cluster has an associated SpacePoint and is associated with a primary pion MCParticle"""
177  mc_particles = cluster.getRelationsTo("MCParticles")
178 
179  space_point = cluster.getRelated("SpacePoints")
180 
181  if space_point is None:
182  return False
183 
184  for mc_particle in mc_particles:
185  if (mc_particle.hasStatus(Belle2.MCParticle.c_PrimaryParticle) and
186  abs(mc_particle.getPDG()) == 211):
187 
188  return True
189 
190  return False
191 
192  def get_tools(self, cluster):
193  """Get the PXD tools for a PXD cluster or the SVD tools for an SVD cluster"""
194  return VXDMomentumEnergyEstimator.do_for_each_hit_type(
195  cluster,
196  lambda cluster: self.svd_toolssvd_tools,
197  lambda cluster: self.pxd_toolspxd_tools)
198 
199  def peel(self, track_cand):
200  """Aggregate the cluster and MCParticle (for primary pion) information associated with the track candidate"""
201  vxd_hit_ids = track_cand.getHitIDs(self.detectordetector)
202 
203  vxd_hits = Belle2.PyStoreArray(self.clustersclusters)
204 
205  for vxd_hit_id in vxd_hit_ids:
206  cluster = vxd_hits[vxd_hit_id]
207 
208  if not self.is_valid_clusteris_valid_cluster(cluster):
209  continue
210 
211  tools = self.get_toolsget_tools(cluster)
212 
213  mc_particles = cluster.getRelationsTo("MCParticles")
214 
215  for mc_particle in mc_particles:
216  if (mc_particle.hasStatus(Belle2.MCParticle.c_PrimaryParticle) and abs(mc_particle.getPDG()) == 211):
217 
218  if self.use_mc_infouse_mc_info:
219  track_momentum = mc_particle.getMomentum()
220  track_position = mc_particle.getProductionVertex()
221  track_charge = mc_particle.getCharge()
222  else:
223  track_momentum = track_cand.getMomSeed()
224  track_position = track_cand.getPosSeed()
225  track_charge = track_cand.getChargeSeed()
226 
227  b_field = Belle2.BFieldManager.getField(XYZVector(track_position)).Z() / Belle2.Unit.T
228  track_helix = Belle2.Helix(track_position, track_momentum, int(track_charge), b_field)
229 
230  cluster_charge = tools.getCalibratedCharge(cluster)
231  path_length = tools.getPathLength(cluster, track_helix)
232  cluster_thickness = tools.getThicknessOfCluster(cluster)
233  cluster_radius = tools.getRadiusOfCluster(cluster)
234  cluster_width = tools.getWidthOfCluster(cluster)
235  cluster_length = tools.getLengthOfCluster(cluster)
236 
237  perp_s_at_cluster_entry = track_helix.getArcLength2DAtCylindricalR(cluster_radius)
238  perp_s_at_cluster_exit = track_helix.getArcLength2DAtCylindricalR(cluster_radius + cluster_thickness)
239 
240  mc_momentum = tools.getEntryMomentumOfMCParticle(cluster)
241  mc_position = tools.getEntryPositionOfMCParticle(cluster)
242 
243  mc_b_field = Belle2.BFieldManager.getField(mc_position).Z() / Belle2.Unit.T
244  mc_helix = Belle2.Helix(mc_position, mc_momentum, int(track_charge), mc_b_field)
245  mc_path_length = tools.getPathLength(cluster, mc_helix)
246 
247  cluster_is_u = VXDMomentumEnergyEstimator.do_for_each_hit_type(
248  cluster,
249  lambda cluster: cluster.isUCluster(),
250  lambda cluster: np.NaN)
251 
252  cluster_is_pxd = VXDMomentumEnergyEstimator.do_for_each_hit_type(
253  cluster,
254  lambda cluster: False,
255  lambda cluster: True)
256 
257  cluster_layer = tools.getLayerOfCluster(cluster)
258  cluster_segment = tools.getSegmentNumberOfCluster(cluster)
259  cluster_ladder = tools.getLadderOfCluster(cluster)
260  cluster_sensor = tools.getSensorNumberOfCluster(cluster)
261 
262  cluster_dict = dict(cluster_charge=cluster_charge,
263  cluster_thickness=cluster_thickness,
264  cluster_radius=cluster_radius,
265  cluster_is_u=cluster_is_u,
266  cluster_is_pxd=cluster_is_pxd,
267  cluster_layer=cluster_layer,
268  cluster_segment=cluster_segment,
269  cluster_ladder=cluster_ladder,
270  cluster_width=cluster_width,
271  cluster_length=cluster_length,
272  cluster_sensor=cluster_sensor)
273 
274  mc_at_hit_dict = dict(mc_helix_perigee_x=mc_helix.getPerigeeX(),
275  mc_helix_perigee_y=mc_helix.getPerigeeY(),
276  mc_helix_perigee_z=mc_helix.getPerigeeZ(),
277  mc_helix_momentum_x=mc_helix.getMomentumX(mc_b_field),
278  mc_helix_momentum_y=mc_helix.getMomentumY(mc_b_field),
279  mc_helix_momentum_z=mc_helix.getMomentumZ(mc_b_field),
280  mc_position=mc_position.Mag(),
281  mc_position_x=mc_position.X(),
282  mc_position_y=mc_position.Y(),
283  mc_position_z=mc_position.Z(),
284  mc_momentum=mc_momentum.Mag(),
285  mc_momentum_x=mc_momentum.X(),
286  mc_momentum_y=mc_momentum.Y(),
287  mc_momentum_z=mc_momentum.Z())
288 
289  dedx_dict = dict(dedx=cluster_charge / path_length,
290  dedx_with_mc=cluster_charge / mc_path_length,
291  dedx_with_thickness=cluster_charge / cluster_thickness,
292  p=mc_momentum.Mag(),
293  perp_s_at_cluster_entry=perp_s_at_cluster_entry,
294  perp_s_at_cluster_exit=perp_s_at_cluster_exit,
295  track_charge=track_charge,
296  path_length=path_length,
297  mc_path_length=mc_path_length,
298  p_origin=mc_particle.getMomentum().R())
299 
300  track_dict = dict(track_helix_perigee_x=track_helix.getPerigeeX(),
301  track_helix_perigee_y=track_helix.getPerigeeY(),
302  track_helix_perigee_z=track_helix.getPerigeeZ(),
303  track_helix_momentum_x=track_helix.getMomentumX(b_field),
304  track_helix_momentum_y=track_helix.getMomentumY(b_field),
305  track_helix_momentum_z=track_helix.getMomentumZ(b_field))
306 
307  result_dict = dict()
308  result_dict.update(cluster_dict)
309  result_dict.update(mc_at_hit_dict)
310  result_dict.update(dedx_dict)
311  result_dict.update(track_dict)
312 
313  yield result_dict
314 
315 
317  save_tree = refiners.SaveTreeRefiner()
318 
319 
321  """Collect dE/dX and the correct p of each PXD hit of the MC particles. """
322 
323  def __init__(self, output_file_name, use_mc_info):
324  """Constructor"""
325  VXDHarvester.__init__(self, clusters="PXDClusters", detector=Belle2.Const.PXD, output_file_name=output_file_name,
326  use_mc_info=use_mc_info)
327 
328 
330  """Collect dE/dX and the correct p of each SVD hit of the MC particles. """
331 
332  def __init__(self, output_file_name, use_mc_info):
333  """Constructor"""
334  VXDHarvester.__init__(self, clusters="SVDClusters", detector=Belle2.Const.SVD, output_file_name=output_file_name,
335  use_mc_info=use_mc_info)
336 
337 
338 class FitHarvester(QueueHarvester):
339  """Collect dE/dX and the correct p of each VXD hit associated with the fitted tracks. """
340 
341  def __init__(self, output_file_name, queue):
342  """Constructor"""
343  QueueHarvester.__init__(self, queue, foreach="TrackFitResults", output_file_name=output_file_name)
344 
346 
347  def pick(self, track_fit_result):
348  """Select a TrackFitResult if it is associated with exactly one MCParticle"""
349  mc_track_cands = track_fit_result.getRelationsFrom("TrackCands")
350  if len(mc_track_cands) != 1:
351  return False
352 
353  mc_track_cand = mc_track_cands[0]
354  mc_particles = self.data_storedata_store.getRelationsFromObj(mc_track_cand, "MCParticles")
355 
356  return len(mc_particles) == 1
357 
358  def peel(self, track_fit_result):
359  """Aggregate the track-fit information associated with a TrackFitResult"""
360  mc_track_cand = track_fit_result.getRelationsFrom("TrackCands")[0]
361  mc_particle = self.data_storedata_store.getRelated(mc_track_cand, "MCParticles")
362 
363  fit_momentum = track_fit_result.getMomentum()
364  true_momentum = mc_particle.getMomentum()
365 
366  related_reco_track = track_fit_result.getRelated("GF2Tracks")
367  cardinal_rep = related_reco_track.getCardinalRep()
368  kalman_fit_state = related_reco_track.getKalmanFitStatus()
369 
370  number_of_measurements_in_total = 0
371  number_of_measurements_with_smaller_weight = 0
372 
373  number_of_momentum_measurements_in_total = 0
374  number_of_momentum_measurements_with_smaller_weight = 0
375 
376  for track_point_ID in range(related_reco_track.getNumPointsWithMeasurement()):
377  track_point = related_reco_track.getPointWithMeasurement(track_point_ID)
378 
379  is_momentum_measurement = track_point.getRawMeasurement().__class__.__name__ == "genfit::PlanarMomentumMeasurement"
380 
381  if is_momentum_measurement:
382  number_of_momentum_measurements_in_total += 1
383 
384  if track_point.hasFitterInfo(cardinal_rep):
385  fitter_info = track_point.getFitterInfo(cardinal_rep)
386  num_measurements = fitter_info.getNumMeasurements()
387 
388  for measurement_id in range(num_measurements):
389  number_of_measurements_in_total += 1
390  weight = fitter_info.getMeasurementOnPlane(measurement_id).getWeight()
391  if weight != 1:
392  number_of_measurements_with_smaller_weight += 1
393 
394  if is_momentum_measurement:
395  number_of_momentum_measurements_with_smaller_weight += 1
396 
397  return dict(fit_momentum_x=fit_momentum.X(),
398  fit_momentum_y=fit_momentum.Y(),
399  fit_momentum_z=fit_momentum.Z(),
400  p_value=kalman_fit_state.getForwardPVal(),
401  backward_p_value=kalman_fit_state.getBackwardPVal(),
402  true_momentum_x=true_momentum.X(),
403  true_momentum_y=true_momentum.Y(),
404  true_momentum_z=true_momentum.Z(),
405  number_of_measurements_in_total=number_of_measurements_in_total,
406  number_of_measurements_with_smaller_weight=number_of_measurements_with_smaller_weight,
407  number_of_momentum_measurements_in_total=number_of_momentum_measurements_in_total,
408  number_of_momentum_measurements_with_smaller_weight=number_of_momentum_measurements_with_smaller_weight)
409 
410 
412  save_tree = refiners.SaveTreeRefiner()
Helix parameter class.
Definition: Helix.h:48
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:54
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
Tools needed for the VXD momentum estimation to, e.g.
data_store
access the DataStore singletion
Definition: harvester.py:345
def peel(self, track_fit_result)
Definition: harvester.py:358
def __init__(self, output_file_name, queue)
Definition: harvester.py:341
def pick(self, track_fit_result)
Definition: harvester.py:347
def pick(self, mc_particle)
Definition: harvester.py:93
def peel(self, mc_particle)
Definition: harvester.py:122
def generate_cluster_dicts(self, charge_list, path_length_list, normalized_charge_list, name)
Definition: harvester.py:101
def peel(self, mc_particle_trajectory)
Definition: harvester.py:75
def __init__(self, output_file_name, use_mc_info)
Definition: harvester.py:323
def __init__(self, output_file_name, use_mc_info)
Definition: harvester.py:332
detector
cached copy of the detector identifier (PXD or SVD)
Definition: harvester.py:170
def peel(self, track_cand)
Definition: harvester.py:199
clusters
cached copy of the name of the cluster StoreArray
Definition: harvester.py:168
def is_valid_cluster(self, cluster)
Definition: harvester.py:175
def __init__(self, clusters, detector, output_file_name, use_mc_info=True)
Definition: harvester.py:158
pxd_tools
cached accessor to the PXD-tools singleton
Definition: harvester.py:165
use_mc_info
if true the MC information is used
Definition: harvester.py:173
svd_tools
cached accessor to the SVD-tools singleton
Definition: harvester.py:163
def get_tools(self, cluster)
Definition: harvester.py:192
def generate_truncated(charge_list)
Definition: harvester.py:58
def do_for_each_hit_type(cluster, svd_function, pxd_function)
Definition: harvester.py:23
def calculate_charges_and_path_lengths_for_one_type(clusters, mc_particle)
Definition: harvester.py:34
static void getField(const double *pos, double *field)
return the magnetic field at a given position.
Definition: BFieldManager.h:91