CMS 3D CMS Logo

CastorRecHitMonitor.cc
Go to the documentation of this file.
1 //***************************************************
2 // CastorRecHitMonitor
3 // Author: Dmytro Volyanskyy
4 // Date : 23.09.2008 (first version)
5 // last modification: Pedro Cipriano 09.07.2013
6 //----------------------------------------------
7 // critical revision 26.06.2014 (Vladimir Popov)
8 //***************************************************
9 
12 #include <string>
13 
14 using namespace std;
15 
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 }
23 
25 
27  const edm::Run &iRun,
28  const edm::EventSetup &iSetup) {
29  char s[60];
30  if (fVerbosity > 0)
31  std::cout << "CastorRecHitMonitor::bookHistograms" << std::endl;
32  ibooker.setCurrentFolder(subsystemname + "/CastorRecHitMonitor");
33 
34  const int N_Sec = 16;
35  const int nySec = 20;
36  float ySec[nySec + 1];
37  float xSec[N_Sec + 1];
38  double E0sec = 1. / 1024.;
39  ySec[0] = 0.;
40  ySec[1] = E0sec;
41  double lnBsec = log(2.);
42  for (int j = 1; j < nySec; j++)
43  ySec[j + 1] = E0sec * exp(j * lnBsec);
44  for (int i = 0; i <= N_Sec; i++)
45  xSec[i] = i;
46 
47  sprintf(s, "Castor Energy by Sectors #Phi");
48  h2RHvsSec = ibooker.book2D(s, s, N_Sec, xSec, nySec, ySec);
49  h2RHvsSec->setAxisTitle("sector #Phi");
50  h2RHvsSec->setAxisTitle("RecHit / GeV", /* axis */ 2);
51  h2RHvsSec->setOption("colz");
52 
53  const int nxCh = 224;
54  const int nyE = 18;
55  float xCh[nxCh + 1];
56  float yErh[nyE + 1];
57  for (int i = 0; i <= nxCh; i++)
58  xCh[i] = i;
59  double E0 = 1. / 1024.;
60  double lnA = log(2.);
61  yErh[0] = 0.;
62  yErh[1] = E0;
63  for (int j = 1; j < nyE; j++)
64  yErh[j + 1] = E0 * exp(j * lnA);
65 
66  string st = "Castor Cell Energy Map (cell-wise)";
67  h2RHchan = ibooker.book2D(st, st + ";moduleZ*16 + sector #Phi;RecHit / GeV", nxCh, xCh, nyE, yErh);
68  h2RHchan->setOption("colz");
69 
70  sprintf(s, "Castor Cell Energy");
71  hallchan = ibooker.book1D(s, s, nyE, yErh);
72  hallchan->setAxisTitle("GeV");
73 
74  st = "Castor cell avr Energy per event Map Z-Phi";
75  h2RHoccmap = ibooker.bookProfile2D(st, st + ";module Z;sector Phi", 14, 0, 14, 16, 0, 16, 0., 1.e10, "");
76  h2RHoccmap->getTProfile2D()->SetOption("colz");
77 
78  sprintf(s, "CastorRecHitEntriesMap");
79  h2RHentriesMap = ibooker.book2D(s, s, 14, 0, 14, 16, 0, 16);
80  h2RHentriesMap->setAxisTitle("moduleZ");
81  h2RHentriesMap->setAxisTitle("sector #Phi", /* axis */ 2);
82  h2RHentriesMap->setOption("colz");
83 
84  sprintf(s, "CastorRecHitTime");
85  hRHtime = ibooker.book1D(s, s, 301, -101., 200.);
86 
87  sprintf(s, "CASTORTowerDepth");
88  hTowerDepth = ibooker.book1D(s, s, 130, -15500., -14200.);
89  hTowerDepth->setAxisTitle("mm");
90 
91  sprintf(s, "CASTORTowerMultiplicity");
92  hTowerMultipl = ibooker.book1D(s, s, 20, 0., 20.);
93 
94  const int NEtow = 20;
95  float EhadTow[NEtow + 1];
96  float EMTow[NEtow + 1];
97  float ETower[NEtow + 2];
98  double E0tow = 1. / 1024.;
99  EMTow[0] = 0.;
100  EMTow[1] = E0tow;
101  EhadTow[0] = 0.;
102  EhadTow[1] = E0tow;
103  ETower[0] = 0.;
104  ETower[1] = E0tow;
105  double lnBtow = log(2.);
106  for (int j = 1; j < NEtow; j++)
107  EMTow[j + 1] = E0tow * exp(j * lnBtow);
108  for (int j = 1; j < NEtow; j++)
109  EhadTow[j + 1] = E0tow * exp(j * lnBtow);
110  for (int j = 1; j <= NEtow; j++)
111  ETower[j + 1] = E0tow * exp(j * lnBtow);
112 
113  sprintf(s, "CASTORTowerEMvsEhad");
114  h2TowerEMhad = ibooker.book2D(s, s, NEtow, EhadTow, NEtow, EMTow);
115  h2TowerEMhad->setAxisTitle("Ehad / GeV");
116  h2TowerEMhad->setAxisTitle("EM / GeV", /* axis */ 2);
117  h2TowerEMhad->setOption("colz");
118 
119  sprintf(s, "CASTORTowerTotalEnergy");
120  hTowerE = ibooker.book1D(s, s, NEtow + 1, ETower);
121  hTowerE->setAxisTitle("GeV");
122 
123  sprintf(s, "CASTORJetsMultiplicity");
124  hJetsMultipl = ibooker.book1D(s, s, 16, 0., 16.);
125 
126  sprintf(s, "CASTORJetEnergy");
127  hJetEnergy = ibooker.book1D(s, s, 5000, 0., 500.);
128 
129  sprintf(s, "CASTORJetEta");
130  hJetEta = ibooker.book1D(s, s, 126, -6.3, 6.3);
131 
132  sprintf(s, "CASTORJetPhi");
133  hJetPhi = ibooker.book1D(s, s, 63, -3.15, 3.15);
134 
135  if (fVerbosity > 0)
136  std::cout << "CastorRecHitMonitor::bookHistograms(end)" << std::endl;
137  return;
138 }
139 
141  if (castorTowers.empty())
142  return;
143  int nTowers = 0;
144 
145  for (reco::CastorTowerCollection::const_iterator iTower = castorTowers.begin(); iTower != castorTowers.end();
146  iTower++) {
147  hTowerE->Fill(iTower->energy() * 0.001);
148  h2TowerEMhad->Fill(iTower->hadEnergy() * 0.001, iTower->emEnergy() * 0.001);
149  hTowerDepth->Fill(iTower->depth());
150  nTowers++;
151  }
152  hTowerMultipl->Fill(nTowers);
153 }
154 
156  if (fVerbosity > 0)
157  std::cout << "CastorRecHitMonitor::processEvent (begin)" << std::endl;
158  ievt_++;
159  for (int z = 0; z < 14; z++)
160  for (int phi = 0; phi < 16; phi++)
161  energyInEachChannel[z][phi] = 0.;
162 
164  // if (showTiming) { cpu_timer.reset(); cpu_timer.start(); }
165 
166  if (castorHits.empty())
167  return;
168 
169  for (CASTORiter = castorHits.begin(); CASTORiter != castorHits.end(); ++CASTORiter) {
170  float energy = CASTORiter->energy();
171  float time = CASTORiter->time();
172  float time2 = time;
173  if (time < -100.)
174  time2 = -100.;
175  hRHtime->Fill(time2);
176 
177  HcalCastorDetId id(CASTORiter->detid().rawId());
178  // float zside = id.zside();
179  int module = (int)id.module(); //-- get module
180  int sector = (int)id.sector(); //-- get sector
181 
182  energyInEachChannel[module - 1][sector - 1] += energy;
183 
184  h2RHentriesMap->Fill(module - 1, sector - 1);
185  } // end for(CASTORiter=castorHits.begin(); CASTORiter!= ...
186 
187  double etot = 0.;
188  for (int phi = 0; phi < 16; phi++) {
189  double es = 0.;
190  for (int z = 0; z < 14; z++) {
191  float rh = energyInEachChannel[z][phi] * 0.001;
192  int ind = z * 16 + phi + 1;
193  // int ind = phi*14 + z +1;
194  h2RHchan->Fill(ind, rh);
195  hallchan->Fill(rh);
196  if (rh < 0.)
197  continue;
198  h2RHoccmap->Fill(z, phi, rh);
199  es += rh;
200  }
201  h2RHvsSec->Fill(phi, es);
202  etot += es;
203  } // end for(int phi=0;
204 
205  if (fVerbosity > 0)
206  std::cout << "CastorRecHitMonitor::processEvent (end)" << std::endl;
207  return;
208 }
209 
211  int nJets = 0;
212  for (reco::BasicJetCollection::const_iterator ibegin = Jets.begin(), iend = Jets.end(), ijet = ibegin; ijet != iend;
213  ++ijet) {
214  nJets++;
215  float energy = ijet->energy() * 0.001;
216  hJetEnergy->Fill(energy);
217  hJetEta->Fill(ijet->eta());
218  hJetPhi->Fill(ijet->phi());
219  }
220  hJetsMultipl->Fill(nJets);
221 }
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX)
Definition: DQMStore.cc:239
T getUntrackedParameter(std::string const &, T const &) const
virtual void setOption(const char *option)
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
std::vector< T >::const_iterator const_iterator
void processEventJets(const reco::BasicJetCollection &Jets)
CastorRecHitMonitor(const edm::ParameterSet &ps)
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")
Definition: DQMStore.cc:381
virtual TProfile2D * getTProfile2D() const
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &)
const_iterator end() const
void processEventTowers(const reco::CastorTowerCollection &castorTowers)
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Definition: DQMStore.cc:266
void processEvent(const CastorRecHitCollection &castorHits)
std::vector< CastorTower > CastorTowerCollection
collection of CastorTower objects
Definition: CastorTower.h:137
Definition: vlib.h:198
const_iterator begin() const
Definition: Run.h:45
std::vector< BasicJet > BasicJetCollection
collection of BasicJet objects
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)