Belle II Software development
MillepedeCalibration Class Reference

Public Member Functions

 __init__ (self, components=None, tracks=None, particles=None, vertices=None, primary_vertices=None, path=None)
 
 set_components (self, components)
 
 set_command (self, cmd_name, command='')
 
 pre_algo (self, algorithm, iteration)
 
 create (self, name, input_files)
 
 get_param (self, param, module=None)
 
 set_param (self, value, param, module=None)
 
 get_module (self, module=None)
 
 __getstate__ (self)
 
 __setstate__ (self, state)
 
 fixGlobalParam (self, uid, element, param, value=0.)
 
 fixCDCLayerX (self, layer)
 
 fixCDCLayerY (self, layer)
 
 fixCDCLayerRot (self, layer)
 
 fixCDCTimeWalk (self, board_id)
 
 fixCDCTimeZero (self, wire_id)
 
 fixVXDid (self, layer, ladder, sensor, segment=0, parameters=[1, 2, 3, 4, 5, 6])
 
 fixPXDYing (self, parameters=UVWABC)
 
 fixPXDYang (self, parameters=UVWABC)
 
 fixSVDPat (self, parameters=UVWABC)
 
 fixSVDMat (self, parameters=UVWABC)
 
 fixEKLMSector (self, endcap, layer, sector, parameters=[1, 2, 6])
 
 fixEKLMSegment (self, endcap, layer, sector, plane, segment, parameters=[1, 2, 6])
 
 fixBKLMModule (self, sector, layer, forward, parameters=UVWABC)
 
 fixBKLM (self, sectors=range(1, 9), layers=range(1, 16), forbackwards=[0, 1])
 

Public Attributes

 algo = Belle2.MillepedeAlgorithm()
 the algorithm
 
 path = path
 pre-collector path
 
 collector = b2.register_module('MillepedeCollector')
 the collector module
 
 parameters = OrderedDict()
 Parameter config at algo level (fixing)
 
 commands = OrderedDict()
 Commands for pede steering.
 
 components = components
 db objects for calibration
 

Detailed Description

 The generic Millepede calibration collector+algorithm 

Definition at line 49 of file __init__.py.

Constructor & Destructor Documentation

◆ __init__()

__init__ ( self,
components = None,
tracks = None,
particles = None,
vertices = None,
primary_vertices = None,
path = None )
components are the names of DB objects to calibrate (BeamSpot etc.)
tracks are collections of RecoTracks of fitted charged tracks (usually cosmic rays, no associated particle)
particles are names of ParticleLists with single charged particles
vertices are names ParticleLists with at least two body decays fitted with vertex constraint
primary_vertices are names of ParticleLists with at least two body decays fitted with beam+vertex constraint

Definition at line 52 of file __init__.py.

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
90 self.algo = Belle2.MillepedeAlgorithm()
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
Class implementing Millepede calibration algorithm.

Member Function Documentation

◆ __getstate__()

__getstate__ ( self)
 serialization 

Definition at line 204 of file __init__.py.

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

◆ __setstate__()

__setstate__ ( self,
state )
 de-serialization 

Definition at line 212 of file __init__.py.

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

◆ create()

create ( self,
name,
input_files )
 Create the CAF Calibration object 

Definition at line 149 of file __init__.py.

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 remember 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

◆ fixBKLM()

fixBKLM ( self,
sectors = range(1, 9),
layers = range(1, 16),
forbackwards = [0, 1] )
 Fix (all by default) BKLM modules 

Definition at line 308 of file __init__.py.

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)

◆ fixBKLMModule()

fixBKLMModule ( self,
sector,
layer,
forward,
parameters = UVWABC )
 Fix a BKLM module 

Definition at line 297 of file __init__.py.

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
static unsigned short getGlobalUniqueID()
Get global unique identifier.

◆ fixCDCLayerRot()

fixCDCLayerRot ( self,
layer )
 fix CDC layer Phi rotation 

Definition at line 236 of file __init__.py.

236 def fixCDCLayerRot(self, layer):
237 """ fix CDC layer Phi rotation """
238 self.fixGlobalParam(Belle2.CDCLayerAlignment.getGlobalUniqueID(), layer, Belle2.CDCLayerAlignment.layerPhi)
239
static unsigned short getGlobalUniqueID()
Get global unique id.

◆ fixCDCLayerX()

fixCDCLayerX ( self,
layer )
 fix CDC layer X-shift 

Definition at line 228 of file __init__.py.

228 def fixCDCLayerX(self, layer):
229 """ fix CDC layer X-shift """
230 self.fixGlobalParam(Belle2.CDCLayerAlignment.getGlobalUniqueID(), layer, Belle2.CDCLayerAlignment.layerX)
231

◆ fixCDCLayerY()

fixCDCLayerY ( self,
layer )
 fix CDC layer Y-shift 

Definition at line 232 of file __init__.py.

232 def fixCDCLayerY(self, layer):
233 """ fix CDC layer Y-shift """
234 self.fixGlobalParam(Belle2.CDCLayerAlignment.getGlobalUniqueID(), layer, Belle2.CDCLayerAlignment.layerY)
235

◆ fixCDCTimeWalk()

fixCDCTimeWalk ( self,
board_id )
 fix CDC time walk parameter for given board 

Definition at line 240 of file __init__.py.

240 def fixCDCTimeWalk(self, board_id):
241 """ fix CDC time walk parameter for given board """
242 self.fixGlobalParam(Belle2.CDCTimeWalks.getGlobalUniqueID(), board_id, 0)
243
static unsigned short getGlobalUniqueID()
Get global unique id.

◆ fixCDCTimeZero()

