#include <HLTriggerOffline/Egamma/interface/EmDQM.h>
Definition at line 20 of file EmDQM.h.
EmDQM::EmDQM | ( | const edm::ParameterSet & | pset | ) | [explicit] |
Constructor.
Definition at line 44 of file EmDQM.cc.
References dbe, dirname_, gencut_, gencutCollection_, genEtaAcc, genEtAcc, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), isoNames, pdgGen, plotBounds, plotiso, reqNum, DQMStore::setCurrentFolder(), DQMStore::setVerbose(), theHLTCollectionLabels, theHLTOutputTypes, theNbins, thePtMax, and thePtMin.
00045 { 00046 00047 00048 00049 dbe = edm::Service < DQMStore > ().operator->(); 00050 dbe->setVerbose(0); 00051 00052 00054 // Read from configuration file // 00056 dirname_="HLT/HLTEgammaValidation/"+pset.getParameter<std::string>("@module_label"); 00057 dbe->setCurrentFolder(dirname_); 00058 00059 // paramters for generator study 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 // plotting paramters (untracked because they don't affect the physics) 00065 thePtMin = pset.getUntrackedParameter<double>("PtMin",0.); 00066 thePtMax = pset.getUntrackedParameter<double>("PtMax",1000.); 00067 theNbins = pset.getUntrackedParameter<unsigned int>("Nbins",40); 00068 00069 //preselction cuts 00070 gencutCollection_= pset.getParameter<edm::InputTag>("cutcollection"); 00071 gencut_ = pset.getParameter<int>("cutnum"); 00072 00074 // Read in the Vector of Parameter Sets. // 00075 // Information for each filter-step // 00077 std::vector<edm::ParameterSet> filters = 00078 pset.getParameter<std::vector<edm::ParameterSet> >("filters"); 00079 00080 for(std::vector<edm::ParameterSet>::iterator filterconf = filters.begin() ; filterconf != filters.end() ; filterconf++) 00081 { 00082 00083 theHLTCollectionLabels.push_back(filterconf->getParameter<edm::InputTag>("HLTCollectionLabels")); 00084 theHLTOutputTypes.push_back(filterconf->getParameter<unsigned int>("theHLTOutputTypes")); 00085 std::vector<double> bounds = filterconf->getParameter<std::vector<double> >("PlotBounds"); 00086 // If the size of plot "bounds" vector != 2, abort 00087 assert(bounds.size() == 2); 00088 plotBounds.push_back(std::pair<double,double>(bounds[0],bounds[1])); 00089 isoNames.push_back(filterconf->getParameter<std::vector<edm::InputTag> >("IsoCollections")); 00090 // If the size of the isoNames vector is not greater than zero, abort 00091 assert(isoNames.back().size()>0); 00092 if (isoNames.back().at(0).label()=="none") 00093 plotiso.push_back(false); 00094 else 00095 plotiso.push_back(true); 00096 00097 } // END of loop over parameter sets 00098 }
void EmDQM::analyze | ( | const edm::Event & | event, | |
const edm::EventSetup & | setup | |||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 175 of file EmDQM.cc.
References funct::abs(), e, eta, etagen, etgen, funct::exp(), MonitorElement::Fill(), gencut_, gencutCollection_, genEtaAcc, genEtAcc, TtGenEvtProducer_cfi::genEvt, int, edm::Handle< T >::isValid(), n, p, pdgGen, reqNum, funct::sin(), theHLTCollectionLabels, theHLTOutputTypes, theta, and total.
00176 { 00177 00179 // Check if there's enough gen particles // 00180 // of interest // 00182 edm::Handle< edm::View<reco::Candidate> > cutCounter; 00183 event.getByLabel(gencutCollection_,cutCounter); 00184 if (cutCounter->size() < (unsigned int)gencut_) { 00185 //edm::LogWarning("EmDQM") << "Less than "<< reqNum <<" gen particles with pdgId=" << pdgGen; 00186 return; 00187 } 00188 00189 00190 // fill L1 and HLT info 00191 // get objects possed by each filter 00192 edm::Handle<trigger::TriggerEventWithRefs> triggerObj; 00193 event.getByLabel("hltTriggerSummaryRAW",triggerObj); 00194 if(!triggerObj.isValid()) { 00195 edm::LogWarning("EmDQM") << "RAW-type HLT results not found, skipping event"; 00196 return; 00197 } 00198 00199 // total event number 00200 total->Fill(theHLTCollectionLabels.size()+0.5); 00201 00202 00204 // Fill generator info // 00206 edm::Handle<edm::HepMCProduct> genEvt; 00207 event.getByLabel("source", genEvt); 00208 00209 std::vector<HepMC::GenParticle> mcparts; 00210 const HepMC::GenEvent * myGenEvent = genEvt->GetEvent(); 00211 unsigned int ncand = 0; 00212 for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end(); ++p ) { 00213 00214 // If the ID number is not what we're looking for or 00215 // it's status is !=1, go to the next particle 00216 if ( !( abs((*p)->pdg_id())==pdgGen && (*p)->status()==1 ) ) continue; 00217 float eta = (*p)->momentum().eta(); 00218 float e = (*p)->momentum().e(); 00219 float theta = 2*atan(exp(-eta)); 00220 float Et = e*sin(theta); 00221 if ( fabs(eta)<genEtaAcc && Et > genEtAcc ) 00222 { 00223 ncand++; 00224 etgen->Fill(Et); 00225 etagen->Fill(eta); 00226 mcparts.push_back(*(*p)); 00227 } 00228 } // END of loop over Generated particles 00229 if (ncand >= reqNum) total->Fill(theHLTCollectionLabels.size()+1.5); 00230 00231 00232 00234 // Loop over filter modules // 00236 for(unsigned int n=0; n < theHLTCollectionLabels.size() ; n++) { 00237 // These numbers are from the Parameter Set, such as: 00238 // theHLTOutputTypes = cms.uint32(100) 00239 switch(theHLTOutputTypes[n]) 00240 { 00241 case 82: // Non-isolated Level 1 00242 fillHistos<l1extra::L1EmParticleCollection>(triggerObj,event,n,mcparts);break; 00243 case 83: // Isolated Level 1 00244 fillHistos<l1extra::L1EmParticleCollection>(triggerObj,event,n,mcparts);break; 00245 case 91: // Photon 00246 fillHistos<reco::RecoEcalCandidateCollection>(triggerObj,event,n,mcparts);break; 00247 case 92: // Electron 00248 fillHistos<reco::ElectronCollection>(triggerObj,event,n,mcparts);break; 00249 case 100: // TriggerCluster 00250 fillHistos<reco::RecoEcalCandidateCollection>(triggerObj,event,n,mcparts);break; 00251 default: 00252 throw(cms::Exception("Release Validation Error") << "HLT output type not implemented: theHLTOutputTypes[n]" ); 00253 } 00254 } // END of loop over filter modules 00255 }
void EmDQM::beginJob | ( | const edm::EventSetup & | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 105 of file EmDQM.cc.
References DQMStore::book1D(), DQMStore::book2D(), dbe, dirname_, etagen, etahist, etahistiso, etahistmatch, etgen, ethist, ethistiso, ethistmatch, i, label, NULL, plotBounds, plotiso, MonitorElement::setBinLabel(), DQMStore::setCurrentFolder(), theHLTCollectionLabels, theNbins, thePtMax, thePtMin, and total.
00106 { 00107 //edm::Service<TFileService> fs; 00108 dbe->setCurrentFolder(dirname_); 00109 00110 std::string histoname="total eff"; 00111 00112 total = dbe->book1D(histoname.c_str(),histoname.c_str(),theHLTCollectionLabels.size()+2,0,theHLTCollectionLabels.size()+2); 00113 total->setBinLabel(theHLTCollectionLabels.size()+1,"Total"); 00114 total->setBinLabel(theHLTCollectionLabels.size()+2,"Gen"); 00115 for (unsigned int u=0; u<theHLTCollectionLabels.size(); u++){total->setBinLabel(u+1,theHLTCollectionLabels[u].label().c_str());} 00116 00117 MonitorElement* tmphisto; 00118 MonitorElement* tmpiso; 00119 00120 histoname = "gen et"; 00121 etgen = dbe->book1D(histoname.c_str(),histoname.c_str(),theNbins,thePtMin,thePtMax); 00122 histoname = "gen eta"; 00123 etagen = dbe->book1D(histoname.c_str(),histoname.c_str(),theNbins,-2.7,2.7); 00124 00125 for(unsigned int i = 0; i< theHLTCollectionLabels.size() ; i++){ 00126 histoname = theHLTCollectionLabels[i].label()+"et"; 00127 tmphisto = dbe->book1D(histoname.c_str(),histoname.c_str(),theNbins,thePtMin,thePtMax); 00128 ethist.push_back(tmphisto); 00129 00130 histoname = theHLTCollectionLabels[i].label()+"eta"; 00131 tmphisto = dbe->book1D(histoname.c_str(),histoname.c_str(),theNbins,-2.7,2.7); 00132 etahist.push_back(tmphisto); 00133 00134 histoname = theHLTCollectionLabels[i].label()+"et MC matched"; 00135 tmphisto = dbe->book1D(histoname.c_str(),histoname.c_str(),theNbins,thePtMin,thePtMax); 00136 ethistmatch.push_back(tmphisto); 00137 00138 histoname = theHLTCollectionLabels[i].label()+"eta MC matched"; 00139 tmphisto = dbe->book1D(histoname.c_str(),histoname.c_str(),theNbins,-2.7,2.7); 00140 etahistmatch.push_back(tmphisto); 00141 00142 if(plotiso[i]) { 00143 histoname = theHLTCollectionLabels[i].label()+"eta isolation"; 00144 tmpiso = dbe->book2D(histoname.c_str(),histoname.c_str(),theNbins,-2.7,2.7,theNbins,plotBounds[i].first,plotBounds[i].second); 00145 } 00146 else { 00147 tmpiso = NULL; 00148 } 00149 etahistiso.push_back(tmpiso); 00150 00151 if(plotiso[i]){ 00152 histoname = theHLTCollectionLabels[i].label()+"et isolation"; 00153 tmpiso = dbe->book2D(histoname.c_str(),histoname.c_str(),theNbins,thePtMin,thePtMax,theNbins,plotBounds[i].first,plotBounds[i].second); 00154 } 00155 else{ 00156 tmpiso = NULL; 00157 } 00158 ethistiso.push_back(tmpiso); 00159 00160 } 00161 }
Reimplemented from edm::EDAnalyzer.
Definition at line 351 of file EmDQM.cc.
00351 { 00352 // Normalize the histograms 00353 // total->Scale(1./total->GetBinContent(1)); 00354 //for(unsigned int n= theHLTCollectionLabels.size()-1 ; n>0;n--){ 00355 // ethist[n]->Divide(ethist[n-1]); 00356 // etahist[n]->Divide(etahist[n-1]); 00357 //} 00358 }
void EmDQM::fillHistos | ( | edm::Handle< trigger::TriggerEventWithRefs > & | triggerObj, | |
const edm::Event & | iEvent, | |||
unsigned int | n, | |||
std::vector< HepMC::GenParticle > & | mcparts | |||
) | [inline, private] |
Definition at line 263 of file EmDQM.cc.
References asciidump::at, e, eta, etahist, etahistiso, etahistmatch, ethist, ethistiso, ethistmatch, funct::exp(), MonitorElement::Fill(), edm::Event::getByLabel(), i, isoNames, edm::Handle< T >::isValid(), j, label, plotiso, reqNum, funct::sin(), theHLTCollectionLabels, theHLTOutputTypes, theta, and total.
00264 { 00265 std::vector<edm::Ref<T> > recoecalcands; 00266 if (!( triggerObj->filterIndex(theHLTCollectionLabels[n])>=triggerObj->size() )){ // only process if available 00267 00269 // Retrieve saved filter objects // 00271 triggerObj->getObjects(triggerObj->filterIndex(theHLTCollectionLabels[n]),theHLTOutputTypes[n],recoecalcands); 00272 //Danger: special case, L1 non-isolated 00273 // needs to be merged with L1 iso 00274 if (theHLTOutputTypes[n] == 82) 00275 { 00276 std::vector<edm::Ref<T> > isocands; 00277 triggerObj->getObjects(triggerObj->filterIndex(theHLTCollectionLabels[n]),83,isocands); 00278 if (isocands.size()>0) 00279 { 00280 for (unsigned int i=0; i < isocands.size(); i++) 00281 recoecalcands.push_back(isocands[i]); 00282 } 00283 } // END of if theHLTOutputTypes == 82 00284 00286 // Fill filter objects into histograms // 00288 if (recoecalcands.size() != 0) 00289 { 00290 if (recoecalcands.size() >= reqNum ) 00291 total->Fill(n+0.5); 00292 for (unsigned int i=0; i<recoecalcands.size(); i++) { 00293 //unmatched 00294 ethist[n]->Fill(recoecalcands[i]->et() ); 00295 etahist[n]->Fill(recoecalcands[i]->eta() ); 00296 //matched 00297 math::XYZVector candDir=recoecalcands[i]->momentum(); 00298 std::vector<HepMC::GenParticle>::iterator closest = mcparts.end(); 00299 double closestDr = 1000. ; 00300 for(std::vector<HepMC::GenParticle>::iterator mc = mcparts.begin(); mc != mcparts.end() ; mc++){ 00301 math::XYZVector mcDir( mc->momentum().px(), 00302 mc->momentum().py(), 00303 mc->momentum().pz()); 00304 double dr = DeltaR(mcDir,candDir); 00305 if (dr < closestDr) { 00306 closestDr = dr; 00307 closest = mc; 00308 } 00309 } 00310 if (closest == mcparts.end()) 00311 edm::LogWarning("EmDQM") << "no MC match, may skew efficieny"; 00312 else { 00313 float eta = closest->momentum().eta(); 00314 float e = closest->momentum().e(); 00315 float theta = 2*atan(exp(-eta)); 00316 float Et = e*sin(theta); 00317 ethistmatch[n]->Fill( Et ); 00318 etahistmatch[n]->Fill( eta ); 00319 } 00320 00322 // Plot isolation variables (show the not-yet-cut // 00323 // isolation, i.e. associated to next filter) // 00325 if (n+1 < theHLTCollectionLabels.size()) // can't plot beyond last 00326 { 00327 if (plotiso[n+1] ){ 00328 for (unsigned int j = 0 ; j < isoNames[n+1].size() ;j++ ){ 00329 edm::Handle<edm::AssociationMap<edm::OneToValue< T , float > > > depMap; 00330 iEvent.getByLabel(isoNames[n+1].at(j).label(),depMap); 00331 if (depMap.isValid()){ //Map may not exist if only one candidate passes a double filter 00332 typename edm::AssociationMap<edm::OneToValue< T , float > >::const_iterator mapi = depMap->find(recoecalcands[i]); 00333 if (mapi!=depMap->end()){ // found candidate in isolation map! 00334 etahistiso[n+1]->Fill(recoecalcands[i]->eta(),mapi->val); 00335 ethistiso[n+1]->Fill(recoecalcands[i]->et(),mapi->val); 00336 break; // to avoid multiple filling we only look until we found the candidate once. 00337 } 00338 } 00339 } 00340 } 00341 } // END of if n+1 < then the number of hlt collections 00342 } 00343 } 00344 } 00345 }
DQMStore* EmDQM::dbe [private] |
std::string EmDQM::dirname_ [private] |
MonitorElement* EmDQM::etagen [private] |
std::vector<MonitorElement*> EmDQM::etahist [private] |
std::vector<MonitorElement*> EmDQM::etahistiso [private] |
std::vector<MonitorElement*> EmDQM::etahistmatch [private] |
MonitorElement* EmDQM::etgen [private] |
std::vector<MonitorElement*> EmDQM::ethist [private] |
std::vector<MonitorElement*> EmDQM::ethistiso [private] |
std::vector<MonitorElement*> EmDQM::ethistmatch [private] |
int EmDQM::gencut_ [private] |
edm::InputTag EmDQM::gencutCollection_ [private] |
double EmDQM::genEtaAcc [private] |
double EmDQM::genEtAcc [private] |
std::vector<std::vector<edm::InputTag> > EmDQM::isoNames [private] |
int EmDQM::pdgGen [private] |
std::vector<std::pair<double,double> > EmDQM::plotBounds [private] |
std::vector<bool> EmDQM::plotiso [private] |
unsigned int EmDQM::reqNum [private] |
std::vector<edm::InputTag> EmDQM::theHLTCollectionLabels [private] |
Definition at line 36 of file EmDQM.h.
Referenced by analyze(), beginJob(), EmDQM(), and fillHistos().
std::string EmDQM::theHltName [private] |
std::vector<int> EmDQM::theHLTOutputTypes [private] |
edm::InputTag EmDQM::theL1Seed [private] |
unsigned int EmDQM::theNbins [private] |
double EmDQM::thePtMax [private] |
double EmDQM::thePtMin [private] |
MonitorElement* EmDQM::total [private] |