CMS 3D CMS Logo

GSFAndHGCalLinker.cc
Go to the documentation of this file.
6 
8 public:
11  _useKDTree(conf.getParameter<bool>("useKDTree")),
12  _debug(conf.getUntrackedParameter<bool>("debug",false)) {}
13 
14  double testLink
15  ( const reco::PFBlockElement*,
16  const reco::PFBlockElement* ) const override;
17 
18 private:
20 };
21 
24  "GSFAndHGCalLinker");
25 
27  ( const reco::PFBlockElement* elem1,
28  const reco::PFBlockElement* elem2) const {
31  const reco::PFBlockElementCluster *hgcalelem(nullptr);
32  const reco::PFBlockElementGsfTrack *gsfelem(nullptr);
33  double dist(-1.0);
34  if( elem1->type() > elem2->type() ) {
35  hgcalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
36  gsfelem = static_cast<const reco::PFBlockElementGsfTrack*>(elem2);
37  } else {
38  hgcalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
39  gsfelem = static_cast<const reco::PFBlockElementGsfTrack*>(elem1);
40  }
41  const reco::PFRecTrack& track = gsfelem->GsftrackPF();
42  const reco::PFClusterRef& clusterref = hgcalelem->clusterRef();
43  const reco::PFTrajectoryPoint& tkAtECAL =
44  track.extrapolatedPoint( ECALShowerMax );
45  if( tkAtECAL.isValid() ) {
46  dist = LinkByRecHit::computeDist( tkAtECAL.positionREP().eta(),
47  tkAtECAL.positionREP().phi(),
48  clusterref->positionREP().Eta(),
49  clusterref->positionREP().Phi() );
50  if( dist > 0.3 ) dist = -1.0;
51  }
52  if ( _debug ) {
53  if ( dist > 0. ) {
54  std::cout << " Here a link has been established"
55  << " between a GSF track an HGCal with dist "
56  << dist << std::endl;
57  } else {
58  if( _debug ) std::cout << " No link found " << std::endl;
59  }
60  }
61 
62  return dist;
63 }
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
Abstract base class for a PFBlock element (track, cluster...)
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:22
static double computeDist(double eta1, double phi1, double eta2, double phi2, bool etaPhi=true)
computes a chisquare
const PFClusterRef & clusterRef() const override
Type type() const
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
const GsfPFRecTrack & GsftrackPF() const
const reco::PFTrajectoryPoint & extrapolatedPoint(unsigned layerid) const
Definition: PFTrack.cc:76
bool isValid() const
is this point valid ?
GSFAndHGCalLinker(const edm::ParameterSet &conf)
#define DEFINE_EDM_PLUGIN(factory, type, name)
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
LayerType
Define the different layers where the track can be propagated.
#define constexpr