CMS 3D CMS Logo

TrackAndHCALLinker.cc
Go to the documentation of this file.
9 
11 public:
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  }
31 
32  double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
33 
34 private:
35  bool useKDTree_;
40  bool debug_;
41  bool checkExit_;
42 };
43 
45 
47  const reco::PFBlockElementCluster* hcalelem(nullptr);
48  const reco::PFBlockElementTrack* tkelem(nullptr);
49  double dist(-1.0);
50  if (elem1->type() < elem2->type()) {
51  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem1);
52  hcalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
53  } else {
54  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem2);
55  hcalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
56  }
57  const reco::PFRecTrackRef& trackref = tkelem->trackRefPF();
58  const reco::PFClusterRef& clusterref = hcalelem->clusterRef();
59  const reco::PFCluster::REPPoint& hcalreppos = clusterref->positionREP();
60  const reco::PFTrajectoryPoint& tkAtHCALEnt = trackref->extrapolatedPoint(trajectoryLayerEntrance_);
61  // Check exit point
62  double dHEta = 0.;
63  double dHPhi = 0.;
64  double dRHCALEx = 0.;
65  if (checkExit_) {
66  const reco::PFTrajectoryPoint& tkAtHCALEx = trackref->extrapolatedPoint(trajectoryLayerExit_);
67  dHEta = (tkAtHCALEx.positionREP().Eta() - tkAtHCALEnt.positionREP().Eta());
68  dHPhi = reco::deltaPhi(tkAtHCALEx.positionREP().Phi(), tkAtHCALEnt.positionREP().Phi());
69  dRHCALEx = tkAtHCALEx.position().R();
70  }
71  if (useKDTree_ && hcalelem->isMultilinksValide()) { //KDTree Algo
72  const reco::PFMultilinksType& multilinks = hcalelem->getMultilinks();
73  const double tracketa = tkAtHCALEnt.positionREP().Eta();
74  const double trackphi = tkAtHCALEnt.positionREP().Phi();
75 
76  // Check if the link Track/Hcal exist
77  reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
78  for (; mlit != multilinks.end(); ++mlit)
79  if ((mlit->first == trackphi) && (mlit->second == tracketa))
80  break;
81 
82  // If the link exist, we fill dist and linktest.
83 
84  if (mlit != multilinks.end()) {
85  // when checkExit_ is false
86  if (!checkExit_) {
87  dist = LinkByRecHit::computeDist(hcalreppos.Eta(), hcalreppos.Phi(), tracketa, trackphi);
88  }
89  // when checkExit_ is true
90  else {
91  //special case ! A looper can exit the barrel inwards and hit the endcap
92  //In this case calculate the distance based on the first crossing since
93  //the looper will probably never make it to the endcap
94  if (dRHCALEx < tkAtHCALEnt.position().R()) {
95  dist = LinkByRecHit::computeDist(hcalreppos.Eta(), hcalreppos.Phi(), tracketa, trackphi);
96  edm::LogWarning("TrackHCALLinker ")
97  << "Special case of linking with track hitting HCAL and looping back in the tracker ";
98  } else {
100  hcalreppos.Eta(), hcalreppos.Phi(), tracketa + 0.1 * dHEta, trackphi + 0.1 * dHPhi);
101  }
102  } // checkExit_
103  } // multilinks
104 
105  } else { // Old algorithm
106  if (tkAtHCALEnt.isValid())
107  dist = LinkByRecHit::testTrackAndClusterByRecHit(*trackref, *clusterref, false, debug_);
108  }
109  return dist;
110 }
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
std::vector< std::pair< double, double > > PFMultilinksType
Abstract This class is used by the KDTree Track / Ecal Cluster linker to store all found links...
Abstract base class for a PFBlock element (track, cluster...)
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
static double computeDist(double eta1, double phi1, double eta2, double phi2, bool etaPhi=true)
computes a chisquare
const PFClusterRef & clusterRef() const override
const math::XYZPoint & position() const
cartesian position (x, y, z)
Type type() const
static LayerType layerTypeByName(const std::string &name)
bool isValid() const
is this point valid ?
TrackAndHCALLinker(const edm::ParameterSet &conf)
std::string trajectoryLayerEntranceString_
bool isMultilinksValide() const
const PFRecTrackRef & trackRefPF() const override
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:52
#define DEFINE_EDM_PLUGIN(factory, type, name)
reco::PFTrajectoryPoint::LayerType trajectoryLayerExit_
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
reco::PFTrajectoryPoint::LayerType trajectoryLayerEntrance_
std::string trajectoryLayerExitString_
LayerType
Define the different layers where the track can be propagated.
static double testTrackAndClusterByRecHit(const reco::PFRecTrack &track, const reco::PFCluster &cluster, bool isBrem=false, bool debug=false)
Definition: LinkByRecHit.cc:16
const PFMultilinksType & getMultilinks() const