File indexing completed on 2024-04-06 12:18:51
0001 #include <iostream>
0002 #include <string>
0003 #include <memory>
0004
0005 #include <fmt/printf.h>
0006
0007
0008
0009
0010 #include <TDirectory.h>
0011 #include <TFile.h>
0012 #include <TH1F.h>
0013 #include <Math/VectorUtil.h>
0014
0015
0016
0017
0018 #include "HLTriggerOffline/Egamma/interface/EmDQMReco.h"
0019
0020
0021
0022
0023 #include "DataFormats/Common/interface/AssociationMap.h"
0024 #include "DataFormats/Common/interface/Handle.h"
0025 #include "DataFormats/Common/interface/RefToBase.h"
0026 #include "DataFormats/Common/interface/TriggerResults.h"
0027 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0028 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0029 #include "DataFormats/HLTReco/interface/TriggerEventWithRefs.h"
0030 #include "DataFormats/HLTReco/interface/TriggerObject.h"
0031 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0032 #include "DataFormats/L1Trigger/interface/L1EmParticle.h"
0033 #include "DataFormats/L1Trigger/interface/L1EmParticleFwd.h"
0034 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0035 #include "FWCore/Framework/interface/Frameworkfwd.h"
0036 #include "FWCore/Framework/interface/MakerMacros.h"
0037 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0038 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0039 #include "FWCore/Utilities/interface/Exception.h"
0040 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0041
0042 using namespace ROOT::Math::VectorUtil;
0043
0044
0045
0046
0047 EmDQMReco::FourVectorMonitorElements::FourVectorMonitorElements(EmDQMReco *_parent,
0048 DQMStore::IBooker &iBooker,
0049 const std::string &histogramNameTemplate,
0050 const std::string &histogramTitleTemplate)
0051 : parent(_parent) {
0052
0053 std::string histName;
0054 std::string histTitle;
0055
0056
0057 histName = fmt::sprintf(histogramNameTemplate, "et");
0058 histTitle = fmt::sprintf(histogramTitleTemplate, "E_{T}");
0059 etMonitorElement =
0060 iBooker.book1D(histName.c_str(), histTitle.c_str(), parent->plotBins, parent->plotPtMin, parent->plotPtMax);
0061
0062
0063 histName = fmt::sprintf(histogramNameTemplate, "eta");
0064 histTitle = fmt::sprintf(histogramTitleTemplate, "#eta");
0065 etaMonitorElement =
0066 iBooker.book1D(histName.c_str(), histTitle.c_str(), parent->plotBins, -parent->plotEtaMax, parent->plotEtaMax);
0067
0068
0069 histName = fmt::sprintf(histogramNameTemplate, "phi");
0070 histTitle = fmt::sprintf(histogramTitleTemplate, "#phi");
0071 phiMonitorElement =
0072 iBooker.book1D(histName.c_str(), histTitle.c_str(), parent->plotBins, -parent->plotPhiMax, parent->plotPhiMax);
0073 }
0074
0075
0076
0077 void EmDQMReco::FourVectorMonitorElements::fill(const math::XYZTLorentzVector &momentum) {
0078 etMonitorElement->Fill(momentum.Et());
0079 etaMonitorElement->Fill(momentum.eta());
0080 phiMonitorElement->Fill(momentum.phi());
0081 }
0082
0083
0084
0085
0086
0087
0088 EmDQMReco::EmDQMReco(const edm::ParameterSet &pset) {
0089
0090
0091
0092 dirname_ = "HLT/HLTEgammaValidation/" + pset.getParameter<std::string>("@module_label");
0093
0094
0095 reqNum = pset.getParameter<unsigned int>("reqNum");
0096 pdgGen = pset.getParameter<int>("pdgGen");
0097 recoEtaAcc = pset.getParameter<double>("genEtaAcc");
0098 recoEtAcc = pset.getParameter<double>("genEtAcc");
0099
0100 plotPtMin = pset.getUntrackedParameter<double>("PtMin", 0.);
0101 plotPtMax = pset.getUntrackedParameter<double>("PtMax", 1000.);
0102 plotEtaMax = pset.getUntrackedParameter<double>("EtaMax", 2.7);
0103 plotPhiMax = pset.getUntrackedParameter<double>("PhiMax", 3.15);
0104 plotBins = pset.getUntrackedParameter<unsigned int>("Nbins", 50);
0105 useHumanReadableHistTitles = pset.getUntrackedParameter<bool>("useHumanReadableHistTitles", false);
0106
0107 triggerNameRecoMonPath = pset.getUntrackedParameter<std::string>("triggerNameRecoMonPath", "HLT_MinBias");
0108 processNameRecoMonPath = pset.getUntrackedParameter<std::string>("processNameRecoMonPath", "HLT");
0109
0110 recoElectronsInput = consumes<reco::GsfElectronCollection>(
0111 pset.getUntrackedParameter<edm::InputTag>("recoElectrons", edm::InputTag("gsfElectrons")));
0112 recoObjectsEBT = consumes<std::vector<reco::SuperCluster>>(edm::InputTag("correctedHybridSuperClusters"));
0113 recoObjectsEET =
0114 consumes<std::vector<reco::SuperCluster>>(edm::InputTag("correctedMulti5x5SuperClustersWithPreshower"));
0115 hltResultsT = consumes<edm::TriggerResults>(edm::InputTag("TriggerResults", "", processNameRecoMonPath));
0116 triggerObjT = consumes<trigger::TriggerEventWithRefs>(edm::InputTag("hltTriggerSummaryRAW"));
0117
0118
0119
0120 recocut_ = pset.getParameter<int>("cutnum");
0121
0122
0123 eventnum = 0;
0124
0125
0126 isHltConfigInitialized_ = false;
0127
0128
0129
0130
0131
0132 std::vector<edm::ParameterSet> filters = pset.getParameter<std::vector<edm::ParameterSet>>("filters");
0133
0134 int i = 0;
0135 for (std::vector<edm::ParameterSet>::iterator filterconf = filters.begin(); filterconf != filters.end();
0136 filterconf++) {
0137 theHLTCollectionLabels.push_back(filterconf->getParameter<edm::InputTag>("HLTCollectionLabels"));
0138 theHLTOutputTypes.push_back(filterconf->getParameter<int>("theHLTOutputTypes"));
0139
0140
0141 theHLTCollectionHumanNames.push_back(
0142 filterconf->getUntrackedParameter<std::string>("HLTCollectionHumanName", theHLTCollectionLabels[i].label()));
0143
0144 std::vector<double> bounds = filterconf->getParameter<std::vector<double>>("PlotBounds");
0145
0146 assert(bounds.size() == 2);
0147 plotBounds.push_back(std::pair<double, double>(bounds[0], bounds[1]));
0148 isoNames.push_back(filterconf->getParameter<std::vector<edm::InputTag>>("IsoCollections"));
0149
0150 for (unsigned int i = 0; i < isoNames.back().size(); i++) {
0151 switch (theHLTOutputTypes.back()) {
0152 case trigger::TriggerL1NoIsoEG:
0153 histoFillerL1NonIso->isoNameTokens_.push_back(
0154 consumes<edm::AssociationMap<edm::OneToValue<l1extra::L1EmParticleCollection, float>>>(
0155 isoNames.back()[i]));
0156 break;
0157 case trigger::TriggerL1IsoEG:
0158 histoFillerL1Iso->isoNameTokens_.push_back(
0159 consumes<edm::AssociationMap<edm::OneToValue<l1extra::L1EmParticleCollection, float>>>(
0160 isoNames.back()[i]));
0161 break;
0162 case trigger::TriggerPhoton:
0163 histoFillerPho->isoNameTokens_.push_back(
0164 consumes<edm::AssociationMap<edm::OneToValue<reco::RecoEcalCandidateCollection, float>>>(
0165 isoNames.back()[i]));
0166 break;
0167 case trigger::TriggerElectron:
0168 histoFillerEle->isoNameTokens_.push_back(
0169 consumes<edm::AssociationMap<edm::OneToValue<reco::ElectronCollection, float>>>(isoNames.back()[i]));
0170 break;
0171 case trigger::TriggerCluster:
0172 histoFillerClu->isoNameTokens_.push_back(
0173 consumes<edm::AssociationMap<edm::OneToValue<reco::RecoEcalCandidateCollection, float>>>(
0174 isoNames.back()[i]));
0175 break;
0176 default:
0177 throw(cms::Exception("Release Validation Error") << "HLT output type not implemented: theHLTOutputTypes[n]");
0178 }
0179 }
0180
0181
0182 assert(!isoNames.back().empty());
0183 if (isoNames.back().at(0).label() == "none") {
0184 plotiso.push_back(false);
0185 } else {
0186 plotiso.push_back(true);
0187 }
0188 i++;
0189 }
0190
0191
0192 numOfHLTCollectionLabels = theHLTCollectionLabels.size();
0193 }
0194
0195
0196
0197
0198 void EmDQMReco::dqmBeginRun(const edm::Run &iRun, const edm::EventSetup &iSetup) {
0199 bool isHltConfigChanged = false;
0200 isHltConfigInitialized_ = hltConfig_.init(iRun, iSetup, "HLT", isHltConfigChanged);
0201 }
0202
0203
0204
0205
0206 void EmDQMReco::bookHistograms(DQMStore::IBooker &iBooker, edm::Run const &iRun, edm::EventSetup const &iSetup) {
0207 iBooker.setCurrentFolder(dirname_);
0208
0209
0210
0211
0212
0213
0214
0215 std::string histName = "total_eff";
0216 std::string histTitle = "total events passing";
0217
0218
0219 totalreco = iBooker.book1D(
0220 histName.c_str(), histTitle.c_str(), numOfHLTCollectionLabels + 2, 0, numOfHLTCollectionLabels + 2);
0221 totalreco->setBinLabel(numOfHLTCollectionLabels + 1, "Total");
0222 totalreco->setBinLabel(numOfHLTCollectionLabels + 2, "Reco");
0223 for (unsigned int u = 0; u < numOfHLTCollectionLabels; u++) {
0224 totalreco->setBinLabel(u + 1, theHLTCollectionLabels[u].label());
0225 }
0226
0227 histName = "total_eff_RECO_matched";
0228 histTitle = "total events passing (Reco matched)";
0229 totalmatchreco = iBooker.book1D(
0230 histName.c_str(), histTitle.c_str(), numOfHLTCollectionLabels + 2, 0, numOfHLTCollectionLabels + 2);
0231 totalmatchreco->setBinLabel(numOfHLTCollectionLabels + 1, "Total");
0232 totalmatchreco->setBinLabel(numOfHLTCollectionLabels + 2, "Reco");
0233 for (unsigned int u = 0; u < numOfHLTCollectionLabels; u++) {
0234 totalmatchreco->setBinLabel(u + 1, theHLTCollectionLabels[u].label());
0235 }
0236
0237
0238 MonitorElement *tmpiso;
0239
0240
0241
0242
0243 std::string pdgIdString;
0244 switch (pdgGen) {
0245 case 11:
0246 pdgIdString = "Electron";
0247 break;
0248 case 22:
0249 pdgIdString = "Photon";
0250 break;
0251 default:
0252 pdgIdString = "Particle";
0253 }
0254
0255
0256
0257
0258
0259 histReco = std::make_unique<FourVectorMonitorElements>(this,
0260 iBooker,
0261 "reco_%s",
0262 "%s of " + pdgIdString + "s");
0263
0264
0265
0266
0267 histRecoMonpath = std::make_unique<FourVectorMonitorElements>(this,
0268 iBooker,
0269 "reco_%s_monpath",
0270 "%s of " + pdgIdString + "s monpath");
0271
0272
0273
0274
0275
0276 histMonpath = std::make_unique<FourVectorMonitorElements>(this,
0277 iBooker,
0278 "final_%s_monpath",
0279 "Final %s Monpath");
0280
0281
0282
0283
0284
0285
0286
0287
0288 std::vector<std::string> HltHistTitle;
0289 if (theHLTCollectionHumanNames.size() == numOfHLTCollectionLabels && useHumanReadableHistTitles) {
0290 HltHistTitle = theHLTCollectionHumanNames;
0291 } else {
0292 for (unsigned int i = 0; i < numOfHLTCollectionLabels; i++) {
0293 HltHistTitle.push_back(theHLTCollectionLabels[i].label());
0294 }
0295 }
0296
0297 for (unsigned int i = 0; i < numOfHLTCollectionLabels; i++) {
0298
0299
0300
0301
0302
0303
0304
0305
0306
0307
0308
0309
0310
0311
0312
0313
0314
0315
0316
0317
0318
0319
0320
0321
0322
0323 standardHist.push_back(std::make_unique<FourVectorMonitorElements>(
0324 this,
0325 iBooker,
0326 theHLTCollectionLabels[i].label() + "%s_all",
0327 HltHistTitle[i] + " %s (ALL)"
0328 ));
0329
0330
0331
0332
0333
0334
0335
0336
0337
0338
0339
0340
0341
0342
0343
0344
0345
0346
0347
0348
0349
0350
0351
0352
0353
0354 histMatchReco.push_back(std::make_unique<FourVectorMonitorElements>(
0355 this,
0356 iBooker,
0357 theHLTCollectionLabels[i].label() + "%s_RECO_matched",
0358 HltHistTitle[i] + " %s (RECO matched)"
0359 ));
0360
0361
0362
0363
0364
0365
0366
0367
0368
0369
0370
0371
0372
0373
0374
0375
0376
0377
0378
0379
0380
0381
0382
0383
0384
0385
0386
0387
0388 histMatchRecoMonPath.push_back(std::make_unique<FourVectorMonitorElements>(
0389 this,
0390 iBooker,
0391 theHLTCollectionLabels[i].label() + "%s_RECO_matched_monpath",
0392 HltHistTitle[i] + " %s (RECO matched, monpath)"
0393 ));
0394
0395
0396
0397
0398
0399
0400
0401
0402
0403
0404
0405
0406
0407
0408
0409
0410
0411
0412
0413
0414
0415
0416
0417
0418
0419
0420 histHltObjMatchToReco.push_back(std::make_unique<FourVectorMonitorElements>(
0421 this,
0422 iBooker,
0423 theHLTCollectionLabels[i].label() + "%s_reco",
0424 HltHistTitle[i] + " %s (reco)"
0425 ));
0426
0427
0428
0429 if (!plotiso[i]) {
0430 tmpiso = nullptr;
0431 etahistiso.push_back(tmpiso);
0432 ethistiso.push_back(tmpiso);
0433 phiHistIso.push_back(tmpiso);
0434
0435 etahistisomatchreco.push_back(tmpiso);
0436 ethistisomatchreco.push_back(tmpiso);
0437 phiHistIsoMatchReco.push_back(tmpiso);
0438
0439 histEtaIsoOfHltObjMatchToReco.push_back(tmpiso);
0440 histEtIsoOfHltObjMatchToReco.push_back(tmpiso);
0441 histPhiIsoOfHltObjMatchToReco.push_back(tmpiso);
0442
0443 } else {
0444
0445
0446
0447
0448
0449 histName = theHLTCollectionLabels[i].label() + "eta_isolation_all";
0450 histTitle = HltHistTitle[i] + " isolation vs #eta (all)";
0451 tmpiso = iBooker.book2D(histName.c_str(),
0452 histTitle.c_str(),
0453 plotBins,
0454 -plotEtaMax,
0455 plotEtaMax,
0456 plotBins,
0457 plotBounds[i].first,
0458 plotBounds[i].second);
0459 etahistiso.push_back(tmpiso);
0460
0461
0462 histName = theHLTCollectionLabels[i].label() + "et_isolation_all";
0463 histTitle = HltHistTitle[i] + " isolation vs Et (all)";
0464 tmpiso = iBooker.book2D(histName.c_str(),
0465 histTitle.c_str(),
0466 plotBins,
0467 plotPtMin,
0468 plotPtMax,
0469 plotBins,
0470 plotBounds[i].first,
0471 plotBounds[i].second);
0472 ethistiso.push_back(tmpiso);
0473
0474
0475 histName = theHLTCollectionLabels[i].label() + "phi_isolation_all";
0476 histTitle = HltHistTitle[i] + " isolation vs #phi (all)";
0477 tmpiso = iBooker.book2D(histName.c_str(),
0478 histTitle.c_str(),
0479 plotBins,
0480 -plotPhiMax,
0481 plotPhiMax,
0482 plotBins,
0483 plotBounds[i].first,
0484 plotBounds[i].second);
0485 phiHistIso.push_back(tmpiso);
0486
0487
0488
0489
0490
0491
0492 histName = theHLTCollectionLabels[i].label() + "eta_isolation_RECO_matched";
0493 histTitle = HltHistTitle[i] + " isolation vs #eta (reco matched)";
0494 tmpiso = iBooker.book2D(histName.c_str(),
0495 histTitle.c_str(),
0496 plotBins,
0497 -plotEtaMax,
0498 plotEtaMax,
0499 plotBins,
0500 plotBounds[i].first,
0501 plotBounds[i].second);
0502 etahistisomatchreco.push_back(tmpiso);
0503
0504
0505 histName = theHLTCollectionLabels[i].label() + "et_isolation_RECO_matched";
0506 histTitle = HltHistTitle[i] + " isolation vs Et (reco matched)";
0507 tmpiso = iBooker.book2D(histName.c_str(),
0508 histTitle.c_str(),
0509 plotBins,
0510 plotPtMin,
0511 plotPtMax,
0512 plotBins,
0513 plotBounds[i].first,
0514 plotBounds[i].second);
0515 ethistisomatchreco.push_back(tmpiso);
0516
0517
0518 histName = theHLTCollectionLabels[i].label() + "phi_isolation_RECO_matched";
0519 histTitle = HltHistTitle[i] + " isolation vs #phi (reco matched)";
0520 tmpiso = iBooker.book2D(histName.c_str(),
0521 histTitle.c_str(),
0522 plotBins,
0523 -plotPhiMax,
0524 plotPhiMax,
0525 plotBins,
0526 plotBounds[i].first,
0527 plotBounds[i].second);
0528 phiHistIsoMatchReco.push_back(tmpiso);
0529
0530
0531
0532
0533
0534
0535
0536 histName = theHLTCollectionLabels[i].label() + "eta_isolation_reco";
0537 histTitle = HltHistTitle[i] + " isolation vs #eta (reco)";
0538 tmpiso = iBooker.book2D(histName.c_str(),
0539 histTitle.c_str(),
0540 plotBins,
0541 -plotEtaMax,
0542 plotEtaMax,
0543 plotBins,
0544 plotBounds[i].first,
0545 plotBounds[i].second);
0546 histEtaIsoOfHltObjMatchToReco.push_back(tmpiso);
0547
0548
0549 histName = theHLTCollectionLabels[i].label() + "et_isolation_reco";
0550 histTitle = HltHistTitle[i] + " isolation vs Et (reco)";
0551 tmpiso = iBooker.book2D(histName.c_str(),
0552 histTitle.c_str(),
0553 plotBins,
0554 plotPtMin,
0555 plotPtMax,
0556 plotBins,
0557 plotBounds[i].first,
0558 plotBounds[i].second);
0559 histEtIsoOfHltObjMatchToReco.push_back(tmpiso);
0560
0561
0562 histName = theHLTCollectionLabels[i].label() + "phi_isolation_reco";
0563 histTitle = HltHistTitle[i] + " isolation vs #phi (reco)";
0564 tmpiso = iBooker.book2D(histName.c_str(),
0565 histTitle.c_str(),
0566 plotBins,
0567 -plotPhiMax,
0568 plotPhiMax,
0569 plotBins,
0570 plotBounds[i].first,
0571 plotBounds[i].second);
0572 histPhiIsoOfHltObjMatchToReco.push_back(tmpiso);
0573
0574
0575 }
0576 }
0577 }
0578
0579
0580
0581
0582 EmDQMReco::~EmDQMReco() {}
0583
0584
0585
0586
0587 void EmDQMReco::analyze(const edm::Event &event, const edm::EventSetup &setup) {
0588
0589 if (!isHltConfigInitialized_)
0590 return;
0591
0592 eventnum++;
0593 bool plotMonpath = false;
0594 bool plotReco = true;
0595
0596 edm::Handle<edm::View<reco::Candidate>> recoObjects;
0597 edm::Handle<std::vector<reco::SuperCluster>> recoObjectsEB;
0598 edm::Handle<std::vector<reco::SuperCluster>> recoObjectsEE;
0599
0600 if (pdgGen == 11) {
0601 event.getByToken(recoElectronsInput, recoObjects);
0602
0603 if (recoObjects->size() < (unsigned int)recocut_) {
0604
0605
0606
0607 return;
0608 }
0609 } else if (pdgGen == 22) {
0610 event.getByToken(recoObjectsEBT, recoObjectsEB);
0611 event.getByToken(recoObjectsEET, recoObjectsEE);
0612
0613 if (recoObjectsEB->size() + recoObjectsEE->size() < (unsigned int)recocut_) {
0614
0615
0616
0617 return;
0618 }
0619 }
0620
0621 edm::Handle<edm::TriggerResults> HLTR;
0622 event.getByToken(hltResultsT, HLTR);
0623
0624
0625
0626
0627
0628
0629
0630
0631
0632
0633
0634
0635
0636
0637
0638
0639 unsigned int triggerIndex;
0640 triggerIndex = hltConfig_.triggerIndex(triggerNameRecoMonPath);
0641
0642
0643 bool isFired = false;
0644 if (triggerIndex < HLTR->size()) {
0645 isFired = HLTR->accept(triggerIndex);
0646 }
0647
0648
0649
0650 edm::Handle<trigger::TriggerEventWithRefs> triggerObj;
0651 event.getByToken(triggerObjT, triggerObj);
0652
0653 if (!triggerObj.isValid()) {
0654 edm::LogWarning("EmDQMReco") << "RAW-type HLT results not found, skipping event";
0655 return;
0656 }
0657
0658
0659
0660
0661
0662 totalreco->Fill(numOfHLTCollectionLabels + 0.5);
0663 totalmatchreco->Fill(numOfHLTCollectionLabels + .5);
0664
0665
0666
0667
0668
0669
0670
0671
0672
0673
0674
0675
0676
0677
0678 std::vector<reco::Particle> sortedReco;
0679 if (plotReco == true) {
0680 if (pdgGen == 11) {
0681 for (edm::View<reco::Candidate>::const_iterator recopart = recoObjects->begin(); recopart != recoObjects->end();
0682 recopart++) {
0683 reco::Particle tmpcand(
0684 recopart->charge(), recopart->p4(), recopart->vertex(), recopart->pdgId(), recopart->status());
0685 sortedReco.push_back(tmpcand);
0686 }
0687 } else if (pdgGen == 22) {
0688 for (std::vector<reco::SuperCluster>::const_iterator recopart2 = recoObjectsEB->begin();
0689 recopart2 != recoObjectsEB->end();
0690 recopart2++) {
0691 float en = recopart2->energy();
0692 float er = sqrt(pow(recopart2->x(), 2) + pow(recopart2->y(), 2) + pow(recopart2->z(), 2));
0693 float px = recopart2->energy() * recopart2->x() / er;
0694 float py = recopart2->energy() * recopart2->y() / er;
0695 float pz = recopart2->energy() * recopart2->z() / er;
0696 reco::Candidate::LorentzVector thisLV(px, py, pz, en);
0697 reco::Particle tmpcand(0, thisLV, math::XYZPoint(0., 0., 0.), 22, 1);
0698 sortedReco.push_back(tmpcand);
0699 }
0700 for (std::vector<reco::SuperCluster>::const_iterator recopart2 = recoObjectsEE->begin();
0701 recopart2 != recoObjectsEE->end();
0702 recopart2++) {
0703 float en = recopart2->energy();
0704 float er = sqrt(pow(recopart2->x(), 2) + pow(recopart2->y(), 2) + pow(recopart2->z(), 2));
0705 float px = recopart2->energy() * recopart2->x() / er;
0706 float py = recopart2->energy() * recopart2->y() / er;
0707 float pz = recopart2->energy() * recopart2->z() / er;
0708 reco::Candidate::LorentzVector thisLV(px, py, pz, en);
0709 reco::Particle tmpcand(0, thisLV, math::XYZPoint(0., 0., 0.), 22, 1);
0710 sortedReco.push_back(tmpcand);
0711 }
0712 }
0713
0714 std::sort(sortedReco.begin(), sortedReco.end(), pTComparator_);
0715
0716
0717
0718
0719
0720 sortedReco.erase(sortedReco.begin() + recocut_, sortedReco.end());
0721
0722 for (unsigned int i = 0; i < recocut_; i++) {
0723
0724 histReco->fill(sortedReco[i].p4());
0725
0726
0727
0728
0729
0730 if (isFired) {
0731 histRecoMonpath->fill(sortedReco[i].p4());
0732 plotMonpath = true;
0733 }
0734
0735 }
0736
0737 if (recocut_ >= reqNum)
0738 totalreco->Fill(numOfHLTCollectionLabels + 1.5);
0739 if (recocut_ >= reqNum)
0740 totalmatchreco->Fill(numOfHLTCollectionLabels + 1.5);
0741 }
0742
0743
0744
0745
0746 for (unsigned int n = 0; n < numOfHLTCollectionLabels; n++) {
0747
0748
0749 switch (theHLTOutputTypes[n]) {
0750 case trigger::TriggerL1NoIsoEG:
0751 histoFillerL1NonIso->fillHistos(triggerObj, event, n, sortedReco, plotReco, plotMonpath);
0752 break;
0753 case trigger::TriggerL1IsoEG:
0754 histoFillerL1Iso->fillHistos(triggerObj, event, n, sortedReco, plotReco, plotMonpath);
0755 break;
0756 case trigger::TriggerPhoton:
0757 histoFillerPho->fillHistos(triggerObj, event, n, sortedReco, plotReco, plotMonpath);
0758 break;
0759 case trigger::TriggerElectron:
0760 histoFillerEle->fillHistos(triggerObj, event, n, sortedReco, plotReco, plotMonpath);
0761 break;
0762 case trigger::TriggerCluster:
0763 histoFillerClu->fillHistos(triggerObj, event, n, sortedReco, plotReco, plotMonpath);
0764 break;
0765 default:
0766 throw(cms::Exception("Release Validation Error") << "HLT output type not implemented: theHLTOutputTypes[n]");
0767 }
0768 }
0769 }
0770
0771
0772
0773
0774
0775 template <class T>
0776 void HistoFillerReco<T>::fillHistos(edm::Handle<trigger::TriggerEventWithRefs> &triggerObj,
0777 const edm::Event &iEvent,
0778 unsigned int n,
0779 std::vector<reco::Particle> &sortedReco,
0780 bool plotReco,
0781 bool plotMonpath) {
0782 std::vector<edm::Ref<T>> recoecalcands;
0783 if ((triggerObj->filterIndex(dqm->theHLTCollectionLabels[n]) >= triggerObj->size())) {
0784 return;
0785 }
0786
0787
0788
0789
0790 triggerObj->getObjects(
0791 triggerObj->filterIndex(dqm->theHLTCollectionLabels[n]), dqm->theHLTOutputTypes[n], recoecalcands);
0792
0793
0794 if (dqm->theHLTOutputTypes[n] == trigger::TriggerL1NoIsoEG) {
0795 std::vector<edm::Ref<T>> isocands;
0796 triggerObj->getObjects(triggerObj->filterIndex(dqm->theHLTCollectionLabels[n]), trigger::TriggerL1IsoEG, isocands);
0797 if (!isocands.empty()) {
0798 for (unsigned int i = 0; i < isocands.size(); i++)
0799 recoecalcands.push_back(isocands[i]);
0800 }
0801 }
0802
0803 if (recoecalcands.empty()) {
0804 return;
0805 }
0806
0807 if (recoecalcands.size() >= dqm->reqNum)
0808 dqm->totalreco->Fill(n + 0.5);
0809
0810
0811
0812
0813
0814 for (unsigned int j = 0; j < recoecalcands.size(); j++) {
0815 if (!(recoecalcands.at(j).isAvailable())) {
0816 edm::LogError("EmDQMReco") << "Event content inconsistent: TriggerEventWithRefs contains "
0817 "invalid Refs"
0818 << std::endl
0819 << "invalid refs for: " << dqm->theHLTCollectionLabels[n].label();
0820 return;
0821 }
0822 }
0823
0824
0825
0826
0827
0828
0829
0830 for (unsigned int i = 0; i < recoecalcands.size(); i++) {
0831 dqm->standardHist[n]->fill(recoecalcands[i]->p4());
0832
0833
0834
0835
0836
0837 if (n + 1 < dqm->numOfHLTCollectionLabels) {
0838 if (dqm->plotiso[n + 1]) {
0839 for (unsigned int j = 0; j < isoNameTokens_.size(); j++) {
0840 edm::Handle<edm::AssociationMap<edm::OneToValue<T, float>>> depMap;
0841 iEvent.getByToken(isoNameTokens_.at(j), depMap);
0842 if (depMap.isValid()) {
0843
0844 typename edm::AssociationMap<edm::OneToValue<T, float>>::const_iterator mapi =
0845 depMap->find(recoecalcands[i]);
0846 if (mapi != depMap->end()) {
0847 dqm->etahistiso[n + 1]->Fill(recoecalcands[i]->eta(), mapi->val);
0848 dqm->ethistiso[n + 1]->Fill(recoecalcands[i]->et(), mapi->val);
0849 dqm->phiHistIso[n + 1]->Fill(recoecalcands[i]->phi(), mapi->val);
0850 }
0851 }
0852 }
0853 }
0854 }
0855 }
0856
0857
0858
0859
0860
0861 if (plotReco == true) {
0862 for (unsigned int i = 0; i < dqm->recocut_; i++) {
0863 math::XYZVector currentRecoParticleMomentum = sortedReco[i].momentum();
0864
0865
0866 float closestRecoDeltaR = 1000.;
0867 int closestRecoEcalCandIndex = -1;
0868 for (unsigned int j = 0; j < recoecalcands.size(); j++) {
0869 float deltaR = DeltaR(recoecalcands[j]->momentum(), currentRecoParticleMomentum);
0870
0871 if (deltaR < closestRecoDeltaR) {
0872 closestRecoDeltaR = deltaR;
0873 closestRecoEcalCandIndex = j;
0874 }
0875 }
0876
0877
0878
0879 if (closestRecoEcalCandIndex >= 0) {
0880
0881
0882
0883
0884
0885
0886
0887 dqm->histHltObjMatchToReco[n]->fill(recoecalcands[closestRecoEcalCandIndex]->p4());
0888
0889
0890 if (n + 1 < dqm->numOfHLTCollectionLabels) {
0891 if (dqm->plotiso[n + 1]) {
0892 for (unsigned int j = 0; j < isoNameTokens_.size(); j++) {
0893 edm::Handle<edm::AssociationMap<edm::OneToValue<T, float>>> depMap;
0894 iEvent.getByToken(isoNameTokens_.at(j), depMap);
0895 if (depMap.isValid()) {
0896
0897 typename edm::AssociationMap<edm::OneToValue<T, float>>::const_iterator mapi =
0898 depMap->find(recoecalcands[closestRecoEcalCandIndex]);
0899 if (mapi != depMap->end()) {
0900 dqm->histEtaIsoOfHltObjMatchToReco[n + 1]->Fill(recoecalcands[closestRecoEcalCandIndex]->eta(),
0901 mapi->val);
0902 dqm->histEtIsoOfHltObjMatchToReco[n + 1]->Fill(recoecalcands[closestRecoEcalCandIndex]->et(),
0903 mapi->val);
0904 dqm->histPhiIsoOfHltObjMatchToReco[n + 1]->Fill(recoecalcands[closestRecoEcalCandIndex]->phi(),
0905 mapi->val);
0906 }
0907 }
0908 }
0909 }
0910 }
0911 }
0912 }
0913
0914
0915
0916
0917 unsigned int mtachedRecoParts = 0;
0918 float minrecodist = 0.3;
0919 if (n == 0)
0920 minrecodist = 0.5;
0921 for (unsigned int i = 0; i < dqm->recocut_; i++) {
0922
0923 bool matchThis = false;
0924 math::XYZVector candDir = sortedReco[i].momentum();
0925 unsigned int closest = 0;
0926 double closestDr = 1000.;
0927 for (unsigned int trigOb = 0; trigOb < recoecalcands.size(); trigOb++) {
0928 double dr = DeltaR(recoecalcands[trigOb]->momentum(), candDir);
0929 if (dr < closestDr) {
0930 closestDr = dr;
0931 closest = trigOb;
0932 }
0933 if (closestDr > minrecodist) {
0934 closest = -1;
0935 } else {
0936 mtachedRecoParts++;
0937 matchThis = true;
0938 }
0939 }
0940 if (!matchThis)
0941 continue;
0942
0943
0944
0945
0946 dqm->histMatchReco[n]->fill(sortedReco[i].p4());
0947
0948 if (plotMonpath) {
0949
0950
0951
0952 dqm->histMatchRecoMonPath[n]->fill(sortedReco[i].p4());
0953 }
0954
0955
0956
0957
0958 if (n + 1 < dqm->numOfHLTCollectionLabels) {
0959 if (dqm->plotiso[n + 1]) {
0960 for (unsigned int j = 0; j < isoNameTokens_.size(); j++) {
0961 edm::Handle<edm::AssociationMap<edm::OneToValue<T, float>>> depMapReco;
0962 iEvent.getByToken(isoNameTokens_.at(j), depMapReco);
0963 if (depMapReco.isValid()) {
0964
0965 typename edm::AssociationMap<edm::OneToValue<T, float>>::const_iterator mapi =
0966 depMapReco->find(recoecalcands[closest]);
0967 if (mapi != depMapReco->end()) {
0968 dqm->etahistisomatchreco[n + 1]->Fill(sortedReco[i].eta(), mapi->val);
0969 dqm->ethistisomatchreco[n + 1]->Fill(sortedReco[i].et(), mapi->val);
0970 dqm->phiHistIsoMatchReco[n + 1]->Fill(sortedReco[i].eta(), mapi->val);
0971 }
0972 }
0973 }
0974 }
0975 }
0976 }
0977
0978 if (mtachedRecoParts >= dqm->reqNum)
0979 dqm->totalmatchreco->Fill(n + 0.5);
0980 }
0981 }
0982
0983 DEFINE_FWK_MODULE(EmDQMReco);