CMS 3D CMS Logo

NoiseRates.cc
Go to the documentation of this file.
1 //
2 // NoiseRates.cc
3 //
4 // description: Calculation for single particle response corrections
5 //
6 // author: J.P. Chou, Brown
7 //
8 //
9 
12 
13 //
14 // constructors and destructor
15 //
16 
18  : outputFile_(iConfig.getUntrackedParameter<std::string>("outputFile", "myfile.root")),
19  rbxCollName_(iConfig.getUntrackedParameter<edm::InputTag>("rbxCollName")),
20  minRBXEnergy_(iConfig.getUntrackedParameter<double>("minRBXEnergy")),
21  minHitEnergy_(iConfig.getUntrackedParameter<double>("minHitEnergy")),
22  useAllHistos_(iConfig.getUntrackedParameter<bool>("useAllHistos", false)),
23  tok_rbx_(consumes<reco::HcalNoiseRBXCollection>(rbxCollName_)),
24  noisetoken_(consumes<HcalNoiseSummary>(iConfig.getParameter<edm::InputTag>("noiselabel"))) {
25  // DQM ROOT output
26  // set parameters
27  // Hcal Noise Summary
28 }
29 
30 //
31 // member functions
32 //
33 
35  ib.setCurrentFolder("NoiseRatesV/NoiseRatesTask");
36 
37  // book histograms
38  Char_t histo[100];
39 
40  // Lumi block is not drawn; the rest are
41  if (useAllHistos_) {
42  sprintf(histo, "hLumiBlockCount");
43  hLumiBlockCount_ = ib.book1D(histo, histo, 1, -0.5, 0.5);
44  }
45 
46  sprintf(histo, "hRBXEnergy");
47  hRBXEnergy_ = ib.book1D(histo, histo, 300, 0, 3000);
48 
49  sprintf(histo, "hRBXEnergyType1");
50  hRBXEnergyType1_ = ib.book1D(histo, histo, 300, 0, 3000);
51 
52  sprintf(histo, "hRBXEnergyType2");
53  hRBXEnergyType2_ = ib.book1D(histo, histo, 300, 0, 3000);
54 
55  sprintf(histo, "hRBXEnergyType3");
56  hRBXEnergyType3_ = ib.book1D(histo, histo, 300, 0, 3000);
57 
58  sprintf(histo, "hRBXNHits");
59  hRBXNHits_ = ib.book1D(histo, histo, 73, -0.5, 72.5);
60 
61  // HcalNoiseSummary
62 
63  sprintf(histo, "nNNumChannels");
64  nNNumChannels_ = ib.book1D(histo, histo, 100, 0, 100);
65  sprintf(histo, "nNSumE");
66  nNSumE_ = ib.book1D(histo, histo, 100, 0, 5000);
67  sprintf(histo, "nNSumEt");
68  nNSumEt_ = ib.book1D(histo, histo, 100, 0, 2000);
69 
70  sprintf(histo, "sNNumChannels");
71  sNNumChannels_ = ib.book1D(histo, histo, 100, 0, 100);
72  sprintf(histo, "sNSumE");
73  sNSumE_ = ib.book1D(histo, histo, 100, 0, 5000);
74  sprintf(histo, "sNSumEt");
75  sNSumEt_ = ib.book1D(histo, histo, 100, 0, 2000);
76 
77  sprintf(histo, "iNNumChannels");
78  iNNumChannels_ = ib.book1D(histo, histo, 100, 0, 100);
79  sprintf(histo, "iNSumE");
80  iNSumE_ = ib.book1D(histo, histo, 100, 0, 5000);
81  sprintf(histo, "iNSumEt");
82  iNSumEt_ = ib.book1D(histo, histo, 100, 0, 2000);
83 
84  sprintf(histo, "hNoise_maxZeros");
85  hNoise_maxZeros_ = ib.book1D(histo, histo, 80, 0, 80);
86  sprintf(histo, "hNoise_maxHPDHits");
87  hNoise_maxHPDHits_ = ib.book1D(histo, histo, 20, 0, 20);
88  sprintf(histo, "hNoise_maxHPDNoOtherHits");
89  hNoise_maxHPDNoOtherHits_ = ib.book1D(histo, histo, 20, 0, 20);
90 }
91 
92 // ------------ method called to for each event ------------
93 void NoiseRates::analyze(const edm::Event &iEvent, const edm::EventSetup &evSetup) {
94  // get the lumi section
95  int lumiSection = iEvent.luminosityBlock();
96  lumiCountMap_[lumiSection]++;
97 
98  // get the RBX Noise collection
100  if (!handle.isValid()) {
102  << " could not find HcalNoiseRBXCollection named " << rbxCollName_ << ".\n";
103  return;
104  }
105 
106  // get the Noise summary object
107  const edm::Handle<HcalNoiseSummary> &summary_h = iEvent.getHandle(noisetoken_);
108  if (!summary_h.isValid()) {
109  throw edm::Exception(edm::errors::ProductNotFound) << " could not find HcalNoiseSummary.\n";
110  return;
111  }
112  const HcalNoiseSummary summary = *summary_h;
113 
114  // Fill the Noise Summary histograms
115  nNNumChannels_->Fill(summary.numNegativeNoiseChannels());
116  nNSumE_->Fill(summary.NegativeNoiseSumE());
117  nNSumEt_->Fill(summary.NegativeNoiseSumEt());
118 
119  sNNumChannels_->Fill(summary.numSpikeNoiseChannels());
120  sNSumE_->Fill(summary.spikeNoiseSumE());
121  sNSumEt_->Fill(summary.spikeNoiseSumEt());
122 
123  iNNumChannels_->Fill(summary.numIsolatedNoiseChannels());
124  iNSumE_->Fill(summary.isolatedNoiseSumE());
125  iNSumEt_->Fill(summary.isolatedNoiseSumEt());
126 
127  hNoise_maxZeros_->Fill(summary.maxZeros());
128  hNoise_maxHPDHits_->Fill(summary.maxHPDHits());
129  hNoise_maxHPDNoOtherHits_->Fill(summary.maxHPDNoOtherHits());
130 
131  // loop over the RBXs and fill the histograms
132  for (reco::HcalNoiseRBXCollection::const_iterator it = handle->begin(); it != handle->end(); ++it) {
133  const reco::HcalNoiseRBX &rbx = (*it);
134 
135  double energy = rbx.recHitEnergy(minHitEnergy_);
136 
137  int nhits = rbx.numRecHits(minHitEnergy_);
138 
139  if (energy < minRBXEnergy_)
140  continue;
141 
143 
144  if (nhits <= 9)
146  else if (nhits <= 18)
148  else
150 
152 
153  } // done looping over RBXs
154 }
155 
156 // define this as a plug-in
MonitorElement * sNNumChannels_
Definition: NoiseRates.h:76
MonitorElement * hRBXEnergy_
Definition: NoiseRates.h:64
std::vector< HcalNoiseRBX > HcalNoiseRBXCollection
Definition: HcalNoiseRBX.h:24
const double minHitEnergy_
Definition: NoiseRates.h:56
const bool useAllHistos_
Definition: NoiseRates.h:57
double recHitEnergy(double theshold=1.5) const
Definition: HcalNoiseRBX.cc:76
MonitorElement * hNoise_maxHPDNoOtherHits_
Definition: NoiseRates.h:86
const edm::EDGetTokenT< reco::HcalNoiseRBXCollection > tok_rbx_
Definition: NoiseRates.h:59
const edm::EDGetTokenT< HcalNoiseSummary > noisetoken_
Definition: NoiseRates.h:61
MonitorElement * hRBXEnergyType3_
Definition: NoiseRates.h:67
int numRecHits(double threshold=1.5) const
MonitorElement * hLumiBlockCount_
Definition: NoiseRates.h:63
MonitorElement * hRBXNHits_
Definition: NoiseRates.h:68
MonitorElement * nNNumChannels_
Definition: NoiseRates.h:72
std::map< int, int > lumiCountMap_
Definition: NoiseRates.h:89
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: NoiseRates.cc:34
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
MonitorElement * hNoise_maxZeros_
Definition: NoiseRates.h:84
MonitorElement * iNNumChannels_
Definition: NoiseRates.h:80
MonitorElement * hRBXEnergyType2_
Definition: NoiseRates.h:66
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
NoiseRates(const edm::ParameterSet &)
Definition: NoiseRates.cc:17
MonitorElement * hNoise_maxHPDHits_
Definition: NoiseRates.h:85
MonitorElement * iNSumEt_
Definition: NoiseRates.h:82
MonitorElement * sNSumEt_
Definition: NoiseRates.h:78
MonitorElement * hRBXEnergyType1_
Definition: NoiseRates.h:65
MonitorElement * nNSumEt_
Definition: NoiseRates.h:74
MonitorElement * iNSumE_
Definition: NoiseRates.h:81
MonitorElement * nNSumE_
Definition: NoiseRates.h:73
MonitorElement * sNSumE_
Definition: NoiseRates.h:77
bool isValid() const
Definition: HandleBase.h:70
const edm::InputTag rbxCollName_
Definition: NoiseRates.h:54
fixed size matrix
HLT enums.
const double minRBXEnergy_
Definition: NoiseRates.h:55
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: NoiseRates.cc:93
Definition: Run.h:45
ib
Definition: cuy.py:661