Belle II Software development
StereoQuadTreePlotter Class Reference
Inheritance diagram for StereoQuadTreePlotter:
QuadTreePlotter QueueStereoQuadTreePlotter

Public Member Functions

def create_trajectory_from_track (self, track)
 
def create_reco_hit3D (self, cdcHit, trajectory3D, rlInfo)
 
def get_plottable_line (self, recoHit3D)
 
def plot_hit_line (self, recoHit3D, color)
 
def event (self)
 

Public Attributes

 draw_quad_tree_content
 by default, draw the QuadTree
 
 range_x_min
 default lower x bound
 
 range_x_max
 default upper x bound
 
 range_y_min
 default lower y bound
 
 range_y_max
 default upper y bound
 

Static Public Attributes

bool draw_mc_hits = False
 by default, do not draw the MC hits
 
bool draw_mc_tracks = False
 by default, do not draw the MC tracks
 
bool draw_track_hits = False
 by default, do not draw the track hits
 
bool draw_last_track = True
 by default, draw the last track
 
bool delete_bad_hits = False
 by default, do not delete the bad track hits
 

Detailed Description

Implementation of a quad tree plotter for StereoHitAssignment

Definition at line 354 of file quadTreePlotter.py.

Member Function Documentation

◆ create_reco_hit3D()

def create_reco_hit3D (   self,
  cdcHit,
  trajectory3D,
  rlInfo 
)
Use a cdc hit and a trajectory to reconstruct a CDCRecoHit3D

params
------
cdcHit: CDCHit
trajectory3D: TrackFindingCDC.CDCTrajectory3D
rlInfo: RightLeftInfo ( = short)

Definition at line 388 of file quadTreePlotter.py.

388 def create_reco_hit3D(self, cdcHit, trajectory3D, rlInfo):
389 """
390 Use a cdc hit and a trajectory to reconstruct a CDCRecoHit3D
391
392 params
393 ------
394 cdcHit: CDCHit
395 trajectory3D: TrackFindingCDC.CDCTrajectory3D
396 rlInfo: RightLeftInfo ( = short)
397 """
398 storedWireHits = Belle2.PyStoreObj('CDCWireHitVector')
399 wireHits = storedWireHits.obj().get()
400
402 wireHit = wireHits.at(bisect.bisect_left(wireHits, cdcHit))
403 rightLeftWireHit = Belle2.TrackFindingCDC.CDCRLWireHit(wireHit, rlInfo)
404 if rightLeftWireHit.getStereoType() != 0:
405 recoHit3D = CDCRecoHit3D.reconstruct(rightLeftWireHit, trajectory3D.getTrajectory2D())
406 return recoHit3D
407 else:
408 return None
409
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67
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

◆ create_trajectory_from_track()

def create_trajectory_from_track (   self,
  track 
)
Convert a genfit::TrackCand into a TrackFindingCDC.CDCTrajectory3D

params
------
track: genfit::TrackCand

Definition at line 371 of file quadTreePlotter.py.

371 def create_trajectory_from_track(self, track):
372 """
373 Convert a genfit::TrackCand into a TrackFindingCDC.CDCTrajectory3D
374
375 params
376 ------
377 track: genfit::TrackCand
378 """
381
382 position = Vector3D(track.getPosSeed())
383 momentum = Vector3D(track.getMomSeed())
384 charge = track.getChargeSeed()
385
386 return Trajectory3D(position, momentum, charge)
387
Particle full three dimensional trajectory.
A three dimensional vector.
Definition: Vector3D.h:33

◆ event()

def event (   self)
Draw the hit content according to the attributes

Attributes
----------
draw_mc_hits
draw_mc_tracks
draw_track_hits
draw_last_track
delete_bad_hits

Reimplemented from QuadTreePlotter.

Definition at line 430 of file quadTreePlotter.py.

