Belle II Software development
SegmentQuadTreePlotter Class Reference
Inheritance diagram for SegmentQuadTreePlotter:
QuadTreePlotter

Public Member Functions

def calculateIntersectionInQuadTreePicture (self, first, second)
 
def calculatePositionInQuadTreePicture (self, position)
 
def forAllAxialSegments (self, f)
 
def convertToQuadTreePicture (self, phi, mag, charge)
 
def event (self)
 

Public Attributes

 range_x_min
 lower x bound for polar angle
 
 range_x_max
 upper x bound for polar angle
 
 range_y_min
 lower y bound
 
 range_y_max
 upper y bound
 

Static Public Attributes

bool draw_segment_intersection = True and False
 by default, do not draw a segment intersection
 
bool draw_segment = True and False
 by default, do not draw a segment
 
bool draw_segment_averaged = True and False
 by default, do not draw an averaged segment
 
bool draw_segment_fitted = True and False
 by default, do not draw a fitted segment
 
bool draw_mc_information = True and False
 by default, do not draw the MC information
 
bool draw_mc_hits = True and False
 by default, do not draw the MC hits
 
bool theta_shifted = False
 by default, polar angles and cuts are in the range (0,pi) rather than (-pi/2,+pi/2)
 
np maximum_theta = np.pi
 an alias for the maximum value of the polar angle
 

Detailed Description

Implementation of a quad tree plotter for SegmentQuadTrees

Definition at line 125 of file quadTreePlotter.py.

Member Function Documentation

◆ calculateIntersectionInQuadTreePicture()

def calculateIntersectionInQuadTreePicture (   self,
  first,
  second 
)
Calculate the point where the two given hits intersect

params
------
first: hit
second: hit

Definition at line 149 of file quadTreePlotter.py.

149 def calculateIntersectionInQuadTreePicture(self, first, second):
150 """
151 Calculate the point where the two given hits intersect
152
153 params
154 ------
155 first: hit
156 second: hit
157 """
158 positionFront = first.getRecoPos2D().conformalTransformed()
159 positionBack = second.getRecoPos2D().conformalTransformed()
160
161 theta_cut = np.arctan2((positionBack - positionFront).x(), (positionFront - positionBack).y())
162
163 if self.theta_shifted:
164 while theta_cut < - self.maximum_theta / 2:
165 theta_cut += self.maximum_theta
166 else:
167 while theta_cut < 0:
168 theta_cut += self.maximum_theta
169
170 r_cut = positionFront.x() * np.cos(theta_cut) + positionFront.y() * np.sin(theta_cut)
171
172 return theta_cut, r_cut
173

◆ calculatePositionInQuadTreePicture()

def calculatePositionInQuadTreePicture (   self,
  position 
)
Transform a given normal coordinate position to a legendre position (conformal transformed)

params
------
position: TrackFindingCDC.Vector2D

Definition at line 174 of file quadTreePlotter.py.

174 def calculatePositionInQuadTreePicture(self, position):
175 """
176 Transform a given normal coordinate position to a legendre position (conformal transformed)
177
178 params
179 ------
180 position: TrackFindingCDC.Vector2D
181 """
182 position = position.conformalTransformed()
183
184 theta = np.linspace(self.range_x_min, self.range_x_max, 100)
185 r = position.x() * np.cos(theta) + position.y() * np.sin(theta)
186
187 return theta, r
188

◆ convertToQuadTreePicture()

def convertToQuadTreePicture (   self,
  phi,
  mag,
  charge 
)
Convert given track parameters into a point in the legendre space

params
------
phi: phi of the track
mag: magnitude of pt
charge: charge of the track

Definition at line 205 of file quadTreePlotter.py.

205 def convertToQuadTreePicture(self, phi, mag, charge):
206 """
207 Convert given track parameters into a point in the legendre space
208
209 params
210 ------
211 phi: phi of the track
212 mag: magnitude of pt
213 charge: charge of the track
214 """
215 theta = phi + np.pi / 2
216 r = 1 / mag * 1.5 * 0.00299792458 * charge
217 if self.theta_shifted:
218 if theta > self.maximum_theta / 2 or theta < -self.maximum_theta / 2:
219 theta = theta % self.maximum_theta - self.maximum_theta / 2
220 else:
221 r *= -1
222 else:
223 if theta > self.maximum_theta or theta < 0:
224 theta = theta % self.maximum_theta
225 else:
226 r *= -1
227 return theta, r
228

