00001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00002 #include "FWCore/ServiceRegistry/interface/Service.h"
00003 #include "FWCore/Framework/interface/MakerMacros.h"
00004
00005 #include "Validation/HcalDigis/interface/HcalDigiTester.h"
00006 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00007 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
00008 #include "DataFormats/HcalDetId/interface/HcalElectronicsId.h"
00009 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00010 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00011 #include "DataFormats/HcalDigi/interface/HcalQIESample.h"
00012
00013 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00014 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00015 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00016 #include "SimDataFormats/CaloHit/interface/PCaloHitContainer.h"
00017
00018 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00019 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
00020 #include "CalibFormats/HcalObjects/interface/HcalCalibrations.h"
00021
00022 #include "DQMServices/Core/interface/DQMStore.h"
00023
00024 #include "DataFormats/HcalDigi/interface/HBHEDataFrame.h"
00025 #include "DataFormats/HcalDigi/interface/HFDataFrame.h"
00026 #include "DataFormats/HcalDigi/interface/HODataFrame.h"
00027 #include <vector>
00028 #include <utility>
00029 #include <ostream>
00030 #include <string>
00031 #include <algorithm>
00032 #include <cmath>
00033
00034 #include "CondFormats/HcalObjects/interface/HcalGain.h"
00035 #include "CondFormats/HcalObjects/interface/HcalGainWidth.h"
00036 #include "CondFormats/HcalObjects/interface/HcalPedestal.h"
00037 #include "CondFormats/HcalObjects/interface/HcalPedestalWidth.h"
00038
00039 template<class Digi>
00040
00041 void HcalDigiTester::reco(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00042
00043
00044 typename edm::Handle<edm::SortedCollection<Digi> > digiCollection;
00045 typename edm::SortedCollection<Digi>::const_iterator digiItr;
00046
00047
00048 HcalCalibrations calibrations;
00049 CaloSamples tool;
00050
00051
00052 iEvent.getByLabel (inputTag_, digiCollection);
00053
00054 int subdet = 0;
00055 if (hcalselector_ == "HB" ) subdet = 1;
00056 if (hcalselector_ == "HE" ) subdet = 2;
00057 if (hcalselector_ == "HO" ) subdet = 3;
00058 if (hcalselector_ == "HF" ) subdet = 4;
00059
00060 if(subdet == 1) nevent1++;
00061 if(subdet == 2) nevent2++;
00062 if(subdet == 3) nevent3++;
00063 if(subdet == 4) nevent4++;
00064
00065 int zsign = 0;
00066 if (zside_ == "+") zsign = 1;
00067 if (zside_ == "-") zsign = -1;
00068
00069 int ndigis = 0;
00070
00071 double ampl1_c = 0.;
00072 double ampl2_c = 0.;
00073 double ampl3_c = 0.;
00074 double ampl4_c = 0.;
00075 double ampl_c = 0.;
00076
00077
00078 int seedSimHit = 0;
00079
00080
00081
00082
00083 int ieta_Sim = 9999;
00084 int iphi_Sim = 9999;
00085 double emax_Sim = -9999.;
00086
00087
00088
00089 if( mc_ == "yes") {
00090 edm::Handle<edm::PCaloHitContainer> hcalHits ;
00091 iEvent.getByLabel("g4SimHits","HcalHits",hcalHits);
00092 const edm::PCaloHitContainer * simhitResult = hcalHits.product () ;
00093
00094
00095 if ( subdet != 0 && noise_ == 0) {
00096
00097 for (std::vector<PCaloHit>::const_iterator simhits = simhitResult->begin (); simhits != simhitResult->end () ; ++simhits) {
00098
00099 HcalDetId cell(simhits->id());
00100 double en = simhits->energy();
00101 int sub = cell.subdet();
00102 int ieta = cell.ieta();
00103 if(ieta > 0) ieta--;
00104 int iphi = cell.iphi()-1;
00105
00106
00107 if(en > emax_Sim && sub == subdet) {
00108 emax_Sim = en;
00109 ieta_Sim = ieta;
00110 iphi_Sim = iphi;
00111
00112 if (mode_ == "multi" &&
00113 ((sub == 4 && en < 100. && en > 1.)
00114 || ((sub !=4) && en < 1. && en > 0.02)))
00115 {
00116 seedSimHit = 1;
00117 break;
00118 }
00119 }
00120
00121 }
00122
00123
00124
00125 if(mode_ != "multi" && emax_Sim > 0.) seedSimHit = 1;
00126 }
00127
00128 }
00129
00130
00131 int Ndig = 0;
00132
00133
00134
00135
00136
00137
00138
00139 for (digiItr=digiCollection->begin();digiItr!=digiCollection->end();digiItr++) {
00140
00141 HcalDetId cell(digiItr->id());
00142 int depth = cell.depth();
00143 int iphi = cell.iphi()-1;
00144 int ieta = cell.ieta();
00145 if(ieta > 0) ieta--;
00146 int sub = cell.subdet();
00147
00148
00149
00150 double ampl = 0.;
00151 double ampl1 = 0.;
00152 double ampl2 = 0.;
00153 double ampl3 = 0.;
00154 double ampl4 = 0.;
00155
00156
00157 if ( ((nevent1 == 1 && subdet == 1) ||
00158 (nevent2 == 1 && subdet == 2) ||
00159 (nevent3 == 1 && subdet == 3) ||
00160 (nevent4 == 1 && subdet == 4)) && noise_ == 1 && sub == subdet) {
00161
00162 HcalGenericDetId hcalGenDetId(digiItr->id());
00163 const HcalPedestal* pedestal = conditions->getPedestal(hcalGenDetId);
00164 const HcalGain* gain = conditions->getGain(hcalGenDetId);
00165 const HcalGainWidth* gainWidth =
00166 conditions->getGainWidth(hcalGenDetId);
00167 const HcalPedestalWidth* pedWidth =
00168 conditions-> getPedestalWidth(hcalGenDetId);
00169
00170 double gainValue0 = gain->getValue(0);
00171 double gainValue1 = gain->getValue(1);
00172 double gainValue2 = gain->getValue(2);
00173 double gainValue3 = gain->getValue(3);
00174
00175 double gainWidthValue0 = gainWidth->getValue(0);
00176 double gainWidthValue1 = gainWidth->getValue(1);
00177 double gainWidthValue2 = gainWidth->getValue(2);
00178 double gainWidthValue3 = gainWidth->getValue(3);
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191 double pedValue0 = pedestal->getValue(0);
00192 double pedValue1 = pedestal->getValue(1);
00193 double pedValue2 = pedestal->getValue(2);
00194 double pedValue3 = pedestal->getValue(3);
00195
00196 double pedWidth0 = pedWidth->getWidth(0);
00197 double pedWidth1 = pedWidth->getWidth(1);
00198 double pedWidth2 = pedWidth->getWidth(2);
00199 double pedWidth3 = pedWidth->getWidth(3);
00200
00201 if (depth == 1) {
00202
00203
00204
00205 monitor()->fillmeGain0Depth1(gainValue0);
00206 monitor()->fillmeGain1Depth1(gainValue1);
00207 monitor()->fillmeGain2Depth1(gainValue2);
00208 monitor()->fillmeGain3Depth1(gainValue3);
00209
00210 monitor()->fillmeGainWidth0Depth1(gainWidthValue0);
00211 monitor()->fillmeGainWidth1Depth1(gainWidthValue1);
00212 monitor()->fillmeGainWidth2Depth1(gainWidthValue2);
00213 monitor()->fillmeGainWidth3Depth1(gainWidthValue3);
00214
00215 monitor()->fillmePed0Depth1(pedValue0);
00216 monitor()->fillmePed1Depth1(pedValue1);
00217 monitor()->fillmePed2Depth1(pedValue2);
00218 monitor()->fillmePed3Depth1(pedValue3);
00219
00220 monitor()->fillmePedWidth0Depth1(pedWidth0);
00221 monitor()->fillmePedWidth1Depth1(pedWidth1);
00222 monitor()->fillmePedWidth2Depth1(pedWidth2);
00223 monitor()->fillmePedWidth3Depth1(pedWidth3);
00224
00225 monitor()->fillmeGainMap1 (double(ieta), double(iphi), gainValue0);
00226 monitor()->fillmePwidthMap1(double(ieta), double(iphi), pedWidth0) ;
00227 }
00228
00229 if (depth == 2) {
00230
00231
00232
00233 monitor()->fillmeGain0Depth2(gainValue0);
00234 monitor()->fillmeGain1Depth2(gainValue1);
00235 monitor()->fillmeGain2Depth2(gainValue2);
00236 monitor()->fillmeGain3Depth2(gainValue3);
00237
00238 monitor()->fillmeGainWidth0Depth2(gainWidthValue0);
00239 monitor()->fillmeGainWidth1Depth2(gainWidthValue1);
00240 monitor()->fillmeGainWidth2Depth2(gainWidthValue2);
00241 monitor()->fillmeGainWidth3Depth2(gainWidthValue3);
00242
00243 monitor()->fillmePed0Depth2(pedValue0);
00244 monitor()->fillmePed1Depth2(pedValue1);
00245 monitor()->fillmePed2Depth2(pedValue2);
00246 monitor()->fillmePed3Depth2(pedValue3);
00247
00248 monitor()->fillmePedWidth0Depth2(pedWidth0);
00249 monitor()->fillmePedWidth1Depth2(pedWidth1);
00250 monitor()->fillmePedWidth2Depth2(pedWidth2);
00251 monitor()->fillmePedWidth3Depth2(pedWidth3);
00252
00253 monitor()->fillmeGainMap2 (double(ieta), double(iphi), gainValue0);
00254 monitor()->fillmePwidthMap2(double(ieta), double(iphi), pedWidth0) ;
00255 }
00256
00257 if (depth == 3) {
00258
00259
00260
00261 monitor()->fillmeGain0Depth3(gainValue0);
00262 monitor()->fillmeGain1Depth3(gainValue1);
00263 monitor()->fillmeGain2Depth3(gainValue2);
00264 monitor()->fillmeGain3Depth3(gainValue3);
00265
00266 monitor()->fillmeGainWidth0Depth3(gainWidthValue0);
00267 monitor()->fillmeGainWidth1Depth3(gainWidthValue1);
00268 monitor()->fillmeGainWidth2Depth3(gainWidthValue2);
00269 monitor()->fillmeGainWidth3Depth3(gainWidthValue3);
00270
00271 monitor()->fillmePed0Depth3(pedValue0);
00272 monitor()->fillmePed1Depth3(pedValue1);
00273 monitor()->fillmePed2Depth3(pedValue2);
00274 monitor()->fillmePed3Depth3(pedValue3);
00275
00276 monitor()->fillmePedWidth0Depth3(pedWidth0);
00277 monitor()->fillmePedWidth1Depth3(pedWidth1);
00278 monitor()->fillmePedWidth2Depth3(pedWidth2);
00279 monitor()->fillmePedWidth3Depth3(pedWidth3);
00280
00281 monitor()->fillmeGainMap3 (double(ieta), double(iphi), gainValue0);
00282 monitor()->fillmePwidthMap3(double(ieta), double(iphi), pedWidth0) ;
00283 }
00284
00285 if (depth == 4) {
00286
00287
00288
00289 monitor()->fillmeGain0Depth4(gainValue0);
00290 monitor()->fillmeGain1Depth4(gainValue1);
00291 monitor()->fillmeGain2Depth4(gainValue2);
00292 monitor()->fillmeGain3Depth4(gainValue3);
00293
00294 monitor()->fillmeGainWidth0Depth4(gainWidthValue0);
00295 monitor()->fillmeGainWidth1Depth4(gainWidthValue1);
00296 monitor()->fillmeGainWidth2Depth4(gainWidthValue2);
00297 monitor()->fillmeGainWidth3Depth4(gainWidthValue3);
00298
00299 monitor()->fillmePed0Depth4(pedValue0);
00300 monitor()->fillmePed1Depth4(pedValue1);
00301 monitor()->fillmePed2Depth4(pedValue2);
00302 monitor()->fillmePed3Depth4(pedValue3);
00303
00304 monitor()->fillmePedWidth0Depth4(pedWidth0);
00305 monitor()->fillmePedWidth1Depth4(pedWidth1);
00306 monitor()->fillmePedWidth2Depth4(pedWidth2);
00307 monitor()->fillmePedWidth3Depth4(pedWidth3);
00308
00309 monitor()->fillmeGainMap4 (double(ieta), double(iphi), gainValue0);
00310 monitor()->fillmePwidthMap4(double(ieta), double(iphi), pedWidth0) ;
00311
00312 }
00313
00314 }
00315
00316
00317
00318
00319 if ( sub == subdet) Ndig++;
00320
00321
00322
00323 if ( sub == subdet && noise_ == 0 ) {
00324
00325
00326 HcalCalibrations calibrations = conditions->getHcalCalibrations(cell);
00327
00328 const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
00329 const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
00330 HcalCoderDb coder (*channelCoder, *shape);
00331 coder.adc2fC(*digiItr,tool);
00332
00333 double noiseADC = (*digiItr)[0].adc();
00334 double noisefC = tool[0];
00335
00336
00337 if(depth == 1) {
00338 monitor()->fillmeADC0_depth1 (noiseADC);
00339 monitor()->fillmeADC0fC_depth1(noisefC);
00340 }
00341 if(depth == 2) {
00342 monitor()->fillmeADC0_depth2 (noiseADC);
00343 monitor()->fillmeADC0fC_depth2(noisefC);
00344 }
00345 if(depth == 3) {
00346 monitor()->fillmeADC0_depth3 (noiseADC);
00347 monitor()->fillmeADC0fC_depth3(noisefC);
00348 }
00349 if(depth == 4) {
00350 monitor()->fillmeADC0_depth4 (noiseADC);
00351 monitor()->fillmeADC0fC_depth4(noisefC);
00352 }
00353
00354
00355 double deta = double(ieta);
00356 double dphi = double(iphi);
00357 if(depth == 1)
00358 monitor()->fillmeOccupancy_map_depth1(deta, dphi);
00359 if(depth == 2)
00360 monitor()->fillmeOccupancy_map_depth2(deta, dphi);
00361 if(depth == 3)
00362 monitor()->fillmeOccupancy_map_depth3(deta, dphi);
00363 if(depth == 4)
00364 monitor()->fillmeOccupancy_map_depth4(deta, dphi);
00365
00366
00367
00368
00369
00370 int closen = 0;
00371 if(ieta == ieta_Sim && iphi == iphi_Sim ) closen = seedSimHit;
00372
00373 for (int ii=0;ii<tool.size();ii++) {
00374 int capid = (*digiItr)[ii].capid();
00375
00376 double val = (tool[ii]-calibrations.pedestal(capid));
00377
00378 if (val > 10.) {
00379 if (depth == 1)
00380 monitor()->fillmeAll10slices_depth1(double(ii), val);
00381 else
00382 monitor()->fillmeAll10slices_depth2(double(ii), val);
00383 }
00384 if (val > 100.) {
00385 if (depth == 1)
00386 monitor()->fillmeAll10slices1D_depth1(double(ii), val);
00387 else
00388 monitor()->fillmeAll10slices1D_depth2(double(ii), val);
00389 }
00390
00391 if( closen == 1 &&( ieta * zsign >= 0 )) {
00392 monitor()->fillmeSignalTimeSlice(double(ii), val);
00393 }
00394
00395
00396
00397 if (subdet != 4 && ii>=4 && ii<=7) {
00398 ampl += val;
00399 if(depth == 1) ampl1 += val;
00400 if(depth == 2) ampl2 += val;
00401 if(depth == 3) ampl3 += val;
00402 if(depth == 4) ampl4 += val;
00403
00404 if( closen == 1 && ( ieta * zsign >= 0 )) {
00405 ampl_c += val;
00406 if(depth == 1) ampl1_c += val;
00407 if(depth == 2) ampl2_c += val;
00408 if(depth == 3) ampl3_c += val;
00409 if(depth == 4) ampl4_c += val;
00410
00411 }
00412 }
00413
00414
00415 if (subdet == 4 && ii==3 ) {
00416 ampl += val;
00417 if(depth == 1) ampl1 += val;
00418 if(depth == 2) ampl2 += val;
00419 if(depth == 3) ampl3 += val;
00420 if(depth == 4) ampl4 += val;
00421 if( closen == 1 && ( ieta * zsign >= 0 )) {
00422 ampl_c += val;
00423 if(depth == 1) ampl1_c += val;
00424 if(depth == 2) ampl2_c += val;
00425 if(depth == 3) ampl3_c += val;
00426 if(depth == 4) ampl4_c += val;
00427
00428 }
00429 }
00430 }
00431
00432
00433 monitor()->fillmeAmplIetaIphi1(double(ieta),double(iphi), ampl1);
00434 monitor()->fillmeAmplIetaIphi2(double(ieta),double(iphi), ampl2);
00435 monitor()->fillmeAmplIetaIphi3(double(ieta),double(iphi), ampl3);
00436 monitor()->fillmeAmplIetaIphi4(double(ieta),double(iphi), ampl4);
00437 monitor()->fillmeSumAmp (ampl);
00438
00439
00440 if(ampl1 > 10. || ampl2 > 10. || ampl3 > 10. || ampl4 > 10. ) ndigis++;
00441
00442
00443 if(ampl1 > 30. && depth == 1 && closen == 1 ) {
00444 double fBin5 = tool[4] - calibrations.pedestal((*digiItr)[4].capid());
00445 double fBin67 = tool[5] + tool[6]
00446 - calibrations.pedestal((*digiItr)[5].capid())
00447 - calibrations.pedestal((*digiItr)[6].capid());
00448 fBin5 /= ampl1;
00449 fBin67 /= ampl1;
00450 monitor()->fillmeBin5Frac (fBin5);
00451 monitor()->fillmeBin67Frac(fBin67);
00452 }
00453
00454 monitor()->fillmeSignalAmp (ampl);
00455 monitor()->fillmeSignalAmp1(ampl1);
00456 monitor()->fillmeSignalAmp2(ampl2);
00457 monitor()->fillmeSignalAmp3(ampl3);
00458 monitor()->fillmeSignalAmp4(ampl4);
00459
00460
00461 }
00462 }
00463
00464
00465 if ( subdet != 0 && noise_ == 0) {
00466
00467 monitor()->fillmenDigis(ndigis);
00468
00469
00470 double eps = 1.e-3;
00471 double ehits = 0.;
00472 double ehits1 = 0.;
00473 double ehits2 = 0.;
00474 double ehits3 = 0.;
00475 double ehits4 = 0.;
00476
00477 if(mc_ == "yes") {
00478 edm::Handle<edm::PCaloHitContainer> hcalHits ;
00479 iEvent.getByLabel("g4SimHits","HcalHits",hcalHits);
00480 const edm::PCaloHitContainer * simhitResult = hcalHits.product () ;
00481 for (std::vector<PCaloHit>::const_iterator simhits = simhitResult->begin (); simhits != simhitResult->end () ; ++simhits) {
00482
00483 HcalDetId cell(simhits->id());
00484 int ieta = cell.ieta();
00485 if(ieta > 0) ieta--;
00486 int iphi = cell.iphi()-1;
00487 int sub = cell.subdet();
00488
00489
00490 if (sub == subdet && ieta == ieta_Sim && iphi == iphi_Sim){
00491 int depth = cell.depth();
00492 double en = simhits->energy();
00493
00494 ehits += en;
00495 if(depth == 1) ehits1 += en;
00496 if(depth == 2) ehits2 += en;
00497 if(depth == 3) ehits3 += en;
00498 if(depth == 4) ehits4 += en;
00499 }
00500 }
00501
00502 if(ehits > eps) monitor()->fillmeDigiSimhit (ehits, ampl_c );
00503 if(ehits1 > eps) monitor()->fillmeDigiSimhit1(ehits1, ampl1_c);
00504 if(ehits2 > eps) monitor()->fillmeDigiSimhit2(ehits2, ampl2_c);
00505 if(ehits3 > eps) monitor()->fillmeDigiSimhit3(ehits3, ampl3_c);
00506 if(ehits4 > eps) monitor()->fillmeDigiSimhit4(ehits4, ampl4_c);
00507
00508 if(ehits > eps) monitor()->fillmeDigiSimhitProfile (ehits, ampl_c );
00509 if(ehits1 > eps) monitor()->fillmeDigiSimhitProfile1(ehits1, ampl1_c);
00510 if(ehits2 > eps) monitor()->fillmeDigiSimhitProfile2(ehits2, ampl2_c);
00511 if(ehits3 > eps) monitor()->fillmeDigiSimhitProfile3(ehits3, ampl3_c);
00512 if(ehits4 > eps) monitor()->fillmeDigiSimhitProfile4(ehits4, ampl4_c);
00513
00514 if(ehits > eps) monitor()->fillmeRatioDigiSimhit (ampl_c / ehits);
00515 if(ehits1 > eps) monitor()->fillmeRatioDigiSimhit1(ampl1_c / ehits1);
00516 if(ehits2 > eps) monitor()->fillmeRatioDigiSimhit2(ampl2_c / ehits2);
00517 if(ehits3 > eps) monitor()->fillmeRatioDigiSimhit3(ampl3_c / ehits3);
00518 if(ehits4 > eps) monitor()->fillmeRatioDigiSimhit4(ampl4_c / ehits4);
00519 }
00520
00521 monitor()->fillmeNdigis(double(Ndig));
00522
00523 }
00524
00525 }
00526
00527
00528 HcalDigiTester::HcalDigiTester(const edm::ParameterSet& iConfig)
00529 : dbe_(edm::Service<DQMStore>().operator->()),
00530 inputTag_(iConfig.getParameter<edm::InputTag>("digiLabel")),
00531 outputFile_(iConfig.getUntrackedParameter<std::string>("outputFile", "")),
00532 hcalselector_(iConfig.getUntrackedParameter<std::string>("hcalselector", "all")),
00533 zside_(iConfig.getUntrackedParameter<std::string>("zside", "*")),
00534 mode_(iConfig.getUntrackedParameter<std::string>("mode", "multi")),
00535 mc_(iConfig.getUntrackedParameter<std::string>("mc", "no")),
00536 monitors_()
00537 {
00538 if ( outputFile_.size() != 0 ) {
00539 edm::LogInfo("OutputInfo") << " Hcal Digi Task histograms will be saved to '" << outputFile_.c_str() << "'";
00540 } else {
00541 edm::LogInfo("OutputInfo") << " Hcal Digi Task histograms will NOT be saved";
00542 }
00543
00544
00545 }
00546
00547
00548 HcalDigiTester::~HcalDigiTester() { }
00549
00550
00551 void HcalDigiTester::endRun() {
00552
00553 if(noise_ != 1) {
00554
00555 if( hcalselector_ == "all") {
00556 hcalselector_ = "HB";
00557 eval_occupancy();
00558 hcalselector_ = "HE";
00559 eval_occupancy();
00560 hcalselector_ = "HO";
00561 eval_occupancy();
00562 hcalselector_ = "HF";
00563 eval_occupancy();
00564 }
00565 else
00566 eval_occupancy();
00567 }
00568
00569 }
00570
00571
00572
00573 void HcalDigiTester::endJob() {
00574
00575 if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_);
00576
00577 }
00578
00579
00580
00581 void HcalDigiTester::eval_occupancy() {
00582
00583 int nx = 82;
00584 int ny = 72;
00585 float cnorm;
00586 float fev = float (nevtot);
00587
00588
00589 float sumphi_1, sumphi_2, sumphi_3, sumphi_4;
00590 float phi_factor;
00591
00592 for (int i = 1; i <= nx; i++) {
00593 sumphi_1 = 0.;
00594 sumphi_2 = 0.;
00595 sumphi_3 = 0.;
00596 sumphi_4 = 0.;
00597
00598 for (int j = 1; j <= ny; j++) {
00599
00600
00601 cnorm = monitor()->getBinContent_depth1(i,j) / fev;
00602 monitor()->setBinContent_depth1(i,j,cnorm);
00603 cnorm = monitor()->getBinContent_depth2(i,j) / fev;
00604 monitor()->setBinContent_depth2(i,j,cnorm);
00605 cnorm = monitor()->getBinContent_depth3(i,j) / fev;
00606 monitor()->setBinContent_depth3(i,j,cnorm);
00607 cnorm = monitor()->getBinContent_depth4(i,j) / fev;
00608 monitor()->setBinContent_depth4(i,j,cnorm);
00609
00610 sumphi_1 += monitor()->getBinContent_depth1(i,j);
00611 sumphi_2 += monitor()->getBinContent_depth2(i,j);
00612 sumphi_3 += monitor()->getBinContent_depth3(i,j);
00613 sumphi_4 += monitor()->getBinContent_depth4(i,j);
00614
00615 }
00616
00617 int ieta = i - 42;
00618 if(ieta >=0 ) ieta +=1;
00619
00620 if(ieta >= -20 && ieta <= 20 )
00621 {phi_factor = 72.;}
00622 else {
00623 if(ieta >= 40 || ieta <= -40 ) {phi_factor = 18.;}
00624 else
00625 phi_factor = 36.;
00626 }
00627
00628
00629 if(ieta >= 0) ieta -= 1;
00630 double deta = double(ieta);
00631
00632 cnorm = sumphi_1 / phi_factor;
00633 monitor() -> fillmeOccupancy_vs_ieta_depth1(deta, cnorm);
00634 cnorm = sumphi_2 / phi_factor;
00635 monitor() -> fillmeOccupancy_vs_ieta_depth2(deta, cnorm);
00636 cnorm = sumphi_3 / phi_factor;
00637 monitor() -> fillmeOccupancy_vs_ieta_depth3(deta, cnorm);
00638 cnorm = sumphi_4 / phi_factor;
00639 monitor() -> fillmeOccupancy_vs_ieta_depth4(deta, cnorm);
00640
00641
00642 }
00643
00644 }
00645
00646 void HcalDigiTester::beginJob() {
00647
00648 nevent1 = 0;
00649 nevent2 = 0;
00650 nevent3 = 0;
00651 nevent4 = 0;
00652
00653 nevtot = 0;
00654
00655 }
00656
00657
00658 HcalSubdetDigiMonitor * HcalDigiTester::monitor()
00659 {
00660 std::map<std::string, HcalSubdetDigiMonitor*>::iterator monitorItr
00661 = monitors_.find(hcalselector_);
00662
00663 if(monitorItr == monitors_.end())
00664 {
00665 HcalSubdetDigiMonitor* m = new HcalSubdetDigiMonitor(dbe_, hcalselector_, noise_);
00666 std::pair<std::string, HcalSubdetDigiMonitor*> mapElement(
00667 hcalselector_, m);
00668 monitorItr = monitors_.insert(mapElement).first;
00669 }
00670 return monitorItr->second;
00671 }
00672
00673 void
00674 HcalDigiTester::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00675 {
00676
00677
00678 iSetup.get<CaloGeometryRecord>().get (geometry);
00679 iSetup.get<HcalDbRecord>().get(conditions);
00680
00681
00682
00683
00684
00685 if( hcalselector_ != "all") {
00686 noise_ = 0;
00687
00688
00689
00690 if (hcalselector_ == "HB" ) reco<HBHEDataFrame>(iEvent,iSetup);
00691 if (hcalselector_ == "HE" ) reco<HBHEDataFrame>(iEvent,iSetup);
00692 if (hcalselector_ == "HO" ) reco<HODataFrame>(iEvent,iSetup);
00693 if (hcalselector_ == "HF" ) reco<HFDataFrame>(iEvent,iSetup);
00694
00695 if (hcalselector_ == "noise") {
00696 noise_ = 1;
00697
00698
00699
00700
00701
00702 hcalselector_ = "HB";
00703 reco<HBHEDataFrame>(iEvent,iSetup);
00704 hcalselector_ = "HE";
00705 reco<HBHEDataFrame>(iEvent,iSetup);
00706 hcalselector_ = "HO";
00707 reco<HODataFrame>(iEvent,iSetup);
00708 hcalselector_ = "HF";
00709 reco<HFDataFrame>(iEvent,iSetup);
00710 hcalselector_ = "noise";
00711 }
00712 }
00713
00714 else {
00715 noise_ = 0;
00716
00717 hcalselector_ = "HB";
00718 reco<HBHEDataFrame>(iEvent,iSetup);
00719 hcalselector_ = "HE";
00720 reco<HBHEDataFrame>(iEvent,iSetup);
00721 hcalselector_ = "HO";
00722 reco<HODataFrame>(iEvent,iSetup);
00723 hcalselector_ = "HF";
00724 reco<HFDataFrame>(iEvent,iSetup);
00725 hcalselector_ = "all";
00726 }
00727
00728 nevtot++;
00729
00730 }
00731
00732 double HcalDigiTester::dR(double eta1, double phi1, double eta2, double phi2) {
00733 double PI = 3.1415926535898;
00734 double deltaphi= phi1 - phi2;
00735 if( phi2 > phi1 ) { deltaphi= phi2 - phi1;}
00736 if(deltaphi > PI) { deltaphi = 2.*PI - deltaphi;}
00737 double deltaeta = eta2 - eta1;
00738 double tmp = sqrt(deltaeta* deltaeta + deltaphi*deltaphi);
00739 return tmp;
00740 }
00741
00742
00743 DEFINE_FWK_MODULE (HcalDigiTester) ;