Belle II Software  release-06-01-15
Alignment.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 # KLM alignment: alignment using the collected data.
13 
14 import sys
15 import basf2
16 import ROOT
17 from ROOT import Belle2
18 from ROOT.Belle2 import KLMChannelIndex, KLMElementNumbers
19 import numpy as np
20 from alignment import MillepedeCalibration
21 
22 basf2.set_log_level(basf2.LogLevel.INFO)
23 
24 basf2.conditions.append_testing_payloads('localdb/database.txt')
25 
26 # Create the algorithm.
27 millepede = MillepedeCalibration(['BKLMAlignment', 'EKLMAlignment', 'EKLMSegmentAlignment'])
28 
29 # Fix module parameters.
30 index = KLMChannelIndex(KLMChannelIndex.c_IndexLevelLayer)
31 index2 = KLMChannelIndex(KLMChannelIndex.c_IndexLevelLayer)
32 while (index != index2.end()):
33  module = index.getKLMModuleNumber()
34  if (index.getSubdetector() == KLMElementNumbers.c_BKLM):
35  for ipar in [1, 2, 3, 4, 5, 6]:
36  # Free parameters are idU, dV, dGamma.
37  if ipar in [1, 2, 6]:
38  continue
39  millepede.fixGlobalParam(Belle2.BKLMAlignment.getGlobalUniqueID(),
40  module, ipar)
41  else:
42  for ipar in [1, 2, 6]:
43  # No parameters are fixed; if necessary, uncomment the following:
44  # millepede.fixGlobalParam(Belle2.EKLMAlignment.getGlobalUniqueID(),
45  # module, ipar)
46  continue
47  index.increment()
48 
49 # Fix EKLM segment parameters.
50 index.setIndexLevel(KLMChannelIndex.c_IndexLevelStrip)
51 index2.setIndexLevel(KLMChannelIndex.c_IndexLevelStrip)
52 index = index2.beginEKLM()
53 index.useEKLMSegments()
54 while (index != index2.endEKLM()):
55  segment = index.getEKLMSegmentNumber()
56  for ipar in [2, 6]:
57  millepede.fixGlobalParam(
59  segment, ipar)
60  index.increment()
61 
62 # Execute the algorithm over all collected data (auto-merged).
63 input_files = sys.argv[1:]
64 calibration = millepede.create('klm_alignment', input_files)
65 millepede.algo.setInputFileNames(input_files)
66 millepede.algo.ignoreUndeterminedParams(True)
67 millepede.algo.invertSign()
68 millepede.algo.execute()
69 millepede.algo.commit()
70 
71 # Get payloads.
72 payloads = list(millepede.algo.getPayloads())
73 bklm_alignment = None
74 eklm_alignment = None
75 eklm_segment_alignment = None
76 for payload in payloads:
77  if payload.name == 'BKLMAlignment':
78  bklm_alignment = payload.object.IsA().DynamicCast(Belle2.BKLMAlignment().IsA(), payload.object, False)
79  elif payload.name == 'EKLMAlignment':
80  eklm_alignment = payload.object.IsA().DynamicCast(Belle2.EKLMAlignment().IsA(), payload.object, False)
81  elif payload.name == 'EKLMAlignment':
82  eklm_segment_alignment = payload.object.IsA().DynamicCast(Belle2.EKLMSegmentAlignment().IsA(), payload.object, False)
83 
84 # Profile plot for all determined parameters
85 profile = ROOT.TH1F(
86  "profile",
87  "correction & errors",
88  millepede.algo.result().getNoDeterminedParameters(),
89  1,
90  millepede.algo.result().getNoDeterminedParameters())
91 
92 # Define some branch variables
93 param = np.zeros(1, dtype=int)
94 value = np.zeros(1, dtype=np.float32)
95 correction = np.zeros(1, dtype=np.float32)
96 error = np.zeros(1, dtype=np.float32)
97 section = np.zeros(1, dtype=int)
98 sector = np.zeros(1, dtype=int)
99 layer = np.zeros(1, dtype=int)
100 sensor = np.zeros(1, dtype=int)
101 plane = np.zeros(1, dtype=int)
102 segment = np.zeros(1, dtype=int)
103 
104 alignment_file = ROOT.TFile('alignment.root', 'recreate')
105 # Tree with BKLM module data.
106 bklm_module_tree = ROOT.TTree('bklm_module', 'BKLM module alignment data')
107 bklm_module_tree.Branch('section', section, 'section/I')
108 bklm_module_tree.Branch('sector', sector, 'sector/I')
109 bklm_module_tree.Branch('layer', layer, 'layer/I')
110 bklm_module_tree.Branch('param', param, 'param/I')
111 bklm_module_tree.Branch('value', value, 'value/F')
112 bklm_module_tree.Branch('correction', correction, 'correction/F')
113 bklm_module_tree.Branch('error', error, 'error/F')
114 # Tree with EKLM module data.
115 eklm_module_tree = ROOT.TTree('eklm_module', 'EKLM module alignment data')
116 eklm_module_tree.Branch('section', section, 'section/I')
117 eklm_module_tree.Branch('sector', sector, 'sector/I')
118 eklm_module_tree.Branch('layer', layer, 'layer/I')
119 eklm_module_tree.Branch('param', param, 'param/I')
120 eklm_module_tree.Branch('value', value, 'value/F')
121 eklm_module_tree.Branch('correction', correction, 'correction/F')
122 eklm_module_tree.Branch('error', error, 'error/F')
123 # Tree with EKLM segment data.
124 eklm_segment_tree = ROOT.TTree('eklm_segment', 'EKLM segment alignment data')
125 eklm_segment_tree.Branch('section', section, 'section/I')
126 eklm_segment_tree.Branch('sector', sector, 'sector/I')
127 eklm_segment_tree.Branch('layer', layer, 'layer/I')
128 eklm_segment_tree.Branch('plane', plane, 'plane/I')
129 eklm_segment_tree.Branch('segment', segment, 'segment/I')
130 eklm_segment_tree.Branch('param', param, 'param/I')
131 eklm_segment_tree.Branch('value', value, 'value/F')
132 eklm_segment_tree.Branch('correction', correction, 'correction/F')
133 eklm_segment_tree.Branch('error', error, 'error/F')
134 
135 # Index of determined parameter.
136 ibin = 0
137 
138 for ipar in range(0, millepede.algo.result().getNoParameters()):
139  label = Belle2.GlobalLabel(millepede.algo.result().getParameterLabel(ipar))
140  param[0] = label.getParameterId()
141  value[0] = 0.0
142 
143  if millepede.algo.result().isParameterDetermined(ipar):
144  correction[0] = millepede.algo.result().getParameterCorrection(ipar)
145  error[0] = millepede.algo.result().getParameterError(ipar)
146  else:
147  correction[0] = 0.0
148  error[0] = -1.0
149 
150  # Module alignment.
151  if (label.getUniqueId() == Belle2.BKLMAlignment.getGlobalUniqueID() or
152  label.getUniqueId() == Belle2.EKLMAlignment.getGlobalUniqueID()):
153  module = label.getElementId()
154  index.setKLMModule(module)
155  section[0] = index.getSection()
156  sector[0] = index.getSector()
157  layer[0] = index.getLayer()
158  if (index.getSubdetector() == KLMElementNumbers.c_BKLM):
159  if (bklm_alignment is not None):
160  value[0] = bklm_alignment.getGlobalParam(module, int(param[0]))
161  bklm_module_tree.Fill()
162  else:
163  if (eklm_alignment is not None):
164  value[0] = eklm_alignment.getGlobalParam(module, int(param[0]))
165  eklm_module_tree.Fill()
166 
167  # EKLM segments alignment.
168  elif (label.getUniqueId() ==
170  segment_global = label.getElementId()
171  index.setEKLMSegment(segment_global)
172  section[0] = index.getSection()
173  sector[0] = index.getSector()
174  layer[0] = index.getLayer()
175  plane[0] = index.getPlane()
176  segment[0] = index.getStrip()
177  if (eklm_segment_alignment is not None):
178  value[0] = eklm_segment_alignment.getGlobalParam(segment, int(param[0]))
179  eklm_segment_tree.Fill()
180 
181  if not millepede.algo.result().isParameterDetermined(ipar):
182  continue
183 
184  ibin = ibin + 1
185  profile.SetBinContent(ibin, value[0])
186  profile.SetBinError(ibin, error[0])
187 
188 alignment_file.cd()
189 profile.Write()
190 bklm_module_tree.Write()
191 eklm_module_tree.Write()
192 eklm_segment_tree.Write()
193 alignment_file.Close()
Class to store BKLM alignment data in the database.
Definition: BKLMAlignment.h:30
static unsigned short getGlobalUniqueID()
Get global unique identifier.
Definition: BKLMAlignment.h:63
Class to store EKLM alignment data in the database.
Definition: EKLMAlignment.h:30
static unsigned short getGlobalUniqueID()
Get global unique identifier.
Definition: EKLMAlignment.h:63
Class to store EKLM alignment data in the database.
static unsigned short getGlobalUniqueID()
Get global unique identifier.
Class to convert to/from global labels for Millepede II to/from detector & parameter identificators.
Definition: GlobalLabel.h:41