◆ event()

def event (   self)
Draw everything according to the given options

Attributes
----------
draw_segment_intersection
draw_segment
draw_segment_averaged
draw_segment_fitted
draw_mc_information
draw_mc_hits

Reimplemented from QuadTreePlotter.

Definition at line 229 of file quadTreePlotter.py.

229 def event(self):
230 """
231 Draw everything according to the given options
232
233 Attributes
234 ----------
235 draw_segment_intersection
236 draw_segment
237 draw_segment_averaged
238 draw_segment_fitted
239 draw_mc_information
240 draw_mc_hits
241 """
242 if self.theta_shifted:
243
244 self.range_x_min = -self.maximum_theta / 2
245
246 self.range_x_max = self.maximum_theta / 2
247 else:
248 self.range_x_min = 0
249 self.range_x_max = self.maximum_theta
250
251
252 self.range_y_min = -0.08
253
254 self.range_y_max = 0.08
255
256 self.init_plotting()
257 # self.plotQuadTreeContent()
258
259 if self.draw_segment_intersection:
260 map = attributemaps.SegmentMCTrackIdColorMap()
261
262 def f(segment):
263 theta, r = self.calculateIntersectionInQuadTreePicture(segment.front(), segment.back())
264 plt.plot(theta, r, color=list(map(0, segment)), marker="o")
265
266 self.forAllAxialSegments(f)
267
268 if self.draw_segment:
269 map = attributemaps.SegmentMCTrackIdColorMap()
270
271 def f(segment):
272 theta, r = self.calculatePositionInQuadTreePicture(segment.front().getRecoPos2D())
273 plt.plot(theta, r, color=list(map(0, segment)), marker="", ls="-")
274
275 self.forAllAxialSegments(f)
276
277 if self.draw_segment_averaged:
278 map = attributemaps.SegmentMCTrackIdColorMap()
279
280 def f(segment):
281 middle_index = int(np.round(segment.size() / 2.0))
282 middle_point = list(segment.items())[middle_index]
283 theta_front, r_front = self.calculateIntersectionInQuadTreePicture(segment.front(), middle_point)
284 theta_back, r_back = self.calculateIntersectionInQuadTreePicture(middle_point, segment.back())
285
286 plt.plot([theta_front, theta_back], [r_front, r_back], color=list(map(0, segment)), marker="o", ls="-")
287
288 self.forAllAxialSegments(f)
289
290 if self.draw_segment_fitted:
291 map = attributemaps.SegmentMCTrackIdColorMap()
293
294 def f(segment):
295 trajectory = fitter.fit(segment)
296 momentum = trajectory.getUnitMom2D(Belle2.TrackFindingCDC.Vector2D(0, 0)).scale(trajectory.getAbsMom2D())
297 theta, r = self.convertToQuadTreePicture(momentum.phi(), momentum.norm(), trajectory.getChargeSign())
298 plt.plot(theta, r, color=list(map(0, segment)), marker="o")
299
300 self.forAllAxialSegments(f)
301
302 if self.draw_hits:
303 cdcHits = Belle2.PyStoreArray("CDCHits")
304 storedWireHits = Belle2.PyStoreObj('CDCWireHitVector')
305 wireHits = storedWireHits.obj().get()
306
307 array = Belle2.PyStoreArray("MCTrackCands")
308 cdc_hits = [cdcHits[i] for track in array for i in track.getHitIDs()]
309
310 for cdcHit in cdcHits:
311 if cdcHit in cdc_hits:
312 continue
313 wireHit = wireHits.at(bisect.bisect_left(wireHits, cdcHit))
314 theta, r = self.calculatePositionInQuadTreePicture(wireHit.getRefPos2D())
315
316 plt.plot(theta, r, marker="", color="black", ls="-", alpha=0.8)
317
318 if self.draw_mc_hits:
319 storedWireHits = Belle2.PyStoreObj('CDCWireHitVector')
320 wireHits = storedWireHits.obj().get()
321
322 map = attributemaps.listColors
323 array = Belle2.PyStoreArray("MCTrackCands")
324 cdcHits = Belle2.PyStoreArray("CDCHits")
325
326 for track in array:
327 mcTrackID = track.getMcTrackId()
328
329 for cdcHitID in track.getHitIDs(Belle2.Const.CDC):
330 cdcHit = cdcHits[cdcHitID]
331 wireHit = wireHits.at(bisect.bisect_left(wireHits, cdcHit))
332
333 theta, r = self.calculatePositionInQuadTreePicture(wireHit.getRefPos2D())
334
335 plt.plot(theta, r, marker="", color=map[mcTrackID % len(map)], ls="-", alpha=0.2)
336
337 if self.draw_mc_information:
338 map = attributemaps.listColors
339 array = Belle2.PyStoreArray("MCTrackCands")
340
341 for track in array:
342 momentum = track.getMomSeed()
343
344 # HARDCODED!!! Temporary solution
345 theta, r = self.convertToQuadTreePicture(momentum.Phi(), momentum.Mag(), track.getChargeSeed())
346 mcTrackID = track.getMcTrackId()
347
348 plt.plot(theta, r, marker="o", color="black", ms=10)
349 plt.plot(theta, r, marker="o", color=map[mcTrackID % len(map)], ms=5)
350
351 self.save_and_show_file()
352
353
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67
static const CDCRiemannFitter & getOriginCircleFitter()
Static getter for an origin circle fitter.
A two dimensional vector which is equipped with functions for correct handeling of orientation relat...
Definition: Vector2D.h:32

