Belle II Software  release-08-01-10
AlignablePXDRecoHit Class Reference

This class is used to transfer PXD information to the track fit. More...

#include <AlignablePXDRecoHit.h>

Inheritance diagram for AlignablePXDRecoHit:
Collaboration diagram for AlignablePXDRecoHit:

Public Member Functions

virtual ~AlignablePXDRecoHit ()
 Destructor.
 
genfit::AbsMeasurementclone () const override
 Creating a deep 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...
 
 PXDRecoHit ()
 Inherit constructors.
 
 PXDRecoHit (const PXDTrueHit *hit, const genfit::TrackCandHit *trackCandHit=NULL, float sigmaU=-1, float sigmaV=-1)
 Inherit constructors.
 
 PXDRecoHit (const PXDCluster *hit, float sigmaU, float sigmaV, float covUV)
 Inherit constructors.
 
 PXDRecoHit (const PXDCluster *hit, const genfit::TrackCandHit *trackCandHit=NULL)
 Inherit constructors.
 
virtual std::vector< genfit::MeasurementOnPlane * > constructMeasurementsOnPlane (const genfit::StateOnPlane &state) const override
 Methods that actually interface to Genfit. More...
 
VxdID getSensorID () const
 Get the compact ID.
 
const PXDTrueHitgetTrueHit () const
 Get pointer to the TrueHit used when creating this RecoHit, can be NULL if created from something else.
 
const PXDClustergetCluster () const
 Get pointer to the Cluster used when creating this RecoHit, can be NULL if created from something else.
 
float getU () const
 Get u coordinate.
 
float getV () const
 Get v coordinate.
 
float getUVariance () const
 Get u coordinate variance.
 
float getVVariance () const
 Get v coordinate variance.
 
float getUVCov () const
 Get u-v error covariance.
 
float getEnergyDep () const
 Get deposited energy.
 
float getShapeLikelyhood (const genfit::StateOnPlane &state) const
 Get deposited energy error. More...
 
virtual const genfit::AbsHMatrixconstructHMatrix (const genfit::AbsTrackRep *) const override
 Construct the hessian matrix.
 
int getPlaneId () const
 
