CMS 3D CMS Logo

HcalNoiseRates.cc
Go to the documentation of this file.
1 //
2 // HcalNoiseRates.cc
3 //
4 // description: Calculation for single particle response corrections
5 //
6 // author: K. Hatakeyama, H. Liu, Baylor
7 //
8 //
9 
12 
13 //
14 // constructors and destructor
15 //
16 
18  // DQM ROOT output
19  outputFile_ = iConfig.getUntrackedParameter<std::string>("outputFile", "myfile.root");
20 
21  // set parameters
22  rbxCollName_ = iConfig.getUntrackedParameter<edm::InputTag>("rbxCollName");
23  tok_rbx_ = consumes<reco::HcalNoiseRBXCollection>(rbxCollName_);
24  minRBXEnergy_ = iConfig.getUntrackedParameter<double>("minRBXEnergy");
25  minHitEnergy_ = iConfig.getUntrackedParameter<double>("minHitEnergy");
26 
27  useAllHistos_ = iConfig.getUntrackedParameter<bool>("useAllHistos", false);
28 
29  // Hcal Noise Summary
30  noisetoken_ = consumes<HcalNoiseSummary>(iConfig.getParameter<edm::InputTag>("noiselabel"));
31 }
32 
34  edm::Run const & /* iRun*/,
35  edm::EventSetup const & /* iSetup */)
36 
37 {
38  ibooker.setCurrentFolder("HcalNoiseRatesD/HcalNoiseRatesTask");
39 
40  Char_t histo[100];
41 
42  // book histograms
43 
44  // Lumi block is not drawn; the rest are
45  if (useAllHistos_) {
46  sprintf(histo, "hLumiBlockCount");
47  hLumiBlockCount_ = ibooker.book1D(histo, histo, 1, -0.5, 0.5);
48  }
49 
50  sprintf(histo, "hRBXEnergy");
51  hRBXEnergy_ = ibooker.book1D(histo, histo, 300, 0, 3000);
52 
53  sprintf(histo, "hRBXEnergyType1");
54  hRBXEnergyType1_ = ibooker.book1D(histo, histo, 300, 0, 3000);
55 
56  sprintf(histo, "hRBXEnergyType2");
57  hRBXEnergyType2_ = ibooker.book1D(histo, histo, 300, 0, 3000);
58 
59  sprintf(histo, "hRBXEnergyType3");
60  hRBXEnergyType3_ = ibooker.book1D(histo, histo, 300, 0, 3000);
61 
62  sprintf(histo, "hRBXNHits");
63  hRBXNHits_ = ibooker.book1D(histo, histo, 73, -0.5, 72.5);
64 
65  // HcalNoiseSummary
66 
67  sprintf(histo, "nNNumChannels");
68  nNNumChannels_ = ibooker.book1D(histo, histo, 100, 0, 100);
69  sprintf(histo, "nNSumE");
70  nNSumE_ = ibooker.book1D(histo, histo, 100, 0, 5000);
71  sprintf(histo, "nNSumEt");
72  nNSumEt_ = ibooker.book1D(histo, histo, 100, 0, 2000);
73 
74  sprintf(histo, "sNNumChannels");
75  sNNumChannels_ = ibooker.book1D(histo, histo, 100, 0, 100);
76  sprintf(histo, "sNSumE");
77  sNSumE_ = ibooker.book1D(histo, histo, 100, 0, 5000);
78  sprintf(histo, "sNSumEt");
79  sNSumEt_ = ibooker.book1D(histo, histo, 100, 0, 2000);
80 
81  sprintf(histo, "iNNumChannels");
82  iNNumChannels_ = ibooker.book1D(histo, histo, 100, 0, 100);
83  sprintf(histo, "iNSumE");
84  iNSumE_ = ibooker.book1D(histo, histo, 100, 0, 5000);
85  sprintf(histo, "iNSumEt");
86  iNSumEt_ = ibooker.book1D(histo, histo, 100, 0, 2000);
87 
88  sprintf(histo, "hNoise_maxZeros");
89  hNoise_maxZeros_ = ibooker.book1D(histo, histo, 80, 0, 80);
90  sprintf(histo, "hNoise_maxHPDHits");
91  hNoise_maxHPDHits_ = ibooker.book1D(histo, histo, 20, 0, 20);
92  sprintf(histo, "hNoise_maxHPDNoOtherHits");
93  hNoise_maxHPDNoOtherHits_ = ibooker.book1D(histo, histo, 20, 0, 20);
94 }
95 
97 
98 //
99 // member functions
100 //
101 
102 // ------------ method called to for each event ------------
104  // get the lumi section
105  int lumiSection = iEvent.luminosityBlock();
106  lumiCountMap_[lumiSection]++;
107 
108  // get the RBX Noise collection
110  iEvent.getByToken(tok_rbx_, handle);
111  if (!handle.isValid()) {
113  << " could not find HcalNoiseRBXCollection named " << rbxCollName_ << ".\n";
114  return;
115  }
116 
117  // get the Noise summary object
119  iEvent.getByToken(noisetoken_, summary_h);
120  if (!summary_h.isValid()) {
121  throw edm::Exception(edm::errors::ProductNotFound) << " could not find HcalNoiseSummary.\n";
122  return;
123  }
124  const HcalNoiseSummary summary = *summary_h;
125 
126  // Fill the Noise Summary histograms
128  nNSumE_->Fill(summary.NegativeNoiseSumE());
129  nNSumEt_->Fill(summary.NegativeNoiseSumEt());
130 
132  sNSumE_->Fill(summary.spikeNoiseSumE());
133  sNSumEt_->Fill(summary.spikeNoiseSumEt());
134 
136  iNSumE_->Fill(summary.isolatedNoiseSumE());
137  iNSumEt_->Fill(summary.isolatedNoiseSumEt());
138 
139  hNoise_maxZeros_->Fill(summary.maxZeros());
140  hNoise_maxHPDHits_->Fill(summary.maxHPDHits());
142 
143  // loop over the RBXs and fill the histograms
144  for (reco::HcalNoiseRBXCollection::const_iterator it = handle->begin(); it != handle->end(); ++it) {
145  const reco::HcalNoiseRBX &rbx = (*it);
146 
147  double energy = rbx.recHitEnergy(minHitEnergy_);
148 
149  int nhits = rbx.numRecHits(minHitEnergy_);
150 
151  if (energy < minRBXEnergy_)
152  continue;
153 
154  hRBXEnergy_->Fill(energy);
155 
156  if (nhits <= 9)
157  hRBXEnergyType1_->Fill(energy);
158  else if (nhits <= 18)
159  hRBXEnergyType2_->Fill(energy);
160  else
161  hRBXEnergyType3_->Fill(energy);
162 
163  hRBXNHits_->Fill(nhits);
164 
165  } // done looping over RBXs
166 }
167 
168 // ------------ method called once each job just before starting event loop
169 // ------------
171 
172 // ------------ method called once each job just after ending the event loop
173 // ------------
175  if (useAllHistos_)
176  hLumiBlockCount_->Fill(0.0, lumiCountMap_.size());
177 }
178 
179 // define this as a plug-in
float NegativeNoiseSumEt(void) const
int numSpikeNoiseChannels(void) const
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * nNNumChannels_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
HcalNoiseRates(const edm::ParameterSet &)
edm::InputTag rbxCollName_
MonitorElement * hRBXEnergyType1_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
int numIsolatedNoiseChannels(void) const
float spikeNoiseSumE(void) const
void beginJob() override
int numRecHits(double threshold=1.5) const
MonitorElement * iNNumChannels_
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:61
MonitorElement * hRBXEnergy_
void analyze(const edm::Event &, const edm::EventSetup &) override
float spikeNoiseSumEt(void) const
MonitorElement * sNNumChannels_
MonitorElement * hRBXEnergyType2_
void Fill(long long x)
int numNegativeNoiseChannels(void) const
float isolatedNoiseSumEt(void) const
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double recHitEnergy(double theshold=1.5) const
Definition: HcalNoiseRBX.cc:99
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
void endJob() override
MonitorElement * nNSumE_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
MonitorElement * hNoise_maxHPDNoOtherHits_
bool isValid() const
Definition: HandleBase.h:74
MonitorElement * sNSumEt_
edm::EDGetTokenT< reco::HcalNoiseRBXCollection > tok_rbx_
edm::EDGetTokenT< HcalNoiseSummary > noisetoken_
MonitorElement * hNoise_maxHPDHits_
std::string outputFile_
float NegativeNoiseSumE(void) const
MonitorElement * hRBXEnergyType3_
double minHitEnergy_
MonitorElement * nNSumEt_
int maxZeros(void) const
~HcalNoiseRates() override
float isolatedNoiseSumE(void) const
int maxHPDNoOtherHits(void) const
MonitorElement * hLumiBlockCount_
MonitorElement * sNSumE_
std::map< int, int > lumiCountMap_
MonitorElement * iNSumEt_
int maxHPDHits(void) const
MonitorElement * iNSumE_
double minRBXEnergy_
MonitorElement * hNoise_maxZeros_
MonitorElement * hRBXNHits_
Definition: Run.h:45