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