File indexing completed on 2024-04-06 12:06:54
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include "DQM/CastorMonitor/interface/CastorRecHitMonitor.h"
0011 #include "DQMServices/Core/interface/DQMStore.h"
0012 #include <string>
0013
0014 using namespace std;
0015
0016 CastorRecHitMonitor::CastorRecHitMonitor(const edm::ParameterSet &ps) {
0017 fVerbosity = ps.getUntrackedParameter<int>("debug", 0);
0018 if (fVerbosity > 0)
0019 std::cout << "CastorRecHitMonitor Constructor: " << this << std::endl;
0020 subsystemname = ps.getUntrackedParameter<std::string>("subSystemFolder", "Castor");
0021 ievt_ = 0;
0022 }
0023
0024 CastorRecHitMonitor::~CastorRecHitMonitor() {}
0025
0026 void CastorRecHitMonitor::bookHistograms(DQMStore::IBooker &ibooker, const edm::Run &iRun) {
0027 char s[60];
0028 if (fVerbosity > 0)
0029 std::cout << "CastorRecHitMonitor::bookHistograms" << std::endl;
0030 ibooker.setCurrentFolder(subsystemname + "/CastorRecHitMonitor");
0031
0032 const int N_Sec = 16;
0033 const int nySec = 20;
0034 float ySec[nySec + 1];
0035 float xSec[N_Sec + 1];
0036 double E0sec = 1. / 1024.;
0037 ySec[0] = 0.;
0038 ySec[1] = E0sec;
0039 double lnBsec = log(2.);
0040 for (int j = 1; j < nySec; j++)
0041 ySec[j + 1] = E0sec * exp(j * lnBsec);
0042 for (int i = 0; i <= N_Sec; i++)
0043 xSec[i] = i;
0044
0045 sprintf(s, "Castor Energy by Sectors #Phi");
0046 h2RHvsSec = ibooker.book2D(s, s, N_Sec, xSec, nySec, ySec);
0047 h2RHvsSec->setAxisTitle("sector #Phi");
0048 h2RHvsSec->setAxisTitle("RecHit / GeV", 2);
0049 h2RHvsSec->setOption("colz");
0050
0051 const int nxCh = 224;
0052 const int nyE = 18;
0053 float xCh[nxCh + 1];
0054 float yErh[nyE + 1];
0055 for (int i = 0; i <= nxCh; i++)
0056 xCh[i] = i;
0057 double E0 = 1. / 1024.;
0058 double lnA = log(2.);
0059 yErh[0] = 0.;
0060 yErh[1] = E0;
0061 for (int j = 1; j < nyE; j++)
0062 yErh[j + 1] = E0 * exp(j * lnA);
0063
0064 string st = "Castor Cell Energy Map (cell-wise)";
0065 h2RHchan = ibooker.book2D(st, st + ";moduleZ*16 + sector #Phi;RecHit / GeV", nxCh, xCh, nyE, yErh);
0066 h2RHchan->setOption("colz");
0067
0068 sprintf(s, "Castor Cell Energy");
0069 hallchan = ibooker.book1D(s, s, nyE, yErh);
0070 hallchan->setAxisTitle("GeV");
0071
0072 st = "Castor cell avr Energy per event Map Z-Phi";
0073 h2RHoccmap = ibooker.bookProfile2D(st, st + ";module Z;sector Phi", 14, 0, 14, 16, 0, 16, 0., 1.e10, "");
0074 h2RHoccmap->getTProfile2D()->SetOption("colz");
0075
0076 sprintf(s, "CastorRecHitEntriesMap");
0077 h2RHentriesMap = ibooker.book2D(s, s, 14, 0, 14, 16, 0, 16);
0078 h2RHentriesMap->setAxisTitle("moduleZ");
0079 h2RHentriesMap->setAxisTitle("sector #Phi", 2);
0080 h2RHentriesMap->setOption("colz");
0081
0082 sprintf(s, "CastorRecHitTime");
0083 hRHtime = ibooker.book1D(s, s, 301, -101., 200.);
0084
0085 sprintf(s, "CASTORTowerDepth");
0086 hTowerDepth = ibooker.book1D(s, s, 130, -15500., -14200.);
0087 hTowerDepth->setAxisTitle("mm");
0088
0089 sprintf(s, "CASTORTowerMultiplicity");
0090 hTowerMultipl = ibooker.book1D(s, s, 20, 0., 20.);
0091
0092 const int NEtow = 20;
0093 float EhadTow[NEtow + 1];
0094 float EMTow[NEtow + 1];
0095 float ETower[NEtow + 2];
0096 double E0tow = 1. / 1024.;
0097 EMTow[0] = 0.;
0098 EMTow[1] = E0tow;
0099 EhadTow[0] = 0.;
0100 EhadTow[1] = E0tow;
0101 ETower[0] = 0.;
0102 ETower[1] = E0tow;
0103 double lnBtow = log(2.);
0104 for (int j = 1; j < NEtow; j++)
0105 EMTow[j + 1] = E0tow * exp(j * lnBtow);
0106 for (int j = 1; j < NEtow; j++)
0107 EhadTow[j + 1] = E0tow * exp(j * lnBtow);
0108 for (int j = 1; j <= NEtow; j++)
0109 ETower[j + 1] = E0tow * exp(j * lnBtow);
0110
0111 sprintf(s, "CASTORTowerEMvsEhad");
0112 h2TowerEMhad = ibooker.book2D(s, s, NEtow, EhadTow, NEtow, EMTow);
0113 h2TowerEMhad->setAxisTitle("Ehad / GeV");
0114 h2TowerEMhad->setAxisTitle("EM / GeV", 2);
0115 h2TowerEMhad->setOption("colz");
0116
0117 sprintf(s, "CASTORTowerTotalEnergy");
0118 hTowerE = ibooker.book1D(s, s, NEtow + 1, ETower);
0119 hTowerE->setAxisTitle("GeV");
0120
0121 sprintf(s, "CASTORJetsMultiplicity");
0122 hJetsMultipl = ibooker.book1D(s, s, 16, 0., 16.);
0123
0124 sprintf(s, "CASTORJetEnergy");
0125 hJetEnergy = ibooker.book1D(s, s, 5000, 0., 500.);
0126
0127 sprintf(s, "CASTORJetEta");
0128 hJetEta = ibooker.book1D(s, s, 126, -6.3, 6.3);
0129
0130 sprintf(s, "CASTORJetPhi");
0131 hJetPhi = ibooker.book1D(s, s, 63, -3.15, 3.15);
0132
0133 if (fVerbosity > 0)
0134 std::cout << "CastorRecHitMonitor::bookHistograms(end)" << std::endl;
0135 return;
0136 }
0137
0138 void CastorRecHitMonitor::processEventTowers(const reco::CastorTowerCollection &castorTowers) {
0139 if (castorTowers.empty())
0140 return;
0141 int nTowers = 0;
0142
0143 for (reco::CastorTowerCollection::const_iterator iTower = castorTowers.begin(); iTower != castorTowers.end();
0144 iTower++) {
0145 hTowerE->Fill(iTower->energy() * 0.001);
0146 h2TowerEMhad->Fill(iTower->hadEnergy() * 0.001, iTower->emEnergy() * 0.001);
0147 hTowerDepth->Fill(iTower->depth());
0148 nTowers++;
0149 }
0150 hTowerMultipl->Fill(nTowers);
0151 }
0152
0153 void CastorRecHitMonitor::processEvent(const CastorRecHitCollection &castorHits) {
0154 if (fVerbosity > 0)
0155 std::cout << "CastorRecHitMonitor::processEvent (begin)" << std::endl;
0156 ievt_++;
0157 for (int z = 0; z < 14; z++)
0158 for (int phi = 0; phi < 16; phi++)
0159 energyInEachChannel[z][phi] = 0.;
0160
0161 CastorRecHitCollection::const_iterator CASTORiter;
0162
0163
0164 if (castorHits.empty())
0165 return;
0166
0167 for (CASTORiter = castorHits.begin(); CASTORiter != castorHits.end(); ++CASTORiter) {
0168 float energy = CASTORiter->energy();
0169 float time = CASTORiter->time();
0170 float time2 = time;
0171 if (time < -100.)
0172 time2 = -100.;
0173 hRHtime->Fill(time2);
0174
0175 HcalCastorDetId id(CASTORiter->detid().rawId());
0176
0177 int module = (int)id.module();
0178 int sector = (int)id.sector();
0179
0180 energyInEachChannel[module - 1][sector - 1] += energy;
0181
0182 h2RHentriesMap->Fill(module - 1, sector - 1);
0183 }
0184
0185 for (int phi = 0; phi < 16; phi++) {
0186 double es = 0.;
0187 for (int z = 0; z < 14; z++) {
0188 float rh = energyInEachChannel[z][phi] * 0.001;
0189 int ind = z * 16 + phi + 1;
0190
0191 h2RHchan->Fill(ind, rh);
0192 hallchan->Fill(rh);
0193 if (rh < 0.)
0194 continue;
0195 h2RHoccmap->Fill(z, phi, rh);
0196 es += rh;
0197 }
0198 h2RHvsSec->Fill(phi, es);
0199 }
0200
0201 if (fVerbosity > 0)
0202 std::cout << "CastorRecHitMonitor::processEvent (end)" << std::endl;
0203 return;
0204 }
0205
0206 void CastorRecHitMonitor::processEventJets(const reco::BasicJetCollection &Jets) {
0207 int nJets = 0;
0208 for (reco::BasicJetCollection::const_iterator ibegin = Jets.begin(), iend = Jets.end(), ijet = ibegin; ijet != iend;
0209 ++ijet) {
0210 nJets++;
0211 float energy = ijet->energy() * 0.001;
0212 hJetEnergy->Fill(energy);
0213 hJetEta->Fill(ijet->eta());
0214 hJetPhi->Fill(ijet->phi());
0215 }
0216 hJetsMultipl->Fill(nJets);
0217 }