CMS 3D CMS Logo

HcalDigiAnalyzer.cc
Go to the documentation of this file.
1 
19 
20 #include <iostream>
21 #include <string>
22 
24 public:
26  int maxBin,
27  float amplitudeThreshold,
28  float expectedPedestal,
29  float binPrevToBinMax,
30  float binNextToBinMax,
31  CaloHitAnalyzer &amplitudeAnalyzer)
32  : maxBin_(maxBin),
34  pedestal_(name + " pedestal", expectedPedestal, 0.),
35  binPrevToBinMax_(name + " binPrevToBinMax", binPrevToBinMax, 0.),
36  binNextToBinMax_(name + " binNextToBinMax", binNextToBinMax, 0.),
37  amplitudeAnalyzer_(amplitudeAnalyzer) {}
38 
39  template <class Digi>
40  void analyze(const Digi &digi);
41 
42 private:
43  int maxBin_;
49 };
50 
51 template <class Digi>
52 void HcalDigiStatistics::analyze(const Digi &digi) {
53  pedestal_.addEntry(digi[0].adc());
54  pedestal_.addEntry(digi[1].adc());
55 
56  double pedestal_fC = 0.5 * (digi[0].nominal_fC() + digi[1].nominal_fC());
57 
58  double maxAmplitude = digi[maxBin_].nominal_fC() - pedestal_fC;
59 
60  if (maxAmplitude > amplitudeThreshold_) {
61  double binPrevToBinMax = (digi[maxBin_ - 1].nominal_fC() - pedestal_fC) / maxAmplitude;
62  binPrevToBinMax_.addEntry(binPrevToBinMax);
63 
64  double binNextToBinMax = (digi[maxBin_ + 1].nominal_fC() - pedestal_fC) / maxAmplitude;
65  binNextToBinMax_.addEntry(binNextToBinMax);
66 
67  double amplitude = digi[maxBin_].nominal_fC() + digi[maxBin_ + 1].nominal_fC() - 2 * pedestal_fC;
68 
69  amplitudeAnalyzer_.analyze(digi.id().rawId(), amplitude);
70  }
71 }
72 
74 public:
75  explicit HcalDigiAnalyzer(edm::ParameterSet const &conf);
76  void analyze(edm::Event const &e, edm::EventSetup const &c) override;
77 
78 private:
93 
102 };
103 
105  : hitReadoutName_("HcalHits"),
106  simParameterMap_(),
107  hbheFilter_(),
108  hoFilter_(),
109  hfFilter_(),
110  hbheHitAnalyzer_("HBHEDigi", 1., &simParameterMap_, &hbheFilter_),
111  hoHitAnalyzer_("HODigi", 1., &simParameterMap_, &hoFilter_),
112  hfHitAnalyzer_("HFDigi", 1., &simParameterMap_, &hfFilter_),
113  zdcHitAnalyzer_("ZDCDigi", 1., &simParameterMap_, &zdcFilter_),
114  hbheDigiStatistics_("HBHEDigi", 4, 10., 6., 0.1, 0.5, hbheHitAnalyzer_),
115  hoDigiStatistics_("HODigi", 4, 10., 6., 0.1, 0.5, hoHitAnalyzer_),
116  hfDigiStatistics_("HFDigi", 3, 10., 6., 0.1, 0.5, hfHitAnalyzer_),
117  zdcDigiStatistics_("ZDCDigi", 3, 10., 6., 0.1, 0.5, zdcHitAnalyzer_),
118  hbheDigiCollectionTag_(conf.getParameter<edm::InputTag>("hbheDigiCollectionTag")),
119  hoDigiCollectionTag_(conf.getParameter<edm::InputTag>("hoDigiCollectionTag")),
120  hfDigiCollectionTag_(conf.getParameter<edm::InputTag>("hfDigiCollectionTag")),
121  hbheDigiCollectionToken_(consumes(hbheDigiCollectionTag_)),
122  hoDigiCollectionToken_(consumes(hoDigiCollectionTag_)),
123  hfDigiCollectionToken_(consumes(hfDigiCollectionTag_)),
124  cfToken_(consumes(edm::InputTag("mix", "HcalHits"))),
125  zdccfToken_(consumes(edm::InputTag("mix", "ZDCHits"))) {}
126 
128  template <class Collection>
130  const edm::Handle<Collection> &digis = e.getHandle(token);
131  for (unsigned i = 0; i < digis->size(); ++i) {
132  edm::LogVerbatim("HcalSim") << (*digis)[i];
133  statistics.analyze((*digis)[i]);
134  }
135  }
136 } // namespace HcalDigiAnalyzerImpl
137 
139  // Step A: Get Inputs
140  const edm::Handle<CrossingFrame<PCaloHit>> &cf = e.getHandle(cfToken_);
141  //const edm::Handle<CrossingFrame<PCaloHit>>& zdccf = e.getHandle(zdccfToken_);
142 
143  // test access to SimHits for HcalHits and ZDC hits
144  std::unique_ptr<MixCollection<PCaloHit>> hits(new MixCollection<PCaloHit>(cf.product()));
145  // std::unique_ptr<MixCollection<PCaloHit> > zdcHits(new MixCollection<PCaloHit>(zdccf.product()));
149  // zdcHitAnalyzer_.fillHits(*zdcHits);
150  HcalDigiAnalyzerImpl::analyze<HBHEDigiCollection>(e, hbheDigiStatistics_, hbheDigiCollectionToken_);
151  HcalDigiAnalyzerImpl::analyze<HODigiCollection>(e, hoDigiStatistics_, hoDigiCollectionToken_);
152  HcalDigiAnalyzerImpl::analyze<HFDigiCollection>(e, hfDigiStatistics_, hfDigiCollectionToken_);
153  // HcalDigiAnalyzerImpl::analyze<ZDCDigiCollection>(e, zdcDigiStatistics_);
154 }
155 
157 
Log< level::Info, true > LogVerbatim
const edm::EDGetTokenT< HFDigiCollection > hfDigiCollectionToken_
void analyze(int detId, double recEnergy)
to be called for each RecHit
CaloValidationStatistics pedestal_
CaloHitAnalyzer hoHitAnalyzer_
constexpr unsigned int maxBin
CaloValidationStatistics binPrevToBinMax_
HcalDigiAnalyzer(edm::ParameterSet const &conf)
HFHitFilter hfFilter_
T const * product() const
Definition: Handle.h:70
HcalDigiStatistics hfDigiStatistics_
HcalDigiStatistics hoDigiStatistics_
CaloHitAnalyzer hfHitAnalyzer_
CaloHitAnalyzer zdcHitAnalyzer_
CaloHitAnalyzer hbheHitAnalyzer_
HOHitFilter hoFilter_
HcalDigiStatistics zdcDigiStatistics_
void fillHits(MixCollection< PCaloHit > &hits)
should be called each event
HBHEHitFilter hbheFilter_
CaloValidationStatistics binNextToBinMax_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< CrossingFrame< PCaloHit > > cfToken_
void analyze(const Digi &digi)
HcalDigiStatistics hbheDigiStatistics_
ZDCHitFilter zdcFilter_
const edm::InputTag hbheDigiCollectionTag_
const edm::InputTag hfDigiCollectionTag_
void addEntry(float value, float weight=1.)
const edm::EDGetTokenT< HBHEDigiCollection > hbheDigiCollectionToken_
CaloHitAnalyzer & amplitudeAnalyzer_
HcalSimParameterMap simParameterMap_
HLT enums.
void analyze(edm::Event const &e, HcalDigiStatistics &statistics, const edm::EDGetTokenT< Collection > &token)
void analyze(edm::Event const &e, edm::EventSetup const &c) override
const edm::EDGetTokenT< CrossingFrame< PCaloHit > > zdccfToken_
const edm::EDGetTokenT< HODigiCollection > hoDigiCollectionToken_
const edm::InputTag hoDigiCollectionTag_
HcalDigiStatistics(std::string name, int maxBin, float amplitudeThreshold, float expectedPedestal, float binPrevToBinMax, float binNextToBinMax, CaloHitAnalyzer &amplitudeAnalyzer)
std::string hitReadoutName_
uint16_t *__restrict__ uint16_t const *__restrict__ adc