430 def event(self):
431 """
432 Draw the hit content according to the attributes
433
434 Attributes
435 ----------
436 draw_mc_hits
437 draw_mc_tracks
438 draw_track_hits
439 draw_last_track
440 delete_bad_hits
441 """
442
443 self.draw_quad_tree_content = True
444
445
446 self.range_x_min = -2 - np.sqrt(3)
447
448 self.range_x_max = 2 + np.sqrt(3)
449
450
451 self.range_y_min = -100
452
453 self.range_y_max = -self.range_y_min
454
455 self.init_plotting()
456 self.plot_quad_tree_content()
457
458 map = attributemaps.listColors
459 cdcHits = Belle2.PyStoreArray("CDCHits")
460
461 items = Belle2.PyStoreObj("CDCTrackVector")
462 wrapped_vector = items.obj()
463 track_vector = wrapped_vector.get()
464
465 mcHitLookUp = Belle2.TrackFindingCDC.CDCMCHitLookUp().getInstance()
466 mcHitLookUp.fill()
467
468 storedWireHits = Belle2.PyStoreObj('CDCWireHitVector')
469 wireHits = storedWireHits.obj().get()
470
471 if self.draw_mc_hits:
472 mc_track_cands = Belle2.PyStoreArray("MCTrackCands")
473
474 for track in mc_track_cands:
475 mcTrackID = track.getMcTrackId()
476 trajectory = self.create_trajectory_from_track(track)
477
478 for cdcHitID in track.getHitIDs(Belle2.Const.CDC):
479 cdcHit = cdcHits[cdcHitID]
480
481 leftRecoHit3D = self.create_reco_hit3D(cdcHit, trajectory, -1)
482 rightRecoHit3D = self.create_reco_hit3D(cdcHit, trajectory, 1)
483
484 self.plot_hit_line(leftRecoHit3D, color=map[mcTrackID % len(map)])
485 self.plot_hit_line(rightRecoHit3D, color=map[mcTrackID % len(map)])
486
487 if self.draw_mc_tracks:
488 mc_track_cands = Belle2.PyStoreArray("MCTrackCands")
489
490 for track in mc_track_cands:
491 mcTrackID = track.getMcTrackId()
492 trajectory = self.create_trajectory_from_track(track)
493 z0 = trajectory.getTrajectorySZ().getZ0()
494
495 for cdcHitID in track.getHitIDs(Belle2.Const.CDC):
496 cdcHit = cdcHits[cdcHitID]
497 recoHit3D = self.create_reco_hit3D(cdcHit, trajectory, mcHitLookUp.getRLInfo(cdcHit))
498
499 if recoHit3D:
500 arr_l = (recoHit3D.getRecoPos3D().z() - z0) / recoHit3D.getArcLength2D()
501 plt.plot(arr_l, z0, marker="o", color=map[mcTrackID % len(map)], ls="", alpha=0.2)
502
503 if self.draw_track_hits:
504 for id, track in enumerate(track_vector):
505 for recoHit3D in list(track.items()):
506 self.plot_hit_line(recoHit3D, color=map[id % len(map)])
507
508 if self.draw_last_track and len(track_vector) != 0:
509
510 last_track = track_vector[-1]
511 trajectory = last_track.getStartTrajectory3D().getTrajectory2D()
512
513 for wireHit in wireHits:
514 for rlInfo in (-1, 1):
515 recoHit3D = Belle2.TrackFindingCDC.CDCRecoHit3D.reconstruct(wireHit, rlInfo, trajectory)
516
517 if (self.delete_bad_hits and
518 (wireHit.getRLInfo() != mcHitLookUp.getRLInfo(wireHit.getWireHit().getHit()) or
519 not recoHit3D.isInCellZBounds())):
520 continue
521
522 if recoHit3D in list(last_track.items()):
523 color = map[len(track_vector) % len(map)]
524 else:
525 if wireHit.getRLInfo() == 1:
526 color = "black"
527 else:
528 color = "gray"
529 self.plot_hit_line(recoHit3D, color)
530
531 plt.xlabel(r"$\tan \ \lambda$")
532 plt.ylabel(r"$z_0$")
533 self.save_and_show_file()
534
535
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
Interface class to the Monte Carlo information for individual hits.
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

◆ get_plottable_line()

def get_plottable_line (   self,
  recoHit3D 
)
Minim the task of the StereoQuadTree by showing the line of quadtree nodes
a hit belongs to

Definition at line 410 of file quadTreePlotter.py.

410 def get_plottable_line(self, recoHit3D):
411 """
412 Minim the task of the StereoQuadTree by showing the line of quadtree nodes
413 a hit belongs to
414 """
415 z0 = [self.range_y_min, self.range_y_max]
416 arr_l = np.array((np.array(recoHit3D.getRecoPos3D().z()) - z0) / recoHit3D.getArcLength2D())
417 return arr_l, z0
418

◆ plot_hit_line()

def plot_hit_line (   self,
  recoHit3D,
  color 
)
Draw one recoHit3D

Definition at line 419 of file quadTreePlotter.py.

419 def plot_hit_line(self, recoHit3D, color):
420 """
421 Draw one recoHit3D
422 """
423 if recoHit3D:
424 if recoHit3D.getStereoType() == 0:
425 return
426
427 arr_l, z0 = self.get_plottable_line(recoHit3D)
428 plt.plot(arr_l, z0, marker="", ls="-", alpha=0.4, color=color)
429

Member Data Documentation

◆ delete_bad_hits

bool delete_bad_hits = False
static

by default, do not delete the bad track hits

Definition at line 369 of file quadTreePlotter.py.

◆ draw_last_track

bool draw_last_track = True
static

by default, draw the last track

Definition at line 367 of file quadTreePlotter.py.

◆ draw_mc_hits

bool draw_mc_hits = False
static

by default, do not draw the MC hits

Definition at line 361 of file quadTreePlotter.py.

◆ draw_mc_tracks

bool draw_mc_tracks = False
static

by default, do not draw the MC tracks

Definition at line 363 of file quadTreePlotter.py.

◆ draw_quad_tree_content

draw_quad_tree_content

by default, draw the QuadTree

Definition at line 443 of file quadTreePlotter.py.

◆ draw_track_hits

bool draw_track_hits = False
static

by default, do not draw the track hits

Definition at line 365 of file quadTreePlotter.py.

◆ range_x_max

range_x_max

default upper x bound

Definition at line 448 of file quadTreePlotter.py.

◆ range_x_min

range_x_min

default lower x bound

Definition at line 446 of file quadTreePlotter.py.

◆ range_y_max

range_y_max

default upper y bound

Definition at line 453 of file quadTreePlotter.py.

◆ range_y_min

range_y_min

default lower y bound

Definition at line 451 of file quadTreePlotter.py.


The documentation for this class was generated from the following file: