Belle II Software prerelease-11-00-00a
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_reconstruction(path, pruneTracks=False, components=components)
72
73 collector = b2.register_module('MillepedeCollector',
74 tracks=['RecoTracks'],
75 primaryVertices=[],
76 calibrateVertex=True, components=db_components,
77 granularity='all', # time dependence needs granularity=all
78 useGblTree=False,
79 absFilePaths=True,
80 timedepConfig=timedep
81 )
82
83 algorithm = MillepedeAlgorithm()
84 algorithm.invertSign(True)
85 algorithm.ignoreUndeterminedParams(False)
86
87 std_components = ROOT.vector('string')()
88 for component in db_components:
89 std_components.push_back(component)
90 algorithm.setComponents(std_components)
91
92 algorithm.steering().command('method diagonalization 3 0.1')
93 algorithm.steering().command('skipemptycons')
94 algorithm.steering().command('threads 4 4')
95 algorithm.steering().command('matiter 1')
96 algorithm.steering().command('scaleerrors 1. 1.')
97 # algorithm.steering().command('entries 10')
98 algorithm.steering().command('printcounts 1')
99 # algorithm.steering().command('histprint')
100 algorithm.steering().command('monitorresiduals')
101 algorithm.steering().command('closeandreopen')
102 algorithm.steering().command('hugecut 50.')
103 # algorithm.steering().command('chiscut 30. 6.')
104 # algorithm.steering().command('outlierdownweighting 3')
105 # algorithm.steering().command('dwfractioncut 0.1')
106 # algorithm.steering().command('presigmas 1.')
107
108 algorithm.steering().command('Parameters')
109
110 def fix_vxd_id(vxd_id, params=None):
111 if params is None:
112 params = [1, 2, 3, 4, 5, 6] # rigid body parameters
113 for ipar in params:
114 label = Belle2.GlobalLabel()
115 label.construct(Belle2.VXDAlignment.getGlobalUniqueID(), vxd_id.getID(), ipar)
116 algorithm.steering().command('{} 0.0 -1. ! {}'.format(str(label.label()), str(vxd_id)))
117
118 # Halfshells
119 # fix_vxd_id(Belle2.VxdID(1, 0, 0, 1)) # ying
120 # fix_vxd_id(Belle2.VxdID(1, 0, 0, 2)) # yang
121 fix_vxd_id(Belle2.VxdID(3, 0, 0, 1)) # pat
122 fix_vxd_id(Belle2.VxdID(3, 0, 0, 2)) # mat
123
124 ladders = [8, 12, 7, 10, 12, 16]
125 sensors = [2, 2, 2, 3, 4, 5]
126
127 for layer in range(1, 7):
128 for ladder in range(1, ladders[layer - 1] + 1):
129 fix_vxd_id(Belle2.VxdID(layer, ladder, 0), params=[1, 2, 3, 4, 5, 6])
130 for sensor in range(1, sensors[layer - 1] + 1):
131 # sensors
132 fix_vxd_id(
134 layer,
135 ladder,
136 sensor),
137 params=[
138 1,
139 2,
140 3,
141 4,
142 5,
143 6,
144 31,
145 32,
146 33,
147 41,
148 42,
149 43,
150 44,
151 51,
152 52,
153 53,
154 54,
155 55,
156 20])
157
158 algorithm.setTimedepConfig(timedep)
159
160 calibration = Calibration('PXDHalfShellsAlignment',
161 collector=collector,
162 algorithms=algorithm,
163 input_files=files,
164 pre_collector_path=path,
165 database_chain=[CentralDatabase(tag) for tag in tags],
166 output_patterns=None,
167 max_files_per_collector_job=1,
168 backend_args=None
169 )
170
171 calibration.strategies = strategies.SingleIOV
172
173 return calibration
174
175
176def generate_test_data(filename):
177 main = b2.create_path()
178
179 main.add_module("EventInfoSetter", evtNumList=evtNumList, runList=runList, expList=expList)
180 # add_beamparameters(main, "Y4S")
181 main.add_module('Gearbox')
182 main.add_module('Geometry')
183
184 add_kkmc_generator(main, 'mu-mu+')
185
186 add_simulation(main)
187
188 main.add_module("RootOutput", outputFileName=filename)
189 main.add_module("Progress")
190
191 b2.process(main, calculateStatistics=True)
192 print(b2.statistics)
193 return os.path.abspath(filename)
194
195
196if __name__ == "__main__":
197 input_files = [os.path.abspath(file) for file in Belle2.Environment.Instance().getInputFilesOverride()]
198
199 if not len(input_files):
200 outfile = Belle2.Environment.Instance().getOutputFileOverride()
201 if not outfile:
202 outfile = 'TimeDepCalibration_TestData.root'
203 print("No input file provided. This will now generate test data into {} and run over it.".format(outfile))
204 print("It will take couple of minutes.")
205 print("You can set this file as input next time using -i option of basf2")
206
207 input_files = [generate_test_data(outfile)]
208
209 tags = b2.conditions.default_globaltags
210 mp2_beamspot = PXDHalfShellsAlignment(input_files, tags)
211 mp2_beamspot.max_iterations = 1
212
213 cal_fw = CAF()
214 cal_fw.add_calibration(mp2_beamspot)
215 cal_fw.backend = backends.LSF()
216
217 # Try to guess if we are at KEKCC and change the backend to Local if not
218 if multiprocessing.cpu_count() < 10:
219 cal_fw.backend = backends.Local(1)
220
221 cal_fw.run()
static Environment & Instance()
Static method to get a reference to the Environment instance.
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.
Class to uniquely identify a any structure of the PXD and SVD.
Definition VxdID.h:32