CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
CastorRecHitMonitor Class Reference

#include <CastorRecHitMonitor.h>

Public Member Functions

void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &)
 
 CastorRecHitMonitor (const edm::ParameterSet &ps)
 
void processEvent (const CastorRecHitCollection &castorHits)
 
void processEventJets (const reco::BasicJetCollection &Jets)
 
void processEventTowers (const reco::CastorTowerCollection &castorTowers)
 
 ~CastorRecHitMonitor ()
 

Private Attributes

float energyInEachChannel [14][16]
 
int fVerbosity = 0
 
TH2F * h2RecHitMap
 
MonitorElementh2RHchan
 
MonitorElementh2RHentriesMap
 
MonitorElementh2RHmap
 
MonitorElementh2RHoccmap
 
MonitorElementh2RHvsSec
 
MonitorElementh2TowerEMhad
 
MonitorElementhallchan
 
MonitorElementhJetEnergy
 
MonitorElementhJetEta
 
MonitorElementhJetPhi
 
MonitorElementhJetsMultipl
 
MonitorElementhRHtime
 
MonitorElementhTowerDepth
 
MonitorElementhTowerE
 
MonitorElementhTowerMultipl
 
int ievt_
 
std::string subsystemname
 

Detailed Description

Definition at line 21 of file CastorRecHitMonitor.h.

Constructor & Destructor Documentation

CastorRecHitMonitor::CastorRecHitMonitor ( const edm::ParameterSet ps)

Definition at line 17 of file CastorRecHitMonitor.cc.

References gather_cfg::cout, edm::ParameterSet::getUntrackedParameter(), and AlCaHLTBitMon_QueryRunRegistry::string.

17  {
18  fVerbosity = ps.getUntrackedParameter<int>("debug", 0);
19  if (fVerbosity > 0)
20  std::cout << "CastorRecHitMonitor Constructor: " << this << std::endl;
21  subsystemname = ps.getUntrackedParameter<std::string>("subSystemFolder", "Castor");
22  ievt_ = 0;
23 }
T getUntrackedParameter(std::string const &, T const &) const
CastorRecHitMonitor::~CastorRecHitMonitor ( )

Definition at line 25 of file CastorRecHitMonitor.cc.

25 {}

Member Function Documentation

void CastorRecHitMonitor::bookHistograms ( DQMStore::IBooker ,
edm::Run const &  ,
edm::EventSetup const &   
)

Definition at line 27 of file CastorRecHitMonitor.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::book2D(), DQMStore::IBooker::bookProfile2D(), gather_cfg::cout, JetChargeProducer_cfi::exp, MonitorElement::getTH1F(), MonitorElement::getTH2F(), MonitorElement::getTProfile2D(), mps_fire::i, cmsBatch::log, alignCSCRings::s, and DQMStore::IBooker::setCurrentFolder().

