Belle II Software  release-05-01-25
AlignableCDCRecoHit Class Reference

This class is used to transfer CDC information to the track fit and Millepede. More...

#include <AlignableCDCRecoHit.h>

Inheritance diagram for AlignableCDCRecoHit:
Collaboration diagram for AlignableCDCRecoHit:

Public Member Functions

 ~AlignableCDCRecoHit ()
 Destructor.
 
AlignableCDCRecoHitclone () const override
 Creating a copy of this hit.
 
virtual std::pair< std::vector< int >, TMatrixD > globalDerivatives (const genfit::StateOnPlane *sop) override
 Labels and derivatives of residuals (local measurement coordinates) w.r.t. More...
 
virtual TMatrixD localDerivatives (const genfit::StateOnPlane *sop) override
 Derivatives for (local) fit parameters. More...
 
 CDCRecoHit ()
 Inherit constructors.
 
 CDCRecoHit (const CDCHit *cdcHit, const genfit::TrackCandHit *trackCandHit)
 Inherit constructors.
 
WireID getWireID () const
 Getter for WireID object.
 
genfit::SharedPlanePtr constructPlane (const genfit::StateOnPlane &state) const override
 Methods that actually interface to Genfit.
 
std::vector< genfit::MeasurementOnPlane * > constructMeasurementsOnPlane (const genfit::StateOnPlane &state) const override
 build MeasurementsOnPlane
 
virtual const genfit::HMatrixUconstructHMatrix (const genfit::AbsTrackRep *) const override
 construct error matrix
 
std::vector< double > timeDerivativesMeasurementsOnPlane (const genfit::StateOnPlane &state) const
 Get the time derivative of the MesuredStateOnPlane (derived from the track fit). More...
 
bool getFlyByDistanceVector (B2Vector3D &pointingVector, B2Vector3D &trackDir, const genfit::AbsTrackRep *rep=NULL, bool usePlaneFromFit=false)
 Get the vector pointing from the wire to the fitted trajectory as well as the direction of the track in the fitted point. More...
 
void setLeftRightResolution (int lr)
 select how to resolve the left/right ambiguity: -1: negative (left) side on vector (wire direction) x (track direction) 0: mirrors enter with same weight, DAF will decide. More...
 
bool isLeftRightMeasurement () const override
 CDC RecoHits always have left-right ambiguity.
 
int getLeftRightResolution () const override
 Getter for left/right passage flag.
 
const CDCHitgetCDCHit () const
 get the pointer to the CDCHit object that was used to create this CDCRecoHit object. More...
 
TrackPoint * getTrackPoint () const
 
void setTrackPoint (TrackPoint *tp)
 
const TVectorD & getRawHitCoords () const
 
TVectorD & getRawHitCoords ()
 
const TMatrixDSym & getRawHitCov () const
 
TMatrixDSym & getRawHitCov ()
 
int getDetId () const
 
int getHitId () const
 
unsigned int getDim () const
 
void setRawHitCoords (const TVectorD &coords)
 
void setRawHitCov (const TMatrixDSym &cov)
 
void setDetId (int detId)
 
void setHitId (int hitId)
 
virtual void Print (const Option_t *="") const
 
virtual std::vector< int > labels ()
 Vector of integer labels for calibration/alignment parameters available (must match #columns of derivatives(...)) More...
 
virtual TMatrixD derivatives (const genfit::StateOnPlane *)
 Derivatives of residuals (local measurement coordinates) w.r.t. More...
 
virtual std::vector< int > localLabels ()
 Vector of integer labels for local calibration parameters available (must match #columns of localDerivatives(...)) More...
 

Static Public Member Functions

static void setTranslators (CDC::ADCCountTranslatorBase *const adcCountTranslator, CDC::CDCGeometryTranslatorBase *const cdcGeometryTranslator, CDC::TDCCountTranslatorBase *const tdcCountTranslator, bool useTrackTime=false, bool cosmics=false)
 Setter for the Translators.
 

Static Public Attributes

static bool s_enableTrackT0LocalDerivative = true
 Static enabling(true) or disabling(false) addition of local derivative for track T0.
 
static bool s_enableWireSaggingGlobalDerivative = false
 Static enabling(true) or disabling(false) addition of global derivative for wire sagging coefficient (per wire)
 
static bool s_enableWireByWireAlignmentGlobalDerivatives = false
 Static enabling(true) or disabling(false) addition of global derivatives for wire-by-wire alignment.
 

Protected Member Functions

 ClassDefOverride (CDCRecoHit, 10)
 ROOT Macro.
 

Protected Attributes

unsigned short m_tdcCount
 TDC Count as out of CDCHit.
 
unsigned short m_adcCount
 ADC Count as out of CDCHit.
 
WireID m_wireID
 Wire Identifier.
 
const CDCHitm_cdcHit
 Pointer to the CDCHit used to created this CDCRecoHit.
 
signed char m_leftRight
 Flag showing left/right passage.
 
TVectorD rawHitCoords_
 
TMatrixDSym rawHitCov_
 
int detId_
 
int hitId_
 
TrackPoint * trackPoint_
 Pointer to TrackPoint where the measurement belongs to.
 

Static Protected Attributes

static std::unique_ptr< CDC::ADCCountTranslatorBases_adcCountTranslator = 0
 Object for ADC Count translation.
 
static std::unique_ptr< CDC::CDCGeometryTranslatorBases_cdcGeometryTranslator = 0
 Object for geometry translation.
 
static std::unique_ptr< CDC::TDCCountTranslatorBases_tdcCountTranslator = 0
 Object for getting drift-length and -resolution.
 
static bool s_useTrackTime = false
 Whether to use the track time or not when building the measurementOnPlane. More...
 
static bool s_cosmics = false
 Switch to use cosmic events, or physics events from IP. More...
 

Private Member Functions

 ClassDefOverride (AlignableCDCRecoHit, 1)
 ROOT Macro.
 

Detailed Description

This class is used to transfer CDC information to the track fit and Millepede.

Definition at line 32 of file AlignableCDCRecoHit.h.

Member Function Documentation

◆ derivatives()

virtual TMatrixD derivatives ( const genfit::StateOnPlane )
inlinevirtualinherited

Derivatives of residuals (local measurement coordinates) w.r.t.

alignment/calibration parameters Matrix "G" of derivatives valid for given prediction of track state:

G(i, j) = d_residual_i/d_parameter_j

For 2D measurement (u,v):

G = ( du/da du/db du/dc ... ) ( dv/da dv/db dv/dc ... )

for calibration parameters a, b, c.

For 1D measurement both forms are allowed:

G = ( 0 0 0 ... ) ( dv/da dv/db dv/dc ... ) for V-strip,

G = ( du/da du/db du/dc ... ) ( 0 0 0 ... ) for U-strip,

or :

G = ( d_sensitive/da d_sensitive/db d_sensitive/dc ... ) as matrix with one row.

A possible algorithm using these derivatives should be able to resolve this based on the measurement HMatrix. Measurements with more dimesions (slopes, curvature) should provide full 4-5Dx(n params) matrix (state as (q/p, u', v', u, v) or (u', v', u, v))

Parameters
sopPredicted state of the track as linearization point around which derivatives of alignment/calibration parameters shall be computed
Returns
TMatrixD Matrix with #rows = dimension of residual, #columns = number of parameters. #columns must match labels().size().

Definition at line 129 of file ICalibrationParametersDerivatives.h.

◆ getCDCHit()

const CDCHit* getCDCHit ( ) const
inlineinherited

get the pointer to the CDCHit object that was used to create this CDCRecoHit object.

Can be NULL if CDCRecoHit was not created with the CDCRecoHit(const CDCHit* cdcHit) constructor

Definition at line 123 of file CDCRecoHit.h.

◆ getFlyByDistanceVector()

bool getFlyByDistanceVector ( B2Vector3D pointingVector,
B2Vector3D trackDir,
const genfit::AbsTrackRep rep = NULL,
bool  usePlaneFromFit = false 
)
inherited

Get the vector pointing from the wire to the fitted trajectory as well as the direction of the track in the fitted point.

Uses the cardinal TrackRep if rep == NULL. Returns false if the track is not fitted in this point for the requested TrackRep. If the hit does not belong to a track, an error is issued and false returned. If usePlaneFromFit is used, constructPlane is not called to evaluate the closest point on the wire to the track. Instead the origin of the plane of the fitted state is used (which will be the same point if the wire is not bent).

Definition at line 370 of file CDCRecoHit.cc.

374 {
375  const genfit::TrackPoint* tp = this->getTrackPoint();
376  if (!tp) {
377  B2ERROR("No genfit::TrackPoint for CDCRecoHit.");
378  return false;
379  }
380  const genfit::AbsFitterInfo* fi = tp->getFitterInfo(rep);
381  if (!fi) {
382  B2DEBUG(200, "No genfit::AbsFitterInfo for this CDCRecoHit.");
383  return false;
384  }
385 
386  const genfit::MeasuredStateOnPlane& mop = fi->getFittedState();
387  B2Vector3D fittedPoca = mop.getPos();
388  // constructPlane places the coordinate center in the POCA to the
389  // wire. Using this is the default behavior. If this should be too
390  // slow, as it has to re-evaluate the POCA, alternatively the user
391  // can set usePlaneFromFit which uses the plane determined by the
392  // track fit.
393  B2Vector3D pocaOnWire = (usePlaneFromFit
394  ? mop.getPlane()->getO()
395  : this->constructPlane(mop)->getO());
396 
397  // The vector from the wire to the track.
398  pointingVector = fittedPoca - pocaOnWire;
399 
400  trackDir = mop.getMom();
401  trackDir.SetMag(1.);
402  return true;
403 }

◆ globalDerivatives()

std::pair< std::vector< int >, TMatrixD > globalDerivatives ( const genfit::StateOnPlane sop)
overridevirtual

Labels and derivatives of residuals (local measurement coordinates) w.r.t.

alignment/calibration parameters Matrix "G" of derivatives valid for given prediction of track state:

G(i, j) = d_residual_i/d_parameter_j

For 2D measurement (u,v):

G = ( du/da du/db du/dc ... ) ( dv/da dv/db dv/dc ... )

for calibration parameters a, b, c.

For 1D measurement:

G = ( 0 0 0 ... ) ( dv/da dv/db dv/dc ... ) for V-strip,

G = ( du/da du/db du/dc ... ) ( 0 0 0 ... ) for U-strip,

Measurements with more dimesions (slopes, curvature) should provide full 4-5Dx(n params) matrix (state as (q/p, u', v', u, v) or (u', v', u, v))

Parameters
sopPredicted state of the track as linearization point around which derivatives of alignment/calibration parameters shall be computed
Returns
pair<vector<int>, TMatrixD> With matrix with #rows = dimension of residual, #columns = number of parameters. #columns must match vector<int>.size().

Reimplemented from ICalibrationParametersDerivatives.

Definition at line 29 of file AlignableCDCRecoHit.cc.

30 {
31  GlobalDerivatives globals;
32  unsigned short LR = (int(m_leftRight) > 0.) ? 1 : 0;
33 
34  const TVector3& mom = sop->getMom();
35  const TVector3& wirePositon = sop->getPlane()->getO();
36  const unsigned short layer = getWireID().getICLayer();
37 
39  const double alpha = cdcgeo.getAlpha(wirePositon, mom);
40  const double theta = cdcgeo.getTheta(mom);
41  const TVectorD& stateOnPlane = sop->getState();
42  const double driftLengthEstimate = std::abs(stateOnPlane(3));
43  const double driftTime = cdcgeo.getDriftTime(driftLengthEstimate, layer, LR, alpha, theta);
44  const double driftVelocity = cdcgeo.getDriftV(driftTime, layer, LR, alpha, theta);
45 
46  // CDC Calibration -------------------------------------------------
47 
48  // T0 calibration (per wire) TODO check sign!!!
49  if (driftTime > 20 && driftTime < 200 && fabs(driftVelocity) < 1.0e-2) {
50  globals.add(
51  GlobalLabel::construct<CDCTimeZeros>(getWireID(), 0),
52  driftVelocity * double(int(m_leftRight))
53  );
54  }
55 
56  // Time-walk calibration (per board) TODO checksign!!!
57  globals.add(
58  GlobalLabel::construct<CDCTimeWalks>(CDCGeometryPar::Instance().getBoardID(getWireID()), 0),
59  driftVelocity * 1. / sqrt(m_adcCount) * double(int(m_leftRight))
60  );
61 
62  // CDC Alignment ---------------------------------------------------
63 
64  auto tdir = sop->getDir();
65  auto ndir = sop->getPlane()->getNormal();
66  auto udir = sop->getPlane()->getU();
67  auto vdir = sop->getPlane()->getV();
68  auto pos = sop->getPos();
69 
70  auto tn = tdir[0] * ndir[0] + tdir[1] * ndir[1] + tdir[2] * ndir [2];
71  // d residual / d measurement
72  auto drdm = TMatrixD(3, 3);
73  drdm.UnitMatrix();
74  for (int i = 0; i < 3; ++i) {
75  for (int j = 0; j < 3; ++j) {
76  drdm(i, j) -= tdir[i] * ndir[j] / tn;
77  }
78  }
79  // d measurement / d global rigid body param
80  auto dmdg = TMatrixD(3, 6);
81  dmdg.Zero();
82  dmdg[0][0] = 1.; dmdg[0][4] = -pos[2]; dmdg[0][5] = pos[1];
83  dmdg[1][1] = 1.; dmdg[1][3] = pos[2]; dmdg[1][5] = -pos[0];
84  dmdg[2][2] = 1.; dmdg[2][3] = -pos[1]; dmdg[2][4] = pos[0];
85  // d local residual / d global residual
86  auto drldrg = TMatrixD(3, 3);
87  for (int i = 0; i < 3; ++i) {
88  drldrg(0, i) = udir[i];
89  drldrg(1, i) = vdir[i];
90  drldrg(2, i) = ndir[i];
91  }
92  // d local residual / d global rigid body param
93  auto drldg = drldrg * (drdm * dmdg);
94 
95  // wire ends
96  const double zWireM = s_cdcGeometryTranslator->getWireBackwardPosition(getWireID(), CDCGeometryPar::c_Aligned)[2];
97  const double zWireP = s_cdcGeometryTranslator->getWireForwardPosition(getWireID(), CDCGeometryPar::c_Aligned)[2];
98  // relative Z position [0..1]
99  const double zRel = std::max(0., std::min(1., (pos[2] - zWireM) / (zWireP - zWireM)));
100 
101  // Layer alignment
102  // wire 511 = no wire (0 is reserved for existing wires) - this has to be compatible with code in CDCGeometryPar::setWirPosAlignParams
103  auto layerID = WireID(getWireID().getICLayer(), 511);
104 
105  // Alignment of layer X (bwd)
106  globals.add(
107  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerX),
108  drldg(0, 0)
109  );
110 
111  // Alignment of layer Y (bwd)
112  globals.add(
113  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerY),
114  drldg(0, 1)
115  );
116 
117  // Alignment of layer rotation (gamma) (bwd)
118  globals.add(
119  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerPhi),
120  drldg(0, 5)
121  );
122 
123  // Difference between wire ends (end plates)
124  // Alignment of layer dX, dX = foward - backward endplate
125  globals.add(
126  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerDx),
127  drldg(0, 0) * zRel
128  );
129 
130  // Alignment of layer dY, dY = foward - backward endplate
131  globals.add(
132  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerDy),
133  drldg(0, 1) * zRel
134  );
135 
136  // Alignment of layer rotation difference d(gamma or phi), dPhi = foward - backward endplate
137  globals.add(
138  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerDPhi),
139  drldg(0, 5) * zRel
140  );
141 
142  //WARNING: experimental (disabled by default)
143  // Wire-by-wire alignment
145  // How much shift (in X or Y) on BWD wire-end will change the residual at estimated track crossing
146  // with the wire (at relative z-position on wire = zRel)
147  double zRelM = fabs(1. - zRel);
148  // Same as above but for FWD wire-end (residual at zRel = zRel * delta(X or Y at FWD wire-end)
149  double zRelP = fabs(zRel - 0.);
150 
151  // Alignment of wires X in global coords at BWD wire-end
152  globals.add(
153  GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireBwdX),
154  drldg(0, 0) * zRelM
155  );
156  // Alignment of wires X in global coords at FWD wire-end
157  globals.add(
158  GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireFwdX),
159  drldg(0, 0) * zRelP
160  );
161 
162  // Alignment of wires Y in global coords at BWD wire-end
163  globals.add(
164  GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireBwdY),
165  drldg(0, 1) * zRelM
166  );
167  // Alignment of wires Y in global coords at FWD wire-end
168  globals.add(
169  GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireFwdY),
170  drldg(0, 1) * zRelP
171  );
172  }
173 
174  //WARNING: experimental (disabled by default)
175  //TODO: missing some factors (we do not align directly the wire-sag coefficient, but
176  // wire tension ... coef = pi * ro * r * r / 8 / tension
177  //TODO: need to get these numbers from CDCGeometryPar!
179  globals.add(
180  GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireTension),
181  drldg(0, 1) * 4.0 * zRel * (1.0 - zRel)
182  );
183  }
184 
185  return globals;
186 }

◆ labels()

virtual std::vector<int> labels ( )
inlinevirtualinherited

Vector of integer labels for calibration/alignment parameters available (must match #columns of derivatives(...))

unique across all sub-detectors in calibration

Returns
std::vector< int > Vector of integer labels

Definition at line 90 of file ICalibrationParametersDerivatives.h.

◆ localDerivatives()

TMatrixD localDerivatives ( const genfit::StateOnPlane sop)
overridevirtual

Derivatives for (local) fit parameters.

Parameters
sopState on virtual plane to calculate derivatives
Returns
TMatrixD of local derivatives, #columns=#params, #row=2 (or measurement dimension if > 2)

Reimplemented from ICalibrationParametersDerivatives.

Definition at line 189 of file AlignableCDCRecoHit.cc.

◆ localLabels()

virtual std::vector<int> localLabels ( )
inlinevirtualinherited

Vector of integer labels for local calibration parameters available (must match #columns of localDerivatives(...))

This will be usually ignored (e.g. does not have to match localDerivatives), but it is a good practice to return vector of zeros of correct size

Returns
std::vector< int > Vector of integer labels

Definition at line 151 of file ICalibrationParametersDerivatives.h.

◆ setLeftRightResolution()

void setLeftRightResolution ( int  lr)
inlineinherited

select how to resolve the left/right ambiguity: -1: negative (left) side on vector (wire direction) x (track direction) 0: mirrors enter with same weight, DAF will decide.

1: positive (right) side on vector (wire direction) x (track direction) where the wire direction is pointing towards +z except for small corrections such as stereo angle, sagging

Definition at line 111 of file CDCRecoHit.h.

◆ timeDerivativesMeasurementsOnPlane()

std::vector< double > timeDerivativesMeasurementsOnPlane ( const genfit::StateOnPlane state) const
inherited

Get the time derivative of the MesuredStateOnPlane (derived from the track fit).


Definition at line 263 of file CDCRecoHit.cc.

Member Data Documentation

◆ s_cosmics

bool s_cosmics = false
staticprotectedinherited

Switch to use cosmic events, or physics events from IP.

Default value is 0, which means "physics event" mode.

Definition at line 147 of file CDCRecoHit.h.

◆ s_useTrackTime

bool s_useTrackTime = false
staticprotectedinherited

Whether to use the track time or not when building the measurementOnPlane.

This needs to be in sync with the TDCCountTranslator.

Definition at line 142 of file CDCRecoHit.h.


The documentation for this class was generated from the following files:
genfit::TrackPoint
Object containing AbsMeasurement and AbsFitterInfo objects.
Definition: TrackPoint.h:46
Belle2::CDCAlignment::wireBwdX
static const baseType wireBwdX
Wire X position w.r.t. nominal on backward endplate.
Definition: CDCAlignment.h:39
Belle2::WireID
Class to identify a wire inside the CDC.
Definition: WireID.h:44
Belle2::B2Vector3::SetMag
void SetMag(DataType mag)
Set magnitude keeping theta and phi constant.
Definition: B2Vector3.h:181
genfit::MeasuredStateOnPlane
#StateOnPlane with additional covariance matrix.
Definition: MeasuredStateOnPlane.h:39
Belle2::CDCAlignment::layerPhi
static const baseType layerPhi
Layer rotation in global X-Y plane (gamma) at backward endplate.
Definition: CDCAlignment.h:64
genfit::AbsFitterInfo
This class collects all information needed and produced by a specific AbsFitter and is specific to on...
Definition: AbsFitterInfo.h:42
Belle2::CDCAlignment::wireFwdX
static const baseType wireFwdX
Wire X position w.r.t. nominal on forward endplate.
Definition: CDCAlignment.h:45
Belle2::CDCRecoHit::m_leftRight
signed char m_leftRight
Flag showing left/right passage.
Definition: CDCRecoHit.h:164
Belle2::CDCAlignment::wireBwdY
static const baseType wireBwdY
Wire Y position w.r.t. nominal on backward endplate.
Definition: CDCAlignment.h:41
Belle2::CDCRecoHit::constructPlane
genfit::SharedPlanePtr constructPlane(const genfit::StateOnPlane &state) const override
Methods that actually interface to Genfit.
Definition: CDCRecoHit.cc:82
Belle2::AlignableCDCRecoHit::s_enableWireByWireAlignmentGlobalDerivatives
static bool s_enableWireByWireAlignmentGlobalDerivatives
Static enabling(true) or disabling(false) addition of global derivatives for wire-by-wire alignment.
Definition: AlignableCDCRecoHit.h:40
Belle2::B2Vector3< double >
genfit::TrackPoint::getFitterInfo
AbsFitterInfo * getFitterInfo(const AbsTrackRep *rep=nullptr) const
Get fitterInfo for rep. Per default, use cardinal rep.
Definition: TrackPoint.cc:170
Belle2::CDC::CDCGeometryPar
The Class for CDC Geometry Parameters.
Definition: CDCGeometryPar.h:75
Belle2::CDCRecoHit::getWireID
WireID getWireID() const
Getter for WireID object.
Definition: CDCRecoHit.h:60
Belle2::CDC::CDCGeometryPar::Instance
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
Definition: CDCGeometryPar.cc:41
Belle2::CDC::CDCGeometryPar::getTheta
double getTheta(const TVector3 &momentum) const
Returns track incident angle (theta in rad.).
Definition: CDCGeometryPar.cc:2674
Belle2::CDCAlignment::wireTension
static const baseType wireTension
Wire tension w.r.t. nominal (=50. ?)
Definition: CDCAlignment.h:53
Belle2::CDC::CDCGeometryPar::getDriftTime
double getDriftTime(double dist, unsigned short layer, unsigned short lr, double alpha, double theta) const
Return the drift time to the sense wire.
Definition: CDCGeometryPar.cc:2457
Belle2::CDC::CDCGeometryPar::getDriftV
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.
Definition: CDCGeometryPar.cc:2006
Belle2::CDCRecoHit::m_adcCount
unsigned short m_adcCount
ADC Count as out of CDCHit.
Definition: CDCRecoHit.h:155
prepareAsicCrosstalkSimDB.fi
fi
open root file to store x-talk probability histogram
Definition: prepareAsicCrosstalkSimDB.py:88
Belle2::AlignableCDCRecoHit::s_enableWireSaggingGlobalDerivative
static bool s_enableWireSaggingGlobalDerivative
Static enabling(true) or disabling(false) addition of global derivative for wire sagging coefficient ...
Definition: AlignableCDCRecoHit.h:38
Belle2::CDCAlignment::layerY
static const baseType layerY
Layer shift in global Y at backward endplate.
Definition: CDCAlignment.h:62
Belle2::CDCAlignment::layerDy
static const baseType layerDy
Layer shift in global Y dY = foward - backward endplate.
Definition: CDCAlignment.h:69
Belle2::CDCAlignment::layerDx
static const baseType layerDx
Layer shift in global X dX = foward - backward endplate.
Definition: CDCAlignment.h:67
Belle2::CDCRecoHit::s_cdcGeometryTranslator
static std::unique_ptr< CDC::CDCGeometryTranslatorBase > s_cdcGeometryTranslator
Object for geometry translation.
Definition: CDCRecoHit.h:134
Belle2::CDCAlignment::wireFwdY
static const baseType wireFwdY
Wire Y position w.r.t. nominal on forward endplate.
Definition: CDCAlignment.h:47
Belle2::CDCAlignment::layerX
static const baseType layerX
Layer shift in global X at backward endplate.
Definition: CDCAlignment.h:60
Belle2::WireID::getICLayer
unsigned short getICLayer() const
Getter for continuous layer numbering.
Definition: WireID.cc:26
Belle2::CDC::CDCGeometryPar::getAlpha
double getAlpha(const TVector3 &posOnWire, const TVector3 &momentum) const
Returns track incident angle in rphi plane (alpha in rad.).
Definition: CDCGeometryPar.cc:2661
Belle2::CDCAlignment::layerDPhi
static const baseType layerDPhi
Layer rotation in global X-Y plane (gamma) dPhi = foward - backward endplate.
Definition: CDCAlignment.h:71