virtual SharedPlanePtr constructPlane (const StateOnPlane &state) const override
 Construct (virtual) detector plane (use state's AbsTrackRep). More...
 
virtual void setPlane (const SharedPlanePtr &physicalPlane, int planeId=-1)
 
void setStripV (bool v=true)
 Use if the coordinate for 1D hits measured in V direction. 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
 
virtual bool isLeftRightMeasurement () const
 If the AbsMeasurement is a wire hit, the left/right resolution will be used.
 
virtual int getLeftRightResolution () 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 TMatrixD localDerivatives (const genfit::StateOnPlane *)
 Derivatives for additional local parameters to be fitted in global calibration algorithms together with with global parameters. More...
 
virtual std::vector< int > localLabels ()
 Vector of integer labels for local calibration parameters available (must match #columns of localDerivatives(...)) More...
 

Static Public Attributes

static bool s_enableLorentzGlobalDerivatives = false
 Static enabling(true) or disabling(false) addition of global derivatives for Lorentz shift.
 

Protected Attributes

SharedPlanePtr physicalPlane_
 
int planeId_
 This is persistent, but '!' makes ROOT shut up.
 
bool stripV_
 
TVectorD rawHitCoords_
 
TMatrixDSym rawHitCov_
 
int detId_
 
int hitId_
 
TrackPoint * trackPoint_
 Pointer to TrackPoint where the measurement belongs to.
 

Private Types

enum  { HIT_DIMENSIONS = 2 }
 

Private Member Functions

 ClassDefOverride (AlignablePXDRecoHit, 3)
 PXD RecoHit extended for alignment/calibration.
 
void setDetectorPlane ()
 Set up Detector plane information.
 
TVectorD applyPlanarDeformation (TVectorD hitCoords, std::vector< double > planarParameters, const genfit::StateOnPlane &state) const
 Apply planar deformation of sensors.
 

Private Attributes

unsigned short m_sensorID
 Unique sensor identifier.
 
const PXDTrueHitm_trueHit
 Pointer to the TrueHit used when creating this object.
 
const PXDClusterm_cluster
 transient member (not written out during streaming) More...
 
float m_energyDep
 transient member (not written out during streaming) More...
 

Friends

class PXDRecoHit
 

Detailed Description

This class is used to transfer PXD information to the track fit.

Definition at line 28 of file AlignablePXDRecoHit.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
privateinherited
Enumerator
HIT_DIMENSIONS 

sensitive Dimensions of the Hit

Definition at line 134 of file PXDRecoHit.h.

134 { HIT_DIMENSIONS = 2 };
@ HIT_DIMENSIONS
sensitive Dimensions of the Hit
Definition: PXDRecoHit.h:134

Member Function Documentation

◆ constructMeasurementsOnPlane()

std::vector< genfit::MeasurementOnPlane * > constructMeasurementsOnPlane ( const genfit::StateOnPlane state) const
overridevirtualinherited

Methods that actually interface to Genfit.


Reimplemented from PlanarMeasurement.

Definition at line 192 of file PXDRecoHit.cc.

193 {
194  // Track-based update only takes place when the RecoHit has an associated cluster
195  if (this->getCluster()) {
196  // Check if we can correct position coordinates based on track info
197  const TVectorD& state5 = state.getState();
198  auto offset = PXD::PXDClusterPositionEstimator::getInstance().getClusterOffset(*this->getCluster(), state5[1], state5[2]);
199 
200  if (offset != nullptr) {
201  // Found a valid offset, lets apply it
202  const Belle2::VxdID& sensorID = (*this->getCluster()).getSensorID();
203  const Belle2::PXD::SensorInfo& Info = dynamic_cast<const Belle2::PXD::SensorInfo&>(VXD::GeoCache::get(sensorID));
204  double posU = Info.getUCellPosition((*this->getCluster()).getUStart());
205  double posV = Info.getVCellPosition((*this->getCluster()).getVStart());
206 
207  TVectorD hitCoords(2);
208  hitCoords(0) = posU + offset->getU();
209  hitCoords(1) = posV + offset->getV();
210  TMatrixDSym hitCov(2);
211  hitCov(0, 0) = offset->getUSigma2();
212  hitCov(0, 1) = offset->getUVCovariance();
213  hitCov(1, 0) = offset->getUVCovariance();
214  hitCov(1, 1) = offset->getVSigma2();
215 
216  // Apply planar deformation
217  TVectorD pos = applyPlanarDeformation(hitCoords, VXD::GeoCache::get(m_sensorID).getSurfaceParameters(), state);
218 
219  return std::vector<genfit::MeasurementOnPlane*>(1, new genfit::MeasurementOnPlane(
220  pos, hitCov, state.getPlane(), state.getRep(), this->constructHMatrix(state.getRep())
221  ));
222  }
223  }
224 
225  // Apply planar deformation
226  TVectorD pos = applyPlanarDeformation(rawHitCoords_, VXD::GeoCache::get(m_sensorID).getSurfaceParameters(), state);
227 
228  // If we reach here, we can do no better than what we have
229  return std::vector<genfit::MeasurementOnPlane*>(1, new genfit::MeasurementOnPlane(
230  pos, rawHitCov_, state.getPlane(), state.getRep(), this->constructHMatrix(state.getRep())
231  ));
232 }
unsigned short m_sensorID
Unique sensor identifier.
Definition: PXDRecoHit.h:136
const PXDCluster * getCluster() const
Get pointer to the Cluster used when creating this RecoHit, can be NULL if created from something els...
Definition: PXDRecoHit.h:106
TVectorD applyPlanarDeformation(TVectorD hitCoords, std::vector< double > planarParameters, const genfit::StateOnPlane &state) const
Apply planar deformation of sensors.
Definition: PXDRecoHit.cc:149
VxdID getSensorID() const
Get the compact ID.
Definition: PXDRecoHit.h:101
static PXDClusterPositionEstimator & getInstance()
Main (and only) way to access the PXDClusterPositionEstimator.
const PXDClusterOffsetPar * getClusterOffset(const PXDCluster &cluster, double tu, double tv) const
Return pointer to cluster offsets, can be nullptr.
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
static const SensorInfoBase & get(Belle2::VxdID id)
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.h:139
double getVCellPosition(int vID) const
Return the position of a specific strip/pixel in v direction.
double getUCellPosition(int uID, int vID=-1) const
Return the position of a specific strip/pixel in u direction.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
Measured coordinates on a plane.

◆ constructPlane()

SharedPlanePtr constructPlane ( const StateOnPlane state) const
overridevirtualinherited

Construct (virtual) detector plane (use state's AbsTrackRep).

It's possible to make corrections to the plane here. The state should be defined somewhere near the measurement. For virtual planes, the state will be extrapolated to the POCA to point (SpacepointMeasurement) or line (WireMeasurement), and from this info the plane will be constructed.

Implements AbsMeasurement.

Definition at line 46 of file PlanarMeasurement.cc.

◆ 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.

◆ getShapeLikelyhood()

float getShapeLikelyhood ( const genfit::StateOnPlane state) const
inherited

Get deposited energy error.

Get the likelyhood that cluster shape is likely to be created from track state.

Definition at line 137 of file PXDRecoHit.cc.

◆ 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 25 of file AlignablePXDRecoHit.cc.

26 {
27  auto alignment = GlobalCalibrationManager::getInstance().getAlignmentHierarchy().getGlobalDerivatives<VXDAlignment>(getPlaneId(),
28  sop);
29 
30  auto globals = GlobalDerivatives(alignment);
31 
33  auto lorentz = GlobalCalibrationManager::getInstance().getLorentzShiftHierarchy().getGlobalDerivatives<VXDAlignment>(getPlaneId(),
34  sop, BFieldManager::getInstance().getField(ROOT::Math::XYZVector(sop->getPos())));
35  globals.add(lorentz);
36  }
37 
38  const PXD::SensorInfo& geometry = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::get(getSensorID()));
39 
40  // Legendre parametrization of deformation
41  auto L1 = [](double x) {return x;};
42  auto L2 = [](double x) {return (3 * x * x - 1) / 2;};
43  auto L3 = [](double x) {return (5 * x * x * x - 3 * x) / 2;};
44  auto L4 = [](double x) {return (35 * x * x * x * x - 30 * x * x + 3) / 8;};
45 
46  double du_dw = sop->getState()[1]; // slope in U direction
47  double dv_dw = sop->getState()[2]; // slope in V direction
48  double u = getU(); // U coordinate of hit
49  double v = getV(); // V coordinate of hit
50  double width = geometry.getWidth(v); // Width of sensor (U side)
51  double length = geometry.getLength(); // Length of sensor (V side)
52  u = u * 2 / width; // Legendre parametrization required U in (-1, 1)
53  v = v * 2 / length; // Legendre parametrization required V in (-1, 1)
54 
55  // Add parameters of surface deformation to alignment
56  // Numbering of VXD alignment parameters:
57  // -> 0-6: Rigid body alignment
58  // -> 31-33: First level of surface deformation
59  // -> 41-44: Second level of surface deformation
60  // -> 51-55: Third level of surface deformation
61  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 31), std::vector<double> {L2(u)*du_dw, L2(u)*dv_dw});
62  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 32), std::vector<double> {L1(u)*L1(v)*du_dw, L1(u)*L1(v)*dv_dw});
63  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 33), std::vector<double> {L2(v)*du_dw, L2(v)*dv_dw});
64 
65  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 41), std::vector<double> {L3(u)*du_dw, L3(u)*dv_dw});
66  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 42), std::vector<double> {L2(u)*L1(v)*du_dw, L2(u)*L1(v)*dv_dw});
67  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 43), std::vector<double> {L1(u)*L2(v)*du_dw, L1(u)*L2(v)*dv_dw});
68  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 44), std::vector<double> {L3(v)*du_dw, L3(v)*dv_dw});
69 
70  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 51), std::vector<double> {L4(u)*du_dw, L4(u)*dv_dw});
71  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 52), std::vector<double> {L3(u)*L1(v)*du_dw, L3(u)*L1(v)*dv_dw});
72  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 53), std::vector<double> {L2(u)*L2(v)*du_dw, L2(u)*L2(v)*dv_dw});
73  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 54), std::vector<double> {L1(u)*L3(v)*du_dw, L1(u)*L3(v)*dv_dw});
74  globals.add(GlobalLabel::construct<VXDAlignment>(getSensorID(), 55), std::vector<double> {L4(v)*du_dw, L4(v)*dv_dw});
75 
76  return GlobalDerivatives(globals);
77 }
static bool s_enableLorentzGlobalDerivatives
Static enabling(true) or disabling(false) addition of global derivatives for Lorentz shift.
static BFieldManager & getInstance()
Return the instance of the magnetic field manager.
virtual double add(baseType id, baseType param, double value, bool subtractInsteadOfAdd=false)
Add correction to already stored (or to 0. if not set yet) constant value (optionaly with minus sign)
float getV() const
Get v coordinate.
Definition: PXDRecoHit.h:111
float getU() const
Get u coordinate.
Definition: PXDRecoHit.h:109
VXD alignment (and maybe some calibration) parameters.
Definition: VXDAlignment.h:19
static void getField(const double *pos, double *field)
return the magnetic field at a given position.
Definition: BFieldManager.h:91

◆ 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()

virtual TMatrixD localDerivatives ( const genfit::StateOnPlane )
inlinevirtualinherited

Derivatives for additional local parameters to be fitted in global calibration algorithms together with with global parameters.

Local parameters are not neccesarily identified by label because their number is proportional to number of measurements included in calibration (possibly very huge number!)

Returns
TMatrixD Matrix in form d_residual_i/d_parameter_j

Reimplemented in AlignableCDCRecoHit.

Definition at line 141 of file ICalibrationParametersDerivatives.h.

◆ 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.

◆ setStripV()

void setStripV ( bool  v = true)
inlineinherited

Use if the coordinate for 1D hits measured in V direction.

Per default for 1D planar hits, the coordinate is measured in U direction. With this function you can set it to be measured in V direction. This affects the outcoe of constructHMatrix().

Definition at line 70 of file PlanarMeasurement.h.

Member Data Documentation

◆ m_cluster

const PXDCluster* m_cluster
privateinherited

transient member (not written out during streaming)

Pointer to the Cluster used when creating this object

Definition at line 140 of file PXDRecoHit.h.

◆ m_energyDep

float m_energyDep
privateinherited

transient member (not written out during streaming)

deposited energy.

Definition at line 141 of file PXDRecoHit.h.


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