CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PreshowerAndECALLinker.cc
Go to the documentation of this file.
5 
7 public:
10  _useKDTree(conf.getParameter<bool>("useKDTree")),
11  _debug(conf.getUntrackedParameter<bool>("debug",false)) {}
12 
14  const reco::PFBlockElement* ) const override;
15 
16  double testLink( const reco::PFBlockElement*,
17  const reco::PFBlockElement* ) const override;
18 
19 private:
21 };
22 
25  "PreshowerAndECALLinker");
26 
29  const reco::PFBlockElement* elem2 ) const {
30  bool result = false;
31  switch( elem1->type() ){
34  result = ( elem1->isMultilinksValide() &&
35  elem1->getMultilinks().size() > 0 );
36  break;
38  result = ( elem2->isMultilinksValide() &&
39  elem2->getMultilinks().size() > 0 );
40  break;
41  default:
42  break;
43  }
44  return (_useKDTree ? result : true);
45 }
46 
49  const reco::PFBlockElement* elem2) const {
50  const reco::PFBlockElementCluster *pselem(NULL), *ecalelem(NULL);
51  double dist(-1.0);
52  if( elem1->type() < elem2->type() ) {
53  pselem = static_cast<const reco::PFBlockElementCluster*>(elem1);
54  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem2);
55  } else {
56  pselem = static_cast<const reco::PFBlockElementCluster*>(elem2);
57  ecalelem = static_cast<const reco::PFBlockElementCluster*>(elem1);
58  }
59  const reco::PFClusterRef& psref = pselem->clusterRef();
60  const reco::PFClusterRef& ecalref = ecalelem->clusterRef();
61  if( psref.isNull() || ecalref.isNull() ) {
62  throw cms::Exception("BadClusterRefs")
63  << "PFBlockElementCluster's refs are null!";
64  }
65  // Check if the linking has been done using the KDTree algo
66  // Glowinski & Gouzevitch
67  if ( _useKDTree && pselem->isMultilinksValide() ) { // KDTree algo
68  const reco::PFMultilinksType& multilinks = pselem->getMultilinks();
69  const reco::PFCluster::REPPoint& ecalreppos = ecalref->positionREP();
70  const math::XYZPoint& ecalxyzpos = ecalref->position();
71  const math::XYZPoint& psxyzpos = psref->position();
72  const double ecalPhi = ecalreppos.Phi();
73  const double ecalEta = ecalreppos.Eta();
74 
75  // Check if the link PS/Ecal exist
76  reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
77  for (; mlit != multilinks.end(); ++mlit)
78  if ((mlit->first == ecalPhi) && (mlit->second == ecalEta))
79  break;
80 
81  // If the link exist, we fill dist and linktest.
82  if (mlit != multilinks.end()){
83  dist =
84  LinkByRecHit::computeDist(ecalxyzpos.X()/1000.,ecalxyzpos.Y()/1000.,
85  psxyzpos.X()/1000. ,psxyzpos.Y()/1000.,
86  false);
87  }
88  } else { //Old algorithm
89  dist = LinkByRecHit::testECALAndPSByRecHit( *ecalref, *psref ,_debug);
90  }
91  return dist;
92 }
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
static double testECALAndPSByRecHit(const reco::PFCluster &clusterECAL, const reco::PFCluster &clusterPS, bool debug=false)
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
bool linkPrefilter(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
bool isNull() const
Checks for null.
Definition: Ref.h:247
tuple result
Definition: query.py:137
tuple conf
Definition: dbtoconf.py:185
PreshowerAndECALLinker(const edm::ParameterSet &conf)
double testLink(const reco::PFBlockElement *, const reco::PFBlockElement *) const override
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
bool isMultilinksValide() const
#define DEFINE_EDM_PLUGIN(factory, type, name)
volatile std::atomic< bool > shutdown_flag false
const PFMultilinksType & getMultilinks() const