CMS 3D CMS Logo

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