Belle II Software development
peelers.py
1
8
9"""
10A set of common purose translators from complex framework objects to flat dictionaries
11"""
12
13import functools
14
15import math
16import numpy as np
17
18import ROOT
19
20from ROOT import Belle2
21
22import basf2
23from tracking.validation.tolerate_missing_key_formatter import TolerateMissingKeyFormatter
24
25from tracking.validation.utilities import getObjectList
26from tracking.validation.utilities import getHelixFromMCParticle
27
28Belle2.RecoTrack.getRightLeftInformation["Belle2::CDCHit"]
29ROOT.gSystem.Load("libtracking")
30
31formatter = TolerateMissingKeyFormatter()
32
33
34
35
36def format_crop_keys(peel_func):
37 @functools.wraps(peel_func)
38 def peel_func_formatted_keys(*args, key="{part_name}", **kwargs):
39 crops = peel_func(*args, **kwargs, key=key)
40 if key and hasattr(crops, 'items'):
41 crops_with_formatted_keys = dict()
42 for part_name, value in list(crops.items()):
43 formatted_key = formatter.format(key, part_name=part_name)
44 crops_with_formatted_keys[formatted_key] = value
45 return crops_with_formatted_keys
46
47 else:
48 return crops
49
50 return peel_func_formatted_keys
51
52
53
54
55@format_crop_keys
56def peel_mc_particle(mc_particle, key="{part_name}"):
57 if mc_particle:
58 helix = getHelixFromMCParticle(mc_particle)
59 momentum = mc_particle.getMomentum()
60 vertex = mc_particle.getVertex()
61 charge = mc_particle.getCharge()
62 pdg_code = mc_particle.getPDG()
63 is_primary = bool(mc_particle.hasStatus(Belle2.MCParticle.c_PrimaryParticle))
64
65 decay_vertex = mc_particle.getDecayVertex()
66 number_of_daughters = mc_particle.getNDaughters()
67 status = mc_particle.getStatus()
68
69 return dict(
70 # At origin assuming perfect magnetic field
71 d0_truth=helix.getD0(),
72 phi0_truth=helix.getPhi0() % (2.0 * math.pi),
73 omega_truth=helix.getOmega(),
74 z0_truth=helix.getZ0(),
75 tan_lambda_truth=helix.getTanLambda(),
76
77 # At the vertex position
78 pt_truth=momentum.Rho(),
79 px_truth=momentum.X(),
80 py_truth=momentum.Y(),
81 pz_truth=momentum.Z(),
82 x_truth=vertex.X(),
83 y_truth=vertex.Y(),
84 z_truth=vertex.Z(),
85
86 decay_vertex_radius_truth=decay_vertex.R(),
87 decay_vertex_x_truth=decay_vertex.X(),
88 decay_vertex_y_truth=decay_vertex.Y(),
89 decay_vertex_z_truth=decay_vertex.Z(),
90 number_of_daughters_truth=number_of_daughters,
91 status_truth=status,
92
93 # MC Particle information
94 charge_truth=charge,
95 pdg_code_truth=pdg_code,
96 is_primary=is_primary,
97 )
98
99 else:
100 nan = float('nan')
101 return dict(
102 # At origin assuming perfect magnetic field
103 d0_truth=nan,
104 phi0_truth=nan,
105 omega_truth=nan,
106 z0_truth=nan,
107 tan_lambda_truth=nan,
108
109 # At the vertex position
110 pt_truth=nan,
111 px_truth=nan,
112 py_truth=nan,
113 pz_truth=nan,
114 x_truth=nan,
115 y_truth=nan,
116 z_truth=nan,
117
118 decay_vertex_radius_truth=nan,
119 decay_vertex_x_truth=nan,
120 decay_vertex_y_truth=nan,
121 decay_vertex_z_truth=nan,
122 number_of_daughters_truth=nan,
123 status_truth=nan,
124
125
126 # MC Particle information
127 charge_truth=nan,
128 pdg_code_truth=0,
129 is_primary=False,
130 )
131
132
133
134
135@format_crop_keys
136def peel_reco_track_hit_content(reco_track, key="{part_name}"):
137 nan = float('nan')
138
139 first_pxd_layer = nan
140 last_pxd_layer = nan
141 first_svd_layer = nan
142 last_svd_layer = nan
143 first_cdc_layer = nan
144 last_cdc_layer = nan
145 sum_of_weights = nan
146 last_fit_layer = nan
147
148 if reco_track:
149 sum_of_weights = 0
150 last_fit_layer = 0
151
152 n_cdc_hits = reco_track.getNumberOfCDCHits()
153 n_svd_hits = reco_track.getNumberOfSVDHits()
154 n_pxd_hits = reco_track.getNumberOfPXDHits()
155 ndf = 2 * n_pxd_hits + n_svd_hits + n_cdc_hits
156
157 pxd_hits = [hit.getSensorID().getLayerNumber() for hit in getObjectList(reco_track.getPXDHitList())]
158 if pxd_hits:
159 first_pxd_layer = min(pxd_hits)
160 last_pxd_layer = max(pxd_hits)
161 svd_hits = [hit.getSensorID().getLayerNumber() for hit in getObjectList(reco_track.getSVDHitList())]
162 if svd_hits:
163 first_svd_layer = min(svd_hits)
164 last_svd_layer = max(svd_hits)
165 cdc_hits = [hit.getICLayer() for hit in getObjectList(reco_track.getCDCHitList())]
166 if cdc_hits:
167 first_cdc_layer = min(cdc_hits)
168 last_cdc_layer = max(cdc_hits)
169 for hit_info in getObjectList(reco_track.getRelationsWith("RecoHitInformations")):
170 track_point = reco_track.getCreatedTrackPoint(hit_info)
171 if track_point:
172 fitted_state = track_point.getFitterInfo()
173 if fitted_state:
174 W = max(fitted_state.getWeights())
175 sum_of_weights += W
176 if W > 0.5 and hit_info.getTrackingDetector() == Belle2.RecoHitInformation.c_CDC:
177 layer = hit_info.getRelated("CDCHits").getICLayer()
178 if layer > last_fit_layer:
179 last_fit_layer = layer
180 return dict(
181 n_pxd_hits=n_pxd_hits,
182 n_svd_hits=n_svd_hits,
183 n_cdc_hits=n_cdc_hits,
184 n_hits=n_pxd_hits + n_svd_hits + n_cdc_hits,
185 ndf_hits=ndf,
186 first_pxd_layer=first_pxd_layer,
187 last_pxd_layer=last_pxd_layer,
188 first_svd_layer=first_svd_layer,
189 last_svd_layer=last_svd_layer,
190 first_cdc_layer=first_cdc_layer,
191 last_cdc_layer=last_cdc_layer,
192 sum_of_weights=sum_of_weights,
193 last_fit_layer=last_fit_layer
194 )
195 else:
196 return dict(
197 n_pxd_hits=nan,
198 n_svd_hits=nan,
199 n_cdc_hits=nan,
200 n_hits=nan,
201 ndf_hits=nan,
202 first_pxd_layer=first_pxd_layer,
203 last_pxd_layer=last_pxd_layer,
204 first_svd_layer=first_svd_layer,
205 last_svd_layer=last_svd_layer,
206 first_cdc_layer=first_cdc_layer,
207 last_cdc_layer=last_cdc_layer,
208 sum_of_weights=sum_of_weights,
209 last_fit_layer=last_fit_layer
210 )
211
212
213
214
215@format_crop_keys
216def peel_reco_track_seed(reco_track, key="{part_name}"):
217 if reco_track:
218 seed_fit_result = get_seed_track_fit_result(reco_track)
219 return peel_track_fit_result(seed_fit_result, key="seed_{part_name}")
220
221 else:
222 return peel_track_fit_result(None, key="seed_{part_name}")
223
224
225
226
227@format_crop_keys
228def peel_event_info(event_meta_data, key="{part_name}"):
229 return dict(
230 event_number=event_meta_data.getEvent(),
231 run_number=event_meta_data.getRun(),
232 experiment_number=event_meta_data.getExperiment(),
233 )
234
235
236
237
238@format_crop_keys
239def peel_store_array_info(item, key="{part_name}"):
240 if item:
241 return dict(
242 store_array_number=item.getArrayIndex()
243 )
244 else:
245 return dict(
246 store_array_number=float("nan")
247 )
248
249
250
251
252@format_crop_keys
253def peel_store_array_size(array_name, key="{part_name}"):
254 array = Belle2.PyStoreArray(array_name)
255 return {str(array_name) + "_size": array.getEntries() if array else 0}
256
257
258
259
260@format_crop_keys
261def peel_event_level_tracking_info(event_level_tracking_info, key="{part_name}"):
262 if not event_level_tracking_info:
263 return dict(
264 has_vxdtf2_failure_flag=False,
265 has_svdckf_failure_flag=False,
266 has_pxdckf_failure_flag=False,
267 has_unspecified_trackfinding_failure=False,
268 )
269 return dict(has_vxdtf2_failure_flag=event_level_tracking_info.hasVXDTF2AbortionFlag(),
270 has_svdckf_failure_flag=event_level_tracking_info.hasSVDCKFAbortionFlag(),
271 has_pxdckf_failure_flag=event_level_tracking_info.hasPXDCKFAbortionFlag(),
272 has_unspecified_trackfinding_failure=event_level_tracking_info.hasUnspecifiedTrackFindingFailure(),
273 )
274
275
276
277
278@format_crop_keys
279def peel_quality_indicators(reco_track, key="{part_name}"):
280 nan = float("nan")
281
282 svd_qi = nan
283 cdc_qi = nan
284 qi = nan
285
286 if reco_track:
287 qi = reco_track.getQualityIndicator()
288 space_point_track_cand = reco_track.getRelated('SPTrackCands')
289
290 # adjust relations if SVD->CDC CKF enabled
291 if not space_point_track_cand:
292 svd_track_cand = reco_track.getRelated('SVDRecoTracks')
293 if not svd_track_cand:
294 svd_cdc_track_cand = reco_track.getRelated('SVDCDCRecoTracks')
295 if svd_cdc_track_cand:
296 svd_track_cand = svd_cdc_track_cand.getRelated('SVDRecoTracks')
297 if not svd_track_cand:
298 temp_svd_track_cand = svd_cdc_track_cand.getRelated('SVDPlusCDCStandaloneRecoTracks')
299 if temp_svd_track_cand:
300 svd_track_cand = temp_svd_track_cand.getRelated('SVDRecoTracks')
301 if svd_track_cand:
302 space_point_track_cand = svd_track_cand.getRelated('SPTrackCands')
303
304 if space_point_track_cand:
305 svd_qi = space_point_track_cand.getQualityIndicator()
306
307 cdc_track_cand = reco_track.getRelated('CDCRecoTracks')
308 if not cdc_track_cand:
309 svd_cdc_track_cand = reco_track.getRelated('SVDCDCRecoTracks')
310 if svd_cdc_track_cand:
311 cdc_track_cand = svd_cdc_track_cand.getRelated('CDCRecoTracks')
312 if not cdc_track_cand:
313 cdc_track_cand = svd_cdc_track_cand.getRelated('CKFCDCRecoTracks')
314 if not cdc_track_cand:
315 temp_cdc_track_cand = svd_cdc_track_cand.getRelated('SVDPlusCDCStandaloneRecoTracks')
316 if temp_cdc_track_cand:
317 cdc_track_cand = temp_cdc_track_cand.getRelated('CDCRecoTracks')
318
319 if cdc_track_cand:
320 cdc_qi = cdc_track_cand.getQualityIndicator()
321
322 crops = dict(
323 quality_indicator=qi,
324 svd_quality_indicator=svd_qi,
325 cdc_qualityindicator=cdc_qi,
326 )
327
328 return crops
329
330
331
332
333@format_crop_keys
334def peel_trackfinder(reco_track, key="{part_name}"):
335 used_CDCTrackFinder = False
336 used_VXDTrackFinder = False
337 used_SVDHough = False
338 used_SVDtoCDCCKF = False
339 used_ECLtoCDCCKF = False
340 used_CDCtoSVDCKF = False
341
342 if reco_track:
343 if reco_track.getNumberOfSVDHits() > 0:
344 info = get_reco_hit_information(reco_track, reco_track.getSVDHitList()[0])
345 svd_tf = info.getFoundByTrackFinder()
346 used_VXDTrackFinder = svd_tf == Belle2.RecoHitInformation.c_VXDTrackFinder
347 used_CDCtoSVDCKF = svd_tf == Belle2.RecoHitInformation.c_CDCtoSVDCKF
348 used_SVDHough = svd_tf == Belle2.RecoHitInformation.c_SVDHough
349
350 if reco_track.getNumberOfCDCHits() > 0:
351 info = get_reco_hit_information(reco_track, reco_track.getCDCHitList()[0])
352 cdc_tf = info.getFoundByTrackFinder()
353 used_CDCTrackFinder = cdc_tf == Belle2.RecoHitInformation.c_CDCTrackFinder
354 used_SVDtoCDCCKF = cdc_tf == Belle2.RecoHitInformation.c_SVDtoCDCCKF
355 used_ECLtoCDCCKF = cdc_tf == Belle2.RecoHitInformation.c_ECLtoCDCCKF
356
357 crops = dict(
358 foundby_CDCTrackFinder=used_CDCTrackFinder,
359 foundby_VXDTrackFinder=used_VXDTrackFinder,
360 foundby_SVDHough=used_SVDHough,
361 foundby_SVDtoCDCCKF=used_SVDtoCDCCKF,
362 foundby_CDCtoSVDCKF=used_CDCtoSVDCKF,
363 foundby_ECLtoCDCCKF=used_ECLtoCDCCKF,
364 )
365
366 return crops
367
368
369
370
371@format_crop_keys
372def peel_fit_status(reco_track, key="{part_name}"):
373 nan = float("nan")
374
375 crops = dict(
376 is_fitted=nan,
377 fit_pion_ok=nan,
378 ndf_pion=nan,
379 fit_muon_ok=nan,
380 ndf_muon=nan,
381 fit_electron_ok=nan,
382 ndf_electron=nan,
383 fit_proton_ok=nan,
384 ndf_proton=nan,
385 fit_kaon_ok=nan,
386 ndf_kaon=nan,
387 )
388
389 if reco_track:
390 crops["is_fitted"] = reco_track.wasFitSuccessful()
391
392 for rep in getObjectList(reco_track.getRepresentations()):
393 was_successful = reco_track.wasFitSuccessful(rep)
394 pdg_code = rep.getPDG()
395
396 for crop in crops.keys():
397 if crop.startswith("fit_"):
398 particle_name = crop.split("_")[1]
399 if abs(getattr(Belle2.Const, particle_name).getPDGCode()) == abs(pdg_code):
400 crops[crop] = was_successful
401
402 if was_successful:
403 crops[f"ndf_{particle_name}"] = reco_track.getTrackFitStatus(rep).getNdf()
404
405 return crops
406
407
408
409
410@format_crop_keys
411def peel_track_fit_result(track_fit_result, key="{part_name}"):
412 nan = float("nan")
413 if track_fit_result:
414 cov6 = track_fit_result.getCovariance6()
415 mom = track_fit_result.getMomentum()
416 pos = track_fit_result.getPosition()
417
418 pt_estimate = mom.Rho()
419
420 pt_variance = np.divide(
421 mom.X() ** 2 * cov6(3, 3) + mom.Y() ** 2 * cov6(4, 4) - 2 * mom.X() * mom.Y() * cov6(3, 4),
422 mom.Perp2()
423 )
424
425 pt_resolution = np.divide(pt_variance, pt_estimate)
426
427 fit_crops = dict(
428 has_trackFitResult=True,
429 d0_estimate=track_fit_result.getD0(),
430 d0_variance=track_fit_result.getCov()[0],
431 phi0_estimate=track_fit_result.getPhi() % (2.0 * math.pi),
432 phi0_variance=track_fit_result.getCov()[5],
433 omega_estimate=track_fit_result.getOmega(),
434 omega_variance=track_fit_result.getCov()[9],
435 z0_estimate=track_fit_result.getZ0(),
436 z0_variance=track_fit_result.getCov()[12],
437 tan_lambda_estimate=track_fit_result.getCotTheta(),
438 tan_lambda_variance=track_fit_result.getCov()[14],
439
440 x_estimate=pos.X(),
441 x_variance=cov6(0, 0),
442 y_estimate=pos.Y(),
443 y_variance=cov6(1, 1),
444 z_estimate=pos.Z(),
445 z_variance=cov6(2, 2),
446
447 pt_estimate=pt_estimate,
448 pt_variance=pt_variance,
449 pt_resolution=pt_resolution,
450
451 track_charge=track_fit_result.getChargeSign(),
452 b_field=Belle2.BFieldManager.getField(pos).Z(),
453
454 px_estimate=mom.X(),
455 px_variance=cov6(3, 3),
456 py_estimate=mom.Y(),
457 py_variance=cov6(4, 4),
458 pz_estimate=mom.Z(),
459 pz_variance=cov6(5, 5),
460
461 p_value=track_fit_result.getPValue(),
462 )
463
464 else:
465 fit_crops = dict(
466 has_trackFitResult=False,
467 d0_estimate=nan,
468 d0_variance=nan,
469 phi0_estimate=nan,
470 phi0_variance=nan,
471 omega_estimate=nan,
472 omega_variance=nan,
473 z0_estimate=nan,
474 z0_variance=nan,
475 tan_lambda_estimate=nan,
476 tan_lambda_variance=nan,
477
478 x_estimate=nan,
479 x_variance=nan,
480 y_estimate=nan,
481 y_variance=nan,
482 z_estimate=nan,
483 z_variance=nan,
484
485 pt_estimate=nan,
486 pt_variance=nan,
487 pt_resolution=nan,
488
489 track_charge=nan,
490 b_field=nan,
491
492 px_estimate=nan,
493 px_variance=nan,
494 py_estimate=nan,
495 py_variance=nan,
496 pz_estimate=nan,
497 pz_variance=nan,
498
499 p_value=nan,
500 )
501
502 return fit_crops
503
504
505
506
507def get_reco_hit_information(reco_track, hit):
508 """Helper function for getting the correct reco hit info"""
509 for info in getObjectList(hit.getRelationsFrom("RecoHitInformations")):
510 if info.getRelatedFrom(reco_track.getArrayName()) == reco_track:
511 return info
512
513
514
515@format_crop_keys
516def peel_subdetector_hit_efficiency(mc_reco_track, reco_track, key="{part_name}"):
517 def get_efficiency(detector_string):
518 if not reco_track or not mc_reco_track:
519 hit_efficiency = float("nan")
520 else:
521 mc_reco_hits = getObjectList(getattr(mc_reco_track, f"get{detector_string.upper()}HitList")())
522 if len(mc_reco_hits) == 0:
523 hit_efficiency = float('nan')
524 else:
525 mc_reco_hit_size = 0
526 hit_efficiency = 0.
527 for mc_reco_hit in mc_reco_hits:
528 info = get_reco_hit_information(mc_reco_track, mc_reco_hit)
529
530 if info.getFoundByTrackFinder() == Belle2.RecoHitInformation.c_MCTrackFinderAuxiliaryHit:
531 continue
532
533 mc_reco_hit_size += 1
534
535 for reco_hit in getObjectList(getattr(reco_track, f"get{detector_string.upper()}HitList")()):
536 if mc_reco_hit.getArrayIndex() == reco_hit.getArrayIndex():
537 hit_efficiency += 1
538 break
539
540 if not mc_reco_hit_size:
541 hit_efficiency = float('nan')
542 else:
543 hit_efficiency /= mc_reco_hit_size
544
545 return {f"{detector_string.lower()}_hit_efficiency": hit_efficiency}
546
547 return dict(**get_efficiency("CDC"), **get_efficiency("SVD"), **get_efficiency("PXD"))
548
549
550
551@format_crop_keys
552def peel_subdetector_hit_purity(reco_track, mc_reco_track, key="{part_name}"):
553 def get_efficiency(detector_string):
554 if not reco_track or not mc_reco_track:
555 hit_purity = float("nan")
556 else:
557 reco_hits = getObjectList(getattr(reco_track, f"get{detector_string.upper()}HitList")())
558 reco_hit_size = len(reco_hits)
559
560 if reco_hit_size == 0:
561 hit_purity = float('nan')
562 else:
563 hit_purity = 0.
564 for reco_hit in reco_hits:
565 for mc_reco_hit in getObjectList(getattr(mc_reco_track, f"get{detector_string.upper()}HitList")()):
566 if mc_reco_hit.getArrayIndex() == reco_hit.getArrayIndex():
567 hit_purity += 1
568 break
569
570 hit_purity /= reco_hit_size
571
572 return {f"{detector_string.lower()}_hit_purity": hit_purity}
573
574 return dict(**get_efficiency("CDC"), **get_efficiency("SVD"), **get_efficiency("PXD"))
575
576
577
578@format_crop_keys
579def peel_hit_information(hit_info, reco_track, key="{part_name}"):
580 nan = float("nan")
581
582 crops = dict(residual=nan,
583 residual_x=nan,
584 residual_y=nan,
585 residuals=nan,
586 weight=nan,
587 tracking_detector=hit_info.getTrackingDetector(),
588 use_in_fit=hit_info.useInFit(),
589 hit_time=nan,
590 layer_number=nan,
591 )
592
593 if hit_info.useInFit() and reco_track.hasTrackFitStatus():
594 track_point = reco_track.getCreatedTrackPoint(hit_info)
595 fitted_state = track_point.getFitterInfo()
596 if fitted_state:
597 try:
598 res_state = fitted_state.getResidual().getState()
599 crops["residual"] = np.sqrt(res_state.Norm2Sqr())
600 if res_state.GetNoElements() == 2:
601 crops["residual_x"] = res_state[0]
602 crops["residual_y"] = res_state[1]
603 weights = fitted_state.getWeights()
604 crops['weight'] = max(weights)
605 except BaseException:
606 pass
607
608 if hit_info.getTrackingDetector() == Belle2.RecoHitInformation.c_SVD:
609 hit = hit_info.getRelated("SVDClusters")
610 crops["hit_time"] = hit.getClsTime()
611 crops["layer_number"] = hit.getSensorID().getLayerNumber()
612 if hit_info.getTrackingDetector() == Belle2.RecoHitInformation.c_PXD:
613 hit = hit_info.getRelated("PXDClusters")
614 crops["layer_number"] = hit.getSensorID().getLayerNumber()
615 if hit_info.getTrackingDetector() == Belle2.RecoHitInformation.c_CDC:
616 hit = hit_info.getRelated("CDCHits")
617 crops["layer_number"] = hit.getICLayer()
618
619 return crops
620
621
622
623@format_crop_keys
624def peel_module_statistics(modules=None, key="{part_name}"):
625 if modules is None:
626 modules = []
627 module_stats = dict()
628
629 for module in basf2.statistics.modules:
630 if module.name in modules:
631 module_stats[str(module.name) + "_mem"] = module.memory_sum(basf2.statistics.EVENT)
632 module_stats[str(module.name) + "_time"] = module.time_sum(basf2.statistics.EVENT)
633 module_stats[str(module.name) + "_calls"] = module.calls(basf2.statistics.EVENT)
634
635 return module_stats
636
637
638
639
640
641def get_seed_track_fit_result(reco_track):
642 position = reco_track.getPositionSeed()
643 momentum = reco_track.getMomentumSeed()
644 cartesian_covariance = reco_track.getSeedCovariance()
645 charge_sign = (-1 if reco_track.getChargeSeed() < 0 else 1)
646 # It does not matter, which particle we put in here, so we just use a pion
647 particle_type = Belle2.Const.pion
648 p_value = float('nan')
649 b_field = Belle2.BFieldManager.getField(ROOT.Math.XYZVector(position)).Z() / Belle2.Unit.T
650 cdc_hit_pattern = 0
651 svd_hit_pattern = 0
652 # the value 0xFFFF will cause the TrackFitResult::getNDF() to return -1
653 ndf = 0xFFFF
654
655 track_fit_result = Belle2.TrackFitResult(
656 position,
657 momentum,
658 cartesian_covariance,
659 charge_sign,
660 particle_type,
661 p_value,
662 b_field,
663 cdc_hit_pattern,
664 svd_hit_pattern,
665 ndf,
666 )
667
668 return track_fit_result
669
670
671
672
673def is_correct_rl_information(cdc_hit, reco_track, hit_lookup):
674 rl_info = reco_track.getRightLeftInformation["const Belle2::CDCHit"](cdc_hit)
675 truth_rl_info = hit_lookup.getRLInfo(cdc_hit)
676
677 if rl_info == Belle2.RecoHitInformation.c_right and truth_rl_info == 1:
678 return True
679 if rl_info == Belle2.RecoHitInformation.c_left and truth_rl_info == 65535: # -1 as short
680 return True
681
682 return False
This class provides a set of constants for the framework.
Definition: Const.h:34
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
Values of the result of a track fit with a given particle hypothesis.
static void getField(const double *pos, double *field)
return the magnetic field at a given position.
Definition: BFieldManager.h:91