CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
LocalTrajectoryParameters Class Reference

#include <LocalTrajectoryParameters.h>

Public Member Functions

float absdz () const
 
TrackCharge charge () const
 Charge (-1, 0 or 1) More...
 
LocalVector direction () const
 Momentum vector unit in the local frame. More...
 
LocalVector directionNotNormalized () const
 Momentum vector unit in the local frame. More...
 
float dxdz () const
 
float dydz () const
 
 LocalTrajectoryParameters ()
 
 LocalTrajectoryParameters (const AlgebraicVector5 &v, float aPzSign, bool charged=true)
 
 LocalTrajectoryParameters (float aQbp, float aDxdz, float aDydz, float aX, float aY, float aPzSign, bool charged=true)
 
 LocalTrajectoryParameters (const LocalPoint &pos, const LocalVector &p, TrackCharge charge)
 Constructor from local position, momentum and charge. More...
 
AlgebraicVector5 mixedFormatVector () const
 
LocalVector momentum () const
 Momentum vector in the local frame. More...
 
LocalPoint position () const
 Local x and y position coordinates. More...
 
float pzSign () const
 Sign of the z-component of the momentum in the local frame. More...
 
float qbp () const
 
float signedInverseMomentum () const
 Signed inverse momentum q/p (zero for neutrals). More...
 
bool updateP (float dP)
 Update of momentum by a scalar dP. More...
 
AlgebraicVector5 vector () const
 

Private Attributes

short theCharge
 charge More...
 
float theDxdz
 tangent of direction in local x vs. z More...
 
float theDydz
 tangent of direction in local y vs. z More...
 
short thePzSign
 sign of local pz More...
 
float theQbp
 q/p (charged) or 1/p (neutral) More...
 
float theX
 local x position More...
 
float theY
 local y position More...
 

Detailed Description

Class providing access to a set of relevant parameters of a trajectory in a local, Cartesian frame. The set consists of the following parameters:

q/p : charged particles: charge (plus or minus one) divided by magnitude of momentum
neutral particles: inverse magnitude of momentum
dxdz : direction tangent in local xz-plane
dydz : direction tangent in local yz-plane
x : local x-coordinate
y : local y-coordinate

In addition, the sign of local p_z is needed to fully define the direction of the track in this local frame.

Definition at line 25 of file LocalTrajectoryParameters.h.

Constructor & Destructor Documentation

◆ LocalTrajectoryParameters() [1/4]

LocalTrajectoryParameters::LocalTrajectoryParameters ( )
inline

Definition at line 29 of file LocalTrajectoryParameters.h.

29 {}

◆ LocalTrajectoryParameters() [2/4]

LocalTrajectoryParameters::LocalTrajectoryParameters ( const AlgebraicVector5 v,
float  aPzSign,
bool  charged = true 
)
inline

Constructor from vector of parameters.

Expects a vector of parameters as defined above, plus the sign of p_z. For charged particles the charge will be determined by the sign of the first element. For neutral particles the last argument should be false, in which case the charge of the first element will be neglected.

Definition at line 38 of file LocalTrajectoryParameters.h.

References theCharge, theDxdz, theDydz, thePzSign, theQbp, theX, theY, and findQualityFiles::v.

38  {
39  theQbp = v[0];
40  theDxdz = v[1];
41  theDydz = v[2];
42  theX = v[3];
43  theY = v[4];
44  thePzSign = aPzSign;
45  if (charged)
46  theCharge = theQbp > 0 ? 1 : -1;
47  else
48  theCharge = 0;
49  }
float theQbp
q/p (charged) or 1/p (neutral)
float theDxdz
tangent of direction in local x vs. z
float theDydz
tangent of direction in local y vs. z

◆ LocalTrajectoryParameters() [3/4]

LocalTrajectoryParameters::LocalTrajectoryParameters ( float  aQbp,
float  aDxdz,
float  aDydz,
float  aX,
float  aY,
float  aPzSign,
bool  charged = true 
)
inline

Constructor from individual parameters.

Expects parameters as defined above, plus the sign of p_z. For charged particles the charge will be determined by the sign of the first argument. For neutral particles the last argument should be false, in which case the charge of the first argument will be neglected.

Definition at line 58 of file LocalTrajectoryParameters.h.

References theCharge, and theQbp.

