Belle II Software development
runTimeDependentVXD_CAF.py
1
8
9import basf2 as b2
10
11import os
12import multiprocessing
13
14import ROOT
15from ROOT import Belle2
16from ROOT.Belle2 import MillepedeAlgorithm
17
18from caf.framework import CAF, Calibration, CentralDatabase
19from caf import backends
20from caf import strategies
21
22import reconstruction as reco
23
24from generators import add_kkmc_generator
25from simulation import add_simulation
26
27b2.set_log_level(b2.LogLevel.INFO)
28# Generate 4 runs, with run number 4, 5, 6, 7
29runList = [4, 5, 6, 7]
30# all runs in experiment 0 (default final Phase3 config)
31expList = [0, 0, 0, 0]
32# let's have runs 5 and 7 longer and let the PXD move each 20 events below
33evtNumList = [20, 60, 20, 40]
34
35# VXD iDs of the 2 PXD half shells
36ying = Belle2.VxdID(1, 0, 0, 1)
37yang = Belle2.VxdID(1, 0, 0, 2)
38
39# Their rigid body parameter IDs
40# x,y, z shift and alpha, beta, gamma rotation
41params = [1, 2, 3, 4, 5, 6]
42
43label = Belle2.GlobalLabel()
44
45pxd_labels = []
46for shell in [ying, yang]:
47 for par in params:
48 label.construct(Belle2.VXDAlignment.getGlobalUniqueID(), shell.getID(), par)
49 pxd_labels.append(label.label())
50
51# Configure the the time dependence for calibration
52timedep = [(pxd_labels, [(0, run, 0) for run in runList] + [(20, 5, 0), (40, 5, 0)] + [(20, 7, 0)])]
53
54db_components = ['VXDAlignment']
55components = ['PXD', 'SVD', 'CDC']
56
57
58def PXDHalfShellsAlignment(files, tags):
59
60 # Set-up re-processing path
61 path = b2.create_path()
62
63 path.add_module('Progress')
64 # Remove all non-raw data to run the full reco again
65 # path.add_module('RootInput', branchNames=input_branches)#, entrySequences=['0:1000'])
66 path.add_module('Gearbox')
67 path.add_module('Geometry')
68 # Not needed for di-muon skim cdst or mdst, but needed to re-run reconstruction
69 # with possibly changed global tags
70 # raw.add_unpackers(path)
71 reco.add_mc_reconstruction(path, pruneTracks=False, components=components)
72 path.add_module('DAFRecoFitter')
73
74 collector = b2.register_module('MillepedeCollector',
75 tracks=['RecoTracks'],
76 primaryVertices=[],
77 calibrateVertex=True, components=db_components,
78 granularity='all', # time dependence needs granularity=all
79 useGblTree=False,
80 absFilePaths=True,
81 timedepConfig=timedep
82 )
83
84 algorithm = MillepedeAlgorithm()
85 algorithm.invertSign(True)
86 algorithm.ignoreUndeterminedParams(False)
87
88 std_components = ROOT.vector('string')()
89 for component in db_components:
90 std_components.push_back(component)
91 algorithm.setComponents(std_components)
92
93 algorithm.steering().command('method diagonalization 3 0.1')
94 algorithm.steering().command('skipemptycons')
95 algorithm.steering().command('threads 4 4')
96 algorithm.steering().command('matiter 1')
97 algorithm.steering().command('scaleerrors 1. 1.')
98 # algorithm.steering().command('entries 10')
99 algorithm.steering().command('printcounts 1')
100 # algorithm.steering().command('histprint')
101 algorithm.steering().command('monitorresiduals')
102 algorithm.steering().command('closeandreopen')
103 algorithm.steering().command('hugecut 50.')
104 # algorithm.steering().command('chiscut 30. 6.')
105 # algorithm.steering().command('outlierdownweighting 3')
106 # algorithm.steering().command('dwfractioncut 0.1')
107 # algorithm.steering().command('presigmas 1.')
108
109 algorithm.steering().command('Parameters')
110
111 def fix_vxd_id(vxd_id, params=None):
112 if params is None:
113 params = [1, 2, 3, 4, 5, 6] # rigid body parameters
114 for ipar in params:
115 label = Belle2.GlobalLabel()
116 label.construct(Belle2.VXDAlignment.getGlobalUniqueID(), vxd_id.getID(), ipar)
117 algorithm.steering().command('{} 0.0 -1. ! {}'.format(str(label.label()), str(vxd_id)))
118
119 # Halfshells
120 # fix_vxd_id(Belle2.VxdID(1, 0, 0, 1)) # ying
121 # fix_vxd_id(Belle2.VxdID(1, 0, 0, 2)) # yang
122 fix_vxd_id(Belle2.VxdID(3, 0, 0, 1)) # pat
123 fix_vxd_id(Belle2.VxdID(3, 0, 0, 2)) # mat
124
125 ladders = [8, 12, 7, 10, 12, 16]
126 sensors = [2, 2, 2, 3, 4, 5]
127
128 for layer in range(1, 7):
129 for ladder in range(1, ladders[layer - 1] + 1):
130 fix_vxd_id(Belle2.VxdID(layer, ladder, 0), params=[1, 2, 3, 4, 5, 6])
131 for sensor in range(1, sensors[layer - 1] + 1):
132 # sensors
133 fix_vxd_id(
135 layer,
136 ladder,
137 sensor),
138 params=[
139 1,
140 2,
141 3,
142 4,
143 5,
144 6,
145 31,
146 32,
147 33,
148 41,
149 42,
150 43,
151 44,
152 51,
153 52,
154 53,
155 54,
156 55,
157 20])
158
159 algorithm.setTimedepConfig(timedep)
160
161 calibration = Calibration('PXDHalfShellsAlignment',
162 collector=collector,
163 algorithms=algorithm,
164 input_files=files,
165 pre_collector_path=path,
166 database_chain=[CentralDatabase(tag) for tag in tags],
167 output_patterns=None,
168 max_files_per_collector_job=1,
169 backend_args=None
170 )
171
172 calibration.strategies = strategies.SingleIOV
173
174 return calibration
175
176
177def generate_test_data(filename):
178 main = b2.create_path()
179
180 main.add_module("EventInfoSetter", evtNumList=evtNumList, runList=runList, expList=expList)
181 # add_beamparameters(main, "Y4S")
182 main.add_module('Gearbox')
183 main.add_module('Geometry')
184
185 add_kkmc_generator(main, 'mu-mu+')
186
187 add_simulation(main)
188
189 main.add_module("RootOutput", outputFileName=filename)
190 main.add_module("Progress")
191
192 b2.process(main)
193 print(b2.statistics)
194 return os.path.abspath(filename)
195
196
197if __name__ == "__main__":
198 input_files = [os.path.abspath(file) for file in Belle2.Environment.Instance().getInputFilesOverride()]
199
200 if not len(input_files):
201 outfile = Belle2.Environment.Instance().getOutputFileOverride()
202 if not outfile:
203 outfile = 'TimeDepCalibration_TestData.root'
204 print("No input file provided. This will now generate test data into {} and run over it.".format(outfile))
205 print("It will take couple of minutes.")
206 print("You can set this file as input next time using -i option of basf2")
207
208 input_files = [generate_test_data(outfile)]
209
210 tags = b2.conditions.default_globaltags
211 mp2_beamspot = PXDHalfShellsAlignment(input_files, tags)
212 mp2_beamspot.max_iterations = 1
213
214 cal_fw = CAF()
215 cal_fw.add_calibration(mp2_beamspot)
216 cal_fw.backend = backends.LSF()
217
218 # Try to guess if we are at KEKCC and change the backend to Local if not
219 if multiprocessing.cpu_count() < 10:
220 cal_fw.backend = backends.Local(1)
221
222 cal_fw.run()
static Environment & Instance()
Static method to get a reference to the Environment instance.
Definition: Environment.cc:28
Class to convert to/from global labels for Millepede II to/from detector & parameter identificators.
Definition: GlobalLabel.h:41
static unsigned short getGlobalUniqueID()
Get global unique id.
Definition: VXDAlignment.h:47
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33