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