CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackAndECALLinker.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 
15  const reco::PFBlockElement* ) const override;
16 
17  double testLink( const reco::PFBlockElement*,
18  const reco::PFBlockElement* ) const override;
19 
20 private:
21  const bool _useKDTree,_debug;
22 };
23 
26  "TrackAndECALLinker");
27 
30  const reco::PFBlockElement* elem2 ) const {
31  bool result = false;
32  switch( elem1->type() ){
34  result = (elem1->isMultilinksValide() && elem1->getMultilinks().size() > 0);
35  break;
37  result = (elem2->isMultilinksValide() && elem2->getMultilinks().size() > 0);
38  default:
39  break;
40  }
41  return (_useKDTree ? result : true);
42 }
43 
46  const reco::PFBlockElement* elem2 ) const {
49  const reco::PFBlockElementCluster *ecalelem(NULL);
50  const reco::PFBlockElementTrack *tkelem(NULL);
51  double dist(-1.0);
52  if( elem1->type() < elem2->type() ) {
53  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem1);
54  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
55  } else {
56  tkelem = static_cast<const reco::PFBlockElementTrack*>(elem2);
57  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
58  }
59  const reco::PFRecTrackRef& trackref = tkelem->trackRefPF();
60  const reco::PFClusterRef& clusterref = ecalelem->clusterRef();
61  const reco::PFCluster::REPPoint& ecalreppos = clusterref->positionREP();
62  const reco::PFTrajectoryPoint& tkAtECAL =
63  trackref->extrapolatedPoint( ECALShowerMax );
64  const reco::PFCluster::REPPoint& tkreppos = tkAtECAL.positionREP();
65 
66  // Check if the linking has been done using the KDTree algo
67  // Glowinski & Gouzevitch
68  if ( _useKDTree && tkelem->isMultilinksValide() ) { //KDTree Algo
69  const reco::PFMultilinksType& multilinks = tkelem->getMultilinks();
70  const double ecalphi = ecalreppos.Phi();
71  const double ecaleta = ecalreppos.Eta();
72 
73  // Check if the link Track/Ecal exist
74  reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
75  for (; mlit != multilinks.end(); ++mlit)
76  if ((mlit->first == ecalphi) && (mlit->second == ecaleta))
77  break;
78 
79  // If the link exist, we fill dist and linktest.
80  if (mlit != multilinks.end()){
81  dist = LinkByRecHit::computeDist(ecaleta, ecalphi,
82  tkreppos.Eta(), tkreppos.Phi());
83  }
84 
85  } else {// Old algorithm
86  if ( tkAtECAL.isValid() )
88  *clusterref,
89  false, _debug );
90  }
91 
92  if ( _debug ) {
93  if( dist > 0. ) {
94  std::cout << " Here a link has been established"
95  << " between a track an Ecal with dist "
96  << dist << std::endl;
97  } else
98  std::cout << " No link found " << std::endl;
99  }
100  return dist;
101 }
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
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...)
static double computeDist(double eta1, double phi1, double eta2, double phi2, bool etaPhi=true)
computes a chisquare
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< Double32_t > > REPPoint
Definition: PFCluster.h:51
Type type() const
#define NULL
Definition: scimark2.h:8
const PFClusterRef & clusterRef() const
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
#define constexpr
tuple result
Definition: query.py:137
tuple conf
Definition: dbtoconf.py:185
bool isValid() const
is this point valid ?
bool isMultilinksValide() const
bool linkPrefilter(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
const PFRecTrackRef & trackRefPF() const
tuple cout
Definition: gather_cfg.py:121
#define DEFINE_EDM_PLUGIN(factory, type, name)
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
volatile std::atomic< bool > shutdown_flag false
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:10
TrackAndECALLinker(const edm::ParameterSet &conf)
const PFMultilinksType & getMultilinks() const