00001
00002
00004 #include "HLTriggerOffline/Egamma/interface/EmDQM.h"
00005
00007
00009 #include "FWCore/Framework/interface/Frameworkfwd.h"
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011 #include "FWCore/Framework/interface/MakerMacros.h"
00012 #include "DataFormats/HLTReco/interface/TriggerEventWithRefs.h"
00013 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
00014 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00015
00016 #include "DataFormats/L1Trigger/interface/L1EmParticle.h"
00017 #include "DataFormats/L1Trigger/interface/L1EmParticleFwd.h"
00018
00019 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021 #include "DataFormats/Common/interface/AssociationMap.h"
00022
00023 #include "DataFormats/Common/interface/Handle.h"
00024 #include "DataFormats/Common/interface/RefToBase.h"
00025 #include "FWCore/ServiceRegistry/interface/Service.h"
00026 #include "FWCore/Utilities/interface/Exception.h"
00027 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
00028 #include <boost/foreach.hpp>
00029
00031
00033 #include "TFile.h"
00034 #include "TDirectory.h"
00035 #include "TH1F.h"
00036 #include <iostream>
00037 #include <string>
00038 #include <Math/VectorUtil.h>
00039 using namespace ROOT::Math::VectorUtil ;
00040
00041
00043
00045 EmDQM::EmDQM(const edm::ParameterSet& pset)
00046 {
00047
00048 dbe = edm::Service < DQMStore > ().operator->();
00049 dbe->setVerbose(0);
00050
00052
00054 dirname_="HLT/HLTEgammaValidation/"+pset.getParameter<std::string>("@module_label");
00055 dbe->setCurrentFolder(dirname_);
00056
00057 triggerobjwithrefs = pset.getParameter<edm::InputTag>("triggerobject");
00058 pathIndex = pset.getUntrackedParameter<unsigned int>("pathIndex", 0);
00059
00060 reqNum = pset.getParameter<unsigned int>("reqNum");
00061 pdgGen = pset.getParameter<int>("pdgGen");
00062 genEtaAcc = pset.getParameter<double>("genEtaAcc");
00063 genEtAcc = pset.getParameter<double>("genEtAcc");
00064
00065 plotEtMin = pset.getUntrackedParameter<double>("genEtMin",0.);
00066 plotPtMin = pset.getUntrackedParameter<double>("PtMin",0.);
00067 plotPtMax = pset.getUntrackedParameter<double>("PtMax",1000.);
00068 plotEtaMax = pset.getUntrackedParameter<double>("EtaMax", 2.7);
00069 plotPhiMax = pset.getUntrackedParameter<double>("PhiMax", 3.15);
00070 plotBins = pset.getUntrackedParameter<unsigned int>("Nbins",40);
00071 plotMinEtForEtaEffPlot = pset.getUntrackedParameter<unsigned int>("minEtForEtaEffPlot", 15);
00072 useHumanReadableHistTitles = pset.getUntrackedParameter<bool>("useHumanReadableHistTitles", false);
00073 mcMatchedOnly = pset.getUntrackedParameter<bool>("mcMatchedOnly", true);
00074 noPhiPlots = pset.getUntrackedParameter<bool>("noPhiPlots", true);
00075 noIsolationPlots = pset.getUntrackedParameter<bool>("noIsolationPlots", true);
00076 verbosity = pset.getUntrackedParameter<unsigned int>("verbosity",0);
00077
00078
00079 gencutCollection_= pset.getParameter<edm::InputTag>("cutcollection");
00080 gencut_ = pset.getParameter<int>("cutnum");
00081
00083
00084
00086 std::vector<edm::ParameterSet> filters =
00087 pset.getParameter<std::vector<edm::ParameterSet> >("filters");
00088
00089 int i = 0;
00090 for(std::vector<edm::ParameterSet>::iterator filterconf = filters.begin() ; filterconf != filters.end() ; filterconf++)
00091 {
00092
00093 theHLTCollectionLabels.push_back(filterconf->getParameter<edm::InputTag>("HLTCollectionLabels"));
00094 theHLTOutputTypes.push_back(filterconf->getParameter<int>("theHLTOutputTypes"));
00095
00096 theHLTCollectionHumanNames.push_back(filterconf->getUntrackedParameter<std::string>("HLTCollectionHumanName",theHLTCollectionLabels[i].label()));
00097
00098 std::vector<double> bounds = filterconf->getParameter<std::vector<double> >("PlotBounds");
00099
00100 assert(bounds.size() == 2);
00101 plotBounds.push_back(std::pair<double,double>(bounds[0],bounds[1]));
00102 isoNames.push_back(filterconf->getParameter<std::vector<edm::InputTag> >("IsoCollections"));
00103
00104 assert(isoNames.back().size()>0);
00105 if (isoNames.back().at(0).label()=="none") {
00106 plotiso.push_back(false);
00107 } else {
00108 if (!noIsolationPlots) plotiso.push_back(true);
00109 else plotiso.push_back(false);
00110 }
00111 nCandCuts.push_back(filterconf->getParameter<int>("ncandcut"));
00112 i++;
00113 }
00114
00115
00116 numOfHLTCollectionLabels = theHLTCollectionLabels.size();
00117
00118 }
00119
00120
00122
00124 void
00125 EmDQM::beginJob()
00126 {
00127
00128 }
00129
00130 void
00131 EmDQM::beginRun(edm::Run const &iRun, edm::EventSetup const &iSetup)
00132 {
00133 bool changed(true);
00134 if (hltConf_.init(iRun, iSetup, triggerobjwithrefs.process(), changed)) {
00135
00136
00137
00138
00139 dbe->setCurrentFolder(dirname_);
00140
00142
00143
00144
00146
00147 std::string histName="total_eff";
00148 std::string histTitle = "total events passing";
00149 if (!mcMatchedOnly) {
00150
00151
00152 total = dbe->book1D(histName.c_str(),histTitle.c_str(),numOfHLTCollectionLabels+2,0,numOfHLTCollectionLabels+2);
00153 total->setBinLabel(numOfHLTCollectionLabels+1,"Total");
00154 total->setBinLabel(numOfHLTCollectionLabels+2,"Gen");
00155 for (unsigned int u=0; u<numOfHLTCollectionLabels; u++){total->setBinLabel(u+1,theHLTCollectionLabels[u].label().c_str());}
00156 }
00157
00158 histName="total_eff_MC_matched";
00159 histTitle="total events passing (mc matched)";
00160 totalmatch = dbe->book1D(histName.c_str(),histTitle.c_str(),numOfHLTCollectionLabels+2,0,numOfHLTCollectionLabels+2);
00161 totalmatch->setBinLabel(numOfHLTCollectionLabels+1,"Total");
00162 totalmatch->setBinLabel(numOfHLTCollectionLabels+2,"Gen");
00163 for (unsigned int u=0; u<numOfHLTCollectionLabels; u++){totalmatch->setBinLabel(u+1,theHLTCollectionLabels[u].label().c_str());}
00164
00165 MonitorElement* tmphisto;
00166 MonitorElement* tmpiso;
00167
00169
00171 std::string pdgIdString;
00172 switch(pdgGen) {
00173 case 11:
00174 pdgIdString="Electron";break;
00175 case 22:
00176 pdgIdString="Photon";break;
00177 default:
00178 pdgIdString="Particle";
00179 }
00180
00181 histName = "gen_et";
00182 histTitle= "E_{T} of " + pdgIdString + "s" ;
00183 etgen = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00184 histName = "gen_eta";
00185 histTitle= "#eta of "+ pdgIdString +"s " ;
00186 etagen = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00187 histName = "gen_phi";
00188 histTitle= "#phi of "+ pdgIdString +"s " ;
00189 if (!noPhiPlots) phigen = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax);
00190
00191
00192
00194
00196
00197
00198 std::vector<std::string> HltHistTitle;
00199 if ( theHLTCollectionHumanNames.size() == numOfHLTCollectionLabels && useHumanReadableHistTitles ) {
00200 HltHistTitle = theHLTCollectionHumanNames;
00201 } else {
00202 for (unsigned int i =0; i < numOfHLTCollectionLabels; i++) {
00203 HltHistTitle.push_back(theHLTCollectionLabels[i].label());
00204 }
00205 }
00206
00207 for(unsigned int i = 0; i< numOfHLTCollectionLabels ; i++){
00208 if (!mcMatchedOnly) {
00209
00210 histName = theHLTCollectionLabels[i].label()+"et_all";
00211 histTitle = HltHistTitle[i]+" Et (ALL)";
00212 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00213 ethist.push_back(tmphisto);
00214
00215
00216 histName = theHLTCollectionLabels[i].label()+"eta_all";
00217 histTitle = HltHistTitle[i]+" #eta (ALL)";
00218 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00219 etahist.push_back(tmphisto);
00220
00221 if (!noPhiPlots) {
00222
00223 histName = theHLTCollectionLabels[i].label()+"phi_all";
00224 histTitle = HltHistTitle[i]+" #phi (ALL)";
00225 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax);
00226 phihist.push_back(tmphisto);
00227 }
00228
00229
00230
00231 histName = theHLTCollectionLabels[i].label()+"et";
00232 histTitle = HltHistTitle[i]+" Et";
00233 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00234 histEtOfHltObjMatchToGen.push_back(tmphisto);
00235
00236
00237 histName = theHLTCollectionLabels[i].label()+"eta";
00238 histTitle = HltHistTitle[i]+" eta";
00239 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00240 histEtaOfHltObjMatchToGen.push_back(tmphisto);
00241
00242 if (!noPhiPlots) {
00243
00244 histName = theHLTCollectionLabels[i].label()+"phi";
00245 histTitle = HltHistTitle[i]+" phi";
00246 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax);
00247 histPhiOfHltObjMatchToGen.push_back(tmphisto);
00248 }
00249 }
00250
00251
00252 histName = theHLTCollectionLabels[i].label()+"et_MC_matched";
00253 histTitle = HltHistTitle[i]+" Et (MC matched)";
00254 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00255 ethistmatch.push_back(tmphisto);
00256
00257
00258 histName = theHLTCollectionLabels[i].label()+"eta_MC_matched";
00259 histTitle = HltHistTitle[i]+" #eta (MC matched)";
00260 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00261 etahistmatch.push_back(tmphisto);
00262
00263 if (!noPhiPlots) {
00264
00265 histName = theHLTCollectionLabels[i].label()+"phi_MC_matched";
00266 histTitle = HltHistTitle[i]+" #phi (MC matched)";
00267 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax);
00268 phihistmatch.push_back(tmphisto);
00269 }
00270
00271
00272 if (!plotiso[i]) {
00273 tmpiso = NULL;
00274 if (!mcMatchedOnly) {
00275 etahistiso.push_back(tmpiso);
00276 phihistiso.push_back(tmpiso);
00277 ethistiso.push_back(tmpiso);
00278 histEtaIsoOfHltObjMatchToGen.push_back(tmpiso);
00279 histPhiIsoOfHltObjMatchToGen.push_back(tmpiso);
00280 histEtIsoOfHltObjMatchToGen.push_back(tmpiso);
00281 }
00282 etahistisomatch.push_back(tmpiso);
00283 phihistisomatch.push_back(tmpiso);
00284 ethistisomatch.push_back(tmpiso);
00285 } else {
00286 if (!mcMatchedOnly) {
00287
00288 histName = theHLTCollectionLabels[i].label()+"eta_isolation_all";
00289 histTitle = HltHistTitle[i]+" isolation vs #eta (all)";
00290 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00291 etahistiso.push_back(tmpiso);
00292
00293
00294 histName = theHLTCollectionLabels[i].label()+"phi_isolation_all";
00295 histTitle = HltHistTitle[i]+" isolation vs #phi (all)";
00296 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00297 phihistiso.push_back(tmpiso);
00298
00299
00300 histName = theHLTCollectionLabels[i].label()+"et_isolation_all";
00301 histTitle = HltHistTitle[i]+" isolation vs Et (all)";
00302 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00303 ethistiso.push_back(tmpiso);
00304
00305
00306
00307 histName = theHLTCollectionLabels[i].label()+"eta_isolation";
00308 histTitle = HltHistTitle[i]+" isolation vs #eta";
00309 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00310 histEtaIsoOfHltObjMatchToGen.push_back(tmpiso);
00311
00312
00313
00314 histName = theHLTCollectionLabels[i].label()+"phi_isolation";
00315 histTitle = HltHistTitle[i]+" isolation vs #phi";
00316 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00317 histPhiIsoOfHltObjMatchToGen.push_back(tmpiso);
00318
00319
00320
00321 histName = theHLTCollectionLabels[i].label()+"et_isolation";
00322 histTitle = HltHistTitle[i]+" isolation vs Et";
00323 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00324 histEtIsoOfHltObjMatchToGen.push_back(tmpiso);
00325 }
00326
00327
00328 histName = theHLTCollectionLabels[i].label()+"eta_isolation_MC_matched";
00329 histTitle = HltHistTitle[i]+" isolation vs #eta (mc matched)";
00330 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00331 etahistisomatch.push_back(tmpiso);
00332
00333
00334 histName = theHLTCollectionLabels[i].label()+"phi_isolation_MC_matched";
00335 histTitle = HltHistTitle[i]+" isolation vs #phi (mc matched)";
00336 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00337 phihistisomatch.push_back(tmpiso);
00338
00339
00340
00341 histName = theHLTCollectionLabels[i].label()+"et_isolation_MC_matched";
00342 histTitle = HltHistTitle[i]+" isolation vs Et (mc matched)";
00343 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00344 ethistisomatch.push_back(tmpiso);
00345
00346 }
00347
00348 }
00349
00350 if (changed) {
00351
00352
00353 }
00354 } else {
00355
00356
00357 if (verbosity >= OUTPUT_ERRORS)
00358 edm::LogError("EmDQM") << " HLT config extraction failure with process name '" << triggerobjwithrefs.process() << "'.";
00359
00360 }
00361 }
00362
00364
00366 EmDQM::~EmDQM(){
00367 }
00369
00370 bool EmDQM::checkGeneratedParticlesRequirement(const edm::Event & event)
00371 {
00373
00374
00375
00377 edm::Handle< edm::View<reco::Candidate> > genParticles;
00378 event.getByLabel("genParticles", genParticles);
00379 if(!genParticles.isValid()) {
00380 if (verbosity >= OUTPUT_WARNINGS)
00381 edm::LogWarning("EmDQM") << "genParticles invalid.";
00382 return false;
00383 }
00384
00385 std::vector<reco::LeafCandidate> allSortedGenParticles;
00386
00387 for(edm::View<reco::Candidate>::const_iterator currentGenParticle = genParticles->begin(); currentGenParticle != genParticles->end(); currentGenParticle++){
00388
00389
00390
00391
00392 if ( !( abs((*currentGenParticle).pdgId())==pdgGen && (*currentGenParticle).status()==1 && (*currentGenParticle).et() > 2.0) ) continue;
00393
00394 reco::LeafCandidate tmpcand( *(currentGenParticle) );
00395
00396 if (tmpcand.et() < plotEtMin) continue;
00397
00398 allSortedGenParticles.push_back(tmpcand);
00399 }
00400
00401 std::sort(allSortedGenParticles.begin(), allSortedGenParticles.end(),pTGenComparator_);
00402
00403
00404 if (allSortedGenParticles.size() < gencut_)
00405 return false;
00406
00407
00408
00409
00410
00411
00412
00413
00414 for (unsigned int i = 0 ; i < gencut_ ; i++ ) {
00415 bool inECALgap = fabs(allSortedGenParticles[i].eta()) > 1.4442 && fabs(allSortedGenParticles[i].eta()) < 1.556;
00416 if ( (fabs(allSortedGenParticles[i].eta()) > genEtaAcc) || inECALgap ) {
00417
00418 return false;
00419 }
00420 }
00421
00422
00423 return true;
00424 }
00426
00427 bool EmDQM::checkRecoParticlesRequirement(const edm::Event & event)
00428 {
00429
00430
00431
00432 edm::Handle< edm::View<reco::Candidate> > referenceParticles;
00433 event.getByLabel(gencutCollection_,referenceParticles);
00434 if(!referenceParticles.isValid()) {
00435 if (verbosity >= OUTPUT_WARNINGS)
00436 edm::LogWarning("EmDQM") << "referenceParticles invalid.";
00437 return false;
00438 }
00439
00440 std::vector<const reco::Candidate *> allSortedReferenceParticles;
00441
00442 for(edm::View<reco::Candidate>::const_iterator currentReferenceParticle = referenceParticles->begin();
00443 currentReferenceParticle != referenceParticles->end();
00444 currentReferenceParticle++)
00445 {
00446 if ( currentReferenceParticle->et() <= 2.0)
00447 continue;
00448
00449
00450
00451
00452
00453
00454 if (currentReferenceParticle->et() < plotEtMin)
00455 continue;
00456
00457
00458 allSortedReferenceParticles.push_back(&(*currentReferenceParticle));
00459 }
00460
00461
00462
00463
00464 return allSortedReferenceParticles.size() >= gencut_;
00465 }
00466
00467
00469
00471 void
00472 EmDQM::analyze(const edm::Event & event , const edm::EventSetup& setup)
00473 {
00475
00476
00478 edm::Handle< edm::View<reco::Candidate> > cutCounter;
00479 event.getByLabel(gencutCollection_,cutCounter);
00480 if (cutCounter->size() < (unsigned int)gencut_) {
00481
00482 return;
00483 }
00484
00485
00486
00487
00488 edm::Handle<trigger::TriggerEventWithRefs> triggerObj;
00489 event.getByLabel(triggerobjwithrefs,triggerObj);
00490 if(!triggerObj.isValid()) {
00491 if (verbosity >= OUTPUT_WARNINGS)
00492 edm::LogWarning("EmDQM") << "parameter triggerobject (" << triggerobjwithrefs << ") does not corresond to a valid TriggerEventWithRefs product. Please check especially the process name (e.g. when running over reprocessed datasets)";
00493 return;
00494 }
00495
00496
00497 if (event.isRealData())
00498 {
00499
00500
00501
00502
00503
00504 if (!checkRecoParticlesRequirement(event))
00505 return;
00506 }
00507 else
00508 {
00509
00510 if (!checkGeneratedParticlesRequirement(event))
00511
00512 return;
00513 }
00514
00515
00516
00517
00519
00520
00522 if (!mcMatchedOnly) total->Fill(numOfHLTCollectionLabels+0.5);
00523 totalmatch->Fill(numOfHLTCollectionLabels+0.5);
00524
00525
00527
00529
00530
00531 std::vector<reco::Particle> sortedGen;
00532 for(edm::View<reco::Candidate>::const_iterator genpart = cutCounter->begin(); genpart != cutCounter->end();genpart++){
00533 reco::Particle tmpcand( genpart->charge(), genpart->p4(), genpart->vertex(),genpart->pdgId(),genpart->status() );
00534 if (tmpcand.et() >= plotEtMin) {
00535 sortedGen.push_back(tmpcand);
00536 }
00537 }
00538 std::sort(sortedGen.begin(),sortedGen.end(),pTComparator_ );
00539
00540
00541
00542
00543 if (sortedGen.size() < gencut_){
00544 return;
00545 }
00546 sortedGen.erase(sortedGen.begin()+gencut_,sortedGen.end());
00547
00548 for (unsigned int i = 0 ; i < gencut_ ; i++ ) {
00549 etgen ->Fill( sortedGen[i].et() );
00550 if (sortedGen[i].et() > plotMinEtForEtaEffPlot) {
00551 etagen->Fill( sortedGen[i].eta() );
00552 if (!noPhiPlots) phigen->Fill( sortedGen[i].phi() );
00553 }
00554 }
00555 if (gencut_ >= reqNum && !mcMatchedOnly) total->Fill(numOfHLTCollectionLabels+1.5);
00556 if (gencut_ >= reqNum) totalmatch->Fill(numOfHLTCollectionLabels+1.5);
00557
00558 bool accepted = true;
00559 edm::Handle<edm::TriggerResults> hltResults;
00560 event.getByLabel(edm::InputTag("TriggerResults","", triggerobjwithrefs.process()), hltResults);
00562
00564 for(unsigned int n=0; n < numOfHLTCollectionLabels ; n++) {
00565
00566 if (sortedGen.size() < nCandCuts.at(n)) {
00567 if (verbosity >= OUTPUT_ERRORS)
00568 edm::LogError("EmDQM") << "There are less generated particles than the module '" << theHLTCollectionLabels[n].label() << "' requires.";
00569 continue;
00570 }
00571 std::vector<reco::Particle> sortedGenForFilter(sortedGen);
00572 sortedGenForFilter.erase(sortedGenForFilter.begin() + nCandCuts.at(n), sortedGenForFilter.end());
00573
00574
00575 if (pathIndex != 0 && hltConf_.moduleIndex(pathIndex, theHLTCollectionLabels[n].label()) > hltResults->index(pathIndex)) break;
00576
00577
00578 switch(theHLTOutputTypes[n])
00579 {
00580 case trigger::TriggerL1NoIsoEG:
00581 fillHistos<l1extra::L1EmParticleCollection>(triggerObj,event,n,sortedGenForFilter,accepted);break;
00582 case trigger::TriggerL1IsoEG:
00583 fillHistos<l1extra::L1EmParticleCollection>(triggerObj,event,n,sortedGenForFilter,accepted);break;
00584 case trigger::TriggerPhoton:
00585 fillHistos<reco::RecoEcalCandidateCollection>(triggerObj,event,n,sortedGenForFilter,accepted);break;
00586 case trigger::TriggerElectron:
00587 fillHistos<reco::ElectronCollection>(triggerObj,event,n,sortedGenForFilter,accepted);break;
00588 case trigger::TriggerCluster:
00589 fillHistos<reco::RecoEcalCandidateCollection>(triggerObj,event,n,sortedGenForFilter,accepted);break;
00590 default:
00591 throw(cms::Exception("Release Validation Error") << "HLT output type not implemented: theHLTOutputTypes[n]" );
00592 }
00593 }
00594 }
00595
00596
00598
00599
00601 template <class T> void EmDQM::fillHistos(edm::Handle<trigger::TriggerEventWithRefs>& triggerObj,const edm::Event& iEvent ,unsigned int n,std::vector<reco::Particle>& sortedGen, bool &accepted)
00602 {
00603 std::vector<edm::Ref<T> > recoecalcands;
00604 if ( ( triggerObj->filterIndex(theHLTCollectionLabels[n])>=triggerObj->size() )){
00605 hltCollectionLabelsMissed.insert(theHLTCollectionLabels[n].encode());
00606 accepted = false;
00607 return;
00608 }
00609
00610 hltCollectionLabelsFound.insert(theHLTCollectionLabels[n].encode());
00611
00613
00615 triggerObj->getObjects(triggerObj->filterIndex(theHLTCollectionLabels[n]),theHLTOutputTypes[n],recoecalcands);
00616
00617
00618 if (theHLTOutputTypes[n] == trigger::TriggerL1NoIsoEG){
00619 std::vector<edm::Ref<T> > isocands;
00620 triggerObj->getObjects(triggerObj->filterIndex(theHLTCollectionLabels[n]),trigger::TriggerL1IsoEG,isocands);
00621 if (isocands.size()>0)
00622 {
00623 for (unsigned int i=0; i < isocands.size(); i++)
00624 recoecalcands.push_back(isocands[i]);
00625 }
00626 }
00627
00628
00629 if (recoecalcands.size() < 1){
00630 accepted = false;
00631 return;
00632 }
00633
00634
00635 if (recoecalcands.size() >= nCandCuts.at(n) && !mcMatchedOnly)
00636 total->Fill(n+0.5);
00637
00639
00640
00642 for (unsigned int j=0; j<recoecalcands.size(); j++){
00643 if(!( recoecalcands.at(j).isAvailable())){
00644 if (verbosity >= OUTPUT_ERRORS)
00645 edm::LogError("EmDQMInvalidRefs") << "Event content inconsistent: TriggerEventWithRefs contains invalid Refs. Invalid refs for: " << theHLTCollectionLabels[n].label() << ". The collection that this module uses may has been dropped in the event.";
00646 return;
00647 }
00648 }
00649
00650 if (!mcMatchedOnly) {
00652
00653
00655
00656 for (unsigned int i=0; i < nCandCuts.at(n); i++) {
00657 math::XYZVector currentGenParticleMomentum = sortedGen[i].momentum();
00658
00659 float closestDeltaR = 0.5;
00660 int closestEcalCandIndex = -1;
00661 for (unsigned int j=0; j<recoecalcands.size(); j++) {
00662 float deltaR = DeltaR(recoecalcands[j]->momentum(),currentGenParticleMomentum);
00663
00664 if (deltaR < closestDeltaR) {
00665 closestDeltaR = deltaR;
00666 closestEcalCandIndex = j;
00667 }
00668 }
00669
00670
00671
00672 if ( closestEcalCandIndex >= 0 ) {
00673 histEtOfHltObjMatchToGen[n] ->Fill( recoecalcands[closestEcalCandIndex]->et() );
00674 histEtaOfHltObjMatchToGen[n]->Fill( recoecalcands[closestEcalCandIndex]->eta() );
00675 if (!noPhiPlots) histPhiOfHltObjMatchToGen[n]->Fill( recoecalcands[closestEcalCandIndex]->phi() );
00676
00677
00678 if (n+1 < numOfHLTCollectionLabels){
00679 if (plotiso[n+1] ){
00680 for (unsigned int j = 0 ; j < isoNames[n+1].size() ;j++ ){
00681 edm::Handle<edm::AssociationMap<edm::OneToValue< T , float > > > depMap;
00682 iEvent.getByLabel(isoNames[n+1].at(j),depMap);
00683 if (depMap.isValid()){
00684 typename edm::AssociationMap<edm::OneToValue< T , float > >::const_iterator mapi = depMap->find(recoecalcands[closestEcalCandIndex]);
00685 if (mapi!=depMap->end()) {
00686 histEtaIsoOfHltObjMatchToGen[n+1]->Fill( recoecalcands[closestEcalCandIndex]->eta(),mapi->val);
00687 histPhiIsoOfHltObjMatchToGen[n+1]->Fill( recoecalcands[closestEcalCandIndex]->phi(),mapi->val);
00688 histEtIsoOfHltObjMatchToGen[n+1] ->Fill( recoecalcands[closestEcalCandIndex]->et(), mapi->val);
00689 }
00690 }
00691 }
00692 }
00693 }
00694 }
00695 }
00696
00698
00699
00701
00702
00703 for (unsigned int i=0; i<recoecalcands.size(); i++) {
00705
00706
00707
00708
00709
00710
00711
00712
00713
00715
00716
00717
00718
00719
00720
00721 ethist[n] ->Fill(recoecalcands[i]->et() );
00722 etahist[n]->Fill(recoecalcands[i]->eta() );
00723 if (!noPhiPlots) phihist[n]->Fill(recoecalcands[i]->phi() );
00724
00726
00727
00729 if ( n+1 < numOfHLTCollectionLabels ) {
00730 if (plotiso[n+1]) {
00731 for (unsigned int j = 0 ; j < isoNames[n+1].size() ;j++ ){
00732 edm::Handle<edm::AssociationMap<edm::OneToValue< T , float > > > depMap;
00733 iEvent.getByLabel(isoNames[n+1].at(j),depMap);
00734 if (depMap.isValid()){
00735 typename edm::AssociationMap<edm::OneToValue< T , float > >::const_iterator mapi = depMap->find(recoecalcands[i]);
00736 if (mapi!=depMap->end()){
00737 etahistiso[n+1]->Fill(recoecalcands[i]->eta(),mapi->val);
00738 phihistiso[n+1]->Fill(recoecalcands[i]->phi(),mapi->val);
00739 ethistiso[n+1]->Fill(recoecalcands[i]->et(),mapi->val);
00740 }
00741 }
00742 }
00743 }
00744 }
00745 }
00746 }
00747
00748
00750
00752 unsigned int mtachedMcParts = 0;
00753 float mindist=0.3;
00754 if(n==0) mindist=0.5;
00755 for(unsigned int i =0; i < nCandCuts.at(n); ++i){
00756
00757 bool matchThis= false;
00758 math::XYZVector candDir=sortedGen[i].momentum();
00759 unsigned int closest = 0;
00760 double closestDr = 1000.;
00761 for(unsigned int trigOb = 0 ; trigOb < recoecalcands.size(); ++trigOb){
00762 double dr = DeltaR(recoecalcands[trigOb]->momentum(),candDir);
00763 if (dr < closestDr) {
00764 closestDr = dr;
00765 closest = trigOb;
00766 }
00767 if (closestDr > mindist) {
00768 closest = -1;
00769 } else {
00770 mtachedMcParts++;
00771 matchThis = true;
00772 }
00773 }
00774 if ( !matchThis ) {
00775 accepted = false;
00776 continue;
00777 }
00778
00779 ethistmatch[n] ->Fill( sortedGen[i].et() );
00780 if (sortedGen[i].et() > plotMinEtForEtaEffPlot) {
00781 etahistmatch[n]->Fill( sortedGen[i].eta() );
00782 if (!noPhiPlots) phihistmatch[n]->Fill( sortedGen[i].phi() );
00783 }
00785
00786
00788 if (n+1 < numOfHLTCollectionLabels){
00789 if (plotiso[n+1] ){
00790 for (unsigned int j = 0 ; j < isoNames[n+1].size() ;j++ ){
00791 edm::Handle<edm::AssociationMap<edm::OneToValue< T , float > > > depMap;
00792 iEvent.getByLabel(isoNames[n+1].at(j),depMap);
00793 if (depMap.isValid()){
00794 typename edm::AssociationMap<edm::OneToValue< T , float > >::const_iterator mapi = depMap->find(recoecalcands[closest]);
00795 if (mapi!=depMap->end()){
00796
00797 etahistisomatch[n+1]->Fill(sortedGen[i].eta(),mapi->val);
00798 phihistisomatch[n+1]->Fill(sortedGen[i].phi(),mapi->val);
00799 ethistisomatch[n+1]->Fill(sortedGen[i].et(),mapi->val);
00800 }
00801 }
00802 }
00803 }
00804 }
00805 }
00806
00807
00808 if (mtachedMcParts >= nCandCuts.at(n) && accepted == true)
00809 totalmatch->Fill(n+0.5);
00810 }
00811
00812 void
00813 EmDQM::endRun(edm::Run const &iRun, edm::EventSetup const &iSetup)
00814 {
00815
00816
00817
00818
00819 std::vector<std::string> labelsNeverFound;
00820
00821
00822
00823 BOOST_FOREACH(const edm::InputTag &tag, hltCollectionLabelsMissed)
00824 {
00825 if (hltCollectionLabelsFound.count(tag.encode()) == 0)
00826
00827 labelsNeverFound.push_back(tag.encode());
00828
00829 }
00830
00831 if (labelsNeverFound.empty())
00832 return;
00833
00834 std::sort(labelsNeverFound.begin(), labelsNeverFound.end());
00835
00836
00837
00838
00839 if (verbosity >= OUTPUT_WARNINGS)
00840 edm::LogWarning("EmDQM") << "There were some HLTCollectionLabels which were never found:";
00841
00842 BOOST_FOREACH(const edm::InputTag &tag, labelsNeverFound)
00843 {
00844 if (verbosity >= OUTPUT_ALL)
00845 edm::LogPrint("EmDQM") << " " << tag;
00846 }
00847 }
00848
00850
00852 void EmDQM::endJob()
00853 {
00854
00855 }
00856
00857
00858
00859 DEFINE_FWK_MODULE(EmDQM);