CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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),
33  amplitudeThreshold_(amplitudeThreshold),
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 
97 };
98 
100  : hitReadoutName_("HcalHits"),
101  simParameterMap_(),
102  hbheFilter_(),
103  hoFilter_(),
104  hfFilter_(),
105  hbheHitAnalyzer_("HBHEDigi", 1., &simParameterMap_, &hbheFilter_),
106  hoHitAnalyzer_("HODigi", 1., &simParameterMap_, &hoFilter_),
107  hfHitAnalyzer_("HFDigi", 1., &simParameterMap_, &hfFilter_),
108  zdcHitAnalyzer_("ZDCDigi", 1., &simParameterMap_, &zdcFilter_),
109  hbheDigiStatistics_("HBHEDigi", 4, 10., 6., 0.1, 0.5, hbheHitAnalyzer_),
110  hoDigiStatistics_("HODigi", 4, 10., 6., 0.1, 0.5, hoHitAnalyzer_),
111  hfDigiStatistics_("HFDigi", 3, 10., 6., 0.1, 0.5, hfHitAnalyzer_),
112  zdcDigiStatistics_("ZDCDigi", 3, 10., 6., 0.1, 0.5, zdcHitAnalyzer_),
113  hbheDigiCollectionTag_(conf.getParameter<edm::InputTag>("hbheDigiCollectionTag")),
114  hoDigiCollectionTag_(conf.getParameter<edm::InputTag>("hoDigiCollectionTag")),
115  hfDigiCollectionTag_(conf.getParameter<edm::InputTag>("hfDigiCollectionTag")) {}
116 
117 namespace HcalDigiAnalyzerImpl {
118  template <class Collection>
119  void analyze(edm::Event const &e, HcalDigiStatistics &statistics, edm::InputTag &tag) {
121  e.getByLabel(tag, digis);
122  for (unsigned i = 0; i < digis->size(); ++i) {
123  std::cout << (*digis)[i] << std::endl;
124  statistics.analyze((*digis)[i]);
125  }
126  }
127 } // namespace HcalDigiAnalyzerImpl
128 
130  // Step A: Get Inputs
132  e.getByLabel("mix", "HcalHits", cf);
133  // e.getByLabel("mix", "ZDCHits", zdccf);
134 
135  // test access to SimHits for HcalHits and ZDC hits
136  std::unique_ptr<MixCollection<PCaloHit>> hits(new MixCollection<PCaloHit>(cf.product()));
137  // std::unique_ptr<MixCollection<PCaloHit> > zdcHits(new
138  // MixCollection<PCaloHit>(zdccf.product()));
139  hbheHitAnalyzer_.fillHits(*hits);
140  hoHitAnalyzer_.fillHits(*hits);
141  hfHitAnalyzer_.fillHits(*hits);
142  // zdcHitAnalyzer_.fillHits(*zdcHits);
143  HcalDigiAnalyzerImpl::analyze<HBHEDigiCollection>(e, hbheDigiStatistics_, hbheDigiCollectionTag_);
144  HcalDigiAnalyzerImpl::analyze<HODigiCollection>(e, hoDigiStatistics_, hoDigiCollectionTag_);
145  HcalDigiAnalyzerImpl::analyze<HFDigiCollection>(e, hfDigiStatistics_, hfDigiCollectionTag_);
146  // HcalDigiAnalyzerImpl::analyze<ZDCDigiCollection>(e, zdcDigiStatistics_);
147 }
148 
150 
const edm::EventSetup & c
edm::InputTag hfDigiCollectionTag_
void analyze(int detId, double recEnergy)
to be called for each RecHit
CaloValidationStatistics pedestal_
CaloHitAnalyzer hoHitAnalyzer_
constexpr unsigned int maxBin
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
CaloValidationStatistics binPrevToBinMax_
HcalDigiAnalyzer(edm::ParameterSet const &conf)
HFHitFilter hfFilter_
HcalDigiStatistics hfDigiStatistics_
HcalDigiStatistics hoDigiStatistics_
edm::InputTag hbheDigiCollectionTag_
CaloHitAnalyzer hfHitAnalyzer_
CaloHitAnalyzer zdcHitAnalyzer_
CaloHitAnalyzer hbheHitAnalyzer_
HOHitFilter hoFilter_
HcalDigiStatistics zdcDigiStatistics_
edm::InputTag hoDigiCollectionTag_
void fillHits(MixCollection< PCaloHit > &hits)
should be called each event
HBHEHitFilter hbheFilter_
CaloValidationStatistics binNextToBinMax_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:500
void analyze(const Digi &digi)
HcalDigiStatistics hbheDigiStatistics_
ZDCHitFilter zdcFilter_
T const * product() const
Definition: Handle.h:70
void addEntry(float value, float weight=1.)
void analyze(edm::Event const &e, HcalDigiStatistics &statistics, edm::InputTag &tag)
CaloHitAnalyzer & amplitudeAnalyzer_
HcalSimParameterMap simParameterMap_
void analyze(edm::Event const &e, edm::EventSetup const &c) override
tuple cout
Definition: gather_cfg.py:144
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