59  : theDxdz(aDxdz), theDydz(aDydz), theX(aX), theY(aY), thePzSign(aPzSign > 0 ? 1 : -1) {
60  if (charged) {
61  theQbp = aQbp;
62  theCharge = theQbp > 0 ? 1 : -1;
63  } else {
64  theQbp = aQbp;
65  theCharge = 0;
66  }
67  }
float theQbp
q/p (charged) or 1/p (neutral)
float theDxdz
tangent of direction in local x vs. z
float theDydz
tangent of direction in local y vs. z

◆ LocalTrajectoryParameters() [4/4]

LocalTrajectoryParameters::LocalTrajectoryParameters ( const LocalPoint pos,
const LocalVector p,
TrackCharge  charge 
)
inline

Constructor from local position, momentum and charge.

Definition at line 70 of file LocalTrajectoryParameters.h.

References charge(), AlCaHLTBitMon_ParallelJobs::p, and theQbp.

71  : theQbp(charge / p.mag()),
72  theDxdz(p.x() / p.z()),
73  theDydz(p.y() / p.z()),
74  theX(pos.x()),
75  theY(pos.y()),
76  thePzSign(p.z() > 0. ? 1 : -1),
77  theCharge(charge) {
78  if (charge == 0)
79  theQbp = 1.f / p.mag();
80  }
float theQbp
q/p (charged) or 1/p (neutral)
TrackCharge charge() const
Charge (-1, 0 or 1)
float theDxdz
tangent of direction in local x vs. z
float theDydz
tangent of direction in local y vs. z

Member Function Documentation

◆ absdz()

float LocalTrajectoryParameters::absdz ( ) const
inline

Definition at line 161 of file LocalTrajectoryParameters.h.

References f, mathSSE::sqrt(), theDxdz, and theDydz.

161 { return 1.f / std::sqrt(1.f + theDxdz * theDxdz + theDydz * theDydz); }
T sqrt(T t)
Definition: SSEVec.h:19
double f[11][100]
float theDxdz
tangent of direction in local x vs. z
float theDydz
tangent of direction in local y vs. z

◆ charge()

TrackCharge LocalTrajectoryParameters::charge ( ) const
inline

◆ direction()

LocalVector LocalTrajectoryParameters::direction ( ) const
inline

Momentum vector unit in the local frame.

Definition at line 99 of file LocalTrajectoryParameters.h.

References PVValHelper::dx, PVValHelper::dy, PVValHelper::dz, f, dqmMemoryStats::float, mathSSE::sqrt(), theDxdz, theDydz, and thePzSign.

Referenced by JacobianCurvilinearToLocal::JacobianCurvilinearToLocal(), and JacobianLocalToCurvilinear::JacobianLocalToCurvilinear().

99  {
100  float dz = float(thePzSign) / std::sqrt(1.f + theDxdz * theDxdz + theDydz * theDydz);
101  float dx = dz * theDxdz;
102  float dy = dz * theDydz;
103  return LocalVector(dx, dy, dz);
104  }
Local3DVector LocalVector
Definition: LocalVector.h:12
T sqrt(T t)
Definition: SSEVec.h:19
double f[11][100]
float theDxdz
tangent of direction in local x vs. z
float theDydz
tangent of direction in local y vs. z

◆ directionNotNormalized()

LocalVector LocalTrajectoryParameters::directionNotNormalized ( ) const
inline

Momentum vector unit in the local frame.

Definition at line 107 of file LocalTrajectoryParameters.h.

References f, theDxdz, and theDydz.

Referenced by StripCPE::getAlgoParam(), and TkClonerImpl::operator()().

107 { return LocalVector(theDxdz, theDydz, 1.f); }
Local3DVector LocalVector
Definition: LocalVector.h:12
double f[11][100]
float theDxdz
tangent of direction in local x vs. z
float theDydz
tangent of direction in local y vs. z

◆ dxdz()

float LocalTrajectoryParameters::dxdz ( ) const
inline

◆ dydz()

float LocalTrajectoryParameters::dydz ( ) const
inline

◆ mixedFormatVector()

AlgebraicVector5 LocalTrajectoryParameters::mixedFormatVector ( ) const
inline

Vector of parameters in internal representation.

Vector of parameters as defined above, with the first element = q/p for charged and = 1/p for neutral.

