Belle II Software  release-08-01-10
constraints.py
1 
8 
9 import basf2
10 from ROOT import Belle2
11 
12 import os
13 import math
14 
15 
16 class 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.valuevalue = value
32 
33  self.commentcomment = comment
34 
35  self.datadata = []
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.datadata.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.datadata]
55  checksum = hash(str(labels_only))
56  return checksum
57 
58 
59 class Constraints():
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.filenamefilename = filename
73 
74  def generate(self):
75  """
76  Should be overriden by the deriving classes and actually
77  fill the dictionary with the constraints (this runs withing 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 overriden be child classes to pass additional configuration to the
86  MillepedeCollector (activated by the use of the constraints)
87  """
88  pass
89 
90 
91 def 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.typetype = type
152 
153  self.pxdpxd = pxd
154 
155  self.svdsvd = svd
156 
157  def configure_collector(self, collector):
158  """
159  Propagate the hierarchy configuration to the collector
160  """
161  collector.param('hierarchyType', self.typetype)
162  collector.param('enablePXDHierarchy', self.pxdpxd)
163  collector.param('enableSVDHierarchy', self.svdsvd)
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.rigidrigid = rigid
213 
214  self.twisttwist = twist
215 
216  self.z_offsetz_offset = z_offset
217 
218  self.r_scaler_scale = r_scale
219 
220  self.z_scalez_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  for par in self.parameterparameter:
240  nTot = 0
241  for sl in self.cdccdc:
242  nlyr = sl[1]
243  for lyr in range(nlyr):
244  nTot += nlyr
245  nLayer = nTot
246 
247  if self.rigidrigid:
248  # all layers
249  for par in self.parameterparameter:
250  if par[0] == 16:
251  continue
252  const = Constraint()
253 
254  nTot = 0
255  for sl in self.cdccdc:
256  nlyr = sl[1]
257  for lyr in range(nlyr):
258  label = cdc_layer_label(nTot + lyr, par[0])
259  const.add(label, 1.0)
260  nTot += nlyr
261 
262  consts.append(const)
263 
264  if self.twisttwist:
265  # all layers
266  for par in self.parameterparameter:
267  if par[0] != 16:
268  continue
269  const = Constraint()
270 
271  nTot = 0
272  for sl in self.cdccdc:
273  nlyr = sl[1]
274  for lyr in range(nlyr):
275  label = cdc_layer_label(nTot + lyr, par[0])
276  const.add(label, 1.0)
277  nTot += nlyr
278 
279  # f.write('\n')
280  consts.append(const)
281 
282  if self.z_offsetz_offset:
283  # stereo layers (Z offset)
284  par = self.parameterparameter[2]
285  const = Constraint()
286 
287  nTot = 0
288  for sl in self.cdccdc:
289  nlyr = sl[1]
290  rInner = sl[3]
291  rOuter = sl[4]
292  stereo = sl[5]
293  if stereo != 0.:
294  dr = (rOuter - rInner) / (nlyr - 1)
295  rWire = rInner
296  for lyr in range(nlyr):
297  label = cdc_layer_label(nTot + lyr, par[0])
298  const.add(label, stereo * rWire)
299  rWire += dr
300  nTot += nlyr
301 
302  # f.write('\n')
303  consts.append(const)
304 
305  if self.r_scaler_scale:
306  # all layers 2nd
307  for par in self.parameterparameter[:2]:
308  const = Constraint()
309 
310  nTot = 0
311  for sl in self.cdccdc:
312  nlyr = sl[1]
313  for lyr in range(nlyr):
314  label = cdc_layer_label(nTot + lyr, par[0])
315  der = cmp(2. * float(nTot + lyr) + 0.5, float(nLayer))
316  const.add(label, der)
317  nTot += nlyr
318 
319  consts.append(const)
320 
321  if self.z_scalez_scale:
322  # stereo layers (Z -scale)
323  par = self.parameterparameter[5]
324  const = Constraint()
325  nTot = 0
326  for sl in self.cdccdc:
327  nlyr = sl[1]
328  stereo = sl[5]
329  if stereo != 0.:
330  for lyr in range(nlyr):
331  label = cdc_layer_label(nTot + lyr, par[0])
332  const.add(label, stereo)
333  nTot += nlyr
334 
335  # f.write('\n')
336  consts.append(const)
337 
338  return consts
339 
340 
342  """
343  Constraint to fix sum of corrections to CDC T0's (per wire) to zero
344  """
345 
346 
347  wires_in_layer = [
348  160, 160, 160, 160, 160, 160, 160, 160,
349  160, 160, 160, 160, 160, 160,
350  192, 192, 192, 192, 192, 192,
351  224, 224, 224, 224, 224, 224,
352  256, 256, 256, 256, 256, 256,
353  288, 288, 288, 288, 288, 288,
354  320, 320, 320, 320, 320, 320,
355  352, 352, 352, 352, 352, 352,
356  384, 384, 384, 384, 384, 384]
357 
358  def __init__(self, filename='cdc-T0-constraints.txt'):
359  """ Initialize
360  filename : str
361  Can use different filename
362  """
363  super().__init__(filename)
364  pass
365 
366  def generate(self):
367  """
368  Generate the constraints
369  """
370  print("Generating constraints for CDC T0's ...")
371  consts = []
372  const = Constraint()
373 
374  for layer in range(0, 56):
375  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
376  label = Belle2.GlobalLabel()
377  label.construct(Belle2.CDCTimeZeros.getGlobalUniqueID(), Belle2.WireID(layer, wire).getEWire(), 0)
378  const.add(label.label(), 1.0)
379 
380  consts.append(const)
381  return consts
382 
383 
385  """
386  Various constraints for CDC wires (w.r.t. layers)
387  """
388 
389 
390  wires_in_layer = [
391  160, 160, 160, 160, 160, 160, 160, 160,
392  160, 160, 160, 160, 160, 160,
393  192, 192, 192, 192, 192, 192,
394  224, 224, 224, 224, 224, 224,
395  256, 256, 256, 256, 256, 256,
396  288, 288, 288, 288, 288, 288,
397  320, 320, 320, 320, 320, 320,
398  352, 352, 352, 352, 352, 352,
399  384, 384, 384, 384, 384, 384]
400 
401  def __init__(
402  self,
403  filename='cdc-wire-constraints.txt',
404  layers=None,
405  layer_rigid=True,
406  layer_radius=None,
407  cdc_radius=False,
408  hemisphere=None):
409  """Initialize constraint
410 
411  Parameters
412  ----------
413  filename : str
414  Can override the default file name
415  layers : list (int)
416  List of layer numbers for which to generate the constraints (default is 0..55)
417  layer_rigid : bool
418  6 constraints - fix sum of shifts of wire ends in x/y/rotation at fwd/bwd end-plate
419  layer_radius : list(int)
420  2 constraints to fix average change of layer radius (wires in layer moving away from origin) for all layers in list
421  cdc_radius : bool
422  1 constraint - fix average change in CDC radius from all wires
423  hemisphere : list(int)
424  Modifies rigid and layer_radius (layer_radius constraint is added if not present for selected layer(s))
425  constraint by splitting the constraints for a given list of layers to L/R up/down half of CDC
426 
427  """
428 
429  super().__init__(filename)
430 
431  if layers is None:
432  layers = list(range(0, 56))
433  self.layerslayers = layers
434 
438  self.layer_rigidlayer_rigid = layer_rigid
439 
441  if layer_radius is None:
442  layer_radius = []
443  self.layer_radiuslayer_radius = layer_radius
444 
446  self.cdc_radiuscdc_radius = cdc_radius
447 
448  if hemisphere is None:
449  hemisphere = []
450  self.hemispherehemisphere = hemisphere
451 
452  def configure_collector(self, collector):
453  """Enables wire-by-wire derivatives in collector
454  """
455  basf2.B2WARNING("Adding CDC wire constraints -> enabling wire-by-wire alignment derivatives")
456  collector.param('enableWireByWireAlignment', True)
457 
458  def get_label(self, layer, wire, parameter):
459  """Return GlobalLabel for wire parameter"""
460 
461  wireid = Belle2.WireID(layer, wire).getEWire()
462  label = Belle2.GlobalLabel()
463  label.construct(Belle2.CDCAlignment.getGlobalUniqueID(), wireid, parameter)
464  return label.label()
465 
466  def generate(self):
467  """Generate constraints """
468  print("Generating constraints for CDC wires...")
469  consts = []
470 
471  layers = self.layerslayers
472 
473  if self.layer_rigidlayer_rigid:
474  for layer in layers:
475  const = Constraint()
476  # sum of wire X (BWD) in layer
477  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
478  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdX), 1.)
479  consts.append(const)
480 
481  for layer in layers:
482  const = Constraint()
483  # sum of wire Y (BWD) in layer
484  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
485  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdY), 1.)
486  consts.append(const)
487 
488  for layer in set(layers) - set(self.hemispherehemisphere):
489  const = Constraint()
490  # sum of wire rotations (BWD) in layer
491  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
492 
493  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
494 
495  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdX), -math.sin(wirePhi))
496  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.cos(wirePhi))
497  consts.append(const)
498 
499  for layer in layers:
500  const = Constraint()
501  # sum of wire X (FWD) in layer
502  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
503  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdX), 1.)
504  consts.append(const)
505 
506  for layer in layers:
507  const = Constraint()
508  # sum of wire Y (FWD) in layer
509  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
510  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdY), 1.)
511  consts.append(const)
512 
513  for layer in set(layers) - set(self.hemispherehemisphere):
514  const = Constraint()
515  # sum of wire rotations (FWD) in layer
516  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
517 
518  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
519 
520  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdX), -math.sin(wirePhi))
521  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.cos(wirePhi))
522  consts.append(const)
523 
524  for layer in self.hemispherehemisphere:
525  const = Constraint()
526  # sum of wire rotations (BWD) in layer ... RIGHT hemisphere
527  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
528 
529  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
530 
531  if math.cos(wirePhi) <= 0.:
532  continue
533 
534  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdX), -math.sin(wirePhi))
535  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.cos(wirePhi))
536  consts.append(const)
537 
538  for layer in self.hemispherehemisphere:
539  const = Constraint()
540  # sum of wire rotations (BWD) in layer ... LEFT hemisphere
541  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
542 
543  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
544 
545  if math.cos(wirePhi) > 0.:
546  continue
547 
548  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdX), -math.sin(wirePhi))
549  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.cos(wirePhi))
550  consts.append(const)
551 
552  for layer in self.hemispherehemisphere:
553  const = Constraint()
554  # sum of wire rotations (FWD) in layer ... RIGHT
555  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
556 
557  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
558 
559  if math.cos(wirePhi) <= 0.:
560  continue
561 
562  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdX), -math.sin(wirePhi))
563  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.cos(wirePhi))
564  consts.append(const)
565 
566  for layer in self.hemispherehemisphere:
567  const = Constraint()
568  # sum of wire rotations (FWD) in layer ... LEFT
569  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
570 
571  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
572 
573  if math.cos(wirePhi) > 0.:
574  continue
575 
576  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdX), -math.sin(wirePhi))
577  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.cos(wirePhi))
578  consts.append(const)
579 
580  # layer_radius:
581  for layer in [lyr for lyr in self.layer_radiuslayer_radius if lyr not in self.hemispherehemisphere]:
582  const = Constraint()
583  # sum of wire rotations (BWD) in layer
584  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
585 
586  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
587 
588  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdX), +math.cos(wirePhi))
589  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.sin(wirePhi))
590  const.add(self.get_labelget_label(layer, 0, 0), 0.)
591  consts.append(const)
592 
593  for layer in [lyr for lyr in self.layer_radiuslayer_radius if lyr not in self.hemispherehemisphere]:
594  const = Constraint()
595  # sum of wire rotations (FWD) in layer
596  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
597 
598  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
599 
600  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdX), +math.cos(wirePhi))
601  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.sin(wirePhi))
602  const.add(self.get_labelget_label(layer, 0, 0), 0.)
603  consts.append(const)
604 
605  # hemisphere:
606  for layer in self.hemispherehemisphere:
607  const = Constraint()
608  # sum of wire rotations (BWD) in layer
609  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
610 
611  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
612 
613  if math.sin(wirePhi) >= 0:
614  continue
615 
616  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdX), +math.cos(wirePhi))
617  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.sin(wirePhi))
618  const.add(self.get_labelget_label(layer, 0, 0), 0.)
619  consts.append(const)
620 
621  for layer in self.hemispherehemisphere:
622  const = Constraint()
623  # sum of wire rotations (BWD) in layer
624  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
625 
626  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
627 
628  if math.sin(wirePhi) < 0:
629  continue
630 
631  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdX), +math.cos(wirePhi))
632  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.sin(wirePhi))
633  const.add(self.get_labelget_label(layer, 0, 0), 0.)
634  consts.append(const)
635 
636  for layer in self.hemispherehemisphere:
637  const = Constraint()
638  # sum of wire rotations (FWD) in layer
639  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
640 
641  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
642 
643  if math.sin(wirePhi) >= 0:
644  continue
645 
646  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdX), +math.cos(wirePhi))
647  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.sin(wirePhi))
648  const.add(self.get_labelget_label(layer, 0, 0), 0.)
649  consts.append(const)
650 
651  for layer in self.hemispherehemisphere:
652  const = Constraint()
653  # sum of wire rotations (FWD) in layer
654  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
655 
656  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
657 
658  if math.sin(wirePhi) < 0:
659  continue
660 
661  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdX), +math.cos(wirePhi))
662  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.sin(wirePhi))
663  const.add(self.get_labelget_label(layer, 0, 0), 0.)
664  consts.append(const)
665 
666  if self.cdc_radiuscdc_radius:
667  const = Constraint()
668  for layer in layers:
669  # sum of wire rotations (BWD) in layer
670  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
671  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getBackwardPos3D().phi()
672  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdX), +math.cos(wirePhi))
673  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireBwdY), +math.sin(wirePhi))
674  consts.append(const)
675 
676  const = Constraint()
677  for layer in layers:
678  # sum of wire rotations (FWD) in layer
679  for wire in range(0, self.wires_in_layerwires_in_layer[layer]):
680  wirePhi = Belle2.TrackFindingCDC.CDCWire.getInstance(Belle2.WireID(layer, wire)).getForwardPos3D().phi()
681  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdX), +math.cos(wirePhi))
682  const.add(self.get_labelget_label(layer, wire, Belle2.CDCAlignment.wireFwdY), +math.sin(wirePhi))
683  consts.append(const)
684 
685  return consts
686 
687 
688 # ------------ Main: Generate some constraint files with default config (no time-dependence, default global tags) ------
689 if __name__ == '__main__':
690 
691  consts6 = CDCLayerConstraints('cdc-layer-constraints-6D.txt', rigid=True, z_offset=False, r_scale=False, z_scale=False)
692 
693  consts7 = CDCLayerConstraints('cdc-layer-constraints-7D.txt', rigid=True, z_offset=True, r_scale=False, z_scale=False)
694 
695  consts10 = CDCLayerConstraints('cdc-layer-constraints-10D.txt', rigid=True, z_offset=True, r_scale=True, z_scale=True)
696 
697  cdcT0 = CDCTimeZerosConstraint()
698 
699  cdcWires = CDCWireConstraints(
700  filename='cdc-wire-constraints-proc12.txt',
701  layers=None,
702  layer_rigid=True,
703  layer_radius=[53],
704  cdc_radius=True,
705  hemisphere=[55])
706 
707  # phase 2
708  # timedep = [([], [(0, 0, 1002)])]
709  # early phase 3
710  # timedep = [([], [(0, 0, 1003)])]
711 
712  # final detector (phase 3)
713 
714  timedep = [] # [([], [(0, 0, 0)])]
715 
716  init_event = (0, 0, 0)
717 
718  files = generate_constraints(
719  [
720  consts6,
721  consts7,
722  consts10,
723  cdcT0,
724  cdcWires,
725  VXDHierarchyConstraints(type=1, pxd=False),
726  Constraints("my_file.txt")],
727 
728  timedep=timedep,
729  global_tags=None,
730  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:148
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:181
def __init__(self, filename='cdc-constraints.txt', rigid=True, twist=True, z_offset=False, r_scale=False, z_scale=False)
Definition: constraints.py:195
def __init__(self, filename='cdc-T0-constraints.txt')
Definition: constraints.py:358
list wires_in_layer
CDC layer/wire configuration.
Definition: constraints.py:347
layer_rigid
6 x 56 (6/layer) constraints.
Definition: constraints.py:438
def __init__(self, filename='cdc-wire-constraints.txt', layers=None, layer_rigid=True, layer_radius=None, cdc_radius=False, hemisphere=None)
Definition: constraints.py:408
layers
List of layers for whose wires to generate the constraints.
Definition: constraints.py:433
def configure_collector(self, collector)
Definition: constraints.py:452
hemisphere
list of layer subject to hemisphere constraints
Definition: constraints.py:450
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:443
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:446
list wires_in_layer
CDC layer/wire configuration.
Definition: constraints.py:390
def get_label(self, layer, wire, parameter)
Definition: constraints.py:458
def __init__(self, comment="", value=0.)
Definition: constraints.py:21
def add(self, label, coeff)
Definition: constraints.py:37
def configure_collector(self, collector)
Definition: constraints.py:83
def __init__(self, filename)
Definition: constraints.py:66
def __init__(self, type=2, pxd=True, svd=True)
Definition: constraints.py:138