CMS 3D CMS Logo

SuperClusterImporter.cc
Go to the documentation of this file.
9 
10 // for single tower H/E
12 
13 //quick pT for superclusters
14 inline double ptFast(const double energy, const math::XYZPoint& position, const math::XYZPoint& origin) {
15  const auto v = position - origin;
16  return energy * std::sqrt(v.perp2() / v.mag2());
17 }
18 
19 #include <unordered_map>
20 
22 public:
24 
25  void updateEventSetup(const edm::EventSetup& es) override;
26 
27  void importToBlock(const edm::Event&, ElementList&) const override;
28 
29 private:
32  const double _maxHoverE, _pTbyPass, _minSCPt;
33  std::unique_ptr<EgammaHadTower> _hadTower;
35  static const math::XYZPoint _zero;
36 };
37 
39 
41 
43  : BlockElementImporterBase(conf, sumes),
44  _srcEB(sumes.consumes<reco::SuperClusterCollection>(conf.getParameter<edm::InputTag>("source_eb"))),
45  _srcEE(sumes.consumes<reco::SuperClusterCollection>(conf.getParameter<edm::InputTag>("source_ee"))),
46  _srcTowers(sumes.consumes<CaloTowerCollection>(conf.getParameter<edm::InputTag>("source_towers"))),
47  _maxHoverE(conf.getParameter<double>("maximumHoverE")),
48  _pTbyPass(conf.getParameter<double>("minPTforBypass")),
49  _minSCPt(conf.getParameter<double>("minSuperClusterPt")),
51  _superClustersArePF(conf.getParameter<bool>("superClustersArePF")) {}
52 
55 }
56 
58  auto eb_scs = e.getHandle(_srcEB);
59  auto ee_scs = e.getHandle(_srcEE);
60  auto towers = e.getHandle(_srcTowers);
61  _hadTower->setTowerCollection(towers.product());
62  elems.reserve(elems.size() + eb_scs->size() + ee_scs->size());
63  // setup our elements so that all the SCs are grouped together
64  auto SCs_end =
65  std::partition(elems.begin(), elems.end(), [](auto const& a) { return a->type() == reco::PFBlockElement::SC; });
66  // add eb superclusters
67  auto bsc = eb_scs->cbegin();
68  auto esc = eb_scs->cend();
69  reco::PFBlockElementSuperCluster* scbe = nullptr;
71  for (auto sc = bsc; sc != esc; ++sc) {
72  scref = reco::SuperClusterRef(eb_scs, std::distance(bsc, sc));
73  PFBlockElementSCEqual myEqual(scref);
74  auto sc_elem = std::find_if(elems.begin(), SCs_end, myEqual);
75  const double scpT = ptFast(sc->energy(), sc->position(), _zero);
76  const double H_tower = (_hadTower->getDepth1HcalESum(*sc) + _hadTower->getDepth2HcalESum(*sc));
77  const double HoverE = H_tower / sc->energy();
78  if (sc_elem == SCs_end && scpT > _minSCPt && (scpT > _pTbyPass || HoverE < _maxHoverE)) {
79  scbe = new reco::PFBlockElementSuperCluster(scref);
81  SCs_end = elems.emplace(SCs_end, scbe);
82  ++SCs_end; // point to element *after* the new one
83  }
84  } // loop on eb superclusters
85  // add ee superclusters
86  bsc = ee_scs->cbegin();
87  esc = ee_scs->cend();
88  for (auto sc = bsc; sc != esc; ++sc) {
89  scref = reco::SuperClusterRef(ee_scs, std::distance(bsc, sc));
90  PFBlockElementSCEqual myEqual(scref);
91  auto sc_elem = std::find_if(elems.begin(), SCs_end, myEqual);
92  const double scpT = ptFast(sc->energy(), sc->position(), _zero);
93  const double H_tower = (_hadTower->getDepth1HcalESum(*sc) + _hadTower->getDepth2HcalESum(*sc));
94  const double HoverE = H_tower / sc->energy();
95  if (sc_elem == SCs_end && scpT > _minSCPt && (scpT > _pTbyPass || HoverE < _maxHoverE)) {
96  scbe = new reco::PFBlockElementSuperCluster(scref);
98  SCs_end = elems.emplace(SCs_end, scbe);
99  ++SCs_end; // point to element *after* the new one
100  }
101  } // loop on ee superclusters
102  elems.shrink_to_fit();
103 }
#define nullptr
edm::EDGetTokenT< reco::SuperClusterCollection > _srcEE
edm::Ref< SuperClusterCollection > SuperClusterRef
reference to an object in a collection of SuperCluster objects
edm::EDGetTokenT< reco::SuperClusterCollection > _srcEB
double ptFast(const double energy, const math::XYZPoint &position, const math::XYZPoint &origin)
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:547
SuperClusterImporter(const edm::ParameterSet &, edm::ConsumesCollector &)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
T sqrt(T t)
Definition: SSEVec.h:19
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:104
std::unique_ptr< EgammaHadTower > _hadTower
edm::EDGetTokenT< CaloTowerCollection > _srcTowers
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:119
static int position[264][3]
Definition: ReadPGInfo.cc:289
void importToBlock(const edm::Event &, ElementList &) const override
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
static const math::XYZPoint _zero
void updateEventSetup(const edm::EventSetup &es) override