◆ forAllAxialSegments()

def forAllAxialSegments (   self,
  f 
)
Loop over all segments and execute a function

params
------
f: function

Definition at line 189 of file quadTreePlotter.py.

189 def forAllAxialSegments(self, f):
190 """
191 Loop over all segments and execute a function
192
193 params
194 ------
195 f: function
196 """
197 items = Belle2.PyStoreObj("CDCSegment2DVector")
198 wrapped_vector = items.obj()
199 vector = wrapped_vector.get()
200
201 for quad_tree_item in vector:
202 if quad_tree_item.getStereoType() == 0:
203 f(quad_tree_item)
204

Member Data Documentation

◆ draw_mc_hits

bool draw_mc_hits = True and False
static

by default, do not draw the MC hits

Definition at line 142 of file quadTreePlotter.py.

◆ draw_mc_information

bool draw_mc_information = True and False
static

by default, do not draw the MC information

Definition at line 140 of file quadTreePlotter.py.

◆ draw_segment

bool draw_segment = True and False
static

by default, do not draw a segment

Definition at line 134 of file quadTreePlotter.py.

◆ draw_segment_averaged

bool draw_segment_averaged = True and False
static

by default, do not draw an averaged segment

Definition at line 136 of file quadTreePlotter.py.

◆ draw_segment_fitted

bool draw_segment_fitted = True and False
static

by default, do not draw a fitted segment

Definition at line 138 of file quadTreePlotter.py.

◆ draw_segment_intersection

bool draw_segment_intersection = True and False
static

by default, do not draw a segment intersection

Definition at line 132 of file quadTreePlotter.py.

◆ maximum_theta

np maximum_theta = np.pi
static

an alias for the maximum value of the polar angle

Definition at line 147 of file quadTreePlotter.py.

◆ range_x_max

range_x_max

upper x bound for polar angle

Definition at line 246 of file quadTreePlotter.py.

◆ range_x_min

range_x_min

lower x bound for polar angle

Definition at line 244 of file quadTreePlotter.py.

◆ range_y_max

range_y_max

upper y bound

Definition at line 254 of file quadTreePlotter.py.

◆ range_y_min

range_y_min

lower y bound

Definition at line 252 of file quadTreePlotter.py.

◆ theta_shifted

bool theta_shifted = False
static

by default, polar angles and cuts are in the range (0,pi) rather than (-pi/2,+pi/2)

Definition at line 145 of file quadTreePlotter.py.


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