![]() |
![]() |
#include <HLTJetMETValidation.h>
Definition at line 52 of file HLTJetMETValidation.h.
HLTJetMETValidation::HLTJetMETValidation | ( | const edm::ParameterSet & | ps | ) | [explicit] |
Definition at line 5 of file HLTJetMETValidation.cc.
References _meGenHT, _meGenHTTrg, _meGenHTTrgLow, _meGenJetEta, _meGenJetEtaTrg, _meGenJetEtaTrgLow, _meGenJetPhi, _meGenJetPhiTrg, _meGenJetPhiTrgLow, _meGenJetPt, _meGenJetPtTrg, _meGenJetPtTrgLow, _meGenMET, _meGenMETTrg, _meGenMETTrgLow, _meRecoHT, _meRecoHTTrg, _meRecoHTTrgLow, _meRecoJetEta, _meRecoJetEtaTrg, _meRecoJetEtaTrgLow, _meRecoJetPhi, _meRecoJetPhiTrg, _meRecoJetPhiTrgLow, _meRecoJetPt, _meRecoJetPtTrg, _meRecoJetPtTrgLow, _meRecoMET, _meRecoMETTrg, _meRecoMETTrgLow, _triggerResults, DQMStore::book1D(), evtCnt, DQMStore::setCurrentFolder(), and triggerTag_.
: triggerEventObject_(ps.getUntrackedParameter<edm::InputTag>("triggerEventObject")), CaloJetAlgorithm( ps.getUntrackedParameter<edm::InputTag>( "CaloJetAlgorithm" ) ), GenJetAlgorithm( ps.getUntrackedParameter<edm::InputTag>( "GenJetAlgorithm" ) ), CaloMETColl( ps.getUntrackedParameter<edm::InputTag>( "CaloMETCollection" ) ), GenMETColl( ps.getUntrackedParameter<edm::InputTag>( "GenMETCollection" ) ), HLTriggerResults( ps.getParameter<edm::InputTag>( "HLTriggerResults" ) ), triggerTag_(ps.getUntrackedParameter<std::string>("DQMFolder","SingleJet")), _HLTPath(ps.getUntrackedParameter<edm::InputTag>("HLTPath")), _HLTLow(ps.getUntrackedParameter<edm::InputTag>("HLTLow")), outFile_(ps.getUntrackedParameter<std::string>("OutputFileName","")), HLTinit_(false), //JL writeFile_(ps.getUntrackedParameter<bool>("WriteFile",false)) { evtCnt=0; //Declare DQM Store DQMStore* store = &*edm::Service<DQMStore>(); if(store) { //Create the histograms store->setCurrentFolder(triggerTag_); _meRecoJetPt= store->book1D("_meRecoJetPt","Single Reconstructed Jet Pt",100,0,500); _meRecoJetPtTrg= store->book1D("_meRecoJetPtTrg","Single Reconstructed Jet Pt -- HLT Triggered",100,0,500); _meRecoJetPtTrgLow= store->book1D("_meRecoJetPtTrgLow","Single Reconstructed Jet Pt -- HLT Triggered Low",100,0,500); _meRecoJetEta= store->book1D("_meRecoJetEta","Single Reconstructed Jet Eta",100,-10,10); _meRecoJetEtaTrg= store->book1D("_meRecoJetEtaTrg","Single Reconstructed Jet Eta -- HLT Triggered",100,-10,10); _meRecoJetEtaTrgLow= store->book1D("_meRecoJetEtaTrgLow","Single Reconstructed Jet Eta -- HLT Triggered Low",100,-10,10); _meRecoJetPhi= store->book1D("_meRecoJetPhi","Single Reconstructed Jet Phi",100,-4.,4.); _meRecoJetPhiTrg= store->book1D("_meRecoJetPhiTrg","Single Reconstructed Jet Phi -- HLT Triggered",100,-4.,4.); _meRecoJetPhiTrgLow= store->book1D("_meRecoJetPhiTrgLow","Single Reconstructed Jet Phi -- HLT Triggered Low",100,-4.,4.); _meGenJetPt= store->book1D("_meGenJetPt","Single Generated Jet Pt",100,0,500); _meGenJetPtTrg= store->book1D("_meGenJetPtTrg","Single Generated Jet Pt -- HLT Triggered",100,0,500); _meGenJetPtTrgLow= store->book1D("_meGenJetPtTrgLow","Single Generated Jet Pt -- HLT Triggered Low",100,0,500); _meGenJetEta= store->book1D("_meGenJetEta","Single Generated Jet Eta",100,-10,10); _meGenJetEtaTrg= store->book1D("_meGenJetEtaTrg","Single Generated Jet Eta -- HLT Triggered",100,-10,10); _meGenJetEtaTrgLow= store->book1D("_meGenJetEtaTrgLow","Single Generated Jet Eta -- HLT Triggered Low",100,-10,10); _meGenJetPhi= store->book1D("_meGenJetPhi","Single Generated Jet Phi",100,-4.,4.); _meGenJetPhiTrg= store->book1D("_meGenJetPhiTrg","Single Generated Jet Phi -- HLT Triggered",100,-4.,4.); _meGenJetPhiTrgLow= store->book1D("_meGenJetPhiTrgLow","Single Generated Jet Phi -- HLT Triggered Low",100,-4.,4.); _meRecoMET= store->book1D("_meRecoMET","Reconstructed Missing ET",100,0,500); _meRecoMETTrg= store->book1D("_meRecoMETTrg","Reconstructed Missing ET -- HLT Triggered",100,0,500); _meRecoMETTrgLow= store->book1D("_meRecoMETTrgLow","Reconstructed Missing ET -- HLT Triggered Low",100,0,500); _meGenMET= store->book1D("_meGenMET","Generated Missing ET",100,0,500); _meGenMETTrg= store->book1D("_meGenMETTrg","Generated Missing ET -- HLT Triggered",100,0,500); _meGenMETTrgLow= store->book1D("_meGenMETTrgLow","Generated Missing ET -- HLT Triggered Low",100,0,500); _meGenHT= store->book1D("_meGenHT","Generated HT",100,0,1000); _meGenHTTrg= store->book1D("_meGenHTTrg","Generated HT -- HLT Triggered",100,0,1000); _meGenHTTrgLow= store->book1D("_meGenHTTrgLow","Generated HT -- HLT Triggered Low",100,0,1000); _meRecoHT= store->book1D("_meRecoHT","Reconstructed HT",100,0,1000); _meRecoHTTrg= store->book1D("_meRecoHTTrg","Reconstructed HT -- HLT Triggered",100,0,1000); _meRecoHTTrgLow= store->book1D("_meRecoHTTrgLow","Reconstructed HT -- HLT Triggered Low",100,0,1000); _triggerResults = store->book1D( "_triggerResults", "HLT Results", 200, 0, 200 ); } //if (writeFile_) printf("Initializing\n"); }
HLTJetMETValidation::~HLTJetMETValidation | ( | ) |
Definition at line 76 of file HLTJetMETValidation.cc.
{ }
void HLTJetMETValidation::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements edm::EDAnalyzer.
Definition at line 97 of file HLTJetMETValidation.cc.
References _HLTLow, _HLTPath, _meGenHT, _meGenHTTrg, _meGenHTTrgLow, _meGenJetEta, _meGenJetEtaTrg, _meGenJetEtaTrgLow, _meGenJetPhi, _meGenJetPhiTrg, _meGenJetPhiTrgLow, _meGenJetPt, _meGenJetPtTrg, _meGenJetPtTrgLow, _meGenMET, _meGenMETTrg, _meGenMETTrgLow, _meRecoHT, _meRecoHTTrg, _meRecoHTTrgLow, _meRecoJetEta, _meRecoJetEtaTrg, _meRecoJetEtaTrgLow, _meRecoJetPhi, _meRecoJetPhiTrg, _meRecoJetPhiTrgLow, _meRecoJetPt, _meRecoJetPtTrg, _meRecoJetPtTrgLow, _meRecoMET, _meRecoMETTrg, _meRecoMETTrgLow, CaloJetAlgorithm, CaloMETColl, evtCnt, MonitorElement::Fill(), HcalObjRepresent::Fill(), relval_steps::gen(), GenJetAlgorithm, HepMCValidationHelper::genMet(), GenMETColl, edm::Event::getByLabel(), getHLTResults(), HLTriggerResults, hltTriggerMap, i, edm::InputTag::label(), dt_dqm_sourceclient_common_cff::reco, trig_iter, triggerEventObject_, and edm::Event::triggerNames().
{ using namespace std; using namespace edm; using namespace reco; using namespace l1extra; using namespace trigger; evtCnt++; //get The triggerEvent Handle<TriggerEventWithRefs> trigEv; iEvent.getByLabel(triggerEventObject_,trigEv); // get TriggerResults object bool gotHLT=true; bool myTrig=false; bool myTrigLow=false; Handle<TriggerResults> hltresults,hltresultsDummy; iEvent.getByLabel(HLTriggerResults,hltresults); if (! hltresults.isValid() ) { //std::cout << " -- No HLTRESULTS"; //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No HLTRESULTS"; gotHLT=false; } if (gotHLT) { const edm::TriggerNames & triggerNames = iEvent.triggerNames(*hltresults); getHLTResults(*hltresults, triggerNames); // trig_iter=hltTriggerMap.find(MyTrigger); //trig_iter=hltTriggerMap.find(_HLTPath.label()); //std::cout << _HLTPath.label() <<std::endl; for( map<std::string,bool>::iterator ii=hltTriggerMap.begin(); ii!=hltTriggerMap.end(); ++ii) { //std::cout << (*ii).first << ": " << (*ii).second << std::endl; //std::cout << (*ii).first << " : " << ((*ii).first).find(_HLTPath.label()) << " : " << string::npos << std::endl; // if _HLTPath.label() is found in the string if ( ((*ii).first).find(_HLTPath.label()) != string::npos) trig_iter=ii; } if (trig_iter==hltTriggerMap.end()){ //std::cout << "Could not find trigger path with name: " << _probefilter.label() << std::endl; //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << "Could not find trigger path with name: " << _probefilter.label(); }else{ myTrig=trig_iter->second; } //trig_iter=hltTriggerMap.find(_HLTLow.label()); for( map<std::string,bool>::iterator ii=hltTriggerMap.begin(); ii!=hltTriggerMap.end(); ++ii) { // if _HLTPath.label() is found in the string if ( ((*ii).first).find(_HLTLow.label()) != string::npos) trig_iter=ii; } if (trig_iter==hltTriggerMap.end()){ //std::cout << "Could not find trigger path with name: " << _probefilter.label() << std::endl; //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << "Could not find trigger path with name: " << _probefilter.label(); }else{ myTrigLow=trig_iter->second; } } Handle<CaloJetCollection> caloJets,caloJetsDummy; iEvent.getByLabel( CaloJetAlgorithm, caloJets ); double calJetPt=-1.; double calJetEta=-999.; double calJetPhi=-999.; double calHT=0; if (caloJets.isValid()) { //Loop over the CaloJets and fill some histograms int jetInd = 0; for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end(); ++ cal ) { // std::cout << "CALO JET #" << jetInd << std::endl << cal->print() << std::endl; // h_ptCal->Fill( cal->pt() ); if (jetInd == 0){ //h_ptCalLeading->Fill( cal->pt() ); calJetPt=cal->pt(); calJetEta=cal->eta(); calJetPhi=cal->phi(); _meRecoJetPt->Fill( calJetPt ); _meRecoJetEta->Fill( calJetEta ); _meRecoJetPhi->Fill( calJetPhi ); if (myTrig) _meRecoJetPtTrg->Fill( calJetPt ); if (myTrig) _meRecoJetEtaTrg->Fill( calJetEta ); if (myTrig) _meRecoJetPhiTrg->Fill( calJetPhi ); if (myTrigLow) _meRecoJetPtTrgLow->Fill( calJetPt ); if (myTrigLow) _meRecoJetEtaTrgLow->Fill( calJetEta ); if (myTrigLow) _meRecoJetPhiTrgLow->Fill( calJetPhi ); //h_etaCalLeading->Fill( cal->eta() ); //h_phiCalLeading->Fill( cal->phi() ); jetInd++; } if (cal->pt()>30) { calHT+=cal->pt(); } } _meRecoHT->Fill( calHT ); if (myTrig) _meRecoHTTrg->Fill( calHT ); if (myTrigLow) _meRecoHTTrgLow->Fill( calHT ); }else{ //std::cout << " -- No CaloJets" << std::endl; //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No CaloJets"; } Handle<GenJetCollection> genJets,genJetsDummy; iEvent.getByLabel( GenJetAlgorithm, genJets ); double genJetPt=-1.; double genJetEta=-999.; double genJetPhi=-999.; double genHT=0; if (genJets.isValid()) { //Loop over the GenJets and fill some histograms int jetInd = 0; for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end(); ++ gen ) { // std::cout << "CALO JET #" << jetInd << std::endl << cal->print() << std::endl; // h_ptCal->Fill( cal->pt() ); if (jetInd == 0){ //h_ptCalLeading->Fill( gen->pt() ); genJetPt=gen->pt(); genJetEta=gen->eta(); genJetPhi=gen->phi(); _meGenJetPt->Fill( genJetPt ); _meGenJetEta->Fill( genJetEta ); _meGenJetPhi->Fill( genJetPhi ); if (myTrig) _meGenJetPtTrg->Fill( genJetPt ); if (myTrig) _meGenJetEtaTrg->Fill( genJetEta ); if (myTrig) _meGenJetPhiTrg->Fill( genJetPhi ); if (myTrigLow) _meGenJetPtTrgLow->Fill( genJetPt ); if (myTrigLow) _meGenJetEtaTrgLow->Fill( genJetEta ); if (myTrigLow) _meGenJetPhiTrgLow->Fill( genJetPhi ); //h_etaCalLeading->Fill( gen->eta() ); //h_phiCalLeading->Fill( gen->phi() ); //if (myTrig) h_ptCalTrig->Fill( gen->pt() ); jetInd++; } if (gen->pt()>30) { genHT+=gen->pt(); } } _meGenHT->Fill( genHT ); if (myTrig) _meGenHTTrg->Fill( genHT ); if (myTrigLow) _meGenHTTrgLow->Fill( genHT ); }else{ //std::cout << " -- No GenJets" << std::endl; //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No GenJets"; } edm::Handle<CaloMETCollection> recmet, recmetDummy; iEvent.getByLabel(CaloMETColl,recmet); double calMet=-1; if (recmet.isValid()) { typedef CaloMETCollection::const_iterator cmiter; //std::cout << "Size of MET collection" << recmet.size() << std::endl; for ( cmiter i=recmet->begin(); i!=recmet->end(); i++) { calMet = i->pt(); //mcalphi = i->phi(); //mcalsum = i->sumEt(); _meRecoMET -> Fill(calMet); if (myTrig) _meRecoMETTrg -> Fill(calMet); if (myTrigLow) _meRecoMETTrgLow -> Fill(calMet); } }else{ //std::cout << " -- No MET Collection with name: " << CaloMETColl << std::endl; //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No MET Collection with name: "<< CaloMETColl; } edm::Handle<GenMETCollection> genmet, genmetDummy; iEvent.getByLabel(GenMETColl,genmet); double genMet=-1; if (genmet.isValid()) { typedef GenMETCollection::const_iterator cmiter; //std::cout << "Size of GenMET collection" << recmet.size() << std::endl; for ( cmiter i=genmet->begin(); i!=genmet->end(); i++) { genMet = i->pt(); //mcalphi = i->phi(); //mcalsum = i->sumEt(); _meGenMET -> Fill(genMet); if (myTrig) _meGenMETTrg -> Fill(genMet); if (myTrigLow) _meGenMETTrgLow -> Fill(genMet); } }else{ //std::cout << " -- No GenMET Collection with name: " << GenMETColl << std::endl; //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No GenMET Collection with name: "<< GenMETColl; } }
void HLTJetMETValidation::endJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 85 of file HLTJetMETValidation.cc.
References outFile_, and writeFile_.
{ //Write DQM thing.. if(outFile_.size()>0) if (&*edm::Service<DQMStore>() && writeFile_) { edm::Service<DQMStore>()->save (outFile_); } }
void HLTJetMETValidation::getHLTResults | ( | const edm::TriggerResults & | hltresults, |
const edm::TriggerNames & | triggerNames | ||
) | [private] |
Definition at line 291 of file HLTJetMETValidation.cc.
References _triggerResults, edm::HLTGlobalStatus::accept(), accept(), MonitorElement::Fill(), HLTinit_, hltTriggerMap, edm::HLTGlobalStatus::size(), trig_iter, and edm::TriggerNames::triggerName().
Referenced by analyze().
{ int ntrigs=hltresults.size(); if (! HLTinit_){ HLTinit_=true; //if (writeFile_) std::cout << "Number of HLT Paths: " << ntrigs << std::endl; for (int itrig = 0; itrig != ntrigs; ++itrig){ std::string trigName = triggerNames.triggerName(itrig); // std::cout << "trigger " << itrig << ": " << trigName << std::endl; } } for (int itrig = 0; itrig != ntrigs; ++itrig){ std::string trigName = triggerNames.triggerName(itrig); bool accept=hltresults.accept(itrig); if (accept) _triggerResults->Fill(float(itrig)); // fill the trigger map typedef std::map<std::string,bool>::value_type valType; trig_iter=hltTriggerMap.find(trigName); if (trig_iter==hltTriggerMap.end()) hltTriggerMap.insert(valType(trigName,accept)); else trig_iter->second=accept; } }
edm::InputTag HLTJetMETValidation::_HLTLow [private] |
Definition at line 73 of file HLTJetMETValidation.h.
Referenced by analyze().
edm::InputTag HLTJetMETValidation::_HLTPath [private] |
Definition at line 72 of file HLTJetMETValidation.h.
Referenced by analyze().
MonitorElement* HLTJetMETValidation::_meGenHT [private] |
Definition at line 85 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meGenHTTrg [private] |
Definition at line 85 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meGenHTTrgLow [private] |
Definition at line 85 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement* HLTJetMETValidation::_meGenJetEta [private] |
Definition at line 81 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meGenJetEtaTrg [private] |
Definition at line 81 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
Definition at line 81 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement* HLTJetMETValidation::_meGenJetPhi [private] |
Definition at line 82 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meGenJetPhiTrg [private] |
Definition at line 82 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
Definition at line 82 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement* HLTJetMETValidation::_meGenJetPt [private] |
Definition at line 80 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meGenJetPtTrg [private] |
Definition at line 80 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meGenJetPtTrgLow [private] |
Definition at line 80 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement* HLTJetMETValidation::_meGenMET [private] |
Definition at line 84 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meGenMETTrg [private] |
Definition at line 84 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meGenMETTrgLow [private] |
Definition at line 84 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement* HLTJetMETValidation::_meRecoHT [private] |
Definition at line 86 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meRecoHTTrg [private] |
Definition at line 86 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meRecoHTTrgLow [private] |
Definition at line 86 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement* HLTJetMETValidation::_meRecoJetEta [private] |
Definition at line 77 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meRecoJetEtaTrg [private] |
Definition at line 77 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
Definition at line 77 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement* HLTJetMETValidation::_meRecoJetPhi [private] |
Definition at line 78 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meRecoJetPhiTrg [private] |
Definition at line 78 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
Definition at line 78 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement* HLTJetMETValidation::_meRecoJetPt [private] |
Definition at line 79 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meRecoJetPtTrg [private] |
Definition at line 79 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
Definition at line 79 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement* HLTJetMETValidation::_meRecoMET [private] |
Definition at line 83 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meRecoMETTrg [private] |
Definition at line 83 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
MonitorElement * HLTJetMETValidation::_meRecoMETTrgLow [private] |
Definition at line 83 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
Definition at line 87 of file HLTJetMETValidation.h.
Referenced by getHLTResults(), and HLTJetMETValidation().
Definition at line 67 of file HLTJetMETValidation.h.
Referenced by analyze().
Definition at line 67 of file HLTJetMETValidation.h.
Referenced by analyze().
int HLTJetMETValidation::evtCnt [private] |
Definition at line 91 of file HLTJetMETValidation.h.
Referenced by analyze(), and HLTJetMETValidation().
Definition at line 67 of file HLTJetMETValidation.h.
Referenced by analyze().
edm::InputTag HLTJetMETValidation::GenMETColl [private] |
Definition at line 67 of file HLTJetMETValidation.h.
Referenced by analyze().
bool HLTJetMETValidation::HLTinit_ [private] |
Definition at line 98 of file HLTJetMETValidation.h.
Referenced by getHLTResults().
Definition at line 67 of file HLTJetMETValidation.h.
Referenced by analyze().
std::map<std::string,bool> HLTJetMETValidation::hltTriggerMap [private] |
Definition at line 95 of file HLTJetMETValidation.h.
Referenced by analyze(), and getHLTResults().
std::vector<bool> HLTJetMETValidation::hlttrigs [private] |
Definition at line 94 of file HLTJetMETValidation.h.
std::string HLTJetMETValidation::MyTrigger [private] |
Definition at line 70 of file HLTJetMETValidation.h.
std::string HLTJetMETValidation::outFile_ [private] |
Definition at line 75 of file HLTJetMETValidation.h.
Referenced by endJob().
std::map<std::string,bool>::iterator HLTJetMETValidation::trig_iter [private] |
Definition at line 96 of file HLTJetMETValidation.h.
Referenced by analyze(), and getHLTResults().
InputTag of TriggerEventWithRefs to analyze.
Definition at line 66 of file HLTJetMETValidation.h.
Referenced by analyze().
std::string HLTJetMETValidation::triggerTag_ [private] |
Definition at line 70 of file HLTJetMETValidation.h.
Referenced by HLTJetMETValidation().
bool HLTJetMETValidation::writeFile_ [private] |
Definition at line 101 of file HLTJetMETValidation.h.
Referenced by endJob().