Belle II Software development
flight_time_mc_compare.py
1
8
10 save_tree,
11 save_pull_analysis,
12)
13from tracking.harvest.harvesting import harvest
14import tracking
15import basf2
16from ROOT import Belle2
17from ROOT import gSystem
18from generators import add_cosmics_generator
19
20gSystem.Load('libcdc')
21gSystem.Load('libtracking')
22gSystem.Load('libtracking_trackFindingCDC')
23
24
26
27
28def main():
29 path = basf2.create_path()
30
31 # Create Event information
32 path.add_module('EventInfoSetter',
33 evtNumList=[50000],
34 runList=[1],
35 )
36
37 # Show progress of processing
38 path.add_module('Progress')
39
40 # Load gearbox parameters
41 gearbox_module = path.add_module('Gearbox')
42
43 # Create geometry
44 geometry_module = path.add_module('Geometry')
45
46 # Adjust gemetry for the test beam #
47 # ################################ #
48 cdc_top_geometry = False
49 if cdc_top_geometry:
50 gearbox_module.param(dict(fileName="geometry/CDCcosmicTests.xml")) # <- does something mysterious to the reconstruction...
51 # The peculiar behaviour seems to originate from the low magnetic field setup close to the earth magentic field.
52 # In the current state of the reconstruction it seems to be better to set the magnetic field to zero
53 # Override of the magnetic field is made below as a temporary work around until we discussed the topic with
54 # the TOP group how composed the geometry/CDCcosmicTests.xml.
55 gearbox_module.param(dict(overrideMultiple=[
56 # Explicitly override magnetic field to zero instead of earth magnetic field as defined in CDC TOP geometry
57 ("/DetectorComponent[@name='MagneticField']//Component/Z", "0", ""),
58 ]
59 ))
60
61 else:
62 # Use the regular geometry but only the cdc, no magnetic field.
63 geometry_module.param(dict(components=["CDC"],))
64
65 gearbox_module.param(dict(
66 override=[
67 # Reset the top volume: must be larger than the generated surface and higher than the detector
68 # It is the users responsibility to ensure a full angular coverage
69 ("/Global/length", "8.", "m"),
70 ("/Global/width", "8.", "m"),
71 ("/Global/height", "1.5", "m"),
72
73 # Adjustments of the CDC setup
74 ("/DetectorComponent[@name='CDC']//t0FileName", "t0.dat", ""),
75 ("/DetectorComponent[@name='CDC']//xtFileName", "xt_noB_v1.dat", ""),
76 # ("/DetectorComponent[@name='CDC']//bwFileName", "badwire_CDCTOP.dat", ""),
77 ("/DetectorComponent[@name='CDC']//GlobalPhiRotation", "1.875", "deg"),
78 ]
79 ))
80
81 # Generator setup #
82 # ############### #
83 # generator_module = "cry"
84 generator_module = "cosmics"
85 # generator_module = "particle_gun"
86
87 if generator_module == "cry":
88 add_cosmics_generator(path, accept_box=[0.5, 0.5, 0.5], keep_box=[0.5, 0.5, 0.5])
89
90 elif generator_module == "cosmics":
91 path.add_module("Cosmics",
92 # Require particles close to the counter
93 ipRequirement=1,
94 ipdr=5,
95 ipdz=15,
96 )
97
98 elif generator_module == "particle_gun":
99 # Generate muons starting in the IP
100 # Yields outgoing arm only
101
102 # PDG code of a muon
103 muon_pdg_code = 13
104 sector = 0
105 phiBounds = [240 + 60.0 * sector % 360.0, 300 + 60.0 * sector % 360.0]
106
107 path.add_module("ParticleGun",
108 pdgCodes=[muon_pdg_code, -muon_pdg_code],
109 nTracks=1,
110 varyNTracks=False,
111 momentumGeneration='uniform',
112 momentumParams=[2.0, 4.0],
113 phiGeneration='uniform',
114 phiParams=phiBounds,
115 thetaGeneration='uniform',
116 thetaParams=[70., 110.])
117
118 else:
119 raise ValueError("Unknown generator module key " + generator_module)
120
121 # Flight time versus trigger setup #
122 # ################################ #
123
124 # Two alternative tof setups:
125 # * Relative to plane sets the Monte Carlo time of flight
126 # such that zero is on first intersection with the trigger counter plane
127 # * Relative to point sets the Monte Carlo time of flight
128 # such that zero is at the perigee relative to central point of the trigger counter plane
129 use_tof_relative_to_trigger_plane = True
130 use_tof_relative_to_trigger_point = not use_tof_relative_to_trigger_plane
131 if use_tof_relative_to_trigger_plane:
132 tof_mode = 1
133 elif use_tof_relative_to_trigger_point:
134 tof_mode = 2
135
136 # Also counter the adjustments the CRY generator makes to the production time in case it is used.
137 use_cry_generator = generator_module == "cry"
138
139 path.add_module('CDCCosmicSelector',
140 lOfCounter=30.,
141 wOfCounter=10.,
142 xOfCounter=0.3744,
143 yOfCounter=0.,
144 zOfCounter=-1.284,
145
146 # xOfCounter = -0.6,
147 # yOfCounter = 5,
148 # zOfCounter = 16,
149
150 TOF=tof_mode,
151 cryGenerator=use_cry_generator,
152 )
153
154 # Run simulation
155 path.add_module('FullSim',
156 ProductionCut=1000000.
157 )
158
159 path.add_module('CDCDigitizer',
160 # **{'2015AprRun' : True} # <- only hits in superlayer4 are digitized
161 # TrigTimeJitter=8., # <- simulate trigger timing jitter
162 AddTimeOfFlight=True,
163 AddInWirePropagationDelay=True,
164 CorrectForWireSag=True,
165 )
166
167 tracking.add_cdc_cr_track_finding(path,
168 trigger_point=(0.3744, 0, -1.284),
169 )
170
171 # path.add_module("TFCDC_TrackFinderCosmics")
172 # tracking.add_cdc_track_finding(path)
173 # path.add_module("TrackFinderMCTruthRecoTracks")
174
175 if use_tof_relative_to_trigger_plane:
176 path.add_module("PlaneTriggerTrackTimeEstimator",
177 triggerPlanePosition=[0.3744, 0, -1.284],
178 triggerPlaneDirection=[0, 1, 0],
179 useFittedInformation=False)
180
181 path.add_module("SetupGenfitExtrapolation")
182 path.add_module("DAFRecoFitter")
183
184 # Create a harvesting module for the time of flight
185 @save_pull_analysis(part_name="time", aux_names=["wire_idx", ], groupby=["layer_sidx"], outlier_z_score=3,)
186 @save_pull_analysis(part_name="time", aux_names=["layer_sidx", "time_truth", ], outlier_z_score=3,)
187 @save_pull_analysis(part_name="drift_length", aux_names=["layer_sidx", ], outlier_z_score=3,)
188 @save_tree(name="tree")
189 @harvest(foreach="RecoTracks", show_results=True)
190 def harvest_cdc_cr(reco_track):
191 fit_status = reco_track.getTrackFitStatus()
192 if not fit_status or not fit_status.isFitConverged():
193 return
194
195 track_points = reco_track.getHitPointsWithMeasurement()
196 for track_point in track_points:
197 kalman_fitter_info = track_point.getKalmanFitterInfo()
198 try:
199 measurement_on_plane = kalman_fitter_info.getFittedState()
200 except BaseException:
201 # Sometimes a throw happens from Genfit, skip track point
202 continue
203
204 abs_measurement = track_point.getRawMeasurement(0)
205 daf_weight = kalman_fitter_info.getWeights().at(0)
206 cdc_hit = abs_measurement.getCDCHit()
207 wire_id = Belle2.WireID(cdc_hit.getID())
208 layer_cidx = wire_id.getICLayer()
209 wire_idx = wire_id.getIWire()
210 cdc_sim_hit = cdc_hit.getRelated("CDCSimHits")
211
212 fit_data = get_fit_data(measurement_on_plane)
213
214 flight_time = fit_data["flight_time_estimate"]
215 drift_time = DriftTimeUtil.getDriftTime(fit_data["drift_length_estimate"],
216 layer_cidx,
217 fit_data["rl"],
218 fit_data["alpha_estimate"],
219 fit_data["theta_estimate"],
220 )
221
222 prop_time = DriftTimeUtil.getPropTime(wire_id, fit_data["z_estimate"])
223 time_walk = DriftTimeUtil.getTimeWalk(wire_id, cdc_hit.getADCCount())
224
225 reco_time = flight_time + drift_time + prop_time + time_walk
226 smear = True
227 measured_time = DriftTimeUtil.getMeasuredTime(wire_id, cdc_hit.getTDCCount(), smear)
228
229 # MC Information
230 flight_time_truth = cdc_sim_hit.getFlightTime()
231 drift_length_truth = cdc_sim_hit.getDriftLength()
232 incoming_arm = -1 if flight_time_truth < 0 else 1
233
234 if drift_time > 1000:
235 # There seems to be a bug in the CDCGeometryPar::getDriftTime function
236 continue
237 if drift_length_truth < 0.1 or drift_length_truth > 0.5:
238 continue
239 if daf_weight < 0.5:
240 continue
241
242 yield dict(drift_length_estimate=fit_data["drift_length_estimate"],
243 drift_length_variance=fit_data["drift_length_variance"],
244 drift_length_truth=drift_length_truth,
245 time_truth=reco_time,
246 time_estimate=measured_time,
247 time_delta=reco_time - measured_time,
248
249 flight_time_estimate=flight_time,
250 drift_time_estimate=drift_time,
251 prop_time_estimate=prop_time,
252 time_walk_estimate=time_walk,
253
254 rl_estimate=fit_data["rl"],
255 alpha_estimate=fit_data["alpha_estimate"],
256 theta_estimate=fit_data["theta_estimate"],
257 daf_weight=daf_weight,
258 flight_time_truth=flight_time_truth,
259 incoming_arm=incoming_arm,
260 layer_sidx=incoming_arm * layer_cidx,
261 wire_idx=wire_idx,
262 tdc=cdc_hit.getTDCCount(),
263 )
264
265 path.add_module(harvest_cdc_cr)
266
267 # Create a harvesting module for the time of flight mc comparision
268 @save_pull_analysis(part_name="flight_time", groupby=[None, "incoming_arm"], outlier_z_score=5,)
269 @save_tree(name="tree")
270 @harvest(foreach="RecoTracks", show_results=True)
271 def harvest_flight_times(reco_track):
272 track_points = reco_track.getHitPointsWithMeasurement()
273 for track_point in track_points:
274 kalman_fitter_info = track_point.getKalmanFitterInfo()
275 try:
276 measurement_on_plane = kalman_fitter_info.getFittedState()
277 except BaseException:
278 # Sometimes a throw happens from Genfit, skip track point
279 continue
280 abs_measurement = track_point.getRawMeasurement(0)
281 cdc_hit = abs_measurement.getCDCHit()
282 cdc_sim_hit = cdc_hit.getRelated("CDCSimHits")
283
284 flight_time_estimate = measurement_on_plane.getTime()
285 flight_time_truth = cdc_sim_hit.getFlightTime()
286 incoming_arm = flight_time_truth < 0
287
288 yield dict(flight_time_estimate=flight_time_estimate,
289 flight_time_truth=flight_time_truth,
290 incoming_arm=incoming_arm)
291
292 # path.add_module(harvest_flight_times)
293
294 # generate events
295 basf2.print_path(path)
296 basf2.process(path)
297
298 # show call statistics
299 print(basf2.statistics)
300
301
302# Helper functions #
303# ################ #
304def get_fit_data(measurement_on_plane):
305 flight_time = measurement_on_plane.getTime()
306 pos = measurement_on_plane.getPos()
307 mom = measurement_on_plane.getMom()
308 wirepos = measurement_on_plane.getPlane().getO()
309 alpha = wirepos.DeltaPhi(mom)
310 theta = mom.Theta()
311
312 state_on_plane = measurement_on_plane.getState() # TVectorD
313 cov_on_plane = measurement_on_plane.getCov() # TMatrixDSym
314
315 # The state at index 3 is equivalent to the drift length
316 drift_length = abs(state_on_plane(3))
317 drift_length_var = abs(cov_on_plane(3, 3))
318 rl_sign = 1 if state_on_plane(3) > 0 else -1
319 rl = 1 if state_on_plane(3) > 0 else 0
320
321 return dict(flight_time_estimate=flight_time,
322 x_estimate=pos.X(),
323 y_estimate=pos.Y(),
324 z_estimate=pos.Z(),
325 px_estimate=mom.X(),
326 py_estimate=mom.Y(),
327 pz_estimate=mom.Z(),
328 wire_x_estimate=wirepos.X(),
329 wire_y_estimate=wirepos.Y(),
330 wire_z_estimate=wirepos.Z(),
331 alpha_estimate=alpha,
332 theta_estimate=theta,
333 drift_length_estimate=drift_length,
334 drift_length_variance=drift_length_var,
335 rl=rl,
336 rl_sign=rl_sign,
337 )
338
339
340if __name__ == "__main__":
341 main()
Class to identify a wire inside the CDC.
Definition: WireID.h:34
Definition: main.py:1
Structure to expose some drift time and length functions from the CDCGeometryPar to Python.
Definition: DriftTimeUtil.h:22