Belle II Software  release-08-01-10
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=nullptr, 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 21 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 112 of file CDCRecoHit.h.

113  {
114  return m_cdcHit;
115  }
const CDCHit * m_cdcHit
Pointer to the CDCHit used to created this CDCRecoHit.
Definition: CDCRecoHit.h:150

◆ getFlyByDistanceVector()

bool getFlyByDistanceVector ( B2Vector3D pointingVector,
B2Vector3D trackDir,
const genfit::AbsTrackRep rep = nullptr,
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 368 of file CDCRecoHit.cc.

372 {
373  const genfit::TrackPoint* tp = this->getTrackPoint();
374  if (!tp) {
375  B2ERROR("No genfit::TrackPoint for CDCRecoHit.");
376  return false;
377  }
378  const genfit::AbsFitterInfo* fi = tp->getFitterInfo(rep);
379  if (!fi) {
380  B2DEBUG(200, "No genfit::AbsFitterInfo for this CDCRecoHit.");
381  return false;
382  }
383 
384  const genfit::MeasuredStateOnPlane& mop = fi->getFittedState();
385  B2Vector3D fittedPoca = mop.getPos();
386  // constructPlane places the coordinate center in the POCA to the
387  // wire. Using this is the default behavior. If this should be too
388  // slow, as it has to re-evaluate the POCA, alternatively the user
389  // can set usePlaneFromFit which uses the plane determined by the
390  // track fit.
391  B2Vector3D pocaOnWire = (usePlaneFromFit
392  ? mop.getPlane()->getO()
393  : this->constructPlane(mop)->getO());
394 
395  // The vector from the wire to the track.
396  pointingVector = fittedPoca - pocaOnWire;
397 
398  trackDir = mop.getMom();
399  trackDir.SetMag(1.);
400  return true;
401 }
void SetMag(DataType mag)
Set magnitude keeping theta and phi constant.
Definition: B2Vector3.h:182
genfit::SharedPlanePtr constructPlane(const genfit::StateOnPlane &state) const override
Methods that actually interface to Genfit.
Definition: CDCRecoHit.cc:80
This class collects all information needed and produced by a specific AbsFitter and is specific to on...
Definition: AbsFitterInfo.h:42
#StateOnPlane with additional covariance matrix.
Object containing AbsMeasurement and AbsFitterInfo objects.
Definition: TrackPoint.h:46
AbsFitterInfo * getFitterInfo(const AbsTrackRep *rep=nullptr) const
Get fitterInfo for rep. Per default, use cardinal rep.
Definition: TrackPoint.cc:170

◆ 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 number of rows = dimension of residual, number of columns = number of parameters. number of columns must match vector<int>.size().

Reimplemented from ICalibrationParametersDerivatives.

Definition at line 28 of file AlignableCDCRecoHit.cc.

29 {
30  GlobalDerivatives globals;
31  unsigned short LR = (int(m_leftRight) > 0.) ? 1 : 0;
32 
33  const B2Vector3D& mom = sop->getMom();
34  const B2Vector3D& wirePositon = sop->getPlane()->getO();
35  const unsigned short layer = getWireID().getICLayer();
36 
38  double alpha = cdcgeo.getAlpha(wirePositon, mom);
39  double theta = cdcgeo.getTheta(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);
44 
45  // CDC Calibration -------------------------------------------------
46 
47  // Time zero calibration (per wire)
48  if (driftTime > 20 && driftTime < 200 && fabs(driftVelocity) < 1.0e-2) {
49  globals.add(
50  GlobalLabel::construct<CDCTimeZeros>(getWireID(), 0),
51  driftVelocity * double(int(m_leftRight))
52  );
53  }
54 
55  // Time walk calibration (per board)
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);
60  globals.add(
61  GlobalLabel::construct<CDCTimeWalks>(board, 0),
62  -driftVelocity * exp(-param[1]* m_adcCount) * double(int(m_leftRight))
63  );
64  globals.add(
65  GlobalLabel::construct<CDCTimeWalks>(board, 1),
66  driftVelocity * m_adcCount * param[0]*exp(-param[1]* m_adcCount) * double(int(m_leftRight))
67  );
68  }
69 
70  // Space time relations calibration
71  if (driftTime > 20 && driftTime < 500 && fabs(driftVelocity) < 1.0e-2) {
72  // TODO: ugly to need to ask XTRelations for something here...
73  // Can't I get this CDCGeometryPar or sth. like this?
74 
75  theta = cdcgeo.getOutgoingTheta(alpha, theta);
76  alpha = cdcgeo.getOutgoingAlpha(alpha);
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) {
82  globals.add(
83  GlobalLabel::construct<CDCXtRelations>(xtid, 0),
84  ROOT::Math::Chebyshev5(driftTime, 1, 0, 0, 0, 0, 0) * double(int(m_leftRight))
85  );
86  globals.add(
87  GlobalLabel::construct<CDCXtRelations>(xtid, 1),
88  ROOT::Math::Chebyshev5(driftTime, 0, 1, 0, 0, 0, 0) * double(int(m_leftRight))
89  );
90  globals.add(
91  GlobalLabel::construct<CDCXtRelations>(xtid, 2),
92  ROOT::Math::Chebyshev5(driftTime, 0, 0, 1, 0, 0, 0) * double(int(m_leftRight))
93  );
94  globals.add(
95  GlobalLabel::construct<CDCXtRelations>(xtid, 3),
96  ROOT::Math::Chebyshev5(driftTime, 0, 0, 0, 1, 0, 0) * double(int(m_leftRight))
97  );
98  globals.add(
99  GlobalLabel::construct<CDCXtRelations>(xtid, 4),
100  ROOT::Math::Chebyshev5(driftTime, 0, 0, 0, 0, 1, 0) * double(int(m_leftRight))
101  );
102  globals.add(
103  GlobalLabel::construct<CDCXtRelations>(xtid, 5),
104  ROOT::Math::Chebyshev5(driftTime, 0, 0, 0, 0, 0, 1) * double(int(m_leftRight))
105  );
106  } else {
107  globals.add(
108  GlobalLabel::construct<CDCXtRelations>(xtid, 7),
109  double(int(m_leftRight))
110  );
111  globals.add(
112  GlobalLabel::construct<CDCXtRelations>(xtid, 8),
113  (driftTime - boundary) * double(int(m_leftRight))
114  );
115  }
116  }
117 
118 
119  // CDC Alignment ---------------------------------------------------
120 
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();
126 
127  auto tn = tdir[0] * ndir[0] + tdir[1] * ndir[1] + tdir[2] * ndir [2];
128  // d residual / d measurement
129  auto drdm = TMatrixD(3, 3);
130  drdm.UnitMatrix();
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;
134  }
135  }
136  // d measurement / d global rigid body param
137  auto dmdg = TMatrixD(3, 6);
138  dmdg.Zero();
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];
142  // d local residual / d global residual
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];
148  }
149  // d local residual / d global rigid body param
150  auto drldg = drldrg * (drdm * dmdg);
151 
152  // wire ends
153  const double zWireM = s_cdcGeometryTranslator->getWireBackwardPosition(getWireID(), CDCGeometryPar::c_Aligned)[2];
154  const double zWireP = s_cdcGeometryTranslator->getWireForwardPosition(getWireID(), CDCGeometryPar::c_Aligned)[2];
155  // relative Z position [0..1]
156  const double zRel = std::max(0., std::min(1., (pos[2] - zWireM) / (zWireP - zWireM)));
157 
158  // Layer alignment
159  // wire 511 = no wire (0 is reserved for existing wires) - this has to be compatible with code in CDCGeometryPar::setWirPosAlignParams
160  auto layerID = WireID(getWireID().getICLayer(), 511);
161 
162  // Alignment of layer X (bwd)
163  globals.add(
164  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerX),
165  drldg(0, 0)
166  );
167 
168  // Alignment of layer Y (bwd)
169  globals.add(
170  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerY),
171  drldg(0, 1)
172  );
173 
174  // Alignment of layer rotation (gamma) (bwd)
175  globals.add(
176  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerPhi),
177  drldg(0, 5)
178  );
179 
180  // Difference between wire ends (end plates)
181  // Alignment of layer dX, dX = foward - backward endplate
182  globals.add(
183  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerDx),
184  drldg(0, 0) * zRel
185  );
186 
187  // Alignment of layer dY, dY = foward - backward endplate
188  globals.add(
189  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerDy),
190  drldg(0, 1) * zRel
191  );
192 
193  // Alignment of layer rotation difference d(gamma or phi), dPhi = foward - backward endplate
194  globals.add(
195  GlobalLabel::construct<CDCAlignment>(layerID, CDCAlignment::layerDPhi),
196  drldg(0, 5) * zRel
197  );
198 
199  //WARNING: experimental (disabled by default)
200  // Wire-by-wire alignment
202  // How much shift (in X or Y) on BWD wire-end will change the residual at estimated track crossing
203  // with the wire (at relative z-position on wire = zRel)
204  double zRelM = fabs(1. - zRel);
205  // Same as above but for FWD wire-end (residual at zRel = zRel * delta(X or Y at FWD wire-end)
206  double zRelP = fabs(zRel - 0.);
207 
208  // Alignment of wires X in global coords at BWD wire-end
209  globals.add(
210  GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireBwdX),
211  drldg(0, 0) * zRelM
212  );
213  // Alignment of wires X in global coords at FWD wire-end
214  globals.add(
215  GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireFwdX),
216  drldg(0, 0) * zRelP
217  );
218 
219  // Alignment of wires Y in global coords at BWD wire-end
220  globals.add(
221  GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireBwdY),
222  drldg(0, 1) * zRelM
223  );
224  // Alignment of wires Y in global coords at FWD wire-end
225  globals.add(
226  GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireFwdY),
227  drldg(0, 1) * zRelP
228  );
229  }
230 
231  //WARNING: experimental (disabled by default)
232  //TODO: missing some factors (we do not align directly the wire-sag coefficient, but
233  // wire tension ... coef = pi * ro * r * r / 8 / tension
234  //TODO: need to get these numbers from CDCGeometryPar!
236  globals.add(
237  GlobalLabel::construct<CDCAlignment>(getWireID().getEWire(), CDCAlignment::wireTension),
238  drldg(0, 1) * 4.0 * zRel * (1.0 - zRel)
239  );
240  }
241 
242 
243 
244  return globals;
245 }
static bool s_enableWireSaggingGlobalDerivative
Static enabling(true) or disabling(false) addition of global derivative for wire sagging coefficient ...
static bool s_enableWireByWireAlignmentGlobalDerivatives
Static enabling(true) or disabling(false) addition of global derivatives for wire-by-wire alignment.
static const baseType layerDPhi
Layer rotation in global X-Y plane (gamma) dPhi = foward - backward endplate.
Definition: CDCAlignment.h:61
static const baseType layerDy
Layer shift in global Y dY = foward - backward endplate.
Definition: CDCAlignment.h:59
static const baseType layerDx
Layer shift in global X dX = foward - backward endplate.
Definition: CDCAlignment.h:57
static const baseType wireBwdY
Wire Y position w.r.t. nominal on backward endplate.
Definition: CDCAlignment.h:31
static const baseType wireFwdY
Wire Y position w.r.t. nominal on forward endplate.
Definition: CDCAlignment.h:37
static const baseType wireFwdX
Wire X position w.r.t. nominal on forward endplate.
Definition: CDCAlignment.h:35
static const baseType wireBwdX
Wire X position w.r.t. nominal on backward endplate.
Definition: CDCAlignment.h:29
static const baseType layerY
Layer shift in global Y at backward endplate.
Definition: CDCAlignment.h:52
static const baseType layerX
Layer shift in global X at backward endplate.
Definition: CDCAlignment.h:50
static const baseType layerPhi
Layer rotation in global X-Y plane (gamma) at backward endplate.
Definition: CDCAlignment.h:54
static const baseType wireTension
Wire tension w.r.t. nominal (=50. ?)
Definition: CDCAlignment.h:43
signed char m_leftRight
Flag showing left/right passage.
Definition: CDCRecoHit.h:153
WireID getWireID() const
Getter for WireID object.
Definition: CDCRecoHit.h:49
unsigned short m_adcCount
ADC Count as out of CDCHit.
Definition: CDCRecoHit.h:144
static std::unique_ptr< CDC::CDCGeometryTranslatorBase > s_cdcGeometryTranslator
Object for geometry translation.
Definition: CDCRecoHit.h:123
The Class for CDC Geometry Parameters.
double getTheta(const B2Vector3D &momentum) const
Returns track incident angle (theta in rad.).
unsigned short getBoardID(const WireID &wID) const
Returns frontend board id. corresponding to the wire id.
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.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Class to identify a wire inside the CDC.
Definition: WireID.h:34
unsigned short getICLayer() const
Getter for continuous layer numbering.
Definition: WireID.cc:24

◆ 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, number of columns = number of params, number of rows = 2 (or measurement dimension if > 2)

Reimplemented from ICalibrationParametersDerivatives.

Definition at line 248 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 100 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 261 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 136 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 131 of file CDCRecoHit.h.


The documentation for this class was generated from the following files: