9 #include <alignment/reconstruction/AlignableCDCRecoHit.h>
11 #include <alignment/GlobalLabel.h>
12 #include <alignment/GlobalDerivatives.h>
13 #include <cdc/dataobjects/WireID.h>
14 #include <cdc/dbobjects/CDCTimeZeros.h>
15 #include <cdc/dbobjects/CDCTimeWalks.h>
16 #include <cdc/geometry/CDCGeometryPar.h>
17 #include "Math/ChebyshevPol.h"
22 using namespace alignment;
24 bool AlignableCDCRecoHit::s_enableTrackT0LocalDerivative =
true;
25 bool AlignableCDCRecoHit::s_enableWireSaggingGlobalDerivative =
false;
26 bool AlignableCDCRecoHit::s_enableWireByWireAlignmentGlobalDerivatives =
false;
28 std::pair<std::vector<int>, TMatrixD> AlignableCDCRecoHit::globalDerivatives(
const genfit::StateOnPlane* sop)
30 GlobalDerivatives globals;
31 unsigned short LR = (int(m_leftRight) > 0.) ? 1 : 0;
34 const B2Vector3D& wirePositon = sop->getPlane()->getO();
35 const unsigned short layer = getWireID().getICLayer();
38 double alpha = cdcgeo.
getAlpha(wirePositon, mom);
40 const TVectorD& stateOnPlane = sop->getState();
41 const double driftLengthEstimate = std::abs(stateOnPlane(3));
42 const double driftTime = cdcgeo.
getDriftTime(driftLengthEstimate, layer, LR, alpha, theta);
43 const double driftVelocity = cdcgeo.
getDriftV(driftTime, layer, LR, alpha, theta);
48 if (driftTime > 20 && driftTime < 200 && fabs(driftVelocity) < 1.0e-2) {
50 GlobalLabel::construct<CDCTimeZeros>(getWireID(), 0),
51 driftVelocity *
double(
int(m_leftRight))
56 if (driftTime > 20 && driftTime < 200 && fabs(driftVelocity) < 1.0e-2 && m_adcCount < 400 && m_adcCount > 2) {
58 unsigned short board = CDCGeometryPar::Instance().getBoardID(getWireID());
59 const std::vector<float> param = tws->getTimeWalkParams(board);
61 GlobalLabel::construct<CDCTimeWalks>(board, 0),
62 -driftVelocity * exp(-param[1]* m_adcCount) *
double(
int(m_leftRight))
65 GlobalLabel::construct<CDCTimeWalks>(board, 1),
66 driftVelocity * m_adcCount * param[0]*exp(-param[1]* m_adcCount) *
double(
int(m_leftRight))
71 if (driftTime > 20 && driftTime < 500 && fabs(driftVelocity) < 1.0e-2) {
78 auto xtid = xts->getXtID(getWireID().getICLayer(), LR, (
float)alpha, (
float)theta);
79 const auto& par = xts->getXtParams(xtid);
80 auto boundary = par.at(6);
81 if (driftTime < boundary) {
83 GlobalLabel::construct<CDCXtRelations>(xtid, 0),
84 ROOT::Math::Chebyshev5(driftTime, 1, 0, 0, 0, 0, 0) *
double(
int(m_leftRight))
87 GlobalLabel::construct<CDCXtRelations>(xtid, 1),
88 ROOT::Math::Chebyshev5(driftTime, 0, 1, 0, 0, 0, 0) *
double(
int(m_leftRight))
91 GlobalLabel::construct<CDCXtRelations>(xtid, 2),
92 ROOT::Math::Chebyshev5(driftTime, 0, 0, 1, 0, 0, 0) *
double(
int(m_leftRight))
95 GlobalLabel::construct<CDCXtRelations>(xtid, 3),
96 ROOT::Math::Chebyshev5(driftTime, 0, 0, 0, 1, 0, 0) *
double(
int(m_leftRight))
99 GlobalLabel::construct<CDCXtRelations>(xtid, 4),
100 ROOT::Math::Chebyshev5(driftTime, 0, 0, 0, 0, 1, 0) *
double(
int(m_leftRight))
103 GlobalLabel::construct<CDCXtRelations>(xtid, 5),
104 ROOT::Math::Chebyshev5(driftTime, 0, 0, 0, 0, 0, 1) *
double(
int(m_leftRight))
108 GlobalLabel::construct<CDCXtRelations>(xtid, 7),
109 double(
int(m_leftRight))
112 GlobalLabel::construct<CDCXtRelations>(xtid, 8),
113 (driftTime - boundary) *
double(
int(m_leftRight))
121 auto tdir = sop->getDir();
122 auto ndir = sop->getPlane()->getNormal();
123 auto udir = sop->getPlane()->getU();
124 auto vdir = sop->getPlane()->getV();
125 auto pos = sop->getPos();
127 auto tn = tdir[0] * ndir[0] + tdir[1] * ndir[1] + tdir[2] * ndir [2];
129 auto drdm = TMatrixD(3, 3);
131 for (
int i = 0; i < 3; ++i) {
132 for (
int j = 0; j < 3; ++j) {
133 drdm(i, j) -= tdir[i] * ndir[j] / tn;
137 auto dmdg = TMatrixD(3, 6);
139 dmdg[0][0] = 1.; dmdg[0][4] = -pos[2]; dmdg[0][5] = pos[1];
140 dmdg[1][1] = 1.; dmdg[1][3] = pos[2]; dmdg[1][5] = -pos[0];
141 dmdg[2][2] = 1.; dmdg[2][3] = -pos[1]; dmdg[2][4] = pos[0];
143 auto drldrg = TMatrixD(3, 3);
144 for (
int i = 0; i < 3; ++i) {
145 drldrg(0, i) = udir[i];
146 drldrg(1, i) = vdir[i];
147 drldrg(2, i) = ndir[i];
150 auto drldg = drldrg * (drdm * dmdg);
153 const double zWireM = s_cdcGeometryTranslator->getWireBackwardPosition(getWireID(), CDCGeometryPar::c_Aligned)[2];
154 const double zWireP = s_cdcGeometryTranslator->getWireForwardPosition(getWireID(), CDCGeometryPar::c_Aligned)[2];
156 const double zRel = std::max(0., std::min(1., (pos[2] - zWireM) / (zWireP - zWireM)));
160 auto layerID =
WireID(getWireID().getICLayer(), 511);
164 GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerX),
170 GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerY),
176 GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerPhi),
183 GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerDx),
189 GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerDy),
195 GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerDPhi),
201 if (s_enableWireByWireAlignmentGlobalDerivatives) {
204 double zRelM = fabs(1. - zRel);
206 double zRelP = fabs(zRel - 0.);
210 GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireBwdX),
215 GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireFwdX),
221 GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireBwdY),
226 GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireFwdY),
235 if (s_enableWireSaggingGlobalDerivative) {
237 GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireTension),
238 drldg(0, 1) * 4.0 * zRel * (1.0 - zRel)
250 if (!s_enableTrackT0LocalDerivative)
253 unsigned short LR = (int(m_leftRight) > 0.) ? 1 : 0;
256 const B2Vector3D& wirePositon = sop->getPlane()->getO();
257 const unsigned short layer = getWireID().getICLayer();
260 const double alpha = cdcgeo.
getAlpha(wirePositon, mom);
261 const double theta = cdcgeo.
getTheta(mom);
262 const TVectorD& stateOnPlane = sop->getState();
263 const double driftLengthEstimate = std::abs(stateOnPlane(3));
264 const double driftTime = cdcgeo.
getDriftTime(driftLengthEstimate, layer, LR, alpha, theta);
265 const double driftVelocity = cdcgeo.
getDriftV(driftTime, layer, LR, alpha, theta);
267 TMatrixD locals(2, 1);
268 if (driftTime > 20 && driftTime < 200 && fabs(driftVelocity) < 1.0e-2) {
269 locals(0, 0) = - double(
int(m_leftRight)) * driftVelocity;
The Class for CDC Geometry Parameters.
double getTheta(const B2Vector3D &momentum) const
Returns track incident angle (theta in rad.).
double getAlpha(const B2Vector3D &posOnWire, const B2Vector3D &momentum) const
Returns track incident angle in rphi plane (alpha in rad.).
double getDriftV(double dt, unsigned short layer, unsigned short lr, double alpha=0., double theta=0.5 *M_PI) const
Get the realistic drift velocity.
double getOutgoingAlpha(const double alpha) const
Converts incoming- to outgoing-alpha.
double getOutgoingTheta(const double alpha, const double theta) const
Converts incoming- to outgoing-theta.
double getDriftTime(double dist, unsigned short layer, unsigned short lr, double alpha, double theta) const
Return the drift time to the sense wire.
Class for accessing objects in the database.
Class to identify a wire inside the CDC.
A state with arbitrary dimension defined in a DetPlane.
Abstract base class for different kinds of events.