CMS 3D CMS Logo

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 
45 };
46 
48  : simParameterMap_(conf),
49  hbheFilter_(),
50  hoFilter_(),
51  hfFilter_(),
52  zdcFilter_(),
53  hbheAnalyzer_("HBHE", 1., &simParameterMap_, &hbheFilter_),
54  hoAnalyzer_("HO", 1., &simParameterMap_, &hoFilter_),
55  hfAnalyzer_("HF", 1., &simParameterMap_, &hfFilter_),
56  zdcAnalyzer_("ZDC", 1., &simParameterMap_, &zdcFilter_),
57  hbheRecHitCollectionTag_(conf.getParameter<edm::InputTag>("hbheRecHitCollectionTag")),
58  hoRecHitCollectionTag_(conf.getParameter<edm::InputTag>("hoRecHitCollectionTag")),
59  hfRecHitCollectionTag_(conf.getParameter<edm::InputTag>("hfRecHitCollectionTag")),
60  hbheRecHitCollectionToken_(consumes(hbheRecHitCollectionTag_)),
61  hoRecHitCollectionToken_(consumes(hoRecHitCollectionTag_)),
62  hfRecHitCollectionToken_(consumes(hfRecHitCollectionTag_)),
63  cfToken_(consumes(edm::InputTag("mix", "HcalHits"))),
64  zdccfToken_(consumes(edm::InputTag("mix", "ZDCHits"))) {}
65 
67  template <class Collection>
69  const edm::Handle<Collection> &recHits = e.getHandle(token);
70  for (unsigned i = 0; i < recHits->size(); ++i) {
71  analyzer.analyze((*recHits)[i].id().rawId(), (*recHits)[i].energy());
72  }
73  }
74 } // namespace HcalHitAnalyzerImpl
75 
77  // Step A: Get Inputs
78  const edm::Handle<CrossingFrame<PCaloHit>> &cf = e.getHandle(cfToken_);
79  //const edm::Handle<CrossingFrame<PCaloHit>>& zdccf = e.getHandle(zdccfToken_);
80 
81  // test access to SimHits for HcalHits and ZDC hits
82  std::unique_ptr<MixCollection<PCaloHit>> hits(new MixCollection<PCaloHit>(cf.product()));
83  // std::unique_ptr<MixCollection<PCaloHit> > zdcHits(new MixCollection<PCaloHit>(zdccf.product()));
85  // hoAnalyzer_.fillHits(*hits);
86  // hfAnalyzer_.fillHits(*hits);
87  // zdcAnalyzer_.fillHits(*hits);
88  HcalHitAnalyzerImpl::analyze<HBHERecHitCollection>(e, hbheAnalyzer_, hbheRecHitCollectionToken_);
89  HcalHitAnalyzerImpl::analyze<HORecHitCollection>(e, hoAnalyzer_, hoRecHitCollectionToken_);
90  HcalHitAnalyzerImpl::analyze<HFRecHitCollection>(e, hfAnalyzer_, hfRecHitCollectionToken_);
91  // HcalHitAnalyzerImpl::analyze<ZDCRecHitCollection>(e, zdcAnalyzer_);
92 }
93 
95 
ZDCHitFilter zdcFilter_
void analyze(edm::Event const &e, edm::EventSetup const &c) override
CaloHitAnalyzer hfAnalyzer_
T const * product() const
Definition: Handle.h:70
HcalSimParameterMap simParameterMap_
HcalHitAnalyzer(edm::ParameterSet const &conf)
const edm::EDGetTokenT< CrossingFrame< PCaloHit > > zdccfToken_
HOHitFilter hoFilter_
const edm::EDGetTokenT< CrossingFrame< PCaloHit > > cfToken_
const edm::EDGetTokenT< HORecHitCollection > hoRecHitCollectionToken_
void analyze(edm::Event const &e, CaloHitAnalyzer &analyzer, edm::EDGetTokenT< Collection > const &token)
const edm::InputTag hfRecHitCollectionTag_
void fillHits(MixCollection< PCaloHit > &hits)
should be called each event
CaloHitAnalyzer hbheAnalyzer_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
CaloHitAnalyzer hoAnalyzer_
const edm::EDGetTokenT< HFRecHitCollection > hfRecHitCollectionToken_
const edm::InputTag hoRecHitCollectionTag_
HBHEHitFilter hbheFilter_
HFHitFilter hfFilter_
HLT enums.
const edm::EDGetTokenT< HBHERecHitCollection > hbheRecHitCollectionToken_
CaloHitAnalyzer zdcAnalyzer_
const edm::InputTag hbheRecHitCollectionTag_