CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HcalHitAnalyzer.cc
Go to the documentation of this file.
1 
17 
18 #include <iostream>
19 #include <string>
20 
22 public:
23  explicit HcalHitAnalyzer(edm::ParameterSet const &conf);
24  void analyze(edm::Event const &e, edm::EventSetup const &c) override;
25 
26 private:
36 
40 };
41 
43  : simParameterMap_(conf),
44  hbheFilter_(),
45  hoFilter_(),
46  hfFilter_(),
47  zdcFilter_(),
48  hbheAnalyzer_("HBHE", 1., &simParameterMap_, &hbheFilter_),
49  hoAnalyzer_("HO", 1., &simParameterMap_, &hoFilter_),
50  hfAnalyzer_("HF", 1., &simParameterMap_, &hfFilter_),
51  zdcAnalyzer_("ZDC", 1., &simParameterMap_, &zdcFilter_),
52  hbheRecHitCollectionTag_(conf.getParameter<edm::InputTag>("hbheRecHitCollectionTag")),
53  hoRecHitCollectionTag_(conf.getParameter<edm::InputTag>("hoRecHitCollectionTag")),
54  hfRecHitCollectionTag_(conf.getParameter<edm::InputTag>("hfRecHitCollectionTag")) {}
55 
56 namespace HcalHitAnalyzerImpl {
57  template <class Collection>
58  void analyze(edm::Event const &e, CaloHitAnalyzer &analyzer, edm::InputTag &tag) {
60  e.getByLabel(tag, recHits);
61  for (unsigned i = 0; i < recHits->size(); ++i) {
62  analyzer.analyze((*recHits)[i].id().rawId(), (*recHits)[i].energy());
63  }
64  }
65 } // namespace HcalHitAnalyzerImpl
66 
68  // Step A: Get Inputs
70  e.getByLabel("mix", "g4SimHitsHcalHits", cf);
71  // e.getByLabel("mix", "ZDCHits", zdccf);
72 
73  // test access to SimHits for HcalHits and ZDC hits
74  std::unique_ptr<MixCollection<PCaloHit>> hits(new MixCollection<PCaloHit>(cf.product()));
75  // std::unique_ptr<MixCollection<PCaloHit> > zdcHits(new
76  // MixCollection<PCaloHit>(zdccf.product()));
77  hbheAnalyzer_.fillHits(*hits);
78  // hoAnalyzer_.fillHits(*hits);
79  // hfAnalyzer_.fillHits(*hits);
80  // zdcAnalyzer_.fillHits(*hits);
81  HcalHitAnalyzerImpl::analyze<HBHERecHitCollection>(e, hbheAnalyzer_, hbheRecHitCollectionTag_);
82  HcalHitAnalyzerImpl::analyze<HORecHitCollection>(e, hoAnalyzer_, hoRecHitCollectionTag_);
83  HcalHitAnalyzerImpl::analyze<HFRecHitCollection>(e, hfAnalyzer_, hfRecHitCollectionTag_);
84  // HcalHitAnalyzerImpl::analyze<ZDCRecHitCollection>(e, zdcAnalyzer_);
85 }
86 
88 
ZDCHitFilter zdcFilter_
void analyze(edm::Event const &e, CaloHitAnalyzer &analyzer, edm::InputTag &tag)
const edm::EventSetup & c
void analyze(int detId, double recEnergy)
to be called for each RecHit
void analyze(edm::Event const &e, edm::EventSetup const &c) override
CaloHitAnalyzer hfAnalyzer_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
HcalSimParameterMap simParameterMap_
HcalHitAnalyzer(edm::ParameterSet const &conf)
HOHitFilter hoFilter_
edm::InputTag hbheRecHitCollectionTag_
edm::InputTag hoRecHitCollectionTag_
void fillHits(MixCollection< PCaloHit > &hits)
should be called each event
CaloHitAnalyzer hbheAnalyzer_
edm::InputTag hfRecHitCollectionTag_
CaloHitAnalyzer hoAnalyzer_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:500
T const * product() const
Definition: Handle.h:70
HBHEHitFilter hbheFilter_
HFHitFilter hfFilter_
CaloHitAnalyzer zdcAnalyzer_