29  {
30  char s[60];
31  if (fVerbosity > 0)
32  std::cout << "CastorRecHitMonitor::bookHistograms" << std::endl;
33  ibooker.setCurrentFolder(subsystemname + "/CastorRecHitMonitor");
34 
35  const int N_Sec = 16;
36  const int nySec = 20;
37  float ySec[nySec + 1];
38  float xSec[N_Sec + 1];
39  double E0sec = 1. / 1024.;
40  ySec[0] = 0.;
41  ySec[1] = E0sec;
42  double lnBsec = log(2.);
43  for (int j = 1; j < nySec; j++)
44  ySec[j + 1] = E0sec * exp(j * lnBsec);
45  for (int i = 0; i <= N_Sec; i++)
46  xSec[i] = i;
47 
48  sprintf(s, "Castor Energy by Sectors #Phi");
49  h2RHvsSec = ibooker.book2D(s, s, N_Sec, xSec, nySec, ySec);
50  h2RHvsSec->getTH2F()->GetXaxis()->SetTitle("sector #Phi");
51  h2RHvsSec->getTH2F()->GetYaxis()->SetTitle("RecHit / GeV");
52  h2RHvsSec->getTH2F()->SetOption("colz");
53 
54  const int nxCh = 224;
55  const int nyE = 18;
56  float xCh[nxCh + 1];
57  float yErh[nyE + 1];
58  for (int i = 0; i <= nxCh; i++)
59  xCh[i] = i;
60  double E0 = 1. / 1024.;
61  double lnA = log(2.);
62  yErh[0] = 0.;
63  yErh[1] = E0;
64  for (int j = 1; j < nyE; j++)
65  yErh[j + 1] = E0 * exp(j * lnA);
66 
67  string st = "Castor Cell Energy Map (cell-wise)";
68  h2RHchan = ibooker.book2D(st, st + ";moduleZ*16 + sector #Phi;RecHit / GeV", nxCh, xCh, nyE, yErh);
69  h2RHchan->getTH2F()->SetOption("colz");
70 
71  sprintf(s, "Castor Cell Energy");
72  hallchan = ibooker.book1D(s, s, nyE, yErh);
73  hallchan->getTH1F()->GetXaxis()->SetTitle("GeV");
74 
75  st = "Castor cell avr Energy per event Map Z-Phi";
76  h2RHoccmap = ibooker.bookProfile2D(st, st + ";module Z;sector Phi", 14, 0, 14, 16, 0, 16, 0., 1.e10, "");
77  h2RHoccmap->getTProfile2D()->SetOption("colz");
78 
79  sprintf(s, "CastorRecHitEntriesMap");
80  h2RHentriesMap = ibooker.book2D(s, s, 14, 0, 14, 16, 0, 16);
81  h2RHentriesMap->getTH2F()->GetXaxis()->SetTitle("moduleZ");
82  h2RHentriesMap->getTH2F()->GetYaxis()->SetTitle("sector #Phi");
83  h2RHentriesMap->getTH2F()->SetOption("colz");
84 
85  sprintf(s, "CastorRecHitTime");
86  hRHtime = ibooker.book1D(s, s, 301, -101., 200.);
87 
88  sprintf(s, "CASTORTowerDepth");
89  hTowerDepth = ibooker.book1D(s, s, 130, -15500., -14200.);
90  hTowerDepth->getTH1F()->GetXaxis()->SetTitle("mm");
91 
92  sprintf(s, "CASTORTowerMultiplicity");
93  hTowerMultipl = ibooker.book1D(s, s, 20, 0., 20.);
94 
95  const int NEtow = 20;
96  float EhadTow[NEtow + 1];
97  float EMTow[NEtow + 1];
98  float ETower[NEtow + 2];
99  double E0tow = 1. / 1024.;
100  EMTow[0] = 0.;
101  EMTow[1] = E0tow;
102  EhadTow[0] = 0.;
103  EhadTow[1] = E0tow;
104  ETower[0] = 0.;
105  ETower[1] = E0tow;
106  double lnBtow = log(2.);
107  for (int j = 1; j < NEtow; j++)
108  EMTow[j + 1] = E0tow * exp(j * lnBtow);
109  for (int j = 1; j < NEtow; j++)
110  EhadTow[j + 1] = E0tow * exp(j * lnBtow);
111  for (int j = 1; j <= NEtow; j++)
112  ETower[j + 1] = E0tow * exp(j * lnBtow);
113 
114  sprintf(s, "CASTORTowerEMvsEhad");
115  h2TowerEMhad = ibooker.book2D(s, s, NEtow, EhadTow, NEtow, EMTow);
116  h2TowerEMhad->getTH2F()->GetXaxis()->SetTitle("Ehad / GeV");
117  h2TowerEMhad->getTH2F()->GetYaxis()->SetTitle("EM / GeV");
118  h2TowerEMhad->getTH2F()->SetOption("colz");
119 
120  sprintf(s, "CASTORTowerTotalEnergy");
121  hTowerE = ibooker.book1D(s, s, NEtow + 1, ETower);
122  hTowerE->getTH1F()->GetXaxis()->SetTitle("GeV");
123 
124  sprintf(s, "CASTORJetsMultiplicity");
125  hJetsMultipl = ibooker.book1D(s, s, 16, 0., 16.);
126 
127  sprintf(s, "CASTORJetEnergy");
128  hJetEnergy = ibooker.book1D(s, s, 5000, 0., 500.);
129 
130  sprintf(s, "CASTORJetEta");
131  hJetEta = ibooker.book1D(s, s, 126, -6.3, 6.3);
132 
133  sprintf(s, "CASTORJetPhi");
134  hJetPhi = ibooker.book1D(s, s, 63, -3.15, 3.15);
135 
136  if (fVerbosity > 0)
137  std::cout << "CastorRecHitMonitor::bookHistograms(end)" << std::endl;
138  return;
139 }
TProfile2D * getTProfile2D() const
TH1F * getTH1F() const
MonitorElement * hJetEta
MonitorElement * h2RHvsSec
MonitorElement * h2RHentriesMap
MonitorElement * hJetEnergy
MonitorElement * hJetsMultipl
MonitorElement * hRHtime
MonitorElement * hTowerDepth
TH2F * getTH2F() const
MonitorElement * hTowerE
MonitorElement * h2RHoccmap
MonitorElement * hJetPhi
MonitorElement * h2TowerEMhad
MonitorElement * hallchan
MonitorElement * h2RHchan
MonitorElement * hTowerMultipl
void CastorRecHitMonitor::processEvent ( const CastorRecHitCollection castorHits)

