CMS 3D CMS Logo

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

#include <CastorRecHitMonitor.h>

Public Types

typedef dqm::legacy::DQMStore DQMStore
 
typedef dqm::legacy::MonitorElement MonitorElement
 

Public Member Functions

void bookHistograms (DQMStore::IBooker &, edm::Run 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 19 of file CastorRecHitMonitor.h.

Member Typedef Documentation

◆ DQMStore

Definition at line 21 of file CastorRecHitMonitor.h.

◆ MonitorElement

Definition at line 22 of file CastorRecHitMonitor.h.

Constructor & Destructor Documentation

◆ CastorRecHitMonitor()

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

Definition at line 16 of file CastorRecHitMonitor.cc.

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

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

◆ ~CastorRecHitMonitor()

CastorRecHitMonitor::~CastorRecHitMonitor ( )

Definition at line 24 of file CastorRecHitMonitor.cc.

24 {}

Member Function Documentation

◆ bookHistograms()

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

Definition at line 26 of file CastorRecHitMonitor.cc.

References dqm::implementation::IBooker::book1D(), dqm::implementation::IBooker::book2D(), dqm::implementation::IBooker::bookProfile2D(), gather_cfg::cout, JetChargeProducer_cfi::exp, dqm::legacy::MonitorElement::getTProfile2D(), mps_fire::i, dqmiolumiharvest::j, dqm-mbProfile::log, alignCSCRings::s, dqm::impl::MonitorElement::setAxisTitle(), dqm::implementation::NavigatorBase::setCurrentFolder(), and dqm::impl::MonitorElement::setOption().

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

◆ processEvent()

void CastorRecHitMonitor::processEvent ( const CastorRecHitCollection castorHits)

Definition at line 153 of file CastorRecHitMonitor.cc.

References edm::SortedCollection< T, SORT >::begin(), gather_cfg::cout, edm::SortedCollection< T, SORT >::empty(), edm::SortedCollection< T, SORT >::end(), hcalRecHitTable_cff::energy, l1ctLayer2EG_cff::id, createfilelist::int, callgraph::module, nano_mu_digi_cff::sector, and hcalRecHitTable_cff::time.

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

◆ processEventJets()

void CastorRecHitMonitor::processEventJets ( const reco::BasicJetCollection Jets)

Definition at line 206 of file CastorRecHitMonitor.cc.

References hcalRecHitTable_cff::energy, METSkim_cff::Jets, and nJets.

206  {
207  int nJets = 0;
208  for (reco::BasicJetCollection::const_iterator ibegin = Jets.begin(), iend = Jets.end(), ijet = ibegin; ijet != iend;
209  ++ijet) {
210  nJets++;
211  float energy = ijet->energy() * 0.001;
213  hJetEta->Fill(ijet->eta());
214  hJetPhi->Fill(ijet->phi());
215  }
217 }
MonitorElement * hJetEta
static constexpr int nJets
void Fill(long long x)
MonitorElement * hJetEnergy
MonitorElement * hJetsMultipl
MonitorElement * hJetPhi

◆ processEventTowers()

void CastorRecHitMonitor::processEventTowers ( const reco::CastorTowerCollection castorTowers)

Definition at line 138 of file CastorRecHitMonitor.cc.

References ecaldqm::nTowers.

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

Member Data Documentation

◆ energyInEachChannel

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

Definition at line 34 of file CastorRecHitMonitor.h.

◆ fVerbosity

int CastorRecHitMonitor::fVerbosity = 0
private

Definition at line 32 of file CastorRecHitMonitor.h.

◆ h2RecHitMap

TH2F* CastorRecHitMonitor::h2RecHitMap
private

Definition at line 46 of file CastorRecHitMonitor.h.

◆ h2RHchan

MonitorElement* CastorRecHitMonitor::h2RHchan
private

Definition at line 47 of file CastorRecHitMonitor.h.

◆ h2RHentriesMap

MonitorElement* CastorRecHitMonitor::h2RHentriesMap
private

Definition at line 51 of file CastorRecHitMonitor.h.

◆ h2RHmap

MonitorElement* CastorRecHitMonitor::h2RHmap
private

Definition at line 49 of file CastorRecHitMonitor.h.

◆ h2RHoccmap

MonitorElement* CastorRecHitMonitor::h2RHoccmap
private

Definition at line 50 of file CastorRecHitMonitor.h.

◆ h2RHvsSec

MonitorElement* CastorRecHitMonitor::h2RHvsSec
private

Definition at line 48 of file CastorRecHitMonitor.h.

◆ h2TowerEMhad

MonitorElement* CastorRecHitMonitor::h2TowerEMhad
private

Definition at line 39 of file CastorRecHitMonitor.h.

◆ hallchan

MonitorElement * CastorRecHitMonitor::hallchan
private

Definition at line 52 of file CastorRecHitMonitor.h.

◆ hJetEnergy

MonitorElement* CastorRecHitMonitor::hJetEnergy
private

Definition at line 42 of file CastorRecHitMonitor.h.

◆ hJetEta

MonitorElement* CastorRecHitMonitor::hJetEta
private

Definition at line 43 of file CastorRecHitMonitor.h.

◆ hJetPhi

MonitorElement* CastorRecHitMonitor::hJetPhi
private

Definition at line 44 of file CastorRecHitMonitor.h.

◆ hJetsMultipl

MonitorElement* CastorRecHitMonitor::hJetsMultipl
private

Definition at line 41 of file CastorRecHitMonitor.h.

◆ hRHtime

MonitorElement* CastorRecHitMonitor::hRHtime
private

Definition at line 52 of file CastorRecHitMonitor.h.

◆ hTowerDepth

MonitorElement* CastorRecHitMonitor::hTowerDepth
private

Definition at line 38 of file CastorRecHitMonitor.h.

◆ hTowerE

MonitorElement* CastorRecHitMonitor::hTowerE
private

Definition at line 37 of file CastorRecHitMonitor.h.

◆ hTowerMultipl

MonitorElement* CastorRecHitMonitor::hTowerMultipl
private

Definition at line 40 of file CastorRecHitMonitor.h.

◆ ievt_

int CastorRecHitMonitor::ievt_
private

Definition at line 33 of file CastorRecHitMonitor.h.

◆ subsystemname

std::string CastorRecHitMonitor::subsystemname
private

Definition at line 35 of file CastorRecHitMonitor.h.