CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Attributes
TrackAndHCALLinker Class Reference
Inheritance diagram for TrackAndHCALLinker:
BlockElementLinkerBase

Public Member Functions

bool linkPrefilter (const reco::PFBlockElement *, const reco::PFBlockElement *) const override
 
double testLink (const reco::PFBlockElement *, const reco::PFBlockElement *) const override
 
 TrackAndHCALLinker (const edm::ParameterSet &conf)
 
- Public Member Functions inherited from BlockElementLinkerBase
 BlockElementLinkerBase (const edm::ParameterSet &conf)
 
 BlockElementLinkerBase (const BlockElementLinkerBase &)=delete
 
const std::string & name () const
 
BlockElementLinkerBaseoperator= (const BlockElementLinkerBase &)=delete
 
virtual ~BlockElementLinkerBase ()=default
 

Private Attributes

bool checkExit_
 
bool debug_
 
reco::PFTrajectoryPoint::LayerType trajectoryLayerEntrance_
 
std::string trajectoryLayerEntranceString_
 
reco::PFTrajectoryPoint::LayerType trajectoryLayerExit_
 
std::string trajectoryLayerExitString_
 
bool useKDTree_
 

Detailed Description

Definition at line 10 of file TrackAndHCALLinker.cc.

Constructor & Destructor Documentation

TrackAndHCALLinker::TrackAndHCALLinker ( const edm::ParameterSet conf)
inline

Definition at line 12 of file TrackAndHCALLinker.cc.

References cms::cuda::assert(), checkExit_, reco::PFTrajectoryPoint::HCALEntrance, reco::PFTrajectoryPoint::HCALExit, reco::PFTrajectoryPoint::layerTypeByName(), trajectoryLayerEntrance_, trajectoryLayerEntranceString_, trajectoryLayerExit_, trajectoryLayerExitString_, reco::PFTrajectoryPoint::Unknown, and reco::PFTrajectoryPoint::VFcalEntrance.

13  : BlockElementLinkerBase(conf),
14  useKDTree_(conf.getParameter<bool>("useKDTree")),
15  trajectoryLayerEntranceString_(conf.getParameter<std::string>("trajectoryLayerEntrance")),
16  trajectoryLayerExitString_(conf.getParameter<std::string>("trajectoryLayerExit")),
17  debug_(conf.getUntrackedParameter<bool>("debug", false)) {
18  // convert TrajectoryLayers info from string to enum
21  // make sure the requested setting is supported
28  // flag if exit layer should be checked or not
30  }
T getUntrackedParameter(std::string const &, T const &) const
assert(be >=bs)
static LayerType layerTypeByName(const std::string &name)
std::string trajectoryLayerEntranceString_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
BlockElementLinkerBase(const edm::ParameterSet &conf)
reco::PFTrajectoryPoint::LayerType trajectoryLayerExit_
reco::PFTrajectoryPoint::LayerType trajectoryLayerEntrance_
std::string trajectoryLayerExitString_

Member Function Documentation

bool TrackAndHCALLinker::linkPrefilter ( const reco::PFBlockElement elem1,
const reco::PFBlockElement elem2 
) const
overridevirtual

Reimplemented from BlockElementLinkerBase.

Definition at line 48 of file TrackAndHCALLinker.cc.

References relativeConstraints::empty, reco::PFBlockElement::getMultilinks(), reco::PFBlockElement::HCAL, reco::PFBlockElement::isMultilinksValide(), mps_fire::result, reco::PFBlockElement::TRACK, reco::PFBlockElement::type(), and useKDTree_.

48  {
49  bool result = false;
50  // Track-HCAL KDTree multilinks are stored to track's elem
51  switch (elem1->type()) {
53  result = (elem1->isMultilinksValide(elem2->type()) && !elem1->getMultilinks(elem2->type()).empty() &&
54  elem2->isMultilinksValide(elem1->type()));
55  break;
57  result = (elem2->isMultilinksValide(elem1->type()) && !elem2->getMultilinks(elem1->type()).empty() &&
58  elem1->isMultilinksValide(elem2->type()));
59  default:
60  break;
61  }
62  return (useKDTree_ ? result : true);
63 }
Type type() const
bool isMultilinksValide(Type type) const
tuple result
Definition: mps_fire.py:311
const PFMultilinksType & getMultilinks(Type type) const
double TrackAndHCALLinker::testLink ( const reco::PFBlockElement elem1,
const reco::PFBlockElement elem2 
) const
overridevirtual

Implements BlockElementLinkerBase.

Definition at line 65 of file TrackAndHCALLinker.cc.

References checkExit_, reco::PFBlockElementCluster::clusterRef(), LinkByRecHit::computeDist(), debug_, reco::deltaPhi(), reco::PFBlockElement::getMultilinks(), reco::PFBlockElement::isMultilinksValide(), reco::PFTrajectoryPoint::isValid(), reco::PFTrajectoryPoint::position(), reco::PFTrajectoryPoint::positionREP(), LinkByRecHit::testTrackAndClusterByRecHit(), reco::PFBlockElementTrack::trackRefPF(), trajectoryLayerEntrance_, trajectoryLayerExit_, reco::PFBlockElement::type(), and useKDTree_.

