Belle II Software  release-06-01-15
quadTreePlotter.py
1 
8 
9 from trackfindingcdc.cdcdisplay.svgdrawing import attributemaps
10 import bisect
11 from datetime import datetime
12 import numpy as np
13 import matplotlib.pyplot as plt
14 import basf2
15 import ROOT
16 from ROOT import Belle2
17 
18 from ROOT import gSystem
19 gSystem.Load('libtracking')
20 gSystem.Load('libtracking_trackFindingCDC')
21 
22 
23 class QuadTreePlotter(basf2.Module):
24  """
25  This Module is able to draw the content coming from a QuadTreeImplementation with debugOutput = True.
26  """
27 
28  def __init__(self, queue):
29  """
30  Do not forget to set the ranges! Otherwise you will end up with an empty plot..
31  """
32  basf2.Module.__init__(self)
33 
34  self.file_name_of_quad_tree_contentfile_name_of_quad_tree_content = "output.root"
35 
36  self.draw_quad_tree_contentdraw_quad_tree_content = True
37 
38  self.range_x_minrange_x_min = 0
39 
40  self.range_x_maxrange_x_max = 0
41 
42  self.range_y_minrange_y_min = 0
43 
44  self.range_y_maxrange_y_max = 0
45 
46 
47  self.queuequeue = queue
48 
49  self.file_namesfile_names = []
50 
52  """
53  Draw the quad tree content coming from the root file if enabled.
54  """
55 
56  import seaborn as sb
57 
58  if not self.draw_quad_tree_contentdraw_quad_tree_content:
59  return
60 
61  input_file = ROOT.TFile(self.file_name_of_quad_tree_contentfile_name_of_quad_tree_content)
62 
63  hist = input_file.Get("histUsed")
64 
65  xAxis = hist.GetXaxis()
66  yAxis = hist.GetYaxis()
67 
68  x_edges = np.array([xAxis.GetBinLowEdge(iX) for iX in range(1, xAxis.GetNbins() + 2)])
69  y_edges = np.array([yAxis.GetBinLowEdge(iY) for iY in range(1, yAxis.GetNbins() + 2)])
70 
71  arr_l = np.array([[hist.GetBinContent(iX, iY) for iY in range(1, yAxis.GetNbins() + 1)]
72  for iX in range(1, xAxis.GetNbins() + 1)])
73 
74  hist = input_file.Get("histUnused")
75 
76  xAxis = hist.GetXaxis()
77  yAxis = hist.GetYaxis()
78 
79  x_edges = np.array([xAxis.GetBinLowEdge(iX) for iX in range(1, xAxis.GetNbins() + 2)])
80  y_edges = np.array([yAxis.GetBinLowEdge(iY) for iY in range(1, yAxis.GetNbins() + 2)])
81 
82  l2 = np.array([[hist.GetBinContent(iX, iY) for iY in range(1, yAxis.GetNbins() + 1)]
83  for iX in range(1, xAxis.GetNbins() + 1)])
84 
85  cmap = sb.cubehelix_palette(8, start=2, rot=0, dark=0, light=1, reverse=False, as_cmap=True)
86 
87  plt.gca().pcolorfast(x_edges, y_edges, (arr_l + l2).T, cmap=cmap)
88 
89  x_labels = ["{1:0.{0}f}".format(int(not float(x).is_integer()), x) if i % 4 == 0 else "" for i, x in enumerate(x_edges)]
90  plt.xticks(x_edges, x_labels)
91  y_labels = ["{1:0.{0}f}".format(int(not float(y).is_integer()), y) if i % 4 == 0 else "" for i, y in enumerate(y_edges)]
92  plt.yticks(y_edges, y_labels)
93 
94  def save_and_show_file(self):
95  """
96  Save the plot to a svg and show it (maybe a png would be better?)
97  """
98  fileName = "/tmp/" + datetime.now().isoformat() + '.svg'
99  plt.savefig(fileName)
100  self.file_namesfile_names.append(fileName)
101 
102  def init_plotting(self):
103  """
104  Initialize the figure with the plot ranges
105  We need to implement axes labels later!
106  """
107  plt.clf()
108  plt.xlim(self.range_x_minrange_x_min, self.range_x_maxrange_x_max)
109  plt.ylim(self.range_y_minrange_y_min, self.range_y_maxrange_y_max)
110 
111  def event(self):
112  """
113  Draw everything
114  """
115  self.init_plottinginit_plotting()
116  self.plot_quad_tree_contentplot_quad_tree_content()
117  self.save_and_show_filesave_and_show_file()
118 
119  def terminate(self):
120  """Termination signal at the end of the event processing"""
121  self.queuequeue.put("quadTree", self.file_namesfile_names)
122 
123 
125 
126  """
127  Implementation of a quad tree plotter for SegmentQuadTrees
128  """
129 
130 
131  draw_segment_intersection = True and False
132 
133  draw_segment = True and False
134 
135  draw_segment_averaged = True and False
136 
137  draw_segment_fitted = True and False
138 
139  draw_mc_information = True and False
140 
141  draw_mc_hits = True and False
142 
143 
144  theta_shifted = False
145 
146  maximum_theta = np.pi
147 
148  def calculateIntersectionInQuadTreePicture(self, first, second):
149  """
150  Calculate the point where the two given hits intersect
151 
152  params
153  ------
154  first: hit
155  second: hit
156  """
157  positionFront = first.getRecoPos2D().conformalTransformed()
158  positionBack = second.getRecoPos2D().conformalTransformed()
159 
160  theta_cut = np.arctan2((positionBack - positionFront).x(), (positionFront - positionBack).y())
161 
162  if self.theta_shiftedtheta_shifted:
163  while theta_cut < - self.maximum_thetamaximum_theta / 2:
164  theta_cut += self.maximum_thetamaximum_theta
165  else:
166  while theta_cut < 0:
167  theta_cut += self.maximum_thetamaximum_theta
168 
169  r_cut = positionFront.x() * np.cos(theta_cut) + positionFront.y() * np.sin(theta_cut)
170 
171  return theta_cut, r_cut
172 
174  """
175  Transform a given normal coordinate position to a legendre position (conformal transformed)
176 
177  params
178  ------
179  position: TrackFindingCDC.Vector2D
180  """
181  position = position.conformalTransformed()
182 
183  theta = np.linspace(self.range_x_minrange_x_minrange_x_min, self.range_x_maxrange_x_maxrange_x_max, 100)
184  r = position.x() * np.cos(theta) + position.y() * np.sin(theta)
185 
186  return theta, r
187 
188  def forAllAxialSegments(self, f):
189  """
190  Loop over all segments and execute a function
191 
192  params
193  ------
194  f: function
195  """
196  items = Belle2.PyStoreObj("CDCSegment2DVector")
197  wrapped_vector = items.obj()
198  vector = wrapped_vector.get()
199 
200  for quad_tree_item in vector:
201  if quad_tree_item.getStereoType() == 0:
202  f(quad_tree_item)
203 
204  def convertToQuadTreePicture(self, phi, mag, charge):
205  """
206  Convert given track parameters into a point in the legendre space
207 
208  params
209  ------
210  phi: phi of the track
211  mag: magnitude of pt
212  charge: charge of the track
213  """
214  theta = phi + np.pi / 2
215  r = 1 / mag * 1.5 * 0.00299792458 * charge
216  if self.theta_shiftedtheta_shifted:
217  if theta > self.maximum_thetamaximum_theta / 2 or theta < -self.maximum_thetamaximum_theta / 2:
218  theta = theta % self.maximum_thetamaximum_theta - self.maximum_thetamaximum_theta / 2
219  else:
220  r *= -1
221  else:
222  if theta > self.maximum_thetamaximum_theta or theta < 0:
223  theta = theta % self.maximum_thetamaximum_theta
224  else:
225  r *= -1
226  return theta, r
227 
228  def event(self):
229  """
230  Draw everything according to the given options
231 
232  Attributes
233  ----------
234  draw_segment_intersection
235  draw_segment
236  draw_segment_averaged
237  draw_segment_fitted
238  draw_mc_information
239  draw_mc_hits
240  """
241  if self.theta_shiftedtheta_shifted:
242 
243  self.range_x_minrange_x_minrange_x_min = -self.maximum_thetamaximum_theta / 2
244 
245  self.range_x_maxrange_x_maxrange_x_max = self.maximum_thetamaximum_theta / 2
246  else:
247  self.range_x_minrange_x_minrange_x_min = 0
248  self.range_x_maxrange_x_maxrange_x_max = self.maximum_thetamaximum_theta
249 
250 
251  self.range_y_minrange_y_minrange_y_min = -0.08
252 
253  self.range_y_maxrange_y_maxrange_y_max = 0.08
254 
255  self.init_plottinginit_plotting()
256  # self.plotQuadTreeContent()
257 
258  if self.draw_segment_intersectiondraw_segment_intersection:
259  map = attributemaps.SegmentMCTrackIdColorMap()
260 
261  def f(segment):
262  theta, r = self.calculateIntersectionInQuadTreePicturecalculateIntersectionInQuadTreePicture(segment.front(), segment.back())
263  plt.plot(theta, r, color=list(map(0, segment)), marker="o")
264 
265  self.forAllAxialSegmentsforAllAxialSegments(f)
266 
267  if self.draw_segmentdraw_segment:
268  map = attributemaps.SegmentMCTrackIdColorMap()
269 
270  def f(segment):
271  theta, r = self.calculatePositionInQuadTreePicturecalculatePositionInQuadTreePicture(segment.front().getRecoPos2D())
272  plt.plot(theta, r, color=list(map(0, segment)), marker="", ls="-")
273 
274  self.forAllAxialSegmentsforAllAxialSegments(f)
275 
276  if self.draw_segment_averageddraw_segment_averaged:
277  map = attributemaps.SegmentMCTrackIdColorMap()
278 
279  def f(segment):
280  middle_index = int(np.round(segment.size() / 2.0))
281  middle_point = list(segment.items())[middle_index]
282  theta_front, r_front = self.calculateIntersectionInQuadTreePicturecalculateIntersectionInQuadTreePicture(segment.front(), middle_point)
283  theta_back, r_back = self.calculateIntersectionInQuadTreePicturecalculateIntersectionInQuadTreePicture(middle_point, segment.back())
284 
285  plt.plot([theta_front, theta_back], [r_front, r_back], color=list(map(0, segment)), marker="o", ls="-")
286 
287  self.forAllAxialSegmentsforAllAxialSegments(f)
288 
289  if self.draw_segment_fitteddraw_segment_fitted:
290  map = attributemaps.SegmentMCTrackIdColorMap()
292 
293  def f(segment):
294  trajectory = fitter.fit(segment)
295  momentum = trajectory.getUnitMom2D(Belle2.TrackFindingCDC.Vector2D(0, 0)).scale(trajectory.getAbsMom2D())
296  theta, r = self.convertToQuadTreePictureconvertToQuadTreePicture(momentum.phi(), momentum.norm(), trajectory.getChargeSign())
297  plt.plot(theta, r, color=list(map(0, segment)), marker="o")
298 
299  self.forAllAxialSegmentsforAllAxialSegments(f)
300 
301  if self.draw_hits:
302  cdcHits = Belle2.PyStoreArray("CDCHits")
303  storedWireHits = Belle2.PyStoreObj('CDCWireHitVector')
304  wireHits = storedWireHits.obj().get()
305 
306  array = Belle2.PyStoreArray("MCTrackCands")
307  cdc_hits = [cdcHits[i] for track in array for i in track.getHitIDs()]
308 
309  for cdcHit in cdcHits:
310  if cdcHit in cdc_hits:
311  continue
312  wireHit = wireHits.at(bisect.bisect_left(wireHits, cdcHit))
313  theta, r = self.calculatePositionInQuadTreePicturecalculatePositionInQuadTreePicture(wireHit.getRefPos2D())
314 
315  plt.plot(theta, r, marker="", color="black", ls="-", alpha=0.8)
316 
317  if self.draw_mc_hitsdraw_mc_hits:
318  storedWireHits = Belle2.PyStoreObj('CDCWireHitVector')
319  wireHits = storedWireHits.obj().get()
320 
321  map = attributemaps.listColors
322  array = Belle2.PyStoreArray("MCTrackCands")
323  cdcHits = Belle2.PyStoreArray("CDCHits")
324 
325  for track in array:
326  mcTrackID = track.getMcTrackId()
327 
328  for cdcHitID in track.getHitIDs(Belle2.Const.CDC):
329  cdcHit = cdcHits[cdcHitID]
330  wireHit = wireHits.at(bisect.bisect_left(wireHits, cdcHit))
331 
332  theta, r = self.calculatePositionInQuadTreePicturecalculatePositionInQuadTreePicture(wireHit.getRefPos2D())
333 
334  plt.plot(theta, r, marker="", color=map[mcTrackID % len(map)], ls="-", alpha=0.2)
335 
336  if self.draw_mc_informationdraw_mc_information:
337  map = attributemaps.listColors
338  array = Belle2.PyStoreArray("MCTrackCands")
339 
340  for track in array:
341  momentum = track.getMomSeed()
342 
343  # HARDCODED!!! Temporary solution
344  theta, r = self.convertToQuadTreePictureconvertToQuadTreePicture(momentum.Phi(), momentum.Mag(), track.getChargeSeed())
345  mcTrackID = track.getMcTrackId()
346 
347  plt.plot(theta, r, marker="o", color="black", ms=10)
348  plt.plot(theta, r, marker="o", color=map[mcTrackID % len(map)], ms=5)
349 
350  self.save_and_show_filesave_and_show_file()
351 
352 
354 
355  """
356  Implementation of a quad tree plotter for StereoHitAssignment
357  """
358 
359 
360  draw_mc_hits = False
361 
362  draw_mc_tracks = False
363 
364  draw_track_hits = False
365 
366  draw_last_track = True
367 
368  delete_bad_hits = False
369 
371  """
372  Convert a genfit::TrackCand into a TrackFindingCDC.CDCTrajectory3D
373 
374  params
375  ------
376  track: genfit::TrackCand
377  """
380 
381  position = Vector3D(track.getPosSeed())
382  momentum = Vector3D(track.getMomSeed())
383  charge = track.getChargeSeed()
384 
385  return Trajectory3D(position, momentum, charge)
386 
387  def create_reco_hit3D(self, cdcHit, trajectory3D, rlInfo):
388  """
389  Use a cdc hit and a trajectory to reconstruct a CDCRecoHit3D
390 
391  params
392  ------
393  cdcHit: CDCHit
394  trajectory3D: TrackFindingCDC.CDCTrajectory3D
395  rlInfo: RightLeftInfo ( = short)
396  """
397  storedWireHits = Belle2.PyStoreObj('CDCWireHitVector')
398  wireHits = storedWireHits.obj().get()
399 
401  wireHit = wireHits.at(bisect.bisect_left(wireHits, cdcHit))
402  rightLeftWireHit = Belle2.TrackFindingCDC.CDCRLWireHit(wireHit, rlInfo)
403  if rightLeftWireHit.getStereoType() != 0:
404  recoHit3D = CDCRecoHit3D.reconstruct(rightLeftWireHit, trajectory3D.getTrajectory2D())
405  return recoHit3D
406  else:
407  return None
408 
409  def get_plottable_line(self, recoHit3D):
410  """
411  Minim the task of the StereoQuadTree by showing the line of quadtree nodes
412  a hit belongs to
413  """
414  z0 = [self.range_y_minrange_y_minrange_y_min, self.range_y_maxrange_y_maxrange_y_max]
415  arr_l = np.array((np.array(recoHit3D.getRecoPos3D().z()) - z0) / recoHit3D.getArcLength2D())
416  return arr_l, z0
417 
418  def plot_hit_line(self, recoHit3D, color):
419  """
420  Draw one recoHit3D
421  """
422  if recoHit3D:
423  if recoHit3D.getStereoType() == 0:
424  return
425 
426  arr_l, z0 = self.get_plottable_lineget_plottable_line(recoHit3D)
427  plt.plot(arr_l, z0, marker="", ls="-", alpha=0.4, color=color)
428 
429  def event(self):
430  """
431  Draw the hit content according to the attributes
432 
433  Attributes
434  ----------
435  draw_mc_hits
436  draw_mc_tracks
437  draw_track_hits
438  draw_last_track
439  delete_bad_hits
440  """
441 
442  self.draw_quad_tree_contentdraw_quad_tree_contentdraw_quad_tree_content = True
443 
444 
445  self.range_x_minrange_x_minrange_x_min = -2 - np.sqrt(3)
446 
447  self.range_x_maxrange_x_maxrange_x_max = 2 + np.sqrt(3)
448 
449 
450  self.range_y_minrange_y_minrange_y_min = -100
451 
452  self.range_y_maxrange_y_maxrange_y_max = -self.range_y_minrange_y_minrange_y_min
453 
454  self.init_plottinginit_plotting()
455  self.plot_quad_tree_contentplot_quad_tree_content()
456 
457  map = attributemaps.listColors
458  cdcHits = Belle2.PyStoreArray("CDCHits")
459 
460  items = Belle2.PyStoreObj("CDCTrackVector")
461  wrapped_vector = items.obj()
462  track_vector = wrapped_vector.get()
463 
464  mcHitLookUp = Belle2.TrackFindingCDC.CDCMCHitLookUp().getInstance()
465  mcHitLookUp.fill()
466 
467  storedWireHits = Belle2.PyStoreObj('CDCWireHitVector')
468  wireHits = storedWireHits.obj().get()
469 
470  if self.draw_mc_hitsdraw_mc_hits:
471  mc_track_cands = Belle2.PyStoreArray("MCTrackCands")
472 
473  for track in mc_track_cands:
474  mcTrackID = track.getMcTrackId()
475  trajectory = self.create_trajectory_from_trackcreate_trajectory_from_track(track)
476 
477  for cdcHitID in track.getHitIDs(Belle2.Const.CDC):
478  cdcHit = cdcHits[cdcHitID]
479 
480  leftRecoHit3D = self.create_reco_hit3Dcreate_reco_hit3D(cdcHit, trajectory, -1)
481  rightRecoHit3D = self.create_reco_hit3Dcreate_reco_hit3D(cdcHit, trajectory, 1)
482 
483  self.plot_hit_lineplot_hit_line(leftRecoHit3D, color=map[mcTrackID % len(map)])
484  self.plot_hit_lineplot_hit_line(rightRecoHit3D, color=map[mcTrackID % len(map)])
485 
486  if self.draw_mc_tracksdraw_mc_tracks:
487  mc_track_cands = Belle2.PyStoreArray("MCTrackCands")
488 
489  for track in mc_track_cands:
490  mcTrackID = track.getMcTrackId()
491  trajectory = self.create_trajectory_from_trackcreate_trajectory_from_track(track)
492  z0 = trajectory.getTrajectorySZ().getZ0()
493 
494  for cdcHitID in track.getHitIDs(Belle2.Const.CDC):
495  cdcHit = cdcHits[cdcHitID]
496  recoHit3D = self.create_reco_hit3Dcreate_reco_hit3D(cdcHit, trajectory, mcHitLookUp.getRLInfo(cdcHit))
497 
498  if recoHit3D:
499  arr_l = (recoHit3D.getRecoPos3D().z() - z0) / recoHit3D.getArcLength2D()
500  plt.plot(arr_l, z0, marker="o", color=map[mcTrackID % len(map)], ls="", alpha=0.2)
501 
502  if self.draw_track_hitsdraw_track_hits:
503  for id, track in enumerate(track_vector):
504  for recoHit3D in list(track.items()):
505  self.plot_hit_lineplot_hit_line(recoHit3D, color=map[id % len(map)])
506 
507  if self.draw_last_trackdraw_last_track and len(track_vector) != 0:
508 
509  last_track = track_vector[-1]
510  trajectory = last_track.getStartTrajectory3D().getTrajectory2D()
511 
512  for wireHit in wireHits:
513  for rlInfo in (-1, 1):
514  recoHit3D = Belle2.TrackFindingCDC.CDCRecoHit3D.reconstruct(wireHit, rlInfo, trajectory)
515 
516  if (self.delete_bad_hitsdelete_bad_hits and
517  (wireHit.getRLInfo() != mcHitLookUp.getRLInfo(wireHit.getWireHit().getHit()) or
518  not recoHit3D.isInCellZBounds())):
519  continue
520 
521  if recoHit3D in list(last_track.items()):
522  color = map[len(track_vector) % len(map)]
523  else:
524  if wireHit.getRLInfo() == 1:
525  color = "black"
526  else:
527  color = "gray"
528  self.plot_hit_lineplot_hit_line(recoHit3D, color)
529 
530  plt.xlabel(r"$\tan \ \lambda$")
531  plt.ylabel(r"$z_0$")
532  self.save_and_show_filesave_and_show_file()
533 
534 
536 
537  """
538  A wrapper around the svg drawer in the tracking package that
539  writes its output files as a list to the queue
540  """
541 
542  def __init__(self, queue, label, *args, **kwargs):
543  """ The same as the base class, except:
544 
545  Arguments
546  ---------
547 
548  queue: The queue to write to
549  label: The key name in the queue
550  """
551 
552  self.queuequeuequeue = queue
553 
554  self.labellabel = label
555  StereoQuadTreePlotter.__init__(self, *args, **kwargs)
556 
557 
558  self.file_listfile_list = []
559 
560  def terminate(self):
561  """ Overwrite the terminate to put the list to the queue"""
562  StereoQuadTreePlotter.terminate(self)
563  self.queuequeuequeue.put(self.labellabel, self.file_listfile_list)
564 
566  """ Overwrite the function to listen for every new filename """
567 
568  from datetime import datetime
569  from matplotlib import pyplot as plt
570 
571  fileName = "/tmp/" + datetime.now().isoformat() + '.svg'
572  plt.savefig(fileName)
573  self.file_listfile_list.append(fileName)
a (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:56
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67
Interface class to the Monte Carlo information for individual hits.
Class representing an oriented hit wire including a hypotheses whether the causing track passes left ...
Definition: CDCRLWireHit.h:41
Class representing a three dimensional reconstructed hit.
Definition: CDCRecoHit3D.h:52
static CDCRecoHit3D reconstruct(const CDCRecoHit2D &recoHit2D, const CDCTrajectory2D &trajectory2D)
Reconstructs the three dimensional hit from the two dimensional and the two dimensional trajectory.
Definition: CDCRecoHit3D.cc:56
static const CDCRiemannFitter & getOriginCircleFitter()
Static getter for an origin circle fitter.
Particle full three dimensional trajectory.
A two dimensional vector which is equipped with functions for correct handeling of orientation relat...
Definition: Vector2D.h:35
A three dimensional vector.
Definition: Vector3D.h:32
range_x_min
cached minimum x value
draw_quad_tree_content
cached flag to draw QuadTree
range_y_max
cached maximum y value
range_x_max
cached maximum x value
file_name_of_quad_tree_content
cached output filename
file_names
cached array of output filenames (one file per image)
range_y_min
cached minimum y value
queue
cached value of the queue input parameter
def __init__(self, queue, label, *args, **kwargs)
label
The label for writing to the queue.
bool draw_segment_averaged
by default, do not draw an averaged segment
bool draw_segment
by default, do not draw a segment
bool draw_segment_fitted
by default, do not draw a fitted segment
range_x_min
lower x bound for polar angle
bool draw_mc_hits
by default, do not draw the MC hits
bool draw_mc_information
by default, do not draw the MC information
def convertToQuadTreePicture(self, phi, mag, charge)
maximum_theta
an alias for the maximum value of the polar angle
def calculatePositionInQuadTreePicture(self, position)
range_x_max
upper x bound for polar angle
bool draw_segment_intersection
by default, do not draw a segment intersection
bool theta_shifted
by default, polar angles and cuts are in the range (0,pi) rather than (-pi/2,+pi/2)
def calculateIntersectionInQuadTreePicture(self, first, second)
def plot_hit_line(self, recoHit3D, color)
bool draw_mc_hits
by default, do not draw the MC hits
bool delete_bad_hits
by default, do not delete the bad track hits
draw_quad_tree_content
by default, draw the QuadTree
bool draw_last_track
by default, draw the last track
def create_reco_hit3D(self, cdcHit, trajectory3D, rlInfo)
bool draw_track_hits
by default, do not draw the track hits
bool draw_mc_tracks
by default, do not draw the MC tracks
HepGeom::Vector3D< double > Vector3D
3D Vector
Definition: Cell.h:34