Belle II Software  release-06-02-00
Fitter_OldVsNew.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 import basf2
13 import numpy as np
14 import os
15 import glob
16 from ROOT import gROOT, Belle2
17 gROOT.ProcessLine("gErrorIgnoreLevel = 4000;") # ignore endless root errors for background files...
18 
19 """
20 Compare the CDCTrigger2DFitterModule output with the old TRGCDCModule,
21 to confirm that the behaviour of the new code is correct.
22 """
23 
24 # ------------ #
25 # user options #
26 # ------------ #
27 
28 # general options
29 seed = 1
30 evtnum = 1000
31 clock = False
32 particlegun_params = {
33  'pdgCodes': [-13, 13],
34  'nTracks': 1,
35  'momentumGeneration': 'inversePt',
36  'momentumParams': [0.3, 10.],
37  'thetaGeneration': 'uniform',
38  'thetaParams': [35, 123],
39  'phiGeneration': 'uniform',
40  'phiParams': [0, 360],
41  'vertexGeneration': 'fixed',
42  'xVertexParams': [0.],
43  'yVertexParams': [0.],
44  'zVertexParams': [0.]}
45 usebkg = False
46 bkgdir = '/sw/belle2/bkg/'
47 
48 # ------------------------- #
49 # create path up to trigger #
50 # ------------------------- #
51 
52 # set random seed
53 basf2.set_random_seed(seed)
54 # suppress messages and warnings during processing:
55 basf2.set_log_level(basf2.LogLevel.ERROR)
56 
57 main = basf2.create_path()
58 
59 main.add_module('EventInfoSetter', evtNumList=evtnum)
60 main.add_module('Progress')
61 main.add_module('Gearbox')
62 main.add_module('Geometry', components=['CDC',
63  'MagneticFieldConstant4LimitedRCDC'])
64 particlegun = basf2.register_module('ParticleGun')
65 particlegun.param(particlegun_params)
66 main.add_module(particlegun)
67 main.add_module('FullSim')
68 if usebkg:
69  bkgmixer = basf2.register_module('BeamBkgMixer')
70  bkgfiles = glob.glob(os.path.join(bkgdir, '*[!(PXD)(ECL)]??.root'))
71  bkgmixer.param('backgroundFiles', bkgfiles)
72  bkgmixer.param('components', ['CDC'])
73  main.add_module(bkgmixer)
74 cdcdigitizer = basf2.register_module('CDCDigitizer')
75 if clock:
76  cdcdigitizer.param('TrigTimeJitter', 32.)
77 main.add_module(cdcdigitizer)
78 
79 # ----------- #
80 # CDC trigger #
81 # ----------- #
82 
83 xtsimple = False
84 fitdrift = True
85 
86 trgcdc = basf2.register_module('TRGCDC')
87 simMode = 0 # 0: full trigger, 1: only TSF
88 if clock:
89  simMode += 2 # 2: full trigger (TSF with clock), 3: only TSF with clock
90 trgcdc_params = {
91  'ConfigFile': Belle2.FileSystem.findFile("data/trg/cdc/TRGCDCConfig_0_20101111.dat"),
92  'InnerTSLUTFile': Belle2.FileSystem.findFile("data/trg/cdc/innerLUT_v3.0.coe"),
93  'OuterTSLUTFile': Belle2.FileSystem.findFile("data/trg/cdc/outerLUT_v3.0.coe"),
94  'SimulationMode': 1, # only fast simulation
95  'FastSimulationMode': simMode,
96  'HoughFinderMappingFileMinus': Belle2.FileSystem.findFile("data/trg/cdc/HoughMappingMinus20160223.dat"),
97  'HoughFinderMappingFilePlus': Belle2.FileSystem.findFile("data/trg/cdc/HoughMappingPlus20160223.dat"),
98  'DebugLevel': 0,
99  'Fitter3Ds2DFitDrift': fitdrift,
100  'Fitter3DsXtSimple': xtsimple,
101  '2DfinderCollection': 'TRGCDC2DFinderTracks'}
102 trgcdc.param(trgcdc_params)
103 if clock:
104  trgcdc.param('inputCollection', 'CDCHits4Trg')
105 main.add_module(trgcdc, logLevel=basf2.LogLevel.INFO)
106 
107 # fitters
108 main.add_module('CDCTriggerETF')
109 main.add_module('CDCTrigger2DFitter', logLevel=basf2.LogLevel.INFO,
110  minHits=2, useDriftTime=fitdrift, xtSimple=xtsimple,
111  outputCollectionName='FitterTracks') # to distinguish from old output
112 main.add_module('CDCTrigger3DFitter', logLevel=basf2.LogLevel.INFO,
113  xtSimple=xtsimple,
114  inputCollectionName='FitterTracks',
115  outputCollectionName='Fitter3DTracks') # to distinguish from old output
116 
117 
118 # ----------- #
119 # test module #
120 # ----------- #
121 
122 class TestModule(basf2.Module):
123  """
124  test module to compare the output of TRGCDC and CDCTrigger2DFitter/CDCTrigger3DFitter
125  """
126 
127  def event(self):
128  """
129  give info for both output lists and warnings in the case of mismatches
130  """
131  phiMC = Belle2.PyStoreArray("MCParticles")[0].getMomentum().Phi()
132  oldTracks = Belle2.PyStoreArray("Trg3DFitterTracks")
133  newTracks = Belle2.PyStoreArray("Fitter3DTracks")
134  if oldTracks.getEntries() == newTracks.getEntries():
135  basf2.B2INFO("%d tracks" % oldTracks.getEntries())
136  else:
137  basf2.B2WARNING("old version: %d, new version: %d" %
138  (oldTracks.getEntries(), newTracks.getEntries()))
139  for i in range(max(oldTracks.getEntries(), newTracks.getEntries())):
140  if i < oldTracks.getEntries():
141  ptfactor = 0.3 * 1.5 / 100 * 222.376063
142  oldString = "phi %.3f pt %.3f charge %d chi2 %.3f z %.3f cot %.3f chi2 %.3f" % \
143  (oldTracks[i].getPhi0() * 180. / np.pi,
144  oldTracks[i].getTransverseMomentum(1.5) / ptfactor,
145  oldTracks[i].getChargeSign(),
146  oldTracks[i].getChi2D(),
147  oldTracks[i].getZ0(),
148  oldTracks[i].getCotTheta(),
149  oldTracks[i].getChi3D())
150  else:
151  oldString = "no track"
152  if i < newTracks.getEntries():
153  newString = "phi %.3f pt %.3f charge %d chi2 %.3f z %.3f cot %.3f chi2 %.3f" % \
154  (newTracks[i].getPhi0() * 180. / np.pi,
155  newTracks[i].getTransverseMomentum(1.5),
156  newTracks[i].getChargeSign(),
157  newTracks[i].getChi2D(),
158  newTracks[i].getZ0(),
159  newTracks[i].getCotTheta(),
160  newTracks[i].getChi3D())
161  else:
162  newString = "no track"
163  if oldString == newString:
164  basf2.B2INFO(oldString)
165  else:
166  basf2.B2WARNING("old: " + oldString)
167  basf2.B2WARNING("new: " + newString)
168 
169 
170 main.add_module(TestModule(), logLevel=basf2.LogLevel.INFO)
171 
172 # Process events
173 basf2.process(main)
174 
175 # Print call statistics
176 print(basf2.statistics)
static std::string findFile(const std::string &path, bool silent=false)
Search for given file or directory in local or central release directory, and return absolute path if...
Definition: FileSystem.cc:145
a (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:56