11 from ROOT
import Belle2
14 from tracking.ipython_tools.wrapper
import QueueHarvester
19 """ The base class with all static methods to use. """
23 """Apply a PXD selection to a PXDCluster or an SVD selection to an SVDCluster"""
24 cluster_type = cluster.__class__.__name__
25 if cluster_type ==
"Belle2::PXDCluster":
26 return pxd_function(cluster)
27 elif cluster_type ==
"Belle2::SVDCluster":
28 return svd_function(cluster)
30 raise TypeError(
"Unknown hit type")
34 """Return lists of charges and path lengths for the clusters associated with an MCParticle"""
40 momentum = mc_particle.getMomentum()
41 position = mc_particle.getProductionVertex()
42 charge = mc_particle.getCharge()
44 helix =
Belle2.Helix(position, momentum, charge, b_field)
46 for cluster
in clusters:
48 calibrated_charge = tools.getCalibratedCharge(cluster)
49 path_length = tools.getPathLength(cluster, helix)
51 charge_list.append(calibrated_charge)
52 path_length_list.append(path_length)
54 return charge_list, path_length_list, list(np.divide(charge_list, path_length_list))
58 """Sort then truncate a list to all but the last 2 entries or the first 4 entries or the first 6 entries"""
59 sorted_list = sorted(charge_list)
60 if len(sorted_list) > 2:
61 return sorted_list[:-2], sorted_list[:4], sorted_list[:6]
63 return sorted_list, sorted_list[:4], sorted_list[:6]
68 """ A harvester to check for the positions of the track points in the MCParticleTrajectories"""
72 HarvestingModule.__init__(self, foreach=
"MCParticleTrajectorys", output_file_name=
"mc_trajectory.root")
74 def peel(self, mc_particle_trajectory):
75 """Aggregate track-point position information for the trajectory"""
76 for track_point
in mc_particle_trajectory:
77 yield {
"x": track_point.x,
"y": track_point.y,
"z": track_point.z,
"index": self.counter}
81 save_tree = refiners.SaveTreeRefiner()
86 """ A harvester to redo parts of the analysis in the Belle II Paper by Robert """
90 HarvestingModule.__init__(self, foreach=
"MCParticles", output_file_name=
"mc_particle.root")
92 def pick(self, mc_particle):
93 """Select the MCParticle if it is a primary pion and has some PXD and/or SVD clusters"""
94 pxd_clusters = mc_particle.getRelationsFrom(
"PXDClusters")
95 svd_clusters = mc_particle.getRelationsFrom(
"SVDClusters")
96 return (mc_particle.hasStatus(Belle2.MCParticle.c_PrimaryParticle)
and
97 abs(mc_particle.getPDG()) == 211
and
98 len(pxd_clusters) + len(svd_clusters) > 0)
101 """Create a dictionary from the lists of charges, normalized charges, and path lengths of the
102 clusters associated with an MCParticle"""
105 truncated, first4, first6 = VXDMomentumEnergyEstimator.generate_truncated(normalized_charge_list)
107 result.update({
"sum_%s_charges" % name: sum(charge_list)})
108 result.update({
"mean_%s_charges" % name: np.mean(charge_list)})
109 result.update({
"sum_%s_normalized_charges" % name: sum(normalized_charge_list)})
110 result.update({
"mean_%s_normalized_charges" % name: np.mean(normalized_charge_list)})
112 result.update({
"sum_%s_normalized_charges_truncated" % name: sum(truncated)})
113 result.update({
"mean_%s_normalized_charges_truncated" % name: np.mean(truncated)})
114 result.update({
"sum_%s_normalized_charges_first4" % name: sum(first4)})
115 result.update({
"mean_%s_normalized_charges_first4" % name: np.mean(first4)})
116 result.update({
"sum_%s_normalized_charges_first6" % name: sum(first6)})
117 result.update({
"mean_%s_normalized_charges_first6" % name: np.mean(first6)})
122 """Aggregate the PXD and SVD cluster information for an MCParticle"""
125 pxd_clusters = mc_particle.getRelationsFrom(
"PXDClusters")
126 svd_clusters = mc_particle.getRelationsFrom(
"SVDClusters")
128 pxd_results = VXDMomentumEnergyEstimator.calculate_charges_and_path_lengths_for_one_type(pxd_clusters, mc_particle)
129 svd_results = VXDMomentumEnergyEstimator.calculate_charges_and_path_lengths_for_one_type(svd_clusters, mc_particle)
132 pxd_charges, pxd_path_length, pxd_normalized_charges = pxd_results
134 svd_charges, svd_path_length, svd_normalized_charges = svd_results
138 pxd_path_length + svd_path_length,
139 pxd_normalized_charges + svd_normalized_charges,
142 result.update(pxd_cluster_dicts)
143 result.update(svd_cluster_dicts)
144 result.update(combined_cluster_dicts)
150 save_tree = refiners.SaveTreeRefiner()
155 """ A base class for the VXD hitwise analysis. Collect dE/dX and the correct p of each hit of the MC particles. """
157 def __init__(self, clusters, detector, output_file_name, use_mc_info=True):
159 HarvestingModule.__init__(self, foreach=
"TrackCands", output_file_name=output_file_name)
175 """Determine if a cluster has an associated SpacePoint and is associated with a primary pion MCParticle"""
176 mc_particles = cluster.getRelationsTo(
"MCParticles")
178 space_point = cluster.getRelated(
"SpacePoints")
180 if space_point
is None:
183 for mc_particle
in mc_particles:
184 if (mc_particle.hasStatus(Belle2.MCParticle.c_PrimaryParticle)
and
185 abs(mc_particle.getPDG()) == 211):
192 """Get the PXD tools for a PXD cluster or the SVD tools for an SVD cluster"""
193 return VXDMomentumEnergyEstimator.do_for_each_hit_type(
199 """Aggregate the cluster and MCParticle (for primary pion) information associated with the track candidate"""
200 vxd_hit_ids = track_cand.getHitIDs(self.
detectordetector)
204 for vxd_hit_id
in vxd_hit_ids:
205 cluster = vxd_hits[vxd_hit_id]
212 mc_particles = cluster.getRelationsTo(
"MCParticles")
214 for mc_particle
in mc_particles:
215 if (mc_particle.hasStatus(Belle2.MCParticle.c_PrimaryParticle)
and abs(mc_particle.getPDG()) == 211):
218 track_momentum = mc_particle.getMomentum()
219 track_position = mc_particle.getProductionVertex()
220 track_charge = mc_particle.getCharge()
222 track_momentum = track_cand.getMomSeed()
223 track_position = track_cand.getPosSeed()
224 track_charge = track_cand.getChargeSeed()
227 track_helix =
Belle2.Helix(track_position, track_momentum, int(track_charge), b_field)
229 cluster_charge = tools.getCalibratedCharge(cluster)
230 path_length = tools.getPathLength(cluster, track_helix)
231 cluster_thickness = tools.getThicknessOfCluster(cluster)
232 cluster_radius = tools.getRadiusOfCluster(cluster)
233 cluster_width = tools.getWidthOfCluster(cluster)
234 cluster_length = tools.getLengthOfCluster(cluster)
236 perp_s_at_cluster_entry = track_helix.getArcLength2DAtCylindricalR(cluster_radius)
237 perp_s_at_cluster_exit = track_helix.getArcLength2DAtCylindricalR(cluster_radius + cluster_thickness)
239 mc_momentum = tools.getEntryMomentumOfMCParticle(cluster)
240 mc_position = tools.getEntryPositionOfMCParticle(cluster)
243 mc_helix =
Belle2.Helix(mc_position, mc_momentum, int(track_charge), mc_b_field)
244 mc_path_length = tools.getPathLength(cluster, mc_helix)
246 cluster_is_u = VXDMomentumEnergyEstimator.do_for_each_hit_type(
248 lambda cluster: cluster.isUCluster(),
249 lambda cluster: np.NaN)
251 cluster_is_pxd = VXDMomentumEnergyEstimator.do_for_each_hit_type(
253 lambda cluster:
False,
254 lambda cluster:
True)
256 cluster_layer = tools.getLayerOfCluster(cluster)
257 cluster_segment = tools.getSegmentNumberOfCluster(cluster)
258 cluster_ladder = tools.getLadderOfCluster(cluster)
259 cluster_sensor = tools.getSensorNumberOfCluster(cluster)
261 cluster_dict = dict(cluster_charge=cluster_charge,
262 cluster_thickness=cluster_thickness,
263 cluster_radius=cluster_radius,
264 cluster_is_u=cluster_is_u,
265 cluster_is_pxd=cluster_is_pxd,
266 cluster_layer=cluster_layer,
267 cluster_segment=cluster_segment,
268 cluster_ladder=cluster_ladder,
269 cluster_width=cluster_width,
270 cluster_length=cluster_length,
271 cluster_sensor=cluster_sensor)
273 mc_at_hit_dict = dict(mc_helix_perigee_x=mc_helix.getPerigeeX(),
274 mc_helix_perigee_y=mc_helix.getPerigeeY(),
275 mc_helix_perigee_z=mc_helix.getPerigeeZ(),
276 mc_helix_momentum_x=mc_helix.getMomentumX(mc_b_field),
277 mc_helix_momentum_y=mc_helix.getMomentumY(mc_b_field),
278 mc_helix_momentum_z=mc_helix.getMomentumZ(mc_b_field),
279 mc_position=mc_position.Mag(),
280 mc_position_x=mc_position.X(),
281 mc_position_y=mc_position.Y(),
282 mc_position_z=mc_position.Z(),
283 mc_momentum=mc_momentum.Mag(),
284 mc_momentum_x=mc_momentum.X(),
285 mc_momentum_y=mc_momentum.Y(),
286 mc_momentum_z=mc_momentum.Z())
288 dedx_dict = dict(dedx=cluster_charge / path_length,
289 dedx_with_mc=cluster_charge / mc_path_length,
290 dedx_with_thickness=cluster_charge / cluster_thickness,
292 perp_s_at_cluster_entry=perp_s_at_cluster_entry,
293 perp_s_at_cluster_exit=perp_s_at_cluster_exit,
294 track_charge=track_charge,
295 path_length=path_length,
296 mc_path_length=mc_path_length,
297 p_origin=mc_particle.getMomentum().Mag())
299 track_dict = dict(track_helix_perigee_x=track_helix.getPerigeeX(),
300 track_helix_perigee_y=track_helix.getPerigeeY(),
301 track_helix_perigee_z=track_helix.getPerigeeZ(),
302 track_helix_momentum_x=track_helix.getMomentumX(b_field),
303 track_helix_momentum_y=track_helix.getMomentumY(b_field),
304 track_helix_momentum_z=track_helix.getMomentumZ(b_field))
307 result_dict.update(cluster_dict)
308 result_dict.update(mc_at_hit_dict)
309 result_dict.update(dedx_dict)
310 result_dict.update(track_dict)
316 save_tree = refiners.SaveTreeRefiner()
320 """Collect dE/dX and the correct p of each PXD hit of the MC particles. """
324 VXDHarvester.__init__(self, clusters=
"PXDClusters", detector=Belle2.Const.PXD, output_file_name=output_file_name,
325 use_mc_info=use_mc_info)
329 """Collect dE/dX and the correct p of each SVD hit of the MC particles. """
333 VXDHarvester.__init__(self, clusters=
"SVDClusters", detector=Belle2.Const.SVD, output_file_name=output_file_name,
334 use_mc_info=use_mc_info)
338 """Collect dE/dX and the correct p of each VXD hit associated with the fitted tracks. """
342 QueueHarvester.__init__(self, queue, foreach=
"TrackFitResults", output_file_name=output_file_name)
346 def pick(self, track_fit_result):
347 """Select a TrackFitResult if it is associated with exactly one MCParticle"""
348 mc_track_cands = track_fit_result.getRelationsFrom(
"TrackCands")
349 if len(mc_track_cands) != 1:
352 mc_track_cand = mc_track_cands[0]
353 mc_particles = self.
data_storedata_store.getRelationsFromObj(mc_track_cand,
"MCParticles")
355 return len(mc_particles) == 1
357 def peel(self, track_fit_result):
358 """Aggregate the track-fit information associated with a TrackFitResult"""
359 mc_track_cand = track_fit_result.getRelationsFrom(
"TrackCands")[0]
360 mc_particle = self.
data_storedata_store.getRelated(mc_track_cand,
"MCParticles")
362 fit_momentum = track_fit_result.getMomentum()
363 true_momentum = mc_particle.getMomentum()
365 related_reco_track = track_fit_result.getRelated(
"GF2Tracks")
366 cardinal_rep = related_reco_track.getCardinalRep()
367 kalman_fit_state = related_reco_track.getKalmanFitStatus()
369 number_of_measurements_in_total = 0
370 number_of_measurements_with_smaller_weight = 0
372 number_of_momentum_measurements_in_total = 0
373 number_of_momentum_measurements_with_smaller_weight = 0
375 for track_point_ID
in range(related_reco_track.getNumPointsWithMeasurement()):
376 track_point = related_reco_track.getPointWithMeasurement(track_point_ID)
378 is_momentum_measurement = track_point.getRawMeasurement().__class__.__name__ ==
"genfit::PlanarMomentumMeasurement"
380 if is_momentum_measurement:
381 number_of_momentum_measurements_in_total += 1
383 if track_point.hasFitterInfo(cardinal_rep):
384 fitter_info = track_point.getFitterInfo(cardinal_rep)
385 num_measurements = fitter_info.getNumMeasurements()
387 for measurement_id
in range(num_measurements):
388 number_of_measurements_in_total += 1
389 weight = fitter_info.getMeasurementOnPlane(measurement_id).getWeight()
391 number_of_measurements_with_smaller_weight += 1
393 if is_momentum_measurement:
394 number_of_momentum_measurements_with_smaller_weight += 1
396 return dict(fit_momentum_x=fit_momentum.X(),
397 fit_momentum_y=fit_momentum.Y(),
398 fit_momentum_z=fit_momentum.Z(),
399 p_value=kalman_fit_state.getForwardPVal(),
400 backward_p_value=kalman_fit_state.getBackwardPVal(),
401 true_momentum_x=true_momentum.X(),
402 true_momentum_y=true_momentum.Y(),
403 true_momentum_z=true_momentum.Z(),
404 number_of_measurements_in_total=number_of_measurements_in_total,
405 number_of_measurements_with_smaller_weight=number_of_measurements_with_smaller_weight,
406 number_of_momentum_measurements_in_total=number_of_momentum_measurements_in_total,
407 number_of_momentum_measurements_with_smaller_weight=number_of_momentum_measurements_with_smaller_weight)
411 save_tree = refiners.SaveTreeRefiner()
static DataStore & Instance()
Instance of singleton Store.
a (simplified) python wrapper for StoreArray.
data_store
access the DataStore singletion
def peel(self, track_fit_result)
def __init__(self, output_file_name, queue)
def pick(self, track_fit_result)
def pick(self, mc_particle)
def peel(self, mc_particle)
def generate_cluster_dicts(self, charge_list, path_length_list, normalized_charge_list, name)
def peel(self, mc_particle_trajectory)
def __init__(self, output_file_name, use_mc_info)
def __init__(self, output_file_name, use_mc_info)
detector
cached copy of the detector identifier (PXD or SVD)
def peel(self, track_cand)
clusters
cached copy of the name of the cluster StoreArray
def is_valid_cluster(self, cluster)
def __init__(self, clusters, detector, output_file_name, use_mc_info=True)
pxd_tools
cached accessor to the PXD-tools singleton
use_mc_info
if true the MC information is used
svd_tools
cached accessor to the SVD-tools singleton
def get_tools(self, cluster)
def generate_truncated(charge_list)
def do_for_each_hit_type(cluster, svd_function, pxd_function)
def calculate_charges_and_path_lengths_for_one_type(clusters, mc_particle)
static void getField(const double *pos, double *field)
return the magnetic field at a given position.