Belle II Software development
constraints.py
1
8
9import basf2
10from ROOT import Belle2
11
12import os
13import math
14
15
16class Constraint():
17 """
18 Class representing a linear constraint for global parameters
19 """
20
21 def __init__(self, comment="", value=0.):
22 """
23 Initialize constraint
24
25 comment : str
26 TODO: Not yet used
27 value : float
28 The constant term of the constraint: sum(c_i * par_i) = value
29 """
30
31 self.value = value
32
33 self.comment = comment
34
35 self.data = []
36
37 def add(self, label, coeff):
38 """
39 Add coefficient for a global parameter
40
41 Parameters
42 ----------
43 label : int
44 global parameter id (GlobalLabel::label())
45 coeff : float
46 coefficient of the parameter in constraint
47 """
48 entry = (label, coeff)
49 self.data.append(entry)
50
51 def get_checksum(self):
52 """Get a checksum to distinguish constraints quickly (computed only from labels)
53 """
54 labels_only = [label for (label, coeff) in self.data]
55 checksum = hash(str(labels_only))
56 return checksum
57
58
60 """
61 Base class representing a "generator" for file with set of constraints
62
63 Can be used directly if you already have a file with constraints.
64 """
65
66 def __init__(self, filename):
67 """
68 filename : str
69 Desired filename for the constraint file
70 """
71
72 self.filename = filename
73
74 def generate(self):
75 """
76 Should be overridden by the deriving classes and actually
77 fill the dictionary with the constraints (this runs within a basf2 module
78 event function - so you have geometry available)
79 """
80 consts = []
81 return consts
82
83 def configure_collector(self, collector):
84 """
85 Can be overridden be child classes to pass additional configuration to the
86 MillepedeCollector (activated by the use of the constraints)
87 """
88 pass
89
90
91def generate_constraints(constraint_sets, timedep, global_tags, init_event):
92 """
93 Run the constraints generation and return file names
94
95 Parameters
96 ----------
97 constraint_sets : list (alignment.Constraints)
98 List of sets of constraints
99 timedep_config : list(tuple(list(int), list(tuple(int, int, int))))
100 Time-dependence configuration.
101 Each list item is 2-tuple with list of parameter numbers (use alignment.parameters to get them) and
102 the (event, run, exp) numbers at which values of these parameters can change.
103 global_tags : list (str)
104 List of global tag names and/or (absolute) file paths to local databases
105 init_event : tuple( int, int, int)
106 Event (event, run, exp) at which to initialize time-INdependent constraints
107 """
108 files = []
109 for filename in [consts.filename for consts in constraint_sets]:
110 files.append(os.path.abspath(filename))
111
112 from alignment.constraints_generator import save_config
113 ccfn = save_config(constraint_sets, timedep, global_tags, init_event)
114 os.system('basf2 {} {}'.format(Belle2.FileSystem.findFile('alignment/scripts/alignment/constraints_generator.py'), ccfn))
115
116 return files
117
118# ----------------------------- Sub-detector specific constraint classes -----------------------------------------------
119
120
122 """
123 Constraints for VXD hierarchy
124
125 They come in 3 types:
126
127 0 : no hierarchy, just sensors aligned in their local system
128 1 : sensors -> half-shells (no ladders)
129 2 : sensors -> ladders ->half-shells (aka full hierarchy)
130
131 and additionally can use them for PXD, SVD or both (default) (setting pxd/svd to False
132 will disable the hierarchy for that sub-detector)
133
134 NOTE: the filename cannot be currently changed, it is hard-coded in the collector module
135 as these constraints are actually generated by C++ code and not by python modules as others.
136 """
137
138 def __init__(self, type=2, pxd=True, svd=True):
139 """ Initialize
140 type : int
141 Type of constraints (0 - no, 1 - half-shells, 2 -full)
142 pxd : bool
143 Use constraints for PXD?
144 svd : bool
145 Use constraints for SVD?
146 """
147
148 # TODO: cannot currently change the filename in collector, so fixed here
149 super().__init__("constraints.txt")
150
151 self.type = type
152
153 self.pxd = pxd
154
155 self.svd = svd
156
157 def configure_collector(self, collector):
158 """
159 Propagate the hierarchy configuration to the collector
160 """
161 collector.param('hierarchyType', self.type)
162 collector.param('enablePXDHierarchy', self.pxd)
163 collector.param('enableSVDHierarchy', self.svd)
164
165 def generate(self):
166 """
167 Generate the constraints - does nothing. The collector will make the job
168 """
169 print("Generating constraints for VXD (no-op, file created by collector) ...")
170 return []
171
172
174 """
175 Various constraints for CDC layers
176
177 (Code from Claus Kleinwort)
178
179 """
180
181 cdc = [['A1', 8, 160, 16.80, 23.80, 0., -35.9, 67.9],
182 ['U2', 6, 160, 25.70, 34.80, 0.068, -51.4, 98.6],
183 ['A3', 6, 192, 36.52, 45.57, 0., -57.5, 132.9],
184 ['V4', 6, 224, 47.69, 56.69, -0.060, -59.6, 144.7],
185 ['A5', 6, 256, 58.41, 67.41, 0., -61.8, 146.8],
186 ['U6', 6, 288, 69.53, 78.53, 0.064, -63.9, 148.9],
187 ['A7', 6, 320, 80.25, 89.25, 0., -66.0, 151.0],
188 ['V8', 6, 352, 91.37, 100.37, -0.072, -68.2, 153.2],
189 ['A9', 6, 384, 102.00, 111.14, 0., -70.2, 155.3]]
190
191
192 parameter = [(1, 'x-offset bwd'), (2, 'y-offset bwd'), (6, 'z-rotation bwd'),
193 (11, 'x-offset fwd-bwd'), (12, 'y-offset fwd-bwd'), (16, 'z-rotation fwd-bwd')]
194
195 def __init__(self, filename='cdc-constraints.txt', rigid=True, twist=True, z_offset=False, r_scale=False, z_scale=False):
196 """
197 filename : str
198 Filename with constraints (probably generated later)
199 rigid : bool
200 5D CDC constraints: x+y+rot + dx+dy
201 twist : bool
202 CDC twist constraint: drot
203 z_offset : bool
204 Constraint for Z-offset
205 r_scale : bool
206 Constraint for R-scale
207 z_scale : bool
208 Constraint for Z-scale
209 """
210 super().__init__(filename)
211
212 self.rigid = rigid
213
214 self.twist = twist
215
216 self.z_offset = z_offset
217
218 self.r_scale = r_scale
219
220 self.z_scale = z_scale
221
222 def generate(self):
223 """Generate constraints from CDC geometry
224 """
225 print("Generating constraints for CDC layers...")
226
227 def cdc_layer_label(layer, param):
228 wire = 511
229 wireid = Belle2.WireID(layer, wire).getEWire()
230 label = Belle2.GlobalLabel()
231 label.construct(Belle2.CDCAlignment.getGlobalUniqueID(), wireid, param)
232 return label.label()
233
234 def cmp(a, b):
235 return (a > b) - (a < b)
236
237 consts = []
238
239 parameters = self.parameter
240
241 for par in parameters:
242 nTot = 0
243 for sl in self.cdc:
244 nlyr = sl[1]
245 for lyr in range(nlyr):
246 nTot += nlyr
247 nLayer = nTot
248
249 if self.rigid:
250 # all layers
251 for par in self.parameter:
252 if par[0] == 16:
253 continue
254 const = Constraint()
255
256 nTot = 0
257 for sl in self.cdc:
258 nlyr = sl[1]
259 for lyr in range(nlyr):
260 label = cdc_layer_label(nTot + lyr, par[0])
261 const.add(label, 1.0)
262 nTot += nlyr
263
264 consts.append(const)
265
266 if self.twist:
267 # all layers
268 for par in self.parameter:
269 if par[0] != 16:
270 continue
271 const = Constraint()
272
273 nTot = 0
274 for sl in self.cdc:
275 nlyr = sl[1]
276 for lyr in range(nlyr):
277 label = cdc_layer_label(nTot + lyr, par[0])
278 const.add(label, 1.0)
279 nTot += nlyr
280
281 # f.write('\n')
282 consts.append(const)
283
284 if self.z_offset:
285 # stereo layers (Z offset)
286 par = self.parameter[2]
287 const = Constraint()
288
289 nTot = 0
290 for sl in self.cdc:
291 nlyr = sl[1]
292 rInner = sl[3]
293 rOuter = sl[4]
294 stereo = sl[5]
295 if stereo != 0.:
296 dr = (rOuter - rInner) / (nlyr - 1)
297 rWire = rInner
298 for lyr in range(nlyr):
299 label = cdc_layer_label(nTot + lyr, par[0])
300 const.add(label, stereo * rWire)
301 rWire += dr
302 nTot += nlyr
303
304 # f.write('\n')
305 consts.append(const)
306
307 if self.r_scale:
308 # all layers 2nd
309 for par in self.parameter[:2]:
310 const = Constraint()
311
312 nTot = 0
313 for sl in self.cdc:
314 nlyr = sl[1]
315 for lyr in range(nlyr):
316 label = cdc_layer_label(nTot + lyr, par[0])
317 der = cmp(2. * float(nTot + lyr) + 0.5, float(nLayer))
318 const.add(label, der)
319 nTot += nlyr
320
321 consts.append(const)
322
323 if self.z_scale:
324 # stereo layers (Z -scale)
325 par = self.parameter[5]
326 const = Constraint()
327 nTot = 0
328 for sl in self.cdc:
329 nlyr = sl[1]
330 stereo = sl[5]
331 if stereo != 0.:
332 for lyr in range(nlyr):
333 label = cdc_layer_label(nTot + lyr, par[0])
334 const.add(label, stereo)
335 nTot += nlyr
336
337 # f.write('\n')
338 consts.append(const)
339
340 return consts
341
342
344 """
345 Constraint to fix sum of corrections to CDC T0's (per wire) to zero
346 """
347
348
349 wires_in_layer = [
350 160, 160, 160, 160, 160, 160, 160, 160,
351 160, 160, 160, 160, 160, 160,
352 192, 192, 192, 192, 192, 192,
353 224, 224, 224, 224, 224, 224,
354 256, 256, 256, 256, 256, 256,
355 288, 288, 288, 288, 288, 288,
356 320, 320, 320, 320, 320, 320,
357 352, 352, 352, 352, 352, 352,
358 384, 384, 384, 384, 384, 384]
359
360 def __init__(self, filename='cdc-T0-constraints.txt'):
361 """ Initialize
362 filename : str
363 Can use different filename
364 """
365 super().__init__(filename)
366 pass
367
368 def generate(self):
369 """
370 Generate the constraints
371 """
372 print("Generating constraints for CDC T0's ...")
373 consts = []
374 const = Constraint()
375
376 for layer in range(0, 56):
377 for wire in range(0, self.wires_in_layer[layer]):
378 label = Belle2.GlobalLabel()
379 label.construct(Belle2.CDCTimeZeros.getGlobalUniqueID(), Belle2.WireID(layer, wire).getEWire(), 0)
380 const.add(label.label(), 1.0)
381
382 consts.append(const)
383 return consts
384
385
387 """
388 Various constraints for CDC wires (w.r.t. layers)
389 """
390
391
392 wires_in_layer = [
393 160, 160, 160, 160, 160, 160, 160, 160,
394 160, 160, 160, 160, 160, 160,
395 192, 192, 192, 192, 192, 192,
396 224, 224, 224, 224, 224, 224,
397 256, 256, 256, 256, 256, 256,
398 288, 288, 288, 288, 288, 288,
399 320, 320, 320, 320, 320, 320,
400 352, 352, 352, 352, 352, 352,
401 384, 384, 384, 384, 384, 384]
402
404 self,
405 filename='cdc-wire-constraints.txt',
406 layers=None,
407 layer_rigid=True,
408 layer_radius=None,
409 cdc_radius=False,
410 hemisphere=None):
411 """Initialize constraint
412
413 Parameters
414 ----------
415 filename : str
416 Can override the default file name
417 layers : list (int)
418 List of layer numbers for which to generate the constraints (default is 0..55)
419 layer_rigid : bool
420 6 constraints - fix sum of shifts of wire ends in x/y/rotation at fwd/bwd end-plate
421 layer_radius : list(int)
422 2 constraints to fix average change of layer radius (wires in layer moving away from origin) for all layers in list
423 cdc_radius : bool
424 1 constraint - fix average change in CDC radius from all wires
425 hemisphere : list(int)
426 Modifies rigid and layer_radius (layer_radius constraint is added if not present for selected layer(s))
427 constraint by splitting the constraints for a given list of layers to L/R up/down half of CDC
428
429 """
430
431 super().__init__(filename)
432
433 if layers is None:
434 layers = list(range(0, 56))
435
436 self.layers = layers
437
441 self.layer_rigid = layer_rigid
442
444 if layer_radius is None:
445 layer_radius = []
446
447 self.layer_radius = layer_radius
448
450 self.cdc_radius = cdc_radius
451
452 if hemisphere is None:
453 hemisphere = []
454
455 self.hemisphere = hemisphere
456
457 def configure_collector(self, collector):
458 """Enables wire-by-wire derivatives in collector
459 """
460 basf2.B2WARNING("Adding CDC wire constraints -> enabling wire-by-wire alignment derivatives")
461 collector.param('enableWireByWireAlignment', True)
462
463 def get_label(self, layer, wire, parameter):
464 """Return GlobalLabel for wire parameter"""
465
466 wireid = Belle2.WireID(layer, wire).getEWire()
467 label = Belle2.GlobalLabel()
468 label.construct(Belle2.CDCAlignment.getGlobalUniqueID(), wireid, parameter)
469 return label.label()
470
471 def generate(self):
472 """Generate constraints """
473 print("Generating constraints for CDC wires...")
474 consts = []
475
476 layers = self.layers
477
478 layer_rigid = self.layer_rigid
479
480 if layer_rigid:
481 for layer in layers:
482 const = Constraint()
483 # sum of wire X (BWD) in layer
484 for wire in range(0, self.wires_in_layer[layer]):
485 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdX), 1.)
486 consts.append(const)
487
488 for layer in layers:
489 const = Constraint()
490 # sum of wire Y (BWD) in layer
491 for wire in range(0, self.wires_in_layer[layer]):
492 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdY), 1.)
493 consts.append(const)
494
495 for layer in set(layers) - set(self.hemisphere):
496 const = Constraint()
497 # sum of wire rotations (BWD) in layer
498 for wire in range(0, self.wires_in_layer[layer]):
499
500 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
501
502 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdX), -math.sin(wirePhi))
503 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.cos(wirePhi))
504 consts.append(const)
505
506 for layer in layers:
507 const = Constraint()
508 # sum of wire X (FWD) in layer
509 for wire in range(0, self.wires_in_layer[layer]):
510 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdX), 1.)
511 consts.append(const)
512
513 for layer in layers:
514 const = Constraint()
515 # sum of wire Y (FWD) in layer
516 for wire in range(0, self.wires_in_layer[layer]):
517 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdY), 1.)
518 consts.append(const)
519
520 for layer in set(layers) - set(self.hemisphere):
521 const = Constraint()
522 # sum of wire rotations (FWD) in layer
523 for wire in range(0, self.wires_in_layer[layer]):
524
525 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
526
527 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdX), -math.sin(wirePhi))
528 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.cos(wirePhi))
529 consts.append(const)
530
531 for layer in self.hemisphere:
532 const = Constraint()
533 # sum of wire rotations (BWD) in layer ... RIGHT hemisphere
534 for wire in range(0, self.wires_in_layer[layer]):
535
536 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
537
538 if math.cos(wirePhi) <= 0.:
539 continue
540
541 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdX), -math.sin(wirePhi))
542 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.cos(wirePhi))
543 consts.append(const)
544
545 for layer in self.hemisphere:
546 const = Constraint()
547 # sum of wire rotations (BWD) in layer ... LEFT hemisphere
548 for wire in range(0, self.wires_in_layer[layer]):
549
550 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
551
552 if math.cos(wirePhi) > 0.:
553 continue
554
555 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdX), -math.sin(wirePhi))
556 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.cos(wirePhi))
557 consts.append(const)
558
559 for layer in self.hemisphere:
560 const = Constraint()
561 # sum of wire rotations (FWD) in layer ... RIGHT
562 for wire in range(0, self.wires_in_layer[layer]):
563
564 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
565
566 if math.cos(wirePhi) <= 0.:
567 continue
568
569 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdX), -math.sin(wirePhi))
570 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.cos(wirePhi))
571 consts.append(const)
572
573 for layer in self.hemisphere:
574 const = Constraint()
575 # sum of wire rotations (FWD) in layer ... LEFT
576 for wire in range(0, self.wires_in_layer[layer]):
577
578 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
579
580 if math.cos(wirePhi) > 0.:
581 continue
582
583 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdX), -math.sin(wirePhi))
584 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.cos(wirePhi))
585 consts.append(const)
586
587 # layer_radius:
588 for layer in [lyr for lyr in self.layer_radius if lyr not in self.hemisphere]:
589 const = Constraint()
590 # sum of wire rotations (BWD) in layer
591 for wire in range(0, self.wires_in_layer[layer]):
592
593 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
594
595 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdX), +math.cos(wirePhi))
596 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.sin(wirePhi))
597 const.add(self.get_label(layer, 0, 0), 0.)
598 consts.append(const)
599
600 for layer in [lyr for lyr in self.layer_radius if lyr not in self.hemisphere]:
601 const = Constraint()
602 # sum of wire rotations (FWD) in layer
603 for wire in range(0, self.wires_in_layer[layer]):
604
605 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
606
607 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdX), +math.cos(wirePhi))
608 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.sin(wirePhi))
609 const.add(self.get_label(layer, 0, 0), 0.)
610 consts.append(const)
611
612 # hemisphere:
613 for layer in self.hemisphere:
614 const = Constraint()
615 # sum of wire rotations (BWD) in layer
616 for wire in range(0, self.wires_in_layer[layer]):
617
618 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
619
620 if math.sin(wirePhi) >= 0:
621 continue
622
623 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdX), +math.cos(wirePhi))
624 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.sin(wirePhi))
625 const.add(self.get_label(layer, 0, 0), 0.)
626 consts.append(const)
627
628 for layer in self.hemisphere:
629 const = Constraint()
630 # sum of wire rotations (BWD) in layer
631 for wire in range(0, self.wires_in_layer[layer]):
632
633 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
634
635 if math.sin(wirePhi) < 0:
636 continue
637
638 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdX), +math.cos(wirePhi))
639 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.sin(wirePhi))
640 const.add(self.get_label(layer, 0, 0), 0.)
641 consts.append(const)
642
643 for layer in self.hemisphere:
644 const = Constraint()
645 # sum of wire rotations (FWD) in layer
646 for wire in range(0, self.wires_in_layer[layer]):
647
648 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
649
650 if math.sin(wirePhi) >= 0:
651 continue
652
653 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdX), +math.cos(wirePhi))
654 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.sin(wirePhi))
655 const.add(self.get_label(layer, 0, 0), 0.)
656 consts.append(const)
657
658 for layer in self.hemisphere:
659 const = Constraint()
660 # sum of wire rotations (FWD) in layer
661 for wire in range(0, self.wires_in_layer[layer]):
662
663 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
664
665 if math.sin(wirePhi) < 0:
666 continue
667
668 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdX), +math.cos(wirePhi))
669 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.sin(wirePhi))
670 const.add(self.get_label(layer, 0, 0), 0.)
671 consts.append(const)
672
673 if self.cdc_radius:
674 const = Constraint()
675 for layer in layers:
676 # sum of wire rotations (BWD) in layer
677 for wire in range(0, self.wires_in_layer[layer]):
678 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
679 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdX), +math.cos(wirePhi))
680 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.sin(wirePhi))
681 consts.append(const)
682
683 const = Constraint()
684 for layer in layers:
685 # sum of wire rotations (FWD) in layer
686 for wire in range(0, self.wires_in_layer[layer]):
687 wirePhi = Belle2.CDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
688 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdX), +math.cos(wirePhi))
689 const.add(self.get_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.sin(wirePhi))
690 consts.append(const)
691
692 return consts
693
694
695# ------------ Main: Generate some constraint files with default config (no time-dependence, default global tags) ------
696if __name__ == '__main__':
697 ## 6D CDC constraints
698 consts6 = CDCLayerConstraints('cdc-layer-constraints-6D.txt', rigid=True, z_offset=False, r_scale=False, z_scale=False)
699 ## 7D CDC constraints
700 consts7 = CDCLayerConstraints('cdc-layer-constraints-7D.txt', rigid=True, z_offset=True, r_scale=False, z_scale=False)
701 ## 10D CD constraints
702 consts10 = CDCLayerConstraints('cdc-layer-constraints-10D.txt', rigid=True, z_offset=True, r_scale=True, z_scale=True)
703 ## CDC EventT0 constraint
704 cdcT0 = CDCTimeZerosConstraint()
705 ## CDC wire constraints
706 cdcWires = CDCWireConstraints(
707 filename='cdc-wire-constraints-proc12.txt',
708 layers=None,
709 layer_rigid=True,
710 layer_radius=[53],
711 cdc_radius=True,
712 hemisphere=[55])
713
714 # phase 2
715 # timedep = [([], [(0, 0, 1002)])]
716 # Run 1 (early phase 3)
717 # timedep = [([], [(0, 0, 1003)])]
718
719 # final detector (phase 3)
720 ## Time dependency
721 timedep = [] # [([], [(0, 0, 0)])]
722 ## Initial event
723 init_event = (0, 0, 0)
724 ## Files
725 files = generate_constraints(
726 [
727 consts6,
728 consts7,
729 consts10,
730 cdcT0,
731 cdcWires,
732 VXDHierarchyConstraints(type=1, pxd=False),
733 Constraints("my_file.txt")],
734
735 timedep=timedep,
736 global_tags=None,
737 init_event=init_event)
static unsigned short getGlobalUniqueID()
Get global unique id.
static unsigned short getGlobalUniqueID()
Get global unique id.
static std::string findFile(const std::string &path, bool silent=false)
Search for given file or directory in local or central release directory, and return absolute path if...
Class to convert to/from global labels for Millepede II to/from detector & parameter identificators.
Definition GlobalLabel.h:41
Class to identify a wire inside the CDC.
Definition WireID.h:34
__init__(self, filename='cdc-constraints.txt', rigid=True, twist=True, z_offset=False, r_scale=False, z_scale=False)
list cdc
CDC geometry hard-code -> TODO.
__init__(self, filename='cdc-T0-constraints.txt')
list wires_in_layer
CDC layer/wire configuration.
layer_rigid
6 x 56 (6/layer) constraints.
get_label(self, layer, wire, parameter)
cdc_radius
2 Constraints: Sum(dr)=0 for all wires in CDC at each end-plate -> "average CDC radius" kept same by ...
__init__(self, filename='cdc-wire-constraints.txt', layers=None, layer_rigid=True, layer_radius=None, cdc_radius=False, hemisphere=None)
__init__(self, comment="", value=0.)
configure_collector(self, collector)
__init__(self, type=2, pxd=True, svd=True)