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