Definition at line 156 of file CastorRecHitMonitor.cc.

References edm::SortedCollection< T, SORT >::begin(), gather_cfg::cout, edm::SortedCollection< T, SORT >::empty(), edm::SortedCollection< T, SORT >::end(), triggerObjects_cff::id, createfilelist::int, and ntuplemaker::time.

156  {
157  if (fVerbosity > 0)
158  std::cout << "CastorRecHitMonitor::processEvent (begin)" << std::endl;
159  ievt_++;
160  for (int z = 0; z < 14; z++)
161  for (int phi = 0; phi < 16; phi++)
162  energyInEachChannel[z][phi] = 0.;
163 
165  // if (showTiming) { cpu_timer.reset(); cpu_timer.start(); }
166 
167  if (castorHits.empty())
168  return;
169 
170  for (CASTORiter = castorHits.begin(); CASTORiter != castorHits.end(); ++CASTORiter) {
171  float energy = CASTORiter->energy();
172  float time = CASTORiter->time();
173  float time2 = time;
174  if (time < -100.)
175  time2 = -100.;
176  hRHtime->Fill(time2);
177 
178  HcalCastorDetId id(CASTORiter->detid().rawId());
179  // float zside = id.zside();
180  int module = (int)id.module(); //-- get module
181  int sector = (int)id.sector(); //-- get sector
182 
183  energyInEachChannel[module - 1][sector - 1] += energy;
184 
185  h2RHentriesMap->Fill(module - 1, sector - 1);
186  } // end for(CASTORiter=castorHits.begin(); CASTORiter!= ...
187 
188  double etot = 0.;
189  for (int phi = 0; phi < 16; phi++) {
190  double es = 0.;
191  for (int z = 0; z < 14; z++) {
192  float rh = energyInEachChannel[z][phi] * 0.001;
193  int ind = z * 16 + phi + 1;
194  // int ind = phi*14 + z +1;
195  h2RHchan->Fill(ind, rh);
196  hallchan->Fill(rh);
197  if (rh < 0.)
198  continue;
199  h2RHoccmap->Fill(z, phi, rh);
200  es += rh;
201  }
202  h2RHvsSec->Fill(phi, es);
203  etot += es;
204  } // end for(int phi=0;
205 
206  if (fVerbosity > 0)
207  std::cout << "CastorRecHitMonitor::processEvent (end)" << std::endl;
208  return;
209 }
float energyInEachChannel[14][16]
std::vector< T >::const_iterator const_iterator
MonitorElement * h2RHvsSec
MonitorElement * h2RHentriesMap
void Fill(long long x)
MonitorElement * hRHtime
const_iterator end() const
MonitorElement * h2RHoccmap
MonitorElement * hallchan
MonitorElement * h2RHchan
Definition: vlib.h:208
const_iterator begin() const
void CastorRecHitMonitor::processEventJets ( const reco::BasicJetCollection Jets)

