Belle II Software  release-08-01-10
TransformData.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 /* Own header. */
10 #include <klm/eklm/geometry/TransformData.h>
11 
12 /* KLM headers. */
13 #include <klm/dbobjects/eklm/EKLMAlignment.h>
14 #include <klm/dbobjects/eklm/EKLMSegmentAlignment.h>
15 #include <klm/eklm/geometry/AlignmentChecker.h>
16 #include <klm/eklm/geometry/GeometryData.h>
17 
18 /* Basf2 headers. */
19 #include <framework/database/DBObjPtr.h>
20 #include <framework/logging/Logger.h>
21 #include <framework/gearbox/Unit.h>
22 
23 using namespace Belle2;
24 
25 EKLM::TransformData::TransformData(bool global, Displacement displacementType)
26 {
27  /* cppcheck-suppress variableScope */
28  int iSection, iLayer, iSector, iPlane, iSegment, iStrip, sector, segment;
29  int nSections, nLayers, nSectors, nPlanes, nStrips, nSegments, nStripsSegment;
30  int nDetectorLayers;
31  AlignmentChecker alignmentChecker(true);
34  nSections = m_GeoDat->getNSections();
35  nSectors = m_GeoDat->getNSectors();
36  nLayers = m_GeoDat->getNLayers();
37  nPlanes = m_GeoDat->getNPlanes();
38  nStrips = m_GeoDat->getNStrips();
39  nSegments = m_GeoDat->getNSegments();
40  nStripsSegment = m_ElementNumbers->getNStripsSegment();
41  m_Section = new HepGeom::Transform3D[nSections];
42  m_Layer = new HepGeom::Transform3D*[nSections];
43  m_Sector = new HepGeom::Transform3D** [nSections];
44  m_Plane = new HepGeom::Transform3D** *[nSections];
45  m_PlaneDisplacement = new HepGeom::Transform3D** *[nSections];
46  m_Segment = new HepGeom::Transform3D**** [nSections];
47  m_Strip = new HepGeom::Transform3D**** [nSections];
48  m_StripInverse = new HepGeom::Transform3D**** [nSections];
49  for (iSection = 0; iSection < nSections; iSection++) {
50  m_GeoDat->getSectionTransform(&m_Section[iSection], iSection);
51  nDetectorLayers = m_GeoDat->getNDetectorLayers(iSection + 1);
52  m_Layer[iSection] = new HepGeom::Transform3D[nLayers];
53  m_Sector[iSection] = new HepGeom::Transform3D*[nLayers];
54  m_Plane[iSection] = new HepGeom::Transform3D** [nLayers];
55  m_PlaneDisplacement[iSection] = new HepGeom::Transform3D** [nLayers];
56  m_Segment[iSection] = new HepGeom::Transform3D** *[nLayers];
57  m_Strip[iSection] = new HepGeom::Transform3D** *[nLayers];
58  m_StripInverse[iSection] = new HepGeom::Transform3D** *[nLayers];
59  for (iLayer = 0; iLayer < nLayers; iLayer++) {
60  m_GeoDat->getLayerTransform(&m_Layer[iSection][iLayer], iLayer);
61  m_Sector[iSection][iLayer] = new HepGeom::Transform3D[nSectors];
62  if (iLayer < nDetectorLayers) {
63  m_Plane[iSection][iLayer] = new HepGeom::Transform3D*[nSectors];
64  m_PlaneDisplacement[iSection][iLayer] =
65  new HepGeom::Transform3D*[nSectors];
66  m_Segment[iSection][iLayer] = new HepGeom::Transform3D** [nSectors];
67  m_Strip[iSection][iLayer] = new HepGeom::Transform3D** [nSectors];
68  m_StripInverse[iSection][iLayer] = new HepGeom::Transform3D** [nSectors];
69  }
70  for (iSector = 0; iSector < nSectors; iSector++) {
71  m_GeoDat->getSectorTransform(&m_Sector[iSection][iLayer][iSector],
72  iSector);
73  if (iLayer >= nDetectorLayers)
74  continue;
75  m_Plane[iSection][iLayer][iSector] = new HepGeom::Transform3D[nPlanes];
76  m_PlaneDisplacement[iSection][iLayer][iSector] =
77  new HepGeom::Transform3D[nPlanes];
78  m_Segment[iSection][iLayer][iSector] =
79  new HepGeom::Transform3D*[nPlanes];
80  m_Strip[iSection][iLayer][iSector] = new HepGeom::Transform3D*[nPlanes];
81  m_StripInverse[iSection][iLayer][iSector] =
82  new HepGeom::Transform3D*[nPlanes];
83  for (iPlane = 0; iPlane < nPlanes; iPlane++) {
85  &m_Plane[iSection][iLayer][iSector][iPlane], iPlane);
86  m_PlaneDisplacement[iSection][iLayer][iSector][iPlane] =
87  HepGeom::Translate3D(0, 0, 0);
88  m_Segment[iSection][iLayer][iSector][iPlane] =
89  new HepGeom::Transform3D[nSegments];
90  for (iSegment = 0; iSegment < nSegments; iSegment++) {
91  m_Segment[iSection][iLayer][iSector][iPlane][iSegment] =
92  HepGeom::Translate3D(0, 0, 0);
93  }
94  m_Strip[iSection][iLayer][iSector][iPlane] =
95  new HepGeom::Transform3D[nStrips];
96  m_StripInverse[iSection][iLayer][iSector][iPlane] =
97  new HepGeom::Transform3D[nStrips];
98  for (iStrip = 0; iStrip < nStrips; iStrip++) {
100  &m_Strip[iSection][iLayer][iSector][iPlane][iStrip], iStrip);
101  }
102  }
103  }
104  }
105  }
106  /* Read alignment data from the database and modify transformations. */
107  if (displacementType != c_None) {
108  std::string payload, segmentPayload;
109  if (displacementType == c_Displacement) {
110  payload = "EKLMDisplacement";
111  segmentPayload = "EKLMSegmentDisplacement";
112  } else {
113  payload = "EKLMAlignment";
114  segmentPayload = "EKLMSegmentAlignment";
115  }
116  DBObjPtr<EKLMAlignment> alignment(payload);
117  DBObjPtr<EKLMSegmentAlignment> segmentAlignment(segmentPayload);
118  if (!alignment.isValid())
119  B2FATAL("No EKLM displacement (alignment) data.");
120  if (displacementType == c_Displacement) {
121  if (!alignmentChecker.checkAlignment(&(*alignment), &(*segmentAlignment)))
122  B2FATAL("EKLM displacement data are incorrect, overlaps exist.");
123  }
124  for (iSection = 1; iSection <= nSections; iSection++) {
125  nDetectorLayers = m_GeoDat->getNDetectorLayers(iSection);
126  for (iLayer = 1; iLayer <= nDetectorLayers; iLayer++) {
127  for (iSector = 1; iSector <= nSectors; iSector++) {
128  sector = m_ElementNumbers->sectorNumber(iSection, iLayer, iSector);
129  const KLMAlignmentData* sectorAlignment =
130  alignment->getModuleAlignment(sector);
131  if (sectorAlignment == nullptr)
132  B2FATAL("Incomplete EKLM displacement (alignment) data.");
133  for (iPlane = 1; iPlane <= nPlanes; iPlane++) {
134  /* First plane is rotated. */
135  if (iPlane == 1) {
136  m_PlaneDisplacement[iSection - 1][iLayer - 1][iSector - 1][iPlane - 1] =
137  HepGeom::Translate3D(
138  sectorAlignment->getDeltaV() * CLHEP::cm / Unit::cm,
139  sectorAlignment->getDeltaU() * CLHEP::cm / Unit::cm, 0) *
140  HepGeom::RotateZ3D(-sectorAlignment->getDeltaGamma() *
141  CLHEP::rad / Unit::rad);
142  } else {
143  m_PlaneDisplacement[iSection - 1][iLayer - 1][iSector - 1][iPlane - 1] =
144  HepGeom::Translate3D(
145  sectorAlignment->getDeltaU() * CLHEP::cm / Unit::cm,
146  sectorAlignment->getDeltaV() * CLHEP::cm / Unit::cm, 0) *
147  HepGeom::RotateZ3D(sectorAlignment->getDeltaGamma() *
148  CLHEP::rad / Unit::rad);
149  }
150  for (iSegment = 1; iSegment <= nSegments; iSegment++) {
151  segment = m_ElementNumbers->segmentNumber(
152  iSection, iLayer, iSector, iPlane, iSegment);
153  const KLMAlignmentData* segmentAlignmentData =
154  segmentAlignment->getSegmentAlignment(segment);
155  if (segmentAlignmentData == nullptr)
156  B2FATAL("Incomplete EKLM displacement (alignment) data.");
157  m_Segment[iSection - 1][iLayer - 1][iSector - 1][iPlane - 1]
158  [iSegment - 1] =
159  HepGeom::Translate3D(
160  segmentAlignmentData->getDeltaU() * CLHEP::cm / Unit::cm,
161  segmentAlignmentData->getDeltaV() * CLHEP::cm / Unit::cm, 0) *
162  m_Segment[iSection - 1][iLayer - 1][iSector - 1][iPlane - 1]
163  [iSegment - 1] *
164  HepGeom::RotateZ3D(segmentAlignmentData->getDeltaGamma() *
165  CLHEP::rad / Unit::rad);
166  for (iStrip = 1; iStrip <= nStripsSegment; iStrip++) {
167  m_Strip[iSection - 1][iLayer - 1][iSector - 1][iPlane - 1]
168  [nStripsSegment * (iSegment - 1) + iStrip - 1] =
169  HepGeom::Translate3D(
170  segmentAlignmentData->getDeltaU() * CLHEP::cm / Unit::cm,
171  segmentAlignmentData->getDeltaV() * CLHEP::cm / Unit::cm,
172  0) *
173  m_Strip[iSection - 1][iLayer - 1][iSector - 1][iPlane - 1]
174  [nStripsSegment * (iSegment - 1) + iStrip - 1] *
175  HepGeom::RotateZ3D(segmentAlignmentData->getDeltaGamma() *
176  CLHEP::rad / Unit::rad);
177  }
178  }
179  }
180  }
181  }
182  }
183  }
184  if (global)
186 }
187 
189 {
190  int iSection, iLayer, iSector, iPlane;
191  /* cppcheck-suppress variableScope */
192  int nSections, nLayers, nDetectorLayers, nSectors, nPlanes;
193  nSections = m_GeoDat->getNSections();
194  nLayers = m_GeoDat->getNLayers();
195  nSectors = m_GeoDat->getNSectors();
196  nPlanes = m_GeoDat->getNPlanes();
197  for (iSection = 0; iSection < nSections; iSection++) {
198  nDetectorLayers = m_GeoDat->getNDetectorLayers(iSection + 1);
199  for (iLayer = 0; iLayer < nLayers; iLayer++) {
200  delete[] m_Sector[iSection][iLayer];
201  if (iLayer >= nDetectorLayers)
202  continue;
203  for (iSector = 0; iSector < nSectors; iSector++) {
204  for (iPlane = 0; iPlane < nPlanes; iPlane++) {
205  delete[] m_Segment[iSection][iLayer][iSector][iPlane];
206  delete[] m_Strip[iSection][iLayer][iSector][iPlane];
207  delete[] m_StripInverse[iSection][iLayer][iSector][iPlane];
208  }
209  delete[] m_Plane[iSection][iLayer][iSector];
210  delete[] m_PlaneDisplacement[iSection][iLayer][iSector];
211  delete[] m_Segment[iSection][iLayer][iSector];
212  delete[] m_Strip[iSection][iLayer][iSector];
213  delete[] m_StripInverse[iSection][iLayer][iSector];
214  }
215  delete[] m_Plane[iSection][iLayer];
216  delete[] m_PlaneDisplacement[iSection][iLayer];
217  delete[] m_Segment[iSection][iLayer];
218  delete[] m_Strip[iSection][iLayer];
219  delete[] m_StripInverse[iSection][iLayer];
220  }
221  delete[] m_Layer[iSection];
222  delete[] m_Sector[iSection];
223  delete[] m_PlaneDisplacement[iSection];
224  delete[] m_Plane[iSection];
225  delete[] m_Segment[iSection];
226  delete[] m_Strip[iSection];
227  delete[] m_StripInverse[iSection];
228  }
229  delete[] m_Section;
230  delete[] m_Layer;
231  delete[] m_Sector;
232  delete[] m_PlaneDisplacement;
233  delete[] m_Plane;
234  delete[] m_Segment;
235  delete[] m_Strip;
236  delete[] m_StripInverse;
237 }
238 
240 {
241  int iSection, iLayer, iSector, iPlane, iSegment, iStrip;
242  /* cppcheck-suppress variableScope */
243  int nSections, nLayers, nDetectorLayers, nSectors, nPlanes, nSegments, nStrips;
244  nSections = m_GeoDat->getNSections();
245  nLayers = m_GeoDat->getNLayers();
246  nSectors = m_GeoDat->getNSectors();
247  nPlanes = m_GeoDat->getNPlanes();
248  nSegments = m_GeoDat->getNSegments();
249  nStrips = m_GeoDat->getNStrips();
250  for (iSection = 0; iSection < nSections; iSection++) {
251  nDetectorLayers = m_GeoDat->getNDetectorLayers(iSection + 1);
252  for (iLayer = 0; iLayer < nLayers; iLayer++) {
253  m_Layer[iSection][iLayer] = m_Section[iSection] * m_Layer[iSection][iLayer];
254  for (iSector = 0; iSector < nSectors; iSector++) {
255  m_Sector[iSection][iLayer][iSector] =
256  m_Layer[iSection][iLayer] * m_Sector[iSection][iLayer][iSector];
257  if (iLayer >= nDetectorLayers)
258  continue;
259  for (iPlane = 0; iPlane < nPlanes; iPlane++) {
260  m_Plane[iSection][iLayer][iSector][iPlane] =
261  m_Sector[iSection][iLayer][iSector] *
262  m_Plane[iSection][iLayer][iSector][iPlane];
263  for (iSegment = 0; iSegment < nSegments; iSegment++) {
264  m_Segment[iSection][iLayer][iSector][iPlane][iSegment] =
265  m_Plane[iSection][iLayer][iSector][iPlane] *
266  m_PlaneDisplacement[iSection][iLayer][iSector][iPlane] *
267  m_Segment[iSection][iLayer][iSector][iPlane][iSegment];
268  }
269  for (iStrip = 0; iStrip < nStrips; iStrip++) {
270  m_Strip[iSection][iLayer][iSector][iPlane][iStrip] =
271  m_Plane[iSection][iLayer][iSector][iPlane] *
272  m_PlaneDisplacement[iSection][iLayer][iSector][iPlane] *
273  m_Strip[iSection][iLayer][iSector][iPlane][iStrip];
274  m_StripInverse[iSection][iLayer][iSector][iPlane][iStrip] =
275  m_Strip[iSection][iLayer][iSector][iPlane][iStrip].inverse();
276  }
277  }
278  }
279  }
280  }
281 }
282 
283 const HepGeom::Transform3D*
285 {
286  return &m_Section[section - 1];
287 }
288 
289 const HepGeom::Transform3D*
290 EKLM::TransformData::getLayerTransform(int section, int layer) const
291 {
292  return &m_Layer[section - 1][layer - 1];
293 }
294 
295 const HepGeom::Transform3D* EKLM::TransformData::
296 getSectorTransform(int section, int layer, int sector) const
297 {
298  return &m_Sector[section - 1][layer - 1][sector - 1];
299 }
300 
301 const HepGeom::Transform3D* EKLM::TransformData::
302 getPlaneTransform(int section, int layer, int sector, int plane) const
303 {
304  return &m_Plane[section - 1][layer - 1][sector - 1][plane - 1];
305 }
306 
307 const HepGeom::Transform3D* EKLM::TransformData::
308 getPlaneDisplacement(int section, int layer, int sector, int plane) const
309 {
310  return &m_PlaneDisplacement[section - 1][layer - 1][sector - 1][plane - 1];
311 }
312 
313 const HepGeom::Transform3D* EKLM::TransformData::
314 getSegmentTransform(int section, int layer, int sector, int plane,
315  int segment) const
316 {
317  return &m_Segment[section - 1][layer - 1][sector - 1][plane - 1][segment - 1];
318 }
319 
320 const HepGeom::Transform3D* EKLM::TransformData::
321 getStripTransform(int section, int layer, int sector, int plane, int strip) const
322 {
323  return &m_Strip[section - 1][layer - 1][sector - 1][plane - 1][strip - 1];
324 }
325 
326 const HepGeom::Transform3D*
328 {
329  return &(m_Strip[hit->getSection() - 1][hit->getLayer() - 1]
330  [hit->getSector() - 1][hit->getPlane() - 1][hit->getStrip() - 1]);
331 }
332 
333 const HepGeom::Transform3D*
335 {
336  return &(m_StripInverse[hit->getSection() - 1][hit->getLayer() - 1]
337  [hit->getSector() - 1][hit->getPlane() - 1][hit->getStrip() - 1]);
338 }
339 
340 const HepGeom::Transform3D*
341 EKLM::TransformData::getStripGlobalToLocal(int section, int layer, int sector,
342  int plane, int strip) const
343 {
344  return &(m_StripInverse[section - 1][layer - 1][sector - 1][plane - 1]
345  [strip - 1]);
346 }
347 
350  double* d1, double* d2, double* sd,
351  bool segments) const
352 {
353  /* Hits must be from the same sector, */
354  if (hit1->getSection() != hit2->getSection())
355  return false;
356  if (hit1->getLayer() != hit2->getLayer())
357  return false;
358  if (hit1->getSector() != hit2->getSector())
359  return false;
360  /* but different planes. */
361  if (hit1->getPlane() == hit2->getPlane())
362  return false;
363  /* Coordinates of strip 1 ends. */
364  double l1 = m_GeoDat->getStripLength(hit1->getStrip());
365  HepGeom::Point3D<double> s1_1(-0.5 * l1, 0.0, 0.0);
366  HepGeom::Point3D<double> s1_2(0.5 * l1, 0.0, 0.0);
367  const HepGeom::Transform3D* tr1 = getStripLocalToGlobal(hit1);
368  HepGeom::Point3D<double> s1_1g = (*tr1) * s1_1;
369  HepGeom::Point3D<double> s1_2g = (*tr1) * s1_2;
370  /* Coordinates of strip 2 ends. */
371  double l2 = m_GeoDat->getStripLength(hit2->getStrip());
372  HepGeom::Point3D<double> s2_1(-0.5 * l2, 0.0, 0.0);
373  HepGeom::Point3D<double> s2_2(0.5 * l2, 0.0, 0.0);
374  const HepGeom::Transform3D* tr2 = getStripLocalToGlobal(hit2);
375  HepGeom::Point3D<double> s2_1g = (*tr2) * s2_1;
376  HepGeom::Point3D<double> s2_2g = (*tr2) * s2_2;
386  HepGeom::Vector3D<double> v1 = s1_2g - s1_1g;
387  HepGeom::Vector3D<double> v2 = s2_2g - s2_1g;
388  HepGeom::Vector3D<double> d = s1_1g - s2_1g;
389  double v1sq = v1.mag2();
390  double v2sq = v2.mag2();
391  double v1dv2 = v1.dot(v2);
392  double ddv1 = d.dot(v1);
393  double ddv2 = d.dot(v2);
394  double den = v1sq * v2sq - v1dv2 * v1dv2;
395  double t1 = (v1dv2 * ddv2 - v2sq * ddv1) / den;
396  double t2 = (- v1dv2 * ddv1 + v1sq * ddv2) / den;
397  /* Segments do not intersect. */
398  if (segments) {
399  if (t1 < 0.0 || t1 > 1.0)
400  return false;
401  if (t2 < 0.0 || t2 > 1.0)
402  return false;
403  }
404  /* Segments intersect, set return values. */
405  HepGeom::Point3D<double> s1_cg = s1_1g + v1 * t1;
406  HepGeom::Point3D<double> s2_cg = s2_1g + v2 * t2;
407  *d1 = s1_2g.distance(s1_cg) / CLHEP::mm * Unit::mm;
408  *d2 = s2_2g.distance(s2_cg) / CLHEP::mm * Unit::mm;
409  *cross = 0.5 * (s1_cg + s2_cg) / CLHEP::mm * Unit::mm;
410  *sd = s1_cg.distance(s2_cg) / CLHEP::mm * Unit::mm;
411  if (s2_cg.mag2() < s1_cg.mag2())
412  *sd = - *sd;
413  return true;
414 }
415 
417  int section, const HepGeom::Point3D<double>& position) const
418 {
419  int sector;
420  double x;
421  if (section == 1) {
422  x = position.x();
423  } else {
424  x = -position.x();
425  }
426  if (position.y() > 0) {
427  if (x > 0)
428  sector = 1;
429  else
430  sector = 2;
431  } else {
432  if (x > 0)
433  sector = 4;
434  else
435  sector = 3;
436  }
437  return sector;
438 }
439 
441  const HepGeom::Point3D<double>& intersection, int* strip1, int* strip2) const
442 {
443  /* cppcheck-suppress variableScope */
444  int section, layer, sector, plane, segment, strip, stripSegment, stripGlobal;
445  /* cppcheck-suppress variableScope */
446  int nLayers, nPlanes, nSegments, nStripsSegment, minDistanceSegment;
447  double solenoidCenter, firstLayerCenter, layerShift;
448  /* cppcheck-suppress variableScope */
449  double x, y, z, l, minY, maxY;
450  double minDistance = 0, minDistanceNew, stripWidth;
451  HepGeom::Point3D<double> intersectionClhep, intersectionLocal;
452  intersectionClhep = intersection * CLHEP::cm / Unit::cm;
453  solenoidCenter = m_GeoDat->getSolenoidZ() / CLHEP::cm * Unit::cm;
454  if (intersection.z() < solenoidCenter)
455  section = 1;
456  else
457  section = 2;
458  firstLayerCenter =
459  (m_GeoDat->getSectionPosition()->getZ()
460  - 0.5 * m_GeoDat->getSectionPosition()->getLength()
461  + m_GeoDat->getLayerShiftZ()
462  - 0.5 * m_GeoDat->getLayerPosition()->getLength()) /
463  CLHEP::cm * Unit::cm;
464  layerShift = m_GeoDat->getLayerShiftZ() / CLHEP::cm * Unit::cm;
465  z = fabs(intersection.z() - solenoidCenter);
466  layer = round((z - firstLayerCenter) / layerShift) + 1;
467  if (layer <= 0)
468  layer = 1;
469  nLayers = m_GeoDat->getNDetectorLayers(section);
470  if (layer > nLayers)
471  layer = nLayers;
472  sector = getSectorByPosition(section, intersection);
473  nPlanes = m_GeoDat->getNPlanes();
474  nSegments = m_GeoDat->getNSegments();
475  nStripsSegment = m_ElementNumbers->getNStripsSegment();
476  stripWidth = m_GeoDat->getStripGeometry()->getWidth() / CLHEP::cm * Unit::cm;
477  minY = -stripWidth / 2;
478  maxY = (double(nStripsSegment) - 0.5) * stripWidth;
479  for (plane = 1; plane <= nPlanes; plane++) {
480  minDistanceSegment = 1;
481  for (segment = 1; segment <= nSegments; segment++) {
482  strip = (segment - 1) * nStripsSegment;
483  intersectionLocal = m_StripInverse[section - 1][layer - 1]
484  [sector - 1][plane - 1][strip] * intersectionClhep;
485  y = intersectionLocal.y() / CLHEP::cm * Unit::cm;
486  if (y < minY) {
487  minDistanceNew = minY - y;
488  } else if (y > maxY) {
489  minDistanceNew = y - maxY;
490  } else {
491  minDistance = 0;
492  minDistanceSegment = segment;
493  break;
494  }
495  if (segment == 1) {
496  minDistance = minDistanceNew;
497  } else if (minDistanceNew < minDistance) {
498  minDistance = minDistanceNew;
499  minDistanceSegment = segment;
500  }
501  }
502  /*
503  * The intersection is required to be strictly within a segment,
504  * this condition might be adjusted later.
505  */
506  if (minDistance > 0)
507  return -1;
508  strip = (minDistanceSegment - 1) * nStripsSegment;
509  intersectionLocal = m_StripInverse[section - 1][layer - 1]
510  [sector - 1][plane - 1][strip] * intersectionClhep;
511  y = intersectionLocal.y() / CLHEP::cm * Unit::cm;
512  stripSegment = ceil((y - 0.5 * stripWidth) / stripWidth);
513  if (stripSegment <= 0)
514  stripSegment = 1;
515  else if (stripSegment > nStripsSegment)
516  stripSegment = nStripsSegment;
517  strip = stripSegment + (minDistanceSegment - 1) * nStripsSegment;
518  intersectionLocal = m_StripInverse[section - 1][layer - 1]
519  [sector - 1][plane - 1][strip - 1] * intersectionClhep;
520  x = intersectionLocal.x();
521  l = m_GeoDat->getStripLength(strip);
522  /*
523  * The intersection is required to be strictly within the strip length,
524  * this condition might be adjusted later.
525  */
526  if (fabs(x) > 0.5 * l)
527  return -1;
528  stripGlobal = m_ElementNumbers->stripNumber(
529  section, layer, sector, plane, strip);
530  if (plane == 1)
531  *strip1 = stripGlobal;
532  else
533  *strip2 = stripGlobal;
534  }
535  return 0;
536 }
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
static const EKLMElementNumbers & Instance()
Instantiation.
int sectorNumber(int section, int layer, int sector) const
Get sector number.
static constexpr int getNStripsSegment()
Get number of strips in a segment.
int segmentNumber(int section, int layer, int sector, int plane, int segment) const
Get segment number.
int getNPlanes() const
Get number of planes.
int getNSections() const
Get number of sections.
int getNLayers() const
Get number of layers.
int getNDetectorLayers(int section) const
Get number of detector layers.
int getNSegments() const
Get number of segments.
int getNStrips() const
Get number of strips.
int getNSectors() const
Get number of sectors.
Class for EKLM alignment checking.
bool checkAlignment(const EKLMAlignment *alignment, const EKLMSegmentAlignment *segmentAlignment) const
Check alignment.
void getSectorTransform(HepGeom::Transform3D *t, int n) const
Get sector transformation.
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Definition: GeometryData.cc:33
void getStripTransform(HepGeom::Transform3D *t, int n) const
Get strip transformation.
void getSectionTransform(HepGeom::Transform3D *t, int n) const
Get section transformation.
void getLayerTransform(HepGeom::Transform3D *t, int n) const
Get layer transformation.
void getPlaneTransform(HepGeom::Transform3D *t, int n) const
Get plane transformation.
const HepGeom::Transform3D * getStripTransform(int section, int layer, int sector, int plane, int strip) const
Get strip transformation.
Displacement
Source of displacement (alignment) data.
Definition: TransformData.h:42
@ c_None
Displacement is not used.
Definition: TransformData.h:43
@ c_Displacement
Use displacement data (for geometry).
Definition: TransformData.h:44
TransformData(bool global, Displacement displacementType)
Constructor.
HepGeom::Transform3D * m_Section
Section transformations.
const HepGeom::Transform3D * getPlaneDisplacement(int section, int layer, int sector, int plane) const
Get additional displacement for plane internal volumes.
HepGeom::Transform3D ***** m_Strip
Strip transformations.
const HepGeom::Transform3D * getPlaneTransform(int section, int layer, int sector, int plane) const
Get plane transformation.
HepGeom::Transform3D ** m_Layer
Layer transformations.
HepGeom::Transform3D *** m_Sector
Sector transformations.
void transformsToGlobal()
Make transformations global from local.
const HepGeom::Transform3D * getSegmentTransform(int section, int layer, int sector, int plane, int segment) const
Get segment transformation.
int getStripsByIntersection(const HepGeom::Point3D< double > &intersection, int *strip1, int *strip2) const
Find strips by intersection.
const HepGeom::Transform3D * getLayerTransform(int section, int layer) const
Get layer transformation.
HepGeom::Transform3D **** m_Plane
Plane transformations.
HepGeom::Transform3D ***** m_StripInverse
Inverse strip transformations.
int getSectorByPosition(int section, const HepGeom::Point3D< double > &position) const
Get sector by position.
const EKLMElementNumbers * m_ElementNumbers
Element numbers.
const HepGeom::Transform3D * getStripLocalToGlobal(KLMDigit *hit) const
Get strip local to global transformation by hit.
const HepGeom::Transform3D * getStripGlobalToLocal(KLMDigit *hit) const
Get strip global to local transformation by hit.
const HepGeom::Transform3D * getSectorTransform(int section, int layer, int sector) const
Get sector transformation.
HepGeom::Transform3D **** m_PlaneDisplacement
Plane internal volumes displacements.
HepGeom::Transform3D ***** m_Segment
Segment transformations.
const GeometryData * m_GeoDat
Geometry data.
bool intersection(KLMDigit *hit1, KLMDigit *hit2, HepGeom::Point3D< double > *cross, double *d1, double *d2, double *sd, bool segments=true) const
Check if strips intersect, and find intersection point if yes.
const HepGeom::Transform3D * getSectionTransform(int section) const
Get section transformation.
KLM Alignment data.
float getDeltaU() const
Get shift in U.
float getDeltaV() const
Get shift in V.
float getDeltaGamma() const
Get rotation in alpha.
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:29
int getLayer() const
Get layer number.
Definition: KLMDigit.h:126
int getSection() const
Get section number.
Definition: KLMDigit.h:90
int getPlane() const
Get plane number.
Definition: KLMDigit.h:144
int getStrip() const
Get strip number.
Definition: KLMDigit.h:162
int getSector() const
Get sector number.
Definition: KLMDigit.h:108
static const double mm
[millimeters]
Definition: Unit.h:70
static const double rad
Standard of [angle].
Definition: Unit.h:50
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47
Abstract base class for different kinds of events.