CMS 3D CMS Logo

GSFAndHCALLinker.cc
Go to the documentation of this file.
6 
8 public:
10  : BlockElementLinkerBase(conf),
11  useKDTree_(conf.getParameter<bool>("useKDTree")),
12  debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
13 
14  double testLink(const reco::PFBlockElement*, const reco::PFBlockElement*) const override;
15 
16 private:
18 };
19 
21 
22 double GSFAndHCALLinker::testLink(const reco::PFBlockElement* elem1, const reco::PFBlockElement* elem2) const {
24  const reco::PFBlockElementCluster* hcalelem(nullptr);
25  const reco::PFBlockElementGsfTrack* gsfelem(nullptr);
26  double dist(-1.0);
27  if (elem1->type() < elem2->type()) {
28  hcalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
29  gsfelem = static_cast<const reco::PFBlockElementGsfTrack*>(elem2);
30  } else {
31  hcalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
32  gsfelem = static_cast<const reco::PFBlockElementGsfTrack*>(elem1);
33  }
34  const reco::PFRecTrack& track = gsfelem->GsftrackPF();
35  const reco::PFClusterRef& clusterref = hcalelem->clusterRef();
36  const reco::PFTrajectoryPoint& tkAtHCAL = track.extrapolatedPoint(HCALEnt);
37  if (tkAtHCAL.isValid()) {
38  dist = LinkByRecHit::testTrackAndClusterByRecHit(track, *clusterref, false, debug_);
39  }
40  if (debug_) {
41  if (dist > 0.) {
42  std::cout << " Here a link has been established"
43  << " between a GSF track an Hcal with dist " << dist << std::endl;
44  } else {
45  if (debug_)
46  std::cout << " No link found " << std::endl;
47  }
48  }
49 
50  return dist;
51 }
electrons_cff.bool
bool
Definition: electrons_cff.py:393
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11713
funct::false
false
Definition: Factorize.h:29
PFBlockElementCluster.h
GSFAndHCALLinker::GSFAndHCALLinker
GSFAndHCALLinker(const edm::ParameterSet &conf)
Definition: GSFAndHCALLinker.cc:9
reco::PFBlockElementGsfTrack::GsftrackPF
const GsfPFRecTrack & GsftrackPF() const
Definition: PFBlockElementGsfTrack.h:52
gather_cfg.cout
cout
Definition: gather_cfg.py:144
LinkByRecHit::testTrackAndClusterByRecHit
static double testTrackAndClusterByRecHit(const reco::PFRecTrack &track, const reco::PFCluster &cluster, bool isBrem=false, bool debug=false)
Definition: LinkByRecHit.cc:18
GSFAndHCALLinker::useKDTree_
bool useKDTree_
Definition: GSFAndHCALLinker.cc:17
BlockElementLinkerBase
Definition: BlockElementLinkerBase.h:10
edm::Ref< PFClusterCollection >
PFCluster.h
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
edm::ParameterSet
Definition: ParameterSet.h:47
edmplugin::PluginFactory
Definition: PluginFactory.h:34
LinkByRecHit.h
reco::PFTrajectoryPoint::HCALEntrance
HCAL front face.
Definition: PFTrajectoryPoint.h:48
GSFAndHCALLinker::testLink
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
Definition: GSFAndHCALLinker.cc:22
reco::PFTrajectoryPoint::LayerType
LayerType
Define the different layers where the track can be propagated.
Definition: PFTrajectoryPoint.h:34
GSFAndHCALLinker::debug_
bool debug_
Definition: GSFAndHCALLinker.cc:17
GSFAndHCALLinker
Definition: GSFAndHCALLinker.cc:7
reco::PFTrajectoryPoint::isValid
bool isValid() const
is this point valid ?
Definition: PFTrajectoryPoint.h:84
reco::PFBlockElement
Abstract base class for a PFBlock element (track, cluster...)
Definition: PFBlockElement.h:26
reco::PFBlockElementGsfTrack
Track Element.
Definition: PFBlockElementGsfTrack.h:18
reco::PFRecTrack
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:20
reco::PFTrajectoryPoint
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
Definition: PFTrajectoryPoint.h:26
reco::PFBlockElementCluster
Cluster Element.
Definition: PFBlockElementCluster.h:16
BlockElementLinkerBase.h
reco::PFBlockElementCluster::clusterRef
const PFClusterRef & clusterRef() const override
Definition: PFBlockElementCluster.h:29
reco::PFBlockElement::type
Type type() const
Definition: PFBlockElement.h:69
PFBlockElementGsfTrack.h