Definition at line 135 of file LocalTrajectoryParameters.h.

References theDxdz, theDydz, theQbp, theX, theY, and findQualityFiles::v.

Referenced by BzeroReferenceTrajectory::BzeroReferenceTrajectory(), SiStripBackplaneCalibration::derivatives(), TwoBowedSurfacesAlignmentParameters::derivatives(), DualBzeroReferenceTrajectory::extractParameters(), DualReferenceTrajectory::extractParameters(), SegmentAlignmentDerivatives4D::operator()(), KarimakiAlignmentDerivatives::operator()(), BowedSurfaceAlignmentDerivatives::operator()(), operator<<(), and ReferenceTrajectory::ReferenceTrajectory().

135  {
137  v[0] = theQbp; // signed in case of charged particles, 1/p for neutrals
138  v[1] = theDxdz;
139  v[2] = theDydz;
140  v[3] = theX;
141  v[4] = theY;
142  return v;
143  }
float theQbp
q/p (charged) or 1/p (neutral)
ROOT::Math::SVector< double, 5 > AlgebraicVector5
float theDxdz
tangent of direction in local x vs. z
float theDydz
tangent of direction in local y vs. z

◆ momentum()

LocalVector LocalTrajectoryParameters::momentum ( ) const
inline

Momentum vector in the local frame.

Definition at line 88 of file LocalTrajectoryParameters.h.

References funct::abs(), MillePedeFileConverter_cfg::e, f, dqmMemoryStats::float, multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, mathSSE::sqrt(), theDxdz, theDydz, thePzSign, and theQbp.

Referenced by SiPixelTrackResidualSource::analyze(), SiPixelErrorEstimation::analyze(), TwoBodyDecayTrajectory::constructSingleTsosWithErrors(), SiPixelHitEfficiencyModule::fill(), L2TauNNProducer::impactParameter(), JacobianCartesianToLocal::JacobianCartesianToLocal(), BasicTrajectoryState::localMomentum(), StripCPEgeometric::localParameters(), PrintRecoObjects::print(), SeedProducerFromSoAT< TrackerTraits >::produce(), PixelTrackProducerFromSoAT< TrackerTraits >::produce(), and TrackInfoProducerAlgorithm::run().

88  {
89  float op = std::abs(theQbp);
90  if (op < 1.e-9f)
91  op = 1.e-9f;
92  float pz = float(thePzSign) / (op * std::sqrt(1.f + theDxdz * theDxdz + theDydz * theDydz));
93  float px = pz * theDxdz;
94  float py = pz * theDydz;
95  return LocalVector(px, py, pz);
96  }
Local3DVector LocalVector
Definition: LocalVector.h:12
float theQbp
q/p (charged) or 1/p (neutral)
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
float theDxdz
tangent of direction in local x vs. z
float theDydz
tangent of direction in local y vs. z

◆ position()

LocalPoint LocalTrajectoryParameters::position ( ) const
inline

◆ pzSign()

float LocalTrajectoryParameters::pzSign ( ) const
inline

◆ qbp()

float LocalTrajectoryParameters::qbp ( ) const
inline

Definition at line 158 of file LocalTrajectoryParameters.h.

References theQbp.

Referenced by KFTrajectoryFitter::fitOne(), and JacobianLocalToCartesian::JacobianLocalToCartesian().

158 { return theQbp; }
float theQbp
q/p (charged) or 1/p (neutral)

◆ signedInverseMomentum()

float LocalTrajectoryParameters::signedInverseMomentum ( ) const
inline

◆ updateP()

bool LocalTrajectoryParameters::updateP ( float  dP)
inline

Update of momentum by a scalar dP.

Definition at line 149 of file LocalTrajectoryParameters.h.

References funct::abs(), f, AlCaHLTBitMon_ParallelJobs::p, and theQbp.

Referenced by GsfMaterialEffectsUpdator::updateState(), and MaterialEffectsUpdator::updateStateInPlace().

149  {
150  float p = 1.f / std::abs(theQbp);
151  if ((p += dP) <= 0.f)
152  return false;
153  float newQbp = theQbp > 0. ? 1.f / p : -1.f / p;
154  theQbp = newQbp;
155  return true;
156  }
float theQbp
q/p (charged) or 1/p (neutral)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]

◆ vector()

