Belle II Software  release-08-01-10
AlignableEKLMRecoHit.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 #include <alignment/reconstruction/AlignableEKLMRecoHit.h>
10 
11 #include <alignment/GlobalDerivatives.h>
12 #include <alignment/GlobalLabel.h>
13 #include <alignment/reconstruction/AlignableEKLMRecoHit.h>
14 #include <klm/dataobjects/eklm/EKLMElementNumbers.h>
15 #include <klm/dataobjects/KLMDigit.h>
16 #include <klm/dataobjects/KLMElementNumbers.h>
17 #include <klm/dataobjects/KLMHit2d.h>
18 #include <klm/dbobjects/eklm/EKLMAlignment.h>
19 #include <klm/eklm/geometry/GeometryData.h>
20 #include <klm/eklm/geometry/TransformDataGlobalAligned.h>
21 
22 using namespace Belle2;
23 using namespace alignment;
24 
26 {
27 }
28 
30  const EKLMAlignmentHit* hit, const genfit::TrackCandHit* trackCandHit) :
31  genfit::PlanarMeasurement(1)
32 {
33  (void)trackCandHit;
34  int digit, plane, segment, strip;
35  const HepGeom::Transform3D* t;
36  CLHEP::HepRotation r;
37  CLHEP::Hep3Vector origin;
38  CLHEP::Hep3Vector u(1, 0, 0);
39  CLHEP::Hep3Vector v(0, 1, 0);
40  B2Vector3D origin2, u2, v2;
42  const KLMElementNumbers* elementNumbers = &(KLMElementNumbers::Instance());
43  const EKLMElementNumbers* eklmElementNumbers =
45  const EKLM::TransformDataGlobalAligned* transformData =
47  RelationVector<KLMHit2d> hit2ds = hit->getRelationsTo<KLMHit2d>();
48  if (hit2ds.size() != 1)
49  B2FATAL("Incorrect number of related KLMHit2ds.");
50  RelationVector<KLMDigit> eklmDigits = hit2ds[0]->getRelationsTo<KLMDigit>();
51  if (eklmDigits.size() != 2)
52  B2FATAL("Incorrect number of related KLMDigits.");
53  digit = hit->getDigitIdentifier();
54  m_Section = eklmDigits[digit]->getSection();
55  m_Layer = eklmDigits[digit]->getLayer();
56  m_Sector = eklmDigits[digit]->getSector();
57  plane = eklmDigits[digit]->getPlane();
58  segment = (eklmDigits[digit]->getStrip() - 1) /
59  eklmElementNumbers->getNStripsSegment() + 1;
60  strip = (segment - 1) * eklmElementNumbers->getNStripsSegment() + 1;
62  m_Segment = eklmElementNumbers->segmentNumber(
63  m_Section, m_Layer, m_Sector, plane, segment);
64  t = transformData->getStripTransform(
65  m_Section, m_Layer, m_Sector, plane, strip);
66  origin = t->getTranslation();
67  origin2.SetX(origin.x() / CLHEP::cm * Unit::cm);
68  origin2.SetY(origin.y() / CLHEP::cm * Unit::cm);
69  origin2.SetZ(origin.z() / CLHEP::cm * Unit::cm);
70  r = t->getRotation();
71  u = r * u;
72  v = r * v;
73  u2.SetX(u.x());
74  u2.SetY(u.y());
75  u2.SetZ(u.z());
76  v2.SetX(v.x());
77  v2.SetY(v.y());
78  v2.SetZ(v.z());
79  t = transformData->getSectorTransform(m_Section, m_Layer, m_Sector);
80  r = t->getRotation().inverse();
81  v = r * v;
82  m_StripV.SetX(v.unit().x());
83  m_StripV.SetY(v.unit().y());
84  genfit::SharedPlanePtr detPlane(new genfit::DetPlane(origin2, u2, v2, 0));
85  setPlane(detPlane, m_Segment);
86  rawHitCoords_[0] = geoDat->getStripGeometry()->getWidth() *
87  ((eklmDigits[digit]->getStrip() - 1) %
88  eklmElementNumbers->getNStripsSegment()) / CLHEP::cm * Unit::cm;
89  rawHitCov_[0][0] = pow(geoDat->getStripGeometry()->getWidth() /
90  CLHEP::cm * Unit::cm, 2) / 12;
91  setStripV();
92 }
93 
95 {
96 }
97 
98 std::pair<std::vector<int>, TMatrixD> AlignableEKLMRecoHit::globalDerivatives(const genfit::StateOnPlane* sop)
99 {
100  std::vector<int> labGlobal;
101  labGlobal.push_back(GlobalLabel::construct<EKLMAlignment>(m_KLMModule, KLMAlignmentData::c_DeltaU)); // dx
102  labGlobal.push_back(GlobalLabel::construct<EKLMAlignment>(m_KLMModule, KLMAlignmentData::c_DeltaV)); // dy
103  labGlobal.push_back(GlobalLabel::construct<EKLMAlignment>(m_KLMModule, KLMAlignmentData::c_DeltaGamma)); // drot
104 
105  /* Local parameters. */
106  const double dalpha = 0;
107  const double dxs = 0;
108  const double dys = 0;
109  const double sinda = sin(dalpha);
110  const double cosda = cos(dalpha);
111  /* Local position in sector coordinates. */
112  HepGeom::Point3D<double> globalPos;
113  HepGeom::Transform3D t;
114  const EKLM::TransformDataGlobalAligned* transformData =
116  t = (*transformData->getSectorTransform(m_Section, m_Layer, m_Sector)).inverse();
117  globalPos.setX(sop->getPos().X() / Unit::cm * CLHEP::cm);
118  globalPos.setY(sop->getPos().Y() / Unit::cm * CLHEP::cm);
119  globalPos.setZ(sop->getPos().Z() / Unit::cm * CLHEP::cm);
120  globalPos = t * globalPos;
121  double x = globalPos.x() / CLHEP::cm * Unit::cm;
122  double y = globalPos.y() / CLHEP::cm * Unit::cm;
123  /*
124  * Matrix of global derivatives of local coordinate (second dimension is
125  * added because of resizing in GblFitterInfo::constructGblPoint()) by
126  * variations of the detector element position.
127  * The sign is inverted to match the convention for other detectors.
128  */
129  TMatrixD derGlobal(2, 3);
130  derGlobal(0, 0) = 0;
131  derGlobal(0, 1) = 0;
132  derGlobal(0, 2) = 0;
133  derGlobal(1, 0) = cosda * m_StripV.X() - sinda * m_StripV.Y();
134  derGlobal(1, 1) = sinda * m_StripV.X() + cosda * m_StripV.Y();
135  derGlobal(1, 2) =
136  -((x - dxs) * (-sinda * m_StripV.X() - cosda * m_StripV.Y()) +
137  (y - dys) * (cosda * m_StripV.X() - sinda * m_StripV.Y()));
138  return alignment::GlobalDerivatives(labGlobal, derGlobal);
139 
140 }
141 
143 {
144  return new AlignableEKLMRecoHit(*this);
145 }
146 
uint16_t m_Segment
Segment number.
virtual std::pair< std::vector< int >, TMatrixD > globalDerivatives(const genfit::StateOnPlane *sop) override
Labels and derivatives of residuals (local measurement coordinates) w.r.t.
genfit::AbsMeasurement * clone() const override
Clone.
B2Vector3D m_StripV
V direction.
uint16_t m_KLMModule
KLM module number.
void SetX(DataType x)
set X/1st-coordinate
Definition: B2Vector3.h:457
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:431
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:433
void SetZ(DataType z)
set Z/3rd-coordinate
Definition: B2Vector3.h:461
void SetY(DataType y)
set Y/2nd-coordinate
Definition: B2Vector3.h:459
This dataobject is used only for EKLM alignment.
EKLM element numbers.
static const EKLMElementNumbers & Instance()
Instantiation.
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.
double getWidth() const
Get width.
const StripGeometry * getStripGeometry() const
Get strip geometry data.
EKLM geometry data.
Definition: GeometryData.h:38
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Definition: GeometryData.cc:33
Transformation data (global, aligned): singleton version.
static const TransformDataGlobalAligned & Instance()
Instantiation.
const HepGeom::Transform3D * getStripTransform(int section, int layer, int sector, int plane, int strip) const
Get strip transformation.
const HepGeom::Transform3D * getSectorTransform(int section, int layer, int sector) const
Get sector transformation.
@ c_DeltaU
Shift in U (EKLM: local X).
@ c_DeltaGamma
Rotation in gamma (EKLM: rotation in local plane).
@ c_DeltaV
Shift in V (EKLM: local Y).
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:29
KLM element numbers.
static const KLMElementNumbers & Instance()
Instantiation.
KLMModuleNumber moduleNumberEKLM(int section, int sector, int layer) const
Get module number for EKLM.
KLM 2d hit.
Definition: KLMHit2d.h:33
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47
Class for easier manipulation with global derivatives (and their labels)
Contains the measurement and covariance in raw detector coordinates.
Detector plane.
Definition: DetPlane.h:59
void setStripV(bool v=true)
Use if the coordinate for 1D hits measured in V direction.
A state with arbitrary dimension defined in a DetPlane.
Definition: StateOnPlane.h:47
Hit object for use in TrackCand.
Definition: TrackCandHit.h:34
Abstract base class for different kinds of events.
Defines for I/O streams used for error and debug printing.
std::shared_ptr< genfit::DetPlane > SharedPlanePtr
Shared Pointer to a DetPlane.