CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Attributes | Friends
reco::PFTrajectoryPoint Class Reference

A PFTrack holds several trajectory points, which basically contain the position and momentum of a track at a given position. More...

#include <PFTrajectoryPoint.h>

Public Types

enum  LayerType {
  ClosestApproach = 0, BeamPipeOrEndVertex = 1, PS1 = 2, PS2 = 3,
  ECALEntrance = 4, ECALShowerMax = 5, HCALEntrance = 6, HCALExit = 7,
  HOLayer = 8, VFcalEntrance = 9, NLayers = 10, Unknown = -1
}
 Define the different layers where the track can be propagated. More...
 
typedef ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
 

Public Member Functions

void calculatePositionREP ()
 calculate posrep_ once and for all More...
 
int detId () const
 measurement detId More...
 
bool isTrackerLayer () const
 is this point corresponding to an intersection with a tracker layer ? More...
 
bool isValid () const
 is this point valid ? More...
 
int layer () const
 trajectory point layer More...
 
const math::XYZTLorentzVectormomentum () const
 4-momenta quadrivector More...
 
bool operator== (const reco::PFTrajectoryPoint &other) const
 
 PFTrajectoryPoint ()
 default constructor. Set variables at default dummy values More...
 
 PFTrajectoryPoint (int detId, int layer, const math::XYZPoint &posxyz, const math::XYZTLorentzVector &momentum)
 constructor from values. set detId to -1 if this point is not from a tracker layer More...
 
const math::XYZPointposition () const
 cartesian position (x, y, z) More...
 
const REPPointpositionREP () const
 trajectory position in (rho, eta, phi) base More...
 

Static Public Member Functions

static LayerType layerTypeByName (const std::string &name)
 

Static Public Attributes

static const std::array< std::string, NLayerslayerTypeNames
 

Private Attributes

int detId_
 detid if measurement is corresponding to a tracker layer More...
 
bool isTrackerLayer_
 Is the measurement corresponding to a tracker layer? or was it obtained by propagating the track to a certain position? More...
 
int layer_
 propagated layer More...
 
math::XYZTLorentzVector momentum_
 momentum quadrivector More...
 
REPPoint posrep_
 position in (rho, eta, phi) base (transient) More...
 
math::XYZPoint posxyz_
 cartesian position (x, y, z) More...
 

Friends

std::ostream & operator<< (std::ostream &out, const reco::PFTrajectoryPoint &trajPoint)
 

Detailed Description

A PFTrack holds several trajectory points, which basically contain the position and momentum of a track at a given position.

Author
Renaud Bruneliere
Date
July 2006

Definition at line 26 of file PFTrajectoryPoint.h.

Member Typedef Documentation

◆ REPPoint

typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<double> > reco::PFTrajectoryPoint::REPPoint

Definition at line 31 of file PFTrajectoryPoint.h.

Member Enumeration Documentation

◆ LayerType

Define the different layers where the track can be propagated.

Enumerator
ClosestApproach 

Point of closest approach from beam axis (initial point in the case of PFSimParticle)

BeamPipeOrEndVertex 
PS1 

Preshower layer 1.

PS2 

Preshower layer 2.

ECALEntrance 

ECAL front face.

ECALShowerMax 

expected maximum of the shower in ECAL, for an e/gamma particle

HCALEntrance 

HCAL front face.

HCALExit 

HCAL exit.

HOLayer 

HO layer.

VFcalEntrance 

VFcal(HF) front face.

NLayers 
Unknown 

Definition at line 34 of file PFTrajectoryPoint.h.

34  {
36  ClosestApproach = 0,
39  PS1 = 2,
41  PS2 = 3,
43  ECALEntrance = 4,
46  ECALShowerMax = 5,
48  HCALEntrance = 6,
50  HCALExit = 7,
52  HOLayer = 8,
54  VFcalEntrance = 9,
55  // Number of valid layers
56  NLayers = 10,
57  // Unknown
58  Unknown = -1
59  };
Point of closest approach from beam axis (initial point in the case of PFSimParticle) ...

Constructor & Destructor Documentation

◆ PFTrajectoryPoint() [1/2]

PFTrajectoryPoint::PFTrajectoryPoint ( )

default constructor. Set variables at default dummy values

Definition at line 21 of file PFTrajectoryPoint.cc.

21 : isTrackerLayer_(false), detId_(-1), layer_(-1) {}
int detId_
detid if measurement is corresponding to a tracker layer
bool isTrackerLayer_
Is the measurement corresponding to a tracker layer? or was it obtained by propagating the track to a...
int layer_
propagated layer

◆ PFTrajectoryPoint() [2/2]

PFTrajectoryPoint::PFTrajectoryPoint ( int  detId,
int  layer,
const math::XYZPoint posxyz,
const math::XYZTLorentzVector momentum 
)

constructor from values. set detId to -1 if this point is not from a tracker layer

Definition at line 23 of file PFTrajectoryPoint.cc.

References detId(), isTrackerLayer_, posrep_, and posxyz_.

28  if (detId)
29  isTrackerLayer_ = true;
30  posrep_.SetCoordinates(posxyz_.Rho(), posxyz_.Eta(), posxyz_.Phi());
31 }
int detId_
detid if measurement is corresponding to a tracker layer
const math::XYZTLorentzVector & momentum() const
4-momenta quadrivector
int layer() const
trajectory point layer
math::XYZPoint posxyz_
cartesian position (x, y, z)
bool isTrackerLayer_
Is the measurement corresponding to a tracker layer? or was it obtained by propagating the track to a...
int layer_
propagated layer
math::XYZTLorentzVector momentum_
momentum quadrivector
int detId() const
measurement detId
REPPoint posrep_
position in (rho, eta, phi) base (transient)