AlgebraicVector5 LocalTrajectoryParameters::vector ( ) const
inline

Vector of parameters with signed inverse momentum.

Vector of parameters as defined above, with the first element = q/p .

Definition at line 120 of file LocalTrajectoryParameters.h.

References signedInverseMomentum(), theDxdz, theDydz, theX, theY, and findQualityFiles::v.

Referenced by GlobalTrackerMuonAlignment::analyzeTrackTrack(), GlobalTrackerMuonAlignment::analyzeTrackTrajectory(), OverlapValidation::analyzeTrajectory(), SiTrackerMultiRecHitUpdator::calcParameters(), CollinearFitAtTM2::CollinearFitAtTM2(), TrajectoryStateCombiner::combine(), GsfTrackProducerBase::computeModeAtTM(), PFGsfHelper::computeQpMode(), SiTrackerMultiRecHitUpdator::ComputeWeight(), GlobalTrackerMuonAlignment::debugTrajectorySOS(), GlobalTrackerMuonAlignment::debugTrajectorySOSv(), MRHChi2MeasurementEstimator::estimate(), CollinearFitAtTM::fit(), HitResol::getPairParameters(), StripCPEgeometric::localParameters(), StripCPEfromTrackAngle::localParameters(), GsfTrackProducerBase::localParametersFromQpMode(), GlobalMuonTrackMatcher::match_Chi2(), TSGForOIFromL2::match_Chi2(), SeedMatcher::matchSimTrack(), MeasurementExtractor::measuredParameters(), GlobalTrackerMuonAlignment::misalignMuonL(), TRecHit5DParamConstraint::parameters(), TrackAssociatorByPositionImpl::quality(), ChargeSignificanceTrajectoryFilter::TBC(), Strip1DMeasurementTransformator::trajectoryParameters(), Tsos2DPhi::Tsos2DPhi(), Tsos2DZed::Tsos2DZed(), and Tsos4D::Tsos4D().

120  {
122  v[0] = signedInverseMomentum();
123  v[1] = theDxdz;
124  v[2] = theDydz;
125  v[3] = theX;
126  v[4] = theY;
127  return v;
128  }
float signedInverseMomentum() const
Signed inverse momentum q/p (zero for neutrals).
ROOT::Math::SVector< double, 5 > AlgebraicVector5
float theDxdz
tangent of direction in local x vs. z
float theDydz
tangent of direction in local y vs. z

Member Data Documentation

◆ theCharge

short LocalTrajectoryParameters::theCharge
private

charge

Definition at line 172 of file LocalTrajectoryParameters.h.

Referenced by charge(), and LocalTrajectoryParameters().

◆ theDxdz

float LocalTrajectoryParameters::theDxdz
private

tangent of direction in local x vs. z

Definition at line 165 of file LocalTrajectoryParameters.h.

Referenced by absdz(), direction(), directionNotNormalized(), dxdz(), LocalTrajectoryParameters(), mixedFormatVector(), momentum(), and vector().

◆ theDydz

float LocalTrajectoryParameters::theDydz
private

tangent of direction in local y vs. z

Definition at line 166 of file LocalTrajectoryParameters.h.

Referenced by absdz(), direction(), directionNotNormalized(), dydz(), LocalTrajectoryParameters(), mixedFormatVector(), momentum(), and vector().

◆ thePzSign

short LocalTrajectoryParameters::thePzSign
private

sign of local pz

Definition at line 170 of file LocalTrajectoryParameters.h.

Referenced by direction(), LocalTrajectoryParameters(), momentum(), and pzSign().

◆ theQbp

float LocalTrajectoryParameters::theQbp
private

q/p (charged) or 1/p (neutral)

Definition at line 164 of file LocalTrajectoryParameters.h.

Referenced by LocalTrajectoryParameters(), mixedFormatVector(), momentum(), qbp(), signedInverseMomentum(), and updateP().

◆ theX

float LocalTrajectoryParameters::theX
private

local x position

Definition at line 167 of file LocalTrajectoryParameters.h.

Referenced by LocalTrajectoryParameters(), mixedFormatVector(), position(), and vector().

◆ theY

float LocalTrajectoryParameters::theY
private

local y position

Definition at line 168 of file LocalTrajectoryParameters.h.

Referenced by LocalTrajectoryParameters(), mixedFormatVector(), position(), and vector().