fixCDCTimeZero ( self,
wire_id )
 fix CDC tie zero parameter for given wire 

Definition at line 244 of file __init__.py.

244 def fixCDCTimeZero(self, wire_id):
245 """ fix CDC tie zero parameter for given wire """
246 self.fixGlobalParam(Belle2.CDCTimeZeros.getGlobalUniqueID(), wire_id, 0)
247
static unsigned short getGlobalUniqueID()
Get global unique id.

◆ fixEKLMSector()

fixEKLMSector ( self,
endcap,
layer,
sector,
parameters = [1, 2, 6] )
 Fix EKLM sector parameters 

Definition at line 280 of file __init__.py.

280 def fixEKLMSector(self, endcap, layer, sector, parameters=[1, 2, 6]):
281 """ Fix EKLM sector parameters """
282
283 for ipar in parameters:
284 self.fixGlobalParam(Belle2.EKLMAlignment.getGlobalUniqueID(),
285 Belle2.EKLMElementID(endcap, layer, sector).getGlobalNumber(),
286 ipar)
287
static unsigned short getGlobalUniqueID()
Get global unique identifier.

◆ fixEKLMSegment()

fixEKLMSegment ( self,
endcap,
layer,
sector,
plane,
segment,
parameters = [1, 2, 6] )
 Fix EKLM segment parameters 

Definition at line 288 of file __init__.py.

288 def fixEKLMSegment(self, endcap, layer, sector, plane, segment, parameters=[1, 2, 6]):
289 """ Fix EKLM segment parameters """
290 for ipar in parameters:
291 self.fixGlobalParam(Belle2.EKLMAlignment.getGlobalUniqueID(),
292 Belle2.EKLMElementID(endcap, layer, sector, plane, segment).getGlobalNumber(),
293 ipar)
294

◆ fixGlobalParam()

fixGlobalParam ( self,
uid,
element,
param,
value = 0. )
 Generic fcn to manipulate (fix to given value) parameters identified by db object id (uid),
    element and parameter number 

Definition at line 219 of file __init__.py.

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
Class to convert to/from global labels for Millepede II to/from detector & parameter identificators.
Definition GlobalLabel.h:41

◆ fixPXDYang()

fixPXDYang ( self,
parameters = UVWABC )
 fix PXD Yang half-shell 

Definition at line 266 of file __init__.py.

266 def fixPXDYang(self, parameters=UVWABC):
267 """ fix PXD Yang half-shell """
268 self.fixVXDid(1, 0, 0, 2, parameters)
269

◆ fixPXDYing()

fixPXDYing ( self,
parameters = UVWABC )
 fix PXD Ying half-shell 

Definition at line 262 of file __init__.py.

262 def fixPXDYing(self, parameters=UVWABC):
263 """ fix PXD Ying half-shell """
264 self.fixVXDid(1, 0, 0, 1, parameters)
265

◆ fixSVDMat()

fixSVDMat ( self,
parameters = UVWABC )
 fix SVD Mat half-shell 

Definition at line 274 of file __init__.py.

274 def fixSVDMat(self, parameters=UVWABC):
275 """ fix SVD Mat half-shell """
276 self.fixVXDid(3, 0, 0, 2, parameters)
277

◆ fixSVDPat()

fixSVDPat ( self,
parameters = UVWABC )
 fix SVD Pat half-shell 

Definition at line 270 of file __init__.py.

270 def fixSVDPat(self, parameters=UVWABC):
271 """ fix SVD Pat half-shell """
272 self.fixVXDid(3, 0, 0, 1, parameters)
273

◆ fixVXDid()

fixVXDid ( self,
layer,
ladder,
sensor,
segment = 0,
parameters = [1, 2, 3, 4, 5, 6] )
 Fix VXD element parameters by layer,ladder,sensor and segment number 

Definition at line 250 of file __init__.py.

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
static unsigned short getGlobalUniqueID()
Get global unique id.
Class to uniquely identify a any structure of the PXD and SVD.
Definition VxdID.h:33

◆ get_module()

get_module ( self,
module = None )
 Get collector module or any other module from path (given its name) 

Definition at line 191 of file __init__.py.

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

◆ get_param()

get_param ( self,
param,
module = None )
 Get parameter of the collector module or any module in path (given its name) 

Definition at line 177 of file __init__.py.

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

◆ pre_algo()

pre_algo ( self,
algorithm,
iteration )
 Fcn to execute before algorithm... 

Definition at line 145 of file __init__.py.

145 def pre_algo(self, algorithm, iteration):
146 """ Fcn to execute before algorithm... """
147 pass
148

◆ set_command()

set_command ( self,
cmd_name,
command = '' )
 Set command for Millepede steering 

Definition at line 131 of file __init__.py.

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

◆ set_components()

set_components ( self,
components )
 Select db objects for calibration from list of their names 

Definition at line 117 of file __init__.py.

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

◆ set_param()

set_param ( self,
value,
param,
module = None )
 Set parameter of the collector module or any module in path (given its name) 

Definition at line 185 of file __init__.py.

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

Member Data Documentation

◆ algo

the algorithm

Definition at line 90 of file __init__.py.

◆ collector

collector = b2.register_module('MillepedeCollector')

the collector module

Definition at line 94 of file __init__.py.

◆ commands

commands = OrderedDict()

Commands for pede steering.

Definition at line 98 of file __init__.py.

◆ components

components = components

db objects for calibration

Definition at line 120 of file __init__.py.

◆ parameters

parameters = OrderedDict()

Parameter config at algo level (fixing)

Definition at line 96 of file __init__.py.

◆ path

path = path

pre-collector path

Definition at line 92 of file __init__.py.


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