Definition at line 211 of file CastorRecHitMonitor.cc.

211  {
212  int nJets = 0;
213  for (reco::BasicJetCollection::const_iterator ibegin = Jets.begin(), iend = Jets.end(), ijet = ibegin; ijet != iend;
214  ++ijet) {
215  nJets++;
216  float energy = ijet->energy() * 0.001;
217  hJetEnergy->Fill(energy);
218  hJetEta->Fill(ijet->eta());
219  hJetPhi->Fill(ijet->phi());
220  }
221  hJetsMultipl->Fill(nJets);
222 }
MonitorElement * hJetEta
void Fill(long long x)
MonitorElement * hJetEnergy
MonitorElement * hJetsMultipl
MonitorElement * hJetPhi
void CastorRecHitMonitor::processEventTowers ( const reco::CastorTowerCollection castorTowers)

Definition at line 141 of file CastorRecHitMonitor.cc.

References ecaldqm::nTowers.

141  {
142  if (castorTowers.empty())
143  return;
144  int nTowers = 0;
145 
146  for (reco::CastorTowerCollection::const_iterator iTower = castorTowers.begin(); iTower != castorTowers.end();
147  iTower++) {
148  hTowerE->Fill(iTower->energy() * 0.001);
149  h2TowerEMhad->Fill(iTower->hadEnergy() * 0.001, iTower->emEnergy() * 0.001);
150  hTowerDepth->Fill(iTower->depth());
151  nTowers++;
152  }
153  hTowerMultipl->Fill(nTowers);
154 }
void Fill(long long x)
MonitorElement * hTowerDepth
MonitorElement * hTowerE
MonitorElement * h2TowerEMhad
MonitorElement * hTowerMultipl

Member Data Documentation

float CastorRecHitMonitor::energyInEachChannel[14][16]
private

Definition at line 34 of file CastorRecHitMonitor.h.

int CastorRecHitMonitor::fVerbosity = 0
private

Definition at line 32 of file CastorRecHitMonitor.h.

TH2F* CastorRecHitMonitor::h2RecHitMap
private

Definition at line 46 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::h2RHchan
private

Definition at line 47 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::h2RHentriesMap
private

Definition at line 51 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::h2RHmap
private

Definition at line 49 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::h2RHoccmap
private

Definition at line 50 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::h2RHvsSec
private

Definition at line 48 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::h2TowerEMhad
private

Definition at line 39 of file CastorRecHitMonitor.h.

MonitorElement * CastorRecHitMonitor::hallchan
private

Definition at line 52 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::hJetEnergy
private

Definition at line 42 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::hJetEta
private

Definition at line 43 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::hJetPhi
private

Definition at line 44 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::hJetsMultipl
private

Definition at line 41 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::hRHtime
private

Definition at line 52 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::hTowerDepth
private

Definition at line 38 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::hTowerE
private

Definition at line 37 of file CastorRecHitMonitor.h.

MonitorElement* CastorRecHitMonitor::hTowerMultipl
private

Definition at line 40 of file CastorRecHitMonitor.h.

int CastorRecHitMonitor::ievt_
private

Definition at line 33 of file CastorRecHitMonitor.h.

std::string CastorRecHitMonitor::subsystemname
private

Definition at line 35 of file CastorRecHitMonitor.h.