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