Belle II Software development
__init__.py
1
8
9import basf2 as b2
10from basf2.pickle_path import deserialize_path, deserialize_module, serialize_path, serialize_module
11
12import ROOT
13from ROOT import Belle2
14from caf.framework import Calibration
15
16from collections import OrderedDict
17
18b2.set_log_level(b2.LogLevel.INFO)
19UVWABC = [1, 2, 3, 4, 5, 6]
20
21
22def collect(calibration, basf2_args=None):
23 import pickle
24 import subprocess
25
26 if basf2_args is None:
27 basf2_args = []
28
29 main = b2.create_path()
30 main.add_module('RootInput', inputFileNames=calibration.input_files)
31 main.add_module('HistoManager', histoFileName='RootOutput.root')
32 main.add_path(calibration.pre_collector_path)
33 main.add_module(calibration.collector)
34
35 path_file_name = calibration.name + '.path'
36 with open(path_file_name, 'bw') as serialized_path_file:
37 pickle.dump(serialize_path(main), serialized_path_file)
38
39 subprocess.call(["basf2", "--execute-path", path_file_name] + basf2_args)
40
41
42def calibrate(calibration, input_file='RootOutput.root'):
43 for algo in calibration.algorithms:
44 algo.algorithm.setInputFileNames([input_file])
45 algo.algorithm.execute()
46 algo.algorithm.commit()
47
48
50 """ The generic Millepede calibration collector+algorithm """
51
52 def __init__(self,
53 components=None,
54 tracks=None,
55 particles=None,
56 vertices=None,
57 primary_vertices=None,
58 path=None):
59 """
60 components are the names of DB objects to calibrate (BeamSpot etc.)
61 tracks are collections of RecoTracks of fitted charged tracks (usually cosmic rays, no associated particle)
62 particles are names of ParticleLists with single charged particles
63 vertices are names ParticleLists with at least two body decays fitted with vertex constraint
64 primary_vertices are names of ParticleLists with at least two body decays fitted with beam+vertex constraint
65 """
66
67 if components is None:
68 components = []
69
70 # By default, at least the raw RecoTracks are fitted (let's specify its name explicitly here)
71 if (tracks is None) and (particles is None) and (vertices is None) and (primary_vertices is None):
72 tracks = ['RecoTracks']
73
74 if particles is None:
75 particles = []
76
77 if vertices is None:
78 vertices = []
79
80 if primary_vertices is None:
81 primary_vertices = []
82
83 if path is None:
84 path = b2.create_path()
85 path.add_module('Gearbox')
86 path.add_module('Geometry')
87 path.add_module('SetupGenfitExtrapolation', noiseBetheBloch=False, noiseCoulomb=False, noiseBrems=False)
88
89
91
92 self.path = path
93
94 self.collector = b2.register_module('MillepedeCollector')
95
96 self.parameters = OrderedDict()
97
98 self.commands = OrderedDict()
99
100 self.set_components(components)
101
102 self.collector.param('tracks', tracks)
103 self.collector.param('particles', particles)
104 self.collector.param('vertices', vertices)
105 self.collector.param('primaryVertices', primary_vertices)
106 self.collector.param('minPValue', 0.)
107 self.collector.param('useGblTree', True)
108
109 self.set_command('method diagonalization 3 0.1')
110 self.set_command('skipemptycons')
111 self.set_command('entries 10')
112 self.set_command('hugecut 50')
113 self.set_command('chiscut 30. 6.')
114 self.set_command('outlierdownweighting 3')
115 self.set_command('dwfractioncut 0.1')
116
117 def set_components(self, components):
118 """ Select db objects for calibration from list of their names """
119
120 self.components = components
121 calibrate_vertex = ((components == []) or ('BeamSpot' in components))
122
123 std_components = ROOT.vector('string')()
124 for component in components:
125 std_components.push_back(component)
126
127 self.algo.setComponents(std_components)
128 self.collector.param('components', components)
129 self.collector.param('calibrateVertex', calibrate_vertex)
130
131 def set_command(self, cmd_name, command=''):
132 """ Set command for Millepede steering """
133 if command is None:
134 self.commands[cmd_name] = ''
135 return
136
137 if command == '':
138 command = cmd_name
139 else:
140 command = cmd_name + ' ' + str(command)
141
142 cmd_words = str(command).split(' ')
143 self.commands[cmd_words[0]] = command
144
145 def pre_algo(self, algorithm, iteration):
146 """ Fcn to execute before algorithm... """
147 pass
148
149 def create(self, name, input_files):
150 """ Create the CAF Calibration object """
151
152 # We keep steering options for pede in python
153 # as long as possible for the user before passing to algo
154 for command in self.commands.values():
155 self.algo.steering().command(command)
156 # Add parameter fixing/setting
157 self.algo.steering().command('Parameters')
158 for command in self.parameters.values():
159 self.algo.steering().command(command)
160
161 # If you call this, you are going to use calibration framework most likely.
162 # This NOW does not anymore move/copy files from collector to algorithm, so
163 # we have to remeber where they are located at the time of creation
164 if not self.get_param('useGblTree'):
165 self.set_param(True, 'absFilePaths')
166
167 cal = Calibration(name, self.collector, self.algo, input_files)
168 cal.output_patterns.append('Belle2FileCatalog.xml')
169 cal.output_patterns.append('constraints.txt')
170 cal.output_patterns.append('*.mille')
171 # Adding in setup basf2 paths and functions for the collector and algorithm
172 cal.pre_collector_path = self.path
173 cal.pre_algorithms = self.pre_algo
174
175 return cal
176
177 def get_param(self, param, module=None):
178 """ Get parameter of the collector module or any module in path (given its name) """
179 module = self.get_module(module)
180 for mpi in module.available_params():
181 if mpi.name == param:
182 return mpi.values
183 return None
184
185 def set_param(self, value, param, module=None):
186 """ Set parameter of the collector module or any module in path (given its name) """
187 module = self.get_module(module)
188 if module is not None:
189 module.param(param, value)
190
191 def get_module(self, module=None):
192 """ Get collector module or any other module from path (given its name) """
193 if module is None or module == self.collector.name() or module == self.collector:
194 module = self.collector
195 else:
196 for mod in self.path.modules():
197 if mod.name() == module:
198 module = mod
199 break
200 if isinstance(module, b2.Module):
201 return module
202 return None
203
204 def __getstate__(self):
205 """ serialization """
206
207 state = self.__dict__.copy()
208 state['path'] = serialize_path(self.path)
209 state['collector'] = serialize_module(self.collector)
210 return state
211
212 def __setstate__(self, state):
213 """ de-serialization """
214
215 self.__dict__.update(state)
216 self.collector = deserialize_module(state['collector'])
217 self.path = deserialize_path(state['path'])
218
219 def fixGlobalParam(self, uid, element, param, value=0.):
220 """ Generic fcn to manipulate (fix to given value) parameters identified by db object id (uid),
221 element and parameter number """
222 label = Belle2.GlobalLabel()
223 label.construct(uid, element, param)
224 self.parameters[str(label.label())] = (str(label.label()) + ' ' + str(value) + ' -1.')
225
226 # CDC helpers --------------------------------------------------------------------------------------------------
227
228 def fixCDCLayerX(self, layer):
229 """ fix CDC layer X-shift """
230 self.fixGlobalParam(Belle2.CDCLayerAlignment.getGlobalUniqueID(), layer, Belle2.CDCLayerAlignment.layerX)
231
232 def fixCDCLayerY(self, layer):
233 """ fix CDC layer Y-shift """
234 self.fixGlobalParam(Belle2.CDCLayerAlignment.getGlobalUniqueID(), layer, Belle2.CDCLayerAlignment.layerY)
235
236 def fixCDCLayerRot(self, layer):
237 """ fix CDC layer Phi rotation """
238 self.fixGlobalParam(Belle2.CDCLayerAlignment.getGlobalUniqueID(), layer, Belle2.CDCLayerAlignment.layerPhi)
239
240 def fixCDCTimeWalk(self, board_id):
241 """ fix CDC time walk parameter for given board """
243
244 def fixCDCTimeZero(self, wire_id):
245 """ fix CDC tie zero parameter for given wire """
247
248 # VXD helpers --------------------------------------------------------------------------------------------------
249
250 def fixVXDid(self, layer, ladder, sensor, segment=0, parameters=[1, 2, 3, 4, 5, 6]):
251 """ Fix VXD element parameters by layer,ladder,sensor and segment number """
252 for param in parameters:
254 layer, ladder, sensor, segment).getID(), param)
255 """
256 {{"PXD.Ying"}, {Belle2::VxdID(1, 0, 0, 1)}},
257 {{"PXD.Yang"}, {Belle2::VxdID(1, 0, 0, 2)}},
258 {{"SVD.Pat"}, {Belle2::VxdID(3, 0, 0, 1)}},
259 {{"SVD.Mat"}, {Belle2::VxdID(3, 0, 0, 2)}}
260 """
261
262 def fixPXDYing(self, parameters=UVWABC):
263 """ fix PXD Ying half-shell """
264 self.fixVXDid(1, 0, 0, 1, parameters)
265
266 def fixPXDYang(self, parameters=UVWABC):
267 """ fix PXD Yang half-shell """
268 self.fixVXDid(1, 0, 0, 2, parameters)
269
270 def fixSVDPat(self, parameters=UVWABC):
271 """ fix SVD Pat half-shell """
272 self.fixVXDid(3, 0, 0, 1, parameters)
273
274 def fixSVDMat(self, parameters=UVWABC):
275 """ fix SVD Mat half-shell """
276 self.fixVXDid(3, 0, 0, 2, parameters)
277
278 # EKLM helpers --------------------------------------------------------------------------------------------------
279
280 def fixEKLMSector(self, endcap, layer, sector, parameters=[1, 2, 6]):
281 """ Fix EKLM sector parameters """
282
283 for ipar in parameters:
285 Belle2.EKLMElementID(endcap, layer, sector).getGlobalNumber(),
286 ipar)
287
288 def fixEKLMSegment(self, endcap, layer, sector, plane, segment, parameters=[1, 2, 6]):
289 """ Fix EKLM segment parameters """
290 for ipar in parameters:
292 Belle2.EKLMElementID(endcap, layer, sector, plane, segment).getGlobalNumber(),
293 ipar)
294
295 # BKLM helpers --------------------------------------------------------------------------------------------------
296
297 def fixBKLMModule(self, sector, layer, forward, parameters=UVWABC):
298 """ Fix a BKLM module """
299 for ipar in parameters:
300 self.fixGlobalParam(
302 Belle2.BKLMElementID(
303 forward,
304 sector,
305 layer).getID(),
306 ipar)
307
308 def fixBKLM(self, sectors=range(1, 9), layers=range(1, 16), forbackwards=[0, 1]):
309 """ Fix (all by default) BKLM modules """
310 for sector in sectors:
311 for layer in layers:
312 for forward in forbackwards:
313 self.fixBKLMModule(forward, sector, layer)
static unsigned short getGlobalUniqueID()
Get global unique identifier.
Definition: BKLMAlignment.h:63
static unsigned short getGlobalUniqueID()
Get global unique id.
static unsigned short getGlobalUniqueID()
Get global unique id.
Definition: CDCTimeWalks.h:155
static unsigned short getGlobalUniqueID()
Get global unique id.
Definition: CDCTimeZeros.h:140
static unsigned short getGlobalUniqueID()
Get global unique identifier.
Definition: EKLMAlignment.h:63
Class to convert to/from global labels for Millepede II to/from detector & parameter identificators.
Definition: GlobalLabel.h:41
Class implementing Millepede calibration algorithm.
static unsigned short getGlobalUniqueID()
Get global unique id.
Definition: VXDAlignment.h:47
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
commands
Commands for pede steering.
Definition: __init__.py:98
def fixEKLMSegment(self, endcap, layer, sector, plane, segment, parameters=[1, 2, 6])
Definition: __init__.py:288
def get_param(self, param, module=None)
Definition: __init__.py:177
def set_command(self, cmd_name, command='')
Definition: __init__.py:131
components
db objects for calibration
Definition: __init__.py:120
def __setstate__(self, state)
Definition: __init__.py:212
def fixCDCLayerX(self, layer)
Definition: __init__.py:228
def fixCDCLayerY(self, layer)
Definition: __init__.py:232
def fixVXDid(self, layer, ladder, sensor, segment=0, parameters=[1, 2, 3, 4, 5, 6])
Definition: __init__.py:250
def fixCDCLayerRot(self, layer)
Definition: __init__.py:236
def fixCDCTimeWalk(self, board_id)
Definition: __init__.py:240
def create(self, name, input_files)
Definition: __init__.py:149
def pre_algo(self, algorithm, iteration)
Definition: __init__.py:145
def fixBKLM(self, sectors=range(1, 9), layers=range(1, 16), forbackwards=[0, 1])
Definition: __init__.py:308
def fixCDCTimeZero(self, wire_id)
Definition: __init__.py:244
path
pre-collector path
Definition: __init__.py:92
def fixPXDYang(self, parameters=UVWABC)
Definition: __init__.py:266
def fixPXDYing(self, parameters=UVWABC)
Definition: __init__.py:262
def get_module(self, module=None)
Definition: __init__.py:191
def fixSVDPat(self, parameters=UVWABC)
Definition: __init__.py:270
def fixEKLMSector(self, endcap, layer, sector, parameters=[1, 2, 6])
Definition: __init__.py:280
parameters
Parameter config at algo level (fixing)
Definition: __init__.py:96
def fixBKLMModule(self, sector, layer, forward, parameters=UVWABC)
Definition: __init__.py:297
def fixGlobalParam(self, uid, element, param, value=0.)
Definition: __init__.py:219
def fixSVDMat(self, parameters=UVWABC)
Definition: __init__.py:274
def __init__(self, components=None, tracks=None, particles=None, vertices=None, primary_vertices=None, path=None)
Definition: __init__.py:58
collector
the collector module
Definition: __init__.py:94
def set_components(self, components)
Definition: __init__.py:117
def set_param(self, value, param, module=None)
Definition: __init__.py:185