Member Function Documentation

◆ calculatePositionREP()

void reco::PFTrajectoryPoint::calculatePositionREP ( )
inline

calculate posrep_ once and for all

Definition at line 100 of file PFTrajectoryPoint.h.

References posrep_, and posxyz_.

100 { posrep_.SetCoordinates(posxyz_.Rho(), posxyz_.Eta(), posxyz_.Phi()); }
math::XYZPoint posxyz_
cartesian position (x, y, z)
REPPoint posrep_
position in (rho, eta, phi) base (transient)

◆ detId()

int reco::PFTrajectoryPoint::detId ( ) const
inline

measurement detId

Definition at line 72 of file PFTrajectoryPoint.h.

References detId_.

Referenced by reco::operator<<(), and PFTrajectoryPoint().

72 { return detId_; }
int detId_
detid if measurement is corresponding to a tracker layer

◆ isTrackerLayer()

bool reco::PFTrajectoryPoint::isTrackerLayer ( ) const
inline

is this point corresponding to an intersection with a tracker layer ?

Definition at line 86 of file PFTrajectoryPoint.h.

References detId_.

Referenced by reco::PFTrack::addPoint().

86  {
87  if (detId_ >= 0)
88  return true;
89  else
90  return false;
91  }
int detId_
detid if measurement is corresponding to a tracker layer

◆ isValid()

bool reco::PFTrajectoryPoint::isValid ( void  ) const
inline

◆ layer()

int reco::PFTrajectoryPoint::layer ( ) const
inline

trajectory point layer

Definition at line 75 of file PFTrajectoryPoint.h.

References layer_.

Referenced by geometryXMLparser.DTAlignable::index(), geometryXMLparser.CSCAlignable::index(), and reco::operator<<().

75 { return layer_; }
int layer_
propagated layer

◆ layerTypeByName()

PFTrajectoryPoint::LayerType PFTrajectoryPoint::layerTypeByName ( const std::string &  name)
static

Definition at line 33 of file PFTrajectoryPoint.cc.

References HLT_2024v13_cff::distance, spr::find(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, layerTypeNames, Skims_PA_cff::name, and Unknown.

Referenced by KDTreeLinkerTrackHcal::KDTreeLinkerTrackHcal(), and TrackAndHCALLinker::TrackAndHCALLinker().

33  {
34  auto it = std::find(layerTypeNames.begin(), layerTypeNames.end(), name);
35  if (it == layerTypeNames.end()) {
36  return Unknown; // better this or throw()?
37  }
38  int index = std::distance(layerTypeNames.begin(), it);
39  return LayerType(index);
40 }
static const std::array< std::string, NLayers > layerTypeNames
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
LayerType
Define the different layers where the track can be propagated.

◆ momentum()

const math::XYZTLorentzVector& reco::PFTrajectoryPoint::momentum ( ) const
inline

4-momenta quadrivector

Definition at line 103 of file PFTrajectoryPoint.h.

References momentum_.

Referenced by reco::operator<<(), ConvBremPFTrackFinder::runConvBremFinder(), KDTreeLinkerTrackEcal::searchLinks(), and LinkByRecHit::testTrackAndClusterByRecHit().

103 { return momentum_; }
math::XYZTLorentzVector momentum_
momentum quadrivector

◆ operator==()

bool reco::PFTrajectoryPoint::operator== ( const reco::PFTrajectoryPoint other) const

◆ position()

const math::XYZPoint& reco::PFTrajectoryPoint::position ( ) const
inline

◆ positionREP()

const REPPoint& reco::PFTrajectoryPoint::positionREP ( ) const
inline

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  out,
const reco::PFTrajectoryPoint trajPoint 
)
friend

Member Data Documentation

◆ detId_

int reco::PFTrajectoryPoint::detId_
private

detid if measurement is corresponding to a tracker layer

Definition at line 115 of file PFTrajectoryPoint.h.

Referenced by detId(), isTrackerLayer(), and isValid().

◆ isTrackerLayer_

bool reco::PFTrajectoryPoint::isTrackerLayer_
private

Is the measurement corresponding to a tracker layer? or was it obtained by propagating the track to a certain position?

Definition at line 112 of file PFTrajectoryPoint.h.

Referenced by PFTrajectoryPoint().

◆ layer_

int reco::PFTrajectoryPoint::layer_
private

propagated layer

Definition at line 118 of file PFTrajectoryPoint.h.

Referenced by isValid(), and layer().

◆ layerTypeNames

std::array< std::string, PFTrajectoryPoint::NLayers > const PFTrajectoryPoint::layerTypeNames
static
Initial value:
{{"ClosestApproach",
"BeamPipeOrEndVertex",
"PS1",
"PS2",
"ECALEntrance",
"ECALShowerMax",
"HCALEntrance",
"HCALExit",
"HOLayer",
"VFcalEntrance"}}

Definition at line 61 of file PFTrajectoryPoint.h.

Referenced by layerTypeByName().

◆ momentum_

math::XYZTLorentzVector reco::PFTrajectoryPoint::momentum_
private

momentum quadrivector

Definition at line 127 of file PFTrajectoryPoint.h.

Referenced by momentum().

◆ posrep_

REPPoint reco::PFTrajectoryPoint::posrep_
private

position in (rho, eta, phi) base (transient)

Definition at line 124 of file PFTrajectoryPoint.h.

Referenced by calculatePositionREP(), PFTrajectoryPoint(), and positionREP().

◆ posxyz_

math::XYZPoint reco::PFTrajectoryPoint::posxyz_
private

cartesian position (x, y, z)

Definition at line 121 of file PFTrajectoryPoint.h.

Referenced by calculatePositionREP(), PFTrajectoryPoint(), and position().