65  {
66  const reco::PFBlockElementCluster* hcalelem(nullptr);
67  const reco::PFBlockElementTrack* tkelem(nullptr);
68  double dist(-1.0);
69  if (elem1->type() < elem2->type()) {
70  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem1);
71  hcalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
72  } else {
73  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem2);
74  hcalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
75  }
76  const reco::PFRecTrackRef& trackref = tkelem->trackRefPF();
77  const reco::PFClusterRef& clusterref = hcalelem->clusterRef();
78  const reco::PFCluster::REPPoint& hcalreppos = clusterref->positionREP();
79  const reco::PFTrajectoryPoint& tkAtHCALEnt = trackref->extrapolatedPoint(trajectoryLayerEntrance_);
80  const reco::PFCluster::REPPoint& tkreppos = tkAtHCALEnt.positionREP();
81  // Check exit point
82  double dHEta = 0.;
83  double dHPhi = 0.;
84  double dRHCALEx = 0.;
85  if (checkExit_) {
86  const reco::PFTrajectoryPoint& tkAtHCALEx = trackref->extrapolatedPoint(trajectoryLayerExit_);
87  dHEta = (tkAtHCALEx.positionREP().Eta() - tkAtHCALEnt.positionREP().Eta());
88  dHPhi = reco::deltaPhi(tkAtHCALEx.positionREP().Phi(), tkAtHCALEnt.positionREP().Phi());
89  dRHCALEx = tkAtHCALEx.position().R();
90  }
91 
92  // Check if the linking has been done using the KDTree algo
93  // Glowinski & Gouzevitch
94  if (useKDTree_ && tkelem->isMultilinksValide(hcalelem->type())) { //KDTree Algo
95  const reco::PFMultilinksType& multilinks = tkelem->getMultilinks(hcalelem->type());
96 
97  // Check if the link Track/Hcal exist
98  reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
99  for (; mlit != multilinks.end(); ++mlit)
100  if (mlit->clusterRef == clusterref)
101  break;
102 
103  // If the link exist, we fill dist and linktest.
104 
105  if (mlit != multilinks.end()) {
106  // when checkExit_ is false
107  if (!checkExit_) {
108  dist = LinkByRecHit::computeDist(hcalreppos.Eta(), hcalreppos.Phi(), tkreppos.Eta(), tkreppos.Phi());
109  }
110  // when checkExit_ is true
111  else {
112  //special case ! A looper can exit the barrel inwards and hit the endcap
113  //In this case calculate the distance based on the first crossing since
114  //the looper will probably never make it to the endcap
115  if (dRHCALEx < tkAtHCALEnt.position().R()) {
116  dist = LinkByRecHit::computeDist(hcalreppos.Eta(), hcalreppos.Phi(), tkreppos.Eta(), tkreppos.Phi());
117  edm::LogWarning("TrackHCALLinker ")
118  << "Special case of linking with track hitting HCAL and looping back in the tracker ";
119  } else {
121  hcalreppos.Eta(), hcalreppos.Phi(), tkreppos.Eta() + 0.1 * dHEta, tkreppos.Phi() + 0.1 * dHPhi);
122  }
123  } // checkExit_
124  } // multilinks
125 
126  } else { // Old algorithm
127  if (tkAtHCALEnt.isValid())
128  dist = LinkByRecHit::testTrackAndClusterByRecHit(*trackref, *clusterref, false, debug_);
129  }
130  return dist;
131 }
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
static double computeDist(double eta1, double phi1, double eta2, double phi2, bool etaPhi=true)
computes a chisquare
const math::XYZPoint & position() const
cartesian position (x, y, z)
Type type() const
std::vector< PFMultilink > PFMultilinksType
collection of PFSuperCluster objects
bool isValid() const
is this point valid ?
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:48
reco::PFTrajectoryPoint::LayerType trajectoryLayerExit_
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
reco::PFTrajectoryPoint::LayerType trajectoryLayerEntrance_
Log< level::Warning, false > LogWarning
static double testTrackAndClusterByRecHit(const reco::PFRecTrack &track, const reco::PFCluster &cluster, bool isBrem=false, bool debug=false)
Definition: LinkByRecHit.cc:18

Member Data Documentation

bool TrackAndHCALLinker::checkExit_
private

Definition at line 43 of file TrackAndHCALLinker.cc.

Referenced by testLink(), and TrackAndHCALLinker().

bool TrackAndHCALLinker::debug_
private

Definition at line 42 of file TrackAndHCALLinker.cc.

Referenced by testLink().

reco::PFTrajectoryPoint::LayerType TrackAndHCALLinker::trajectoryLayerEntrance_
private

Definition at line 40 of file TrackAndHCALLinker.cc.

Referenced by testLink(), and TrackAndHCALLinker().

std::string TrackAndHCALLinker::trajectoryLayerEntranceString_
private

Definition at line 38 of file TrackAndHCALLinker.cc.

Referenced by TrackAndHCALLinker().

reco::PFTrajectoryPoint::LayerType TrackAndHCALLinker::trajectoryLayerExit_
private

Definition at line 41 of file TrackAndHCALLinker.cc.

Referenced by testLink(), and TrackAndHCALLinker().

std::string TrackAndHCALLinker::trajectoryLayerExitString_
private

Definition at line 39 of file TrackAndHCALLinker.cc.

Referenced by TrackAndHCALLinker().

bool TrackAndHCALLinker::useKDTree_
private

Definition at line 37 of file TrackAndHCALLinker.cc.

Referenced by linkPrefilter(), and testLink().