10from basf2.pickle_path import deserialize_path, deserialize_module, serialize_path, serialize_module
13from ROOT
import Belle2
14from caf.framework
import Calibration
16from collections
import OrderedDict
18b2.set_log_level(b2.LogLevel.INFO)
19UVWABC = [1, 2, 3, 4, 5, 6]
22def collect(calibration, basf2_args=None):
26 if basf2_args
is None:
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)
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)
39 subprocess.call([
"basf2",
"--execute-path", path_file_name] + basf2_args)
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()
50 """ The generic Millepede calibration collector+algorithm """
57 primary_vertices=None,
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
67 if components
is None:
71 if (tracks
is None)
and (particles
is None)
and (vertices
is None)
and (primary_vertices
is None):
72 tracks = [
'RecoTracks']
80 if primary_vertices
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)
94 self.
collector = b2.register_module(
'MillepedeCollector')
103 self.
collector.param(
'particles', particles)
104 self.
collector.param(
'vertices', vertices)
105 self.
collector.param(
'primaryVertices', primary_vertices)
118 """ Select db objects for calibration from list of their names """
121 calibrate_vertex = ((components == [])
or (
'BeamSpot' in components))
123 std_components = ROOT.vector(
'string')()
124 for component
in components:
125 std_components.push_back(component)
127 self.
algo.setComponents(std_components)
128 self.
collector.param(
'components', components)
129 self.
collector.param(
'calibrateVertex', calibrate_vertex)
132 """ Set command for Millepede steering """
140 command = cmd_name +
' ' + str(command)
142 cmd_words = str(command).split(
' ')
143 self.
commands[cmd_words[0]] = command
146 """ Fcn to execute before algorithm... """
150 """ Create the CAF Calibration object """
154 for command
in self.commands.values():
155 self.algo.
steering().command(command)
157 self.algo.
steering().command(
'Parameters')
158 for command
in self.parameters.values():
159 self.algo.
steering().command(command)
164 if not self.get_param(
'useGblTree'):
165 self.set_param(
True,
'absFilePaths')
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')
172 cal.pre_collector_path = self.path
173 cal.pre_algorithms = self.pre_algo
178 """ Get parameter of the collector module or any module in path (given its name) """
180 for mpi
in module.available_params():
181 if mpi.name == param:
186 """ Set parameter of the collector module or any module in path (given its name) """
188 if module
is not None:
189 module.param(param, value)
192 """ Get collector module or any other module from path (given its name) """
197 if mod.name() == module:
200 if isinstance(module, b2.Module):
205 """ serialization """
207 state = self.__dict__.copy()
208 state[
'path'] = serialize_path(self.
path)
209 state[
'collector'] = serialize_module(self.
collector)
213 """ de-serialization """
215 self.__dict__.update(state)
216 self.
collector = deserialize_module(state[
'collector'])
217 self.
path = deserialize_path(state[
'path'])
220 """ Generic fcn to manipulate (fix to given value) parameters identified by db object id (uid),
221 element and parameter number
"""
223 label.construct(uid, element, param)
224 self.parameters[str(label.label())] = (str(label.label()) + ' ' + str(value) +
' -1.')
229 """ fix CDC layer X-shift """
233 """ fix CDC layer Y-shift """
237 """ fix CDC layer Phi rotation """
241 """ fix CDC time walk parameter for given board """
245 """ fix CDC tie zero parameter for given wire """
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)
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)}}
263 """ fix PXD Ying half-shell """
264 self.
fixVXDid(1, 0, 0, 1, parameters)
267 """ fix PXD Yang half-shell """
268 self.
fixVXDid(1, 0, 0, 2, parameters)
271 """ fix SVD Pat half-shell """
272 self.
fixVXDid(3, 0, 0, 1, parameters)
275 """ fix SVD Mat half-shell """
276 self.
fixVXDid(3, 0, 0, 2, parameters)
281 """ Fix EKLM sector parameters """
283 for ipar
in parameters:
285 Belle2.EKLMElementID(endcap, layer, sector).getGlobalNumber(),
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(),
298 """ Fix a BKLM module """
299 for ipar
in parameters:
302 Belle2.BKLMElementID(
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:
312 for forward
in forbackwards:
static unsigned short getGlobalUniqueID()
Get global unique identifier.
static unsigned short getGlobalUniqueID()
Get global unique id.
static unsigned short getGlobalUniqueID()
Get global unique id.
static unsigned short getGlobalUniqueID()
Get global unique id.
static unsigned short getGlobalUniqueID()
Get global unique identifier.
Class to convert to/from global labels for Millepede II to/from detector & parameter identificators.
Class implementing Millepede calibration algorithm.
static unsigned short getGlobalUniqueID()
Get global unique id.
Class to uniquely identify a any structure of the PXD and SVD.
commands
Commands for pede steering.
def fixEKLMSegment(self, endcap, layer, sector, plane, segment, parameters=[1, 2, 6])
def get_param(self, param, module=None)
def set_command(self, cmd_name, command='')
components
db objects for calibration
def __setstate__(self, state)
def fixCDCLayerX(self, layer)
def fixCDCLayerY(self, layer)
def fixVXDid(self, layer, ladder, sensor, segment=0, parameters=[1, 2, 3, 4, 5, 6])
def fixCDCLayerRot(self, layer)
def fixCDCTimeWalk(self, board_id)
def create(self, name, input_files)
def pre_algo(self, algorithm, iteration)
def fixBKLM(self, sectors=range(1, 9), layers=range(1, 16), forbackwards=[0, 1])
def fixCDCTimeZero(self, wire_id)
def fixPXDYang(self, parameters=UVWABC)
def fixPXDYing(self, parameters=UVWABC)
def get_module(self, module=None)
def fixSVDPat(self, parameters=UVWABC)
def fixEKLMSector(self, endcap, layer, sector, parameters=[1, 2, 6])
parameters
Parameter config at algo level (fixing)
def fixBKLMModule(self, sector, layer, forward, parameters=UVWABC)
def fixGlobalParam(self, uid, element, param, value=0.)
def fixSVDMat(self, parameters=UVWABC)
def __init__(self, components=None, tracks=None, particles=None, vertices=None, primary_vertices=None, path=None)
collector
the collector module
def set_components(self, components)
def set_param(self, value, param, module=None)