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,
15  const math::XYZPoint& position,
16  const math::XYZPoint& origin ) {
17  const auto v = position - origin;
18  return energy*std::sqrt(v.perp2()/v.mag2());
19 }
20 
21 #include <unordered_map>
22 
24 public:
26 
27  void updateEventSetup( const edm::EventSetup& es ) override;
28 
29  void importToBlock( const edm::Event& ,
30  ElementList& ) const override;
31 
32 private:
35  const double _maxHoverE, _pTbyPass, _minSCPt;
36  std::unique_ptr<EgammaHadTower> _hadTower;
38  static const math::XYZPoint _zero;
39 };
40 
42 
45  "SuperClusterImporter");
46 
48  edm::ConsumesCollector& sumes) :
49  BlockElementImporterBase(conf,sumes),
50  _srcEB(sumes.consumes<reco::SuperClusterCollection>(conf.getParameter<edm::InputTag>("source_eb"))),
51  _srcEE(sumes.consumes<reco::SuperClusterCollection>(conf.getParameter<edm::InputTag>("source_ee"))),
52  _srcTowers(sumes.consumes<CaloTowerCollection>(conf.getParameter<edm::InputTag>("source_towers"))),
53  _maxHoverE(conf.getParameter<double>("maximumHoverE")),
54  _pTbyPass(conf.getParameter<double>("minPTforBypass")),
55  _minSCPt(conf.getParameter<double>("minSuperClusterPt")),
57  _superClustersArePF(conf.getParameter<bool>("superClustersArePF")) {
58 }
59 
63 }
64 
68  auto eb_scs = e.getHandle(_srcEB);
69  auto ee_scs = e.getHandle(_srcEE);
70  auto towers = e.getHandle(_srcTowers);
71  _hadTower->setTowerCollection(towers.product());
72  elems.reserve(elems.size()+eb_scs->size()+ee_scs->size());
73  // setup our elements so that all the SCs are grouped together
74  auto SCs_end = std::partition(elems.begin(),elems.end(),
75  [](auto const& a){
76  return a->type() == reco::PFBlockElement::SC;
77  });
78  // add eb superclusters
79  auto bsc = eb_scs->cbegin();
80  auto esc = eb_scs->cend();
81  reco::PFBlockElementSuperCluster* scbe = nullptr;
83  for( auto sc = bsc; sc != esc; ++sc ) {
84  scref = reco::SuperClusterRef(eb_scs,std::distance(bsc,sc));
85  PFBlockElementSCEqual myEqual(scref);
86  auto sc_elem = std::find_if(elems.begin(),SCs_end,myEqual);
87  const double scpT = ptFast(sc->energy(),sc->position(),_zero);
88  const double H_tower = ( _hadTower->getDepth1HcalESum(*sc) +
89  _hadTower->getDepth2HcalESum(*sc) );
90  const double HoverE = H_tower/sc->energy();
91  if( sc_elem == SCs_end && scpT > _minSCPt &&
92  (scpT > _pTbyPass || HoverE < _maxHoverE) ) {
93  scbe = new reco::PFBlockElementSuperCluster(scref);
95  SCs_end = elems.emplace(SCs_end, scbe);
96  ++SCs_end; // point to element *after* the new one
97  }
98  }// loop on eb superclusters
99  // add ee superclusters
100  bsc = ee_scs->cbegin();
101  esc = ee_scs->cend();
102  for( auto sc = bsc; sc != esc; ++sc ) {
103  scref = reco::SuperClusterRef(ee_scs,std::distance(bsc,sc));
104  PFBlockElementSCEqual myEqual(scref);
105  auto sc_elem = std::find_if(elems.begin(),SCs_end,myEqual);
106  const double scpT = ptFast(sc->energy(),sc->position(),_zero);
107  const double H_tower = ( _hadTower->getDepth1HcalESum(*sc) +
108  _hadTower->getDepth2HcalESum(*sc) );
109  const double HoverE = H_tower/sc->energy();
110  if( sc_elem == SCs_end && scpT > _minSCPt &&
111  (scpT > _pTbyPass || HoverE < _maxHoverE)) {
112  scbe = new reco::PFBlockElementSuperCluster(scref);
114  SCs_end = elems.emplace(SCs_end, scbe);
115  ++SCs_end; // point to element *after* the new one
116  }
117  }// loop on ee superclusters
118  elems.shrink_to_fit();
119 }
#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:539
SuperClusterImporter(const edm::ParameterSet &, edm::ConsumesCollector &)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
T sqrt(T t)
Definition: SSEVec.h:18
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:121
static int position[264][3]
Definition: ReadPGInfo.cc:509
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