CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/HLTriggerOffline/Egamma/src/EmDQM.cc

Go to the documentation of this file.
00001 
00002 //                    Header file for this                                    //
00004 #include "HLTriggerOffline/Egamma/interface/EmDQM.h"
00005 
00007 //                    Collaborating Class Header                              //
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 //#include "DataFormats/EgammaCandidates/interface/Photon.h"
00016 #include "DataFormats/L1Trigger/interface/L1EmParticle.h"
00017 #include "DataFormats/L1Trigger/interface/L1EmParticleFwd.h"
00018 //#include "SimDataFormats/HepMCProduct/interface/HepMCProduct.h"
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"
00029 //                           Root include files                               //
00031 #include "TFile.h"
00032 #include "TDirectory.h"
00033 #include "TH1F.h"
00034 #include <iostream>
00035 #include <string>
00036 #include <Math/VectorUtil.h>
00037 using namespace ROOT::Math::VectorUtil ;
00038 
00039 
00041 //                             Constructor                                    //
00043 EmDQM::EmDQM(const edm::ParameterSet& pset)  
00044 {
00045 
00046   dbe = edm::Service < DQMStore > ().operator->();
00047   dbe->setVerbose(0);
00048 
00049 
00051   //          Read from configuration file                  //
00053   dirname_="HLT/HLTEgammaValidation/"+pset.getParameter<std::string>("@module_label");
00054   dbe->setCurrentFolder(dirname_);
00055 
00056   triggerobjwithrefs = pset.getParameter<edm::InputTag>("triggerobject");
00057   // paramters for generator study 
00058   reqNum    = pset.getParameter<unsigned int>("reqNum");
00059   pdgGen    = pset.getParameter<int>("pdgGen");
00060   genEtaAcc = pset.getParameter<double>("genEtaAcc");
00061   genEtAcc  = pset.getParameter<double>("genEtAcc");
00062   // plotting paramters (untracked because they don't affect the physics)
00063   plotPtMin  = pset.getUntrackedParameter<double>("PtMin",0.);
00064   plotPtMax  = pset.getUntrackedParameter<double>("PtMax",1000.);
00065   plotEtaMax = pset.getUntrackedParameter<double>("EtaMax", 2.7);
00066   plotPhiMax = pset.getUntrackedParameter<double>("EtaMax", 3.15);
00067   plotBins   = pset.getUntrackedParameter<unsigned int>("Nbins",40);
00068   useHumanReadableHistTitles = pset.getUntrackedParameter<bool>("useHumanReadableHistTitles", false);
00069 
00070   //preselction cuts 
00071   gencutCollection_= pset.getParameter<edm::InputTag>("cutcollection");
00072   gencut_          = pset.getParameter<int>("cutnum");
00073 
00075   //         Read in the Vector of Parameter Sets.          //
00076   //           Information for each filter-step             //
00078   std::vector<edm::ParameterSet> filters = 
00079        pset.getParameter<std::vector<edm::ParameterSet> >("filters");
00080 
00081   int i = 0;
00082   for(std::vector<edm::ParameterSet>::iterator filterconf = filters.begin() ; filterconf != filters.end() ; filterconf++)
00083   {
00084 
00085     theHLTCollectionLabels.push_back(filterconf->getParameter<edm::InputTag>("HLTCollectionLabels"));
00086     theHLTOutputTypes.push_back(filterconf->getParameter<int>("theHLTOutputTypes"));
00087     // Grab the human-readable name, if it is not specified, use the Collection Label
00088     theHLTCollectionHumanNames.push_back(filterconf->getUntrackedParameter<std::string>("HLTCollectionHumanName",theHLTCollectionLabels[i].label()));
00089 
00090     std::vector<double> bounds = filterconf->getParameter<std::vector<double> >("PlotBounds");
00091     // If the size of plot "bounds" vector != 2, abort
00092     assert(bounds.size() == 2);
00093     plotBounds.push_back(std::pair<double,double>(bounds[0],bounds[1]));
00094     isoNames.push_back(filterconf->getParameter<std::vector<edm::InputTag> >("IsoCollections"));
00095     // If the size of the isoNames vector is not greater than zero, abort
00096     assert(isoNames.back().size()>0);
00097     if (isoNames.back().at(0).label()=="none") {
00098       plotiso.push_back(false);
00099     } else {
00100       plotiso.push_back(true);
00101     }
00102     i++;
00103   } // END of loop over parameter sets
00104 
00105   // Record number of HLTCollectionLabels
00106   numOfHLTCollectionLabels = theHLTCollectionLabels.size();
00107   
00108 }
00109 
00110 
00112 //       method called once each job just before starting event loop          //
00114 void 
00115 EmDQM::beginJob()
00116 {
00117   //edm::Service<TFileService> fs;
00118   dbe->setCurrentFolder(dirname_);
00119 
00121   //  Set up Histogram of Effiency vs Step.                 //
00122   //   theHLTCollectionLabels is a vector of InputTags      //
00123   //    from the configuration file.                        //
00125 
00126   std::string histName="total_eff";
00127   std::string histTitle = "total events passing";
00128   // This plot will have bins equal to 2+(number of
00129   //        HLTCollectionLabels in the config file)
00130   total = dbe->book1D(histName.c_str(),histTitle.c_str(),numOfHLTCollectionLabels+2,0,numOfHLTCollectionLabels+2);
00131   total->setBinLabel(numOfHLTCollectionLabels+1,"Total");
00132   total->setBinLabel(numOfHLTCollectionLabels+2,"Gen");
00133   for (unsigned int u=0; u<numOfHLTCollectionLabels; u++){total->setBinLabel(u+1,theHLTCollectionLabels[u].label().c_str());}
00134 
00135   histName="total_eff_MC_matched";
00136   histTitle="total events passing (mc matched)";
00137   totalmatch = dbe->book1D(histName.c_str(),histTitle.c_str(),numOfHLTCollectionLabels+2,0,numOfHLTCollectionLabels+2);
00138   totalmatch->setBinLabel(numOfHLTCollectionLabels+1,"Total");
00139   totalmatch->setBinLabel(numOfHLTCollectionLabels+2,"Gen");
00140   for (unsigned int u=0; u<numOfHLTCollectionLabels; u++){totalmatch->setBinLabel(u+1,theHLTCollectionLabels[u].label().c_str());}
00141 
00142   MonitorElement* tmphisto;
00143   MonitorElement* tmpiso;
00144 
00146   // Set up generator-level histograms                      //
00148   std::string pdgIdString;
00149   switch(pdgGen) {
00150   case 11:
00151     pdgIdString="Electron";break;
00152   case 22:
00153     pdgIdString="Photon";break;
00154   default:
00155     pdgIdString="Particle";
00156   }
00157 
00158   histName = "gen_et";
00159   histTitle= "E_{T} of " + pdgIdString + "s" ;
00160   etgen =  dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00161   histName = "gen_eta";
00162   histTitle= "#eta of "+ pdgIdString +"s " ;
00163   etagen = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00164   histName = "gen_phi";
00165   histTitle= "#phi of "+ pdgIdString +"s " ;
00166   phigen = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax);
00167 
00168   
00169 
00171   //  Set up histograms of HLT objects                      //
00173 
00174   // Determine what strings to use for histogram titles
00175   std::vector<std::string> HltHistTitle;
00176   if ( theHLTCollectionHumanNames.size() == numOfHLTCollectionLabels && useHumanReadableHistTitles ) {
00177     HltHistTitle = theHLTCollectionHumanNames;
00178   } else {
00179     for (unsigned int i =0; i < numOfHLTCollectionLabels; i++) {
00180       HltHistTitle.push_back(theHLTCollectionLabels[i].label());
00181     }
00182   }
00183  
00184   for(unsigned int i = 0; i< numOfHLTCollectionLabels ; i++){
00185     // Et distribution of HLT objects passing filter i
00186     histName = theHLTCollectionLabels[i].label()+"et_all";
00187     histTitle = HltHistTitle[i]+" Et (ALL)";
00188     tmphisto =  dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00189     ethist.push_back(tmphisto);
00190     
00191     // Eta distribution of HLT objects passing filter i
00192     histName = theHLTCollectionLabels[i].label()+"eta_all";
00193     histTitle = HltHistTitle[i]+" #eta (ALL)";
00194     tmphisto =  dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00195     etahist.push_back(tmphisto);          
00196 
00197     // Phi distribution of HLT objects passing filter i
00198     histName = theHLTCollectionLabels[i].label()+"phi_all";
00199     histTitle = HltHistTitle[i]+" #phi (ALL)";
00200     tmphisto =  dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax);
00201     phihist.push_back(tmphisto);
00202 
00203 
00204     // Et distribution of gen object matching HLT object passing filter i
00205     histName = theHLTCollectionLabels[i].label()+"et_MC_matched";
00206     histTitle = HltHistTitle[i]+" Et (MC matched)";
00207     tmphisto =  dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00208     ethistmatch.push_back(tmphisto);
00209     
00210     // Eta distribution of gen object matching HLT object passing filter i
00211     histName = theHLTCollectionLabels[i].label()+"eta_MC_matched";
00212     histTitle = HltHistTitle[i]+" #eta (MC matched)";
00213     tmphisto =  dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00214     etahistmatch.push_back(tmphisto);
00215 
00216     // Phi distribution of gen object matching HLT object passing filter i
00217     histName = theHLTCollectionLabels[i].label()+"phi_MC_matched";
00218     histTitle = HltHistTitle[i]+" #phi (MC matched)";
00219     tmphisto =  dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax);
00220     phihistmatch.push_back(tmphisto);
00221 
00222 
00223 
00224     // Et distribution of HLT object that is closest delta-R match to sorted gen particle(s)
00225     histName  = theHLTCollectionLabels[i].label()+"et";
00226     histTitle = HltHistTitle[i]+" Et";
00227     tmphisto  = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00228     histEtOfHltObjMatchToGen.push_back(tmphisto);
00229 
00230     // eta distribution of HLT object that is closest delta-R match to sorted gen particle(s)
00231     histName  = theHLTCollectionLabels[i].label()+"eta";
00232     histTitle = HltHistTitle[i]+" eta";
00233     tmphisto  = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00234     histEtaOfHltObjMatchToGen.push_back(tmphisto);
00235 
00236     // phi distribution of HLT object that is closest delta-R match to sorted gen particle(s)
00237     histName  = theHLTCollectionLabels[i].label()+"phi";
00238     histTitle = HltHistTitle[i]+" phi";
00239     tmphisto  = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax);
00240     histPhiOfHltObjMatchToGen.push_back(tmphisto);
00241 
00242 
00243 
00244     if (!plotiso[i]) {
00245       tmpiso = NULL;
00246       etahistiso.push_back(tmpiso);
00247       phihistiso.push_back(tmpiso);
00248       ethistiso.push_back(tmpiso);
00249       etahistisomatch.push_back(tmpiso);
00250       phihistisomatch.push_back(tmpiso);
00251       ethistisomatch.push_back(tmpiso);
00252       histEtaIsoOfHltObjMatchToGen.push_back(tmpiso);
00253       histPhiIsoOfHltObjMatchToGen.push_back(tmpiso);
00254       histEtIsoOfHltObjMatchToGen.push_back(tmpiso);
00255     } else {
00256       // 2D plot: Isolation values vs eta for all objects
00257       histName  = theHLTCollectionLabels[i].label()+"eta_isolation_all";
00258       histTitle = HltHistTitle[i]+" isolation vs #eta (all)";
00259       tmpiso    = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00260       etahistiso.push_back(tmpiso);
00261 
00262       // 2D plot: Isolation values vs phi for all objects
00263       histName  = theHLTCollectionLabels[i].label()+"phi_isolation_all";
00264       histTitle = HltHistTitle[i]+" isolation vs #phi (all)";
00265       tmpiso    = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00266       phihistiso.push_back(tmpiso);
00267 
00268 
00269       // 2D plot: Isolation values vs et for all objects
00270       histName  = theHLTCollectionLabels[i].label()+"et_isolation_all";
00271       histTitle = HltHistTitle[i]+" isolation vs Et (all)";
00272       tmpiso    = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00273       ethistiso.push_back(tmpiso);
00274 
00275       // 2D plot: Isolation values vs eta for matched objects
00276       histName  = theHLTCollectionLabels[i].label()+"eta_isolation_MC_matched";
00277       histTitle = HltHistTitle[i]+" isolation vs #eta (mc matched)";
00278       tmpiso    = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00279       etahistisomatch.push_back(tmpiso);
00280 
00281       // 2D plot: Isolation values vs phi for matched objects
00282       histName  = theHLTCollectionLabels[i].label()+"phi_isolation_MC_matched";
00283       histTitle = HltHistTitle[i]+" isolation vs #phi (mc matched)";
00284       tmpiso    = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00285       phihistisomatch.push_back(tmpiso);
00286 
00287 
00288       // 2D plot: Isolation values vs et for matched objects
00289       histName  = theHLTCollectionLabels[i].label()+"et_isolation_MC_matched";
00290       histTitle = HltHistTitle[i]+" isolation vs Et (mc matched)";
00291       tmpiso    = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00292       ethistisomatch.push_back(tmpiso);
00293 
00294 
00295       // 2D plot: Isolation values vs eta for HLT object that 
00296       // is closest delta-R match to sorted gen particle(s)
00297       histName  = theHLTCollectionLabels[i].label()+"eta_isolation";
00298       histTitle = HltHistTitle[i]+" isolation vs #eta";
00299       tmpiso    = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00300       histEtaIsoOfHltObjMatchToGen.push_back(tmpiso);
00301 
00302       // 2D plot: Isolation values vs phi for HLT object that
00303       // is closest delta-R match to sorted gen particle(s)
00304       histName  = theHLTCollectionLabels[i].label()+"phi_isolation";
00305       histTitle = HltHistTitle[i]+" isolation vs #phi";
00306       tmpiso    = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotPhiMax,plotPhiMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00307       histPhiIsoOfHltObjMatchToGen.push_back(tmpiso);
00308 
00309       // 2D plot: Isolation values vs et for HLT object that 
00310       // is closest delta-R match to sorted gen particle(s)
00311       histName  = theHLTCollectionLabels[i].label()+"et_isolation";
00312       histTitle = HltHistTitle[i]+" isolation vs Et";
00313       tmpiso    = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00314       histEtIsoOfHltObjMatchToGen.push_back(tmpiso);
00315     } // END of HLT histograms
00316 
00317   }
00318 }
00319 
00320 
00322 //                                Destructor                                  //
00324 EmDQM::~EmDQM(){
00325 }
00326 
00327 
00329 //                     method called to for each event                        //
00331 void 
00332 EmDQM::analyze(const edm::Event & event , const edm::EventSetup& setup)
00333 {
00334   
00336   //           Check if there's enough gen particles        //
00337   //             of interest                                //
00339   edm::Handle< edm::View<reco::Candidate> > cutCounter;
00340   event.getByLabel(gencutCollection_,cutCounter);
00341   if (cutCounter->size() < (unsigned int)gencut_) {
00342     //edm::LogWarning("EmDQM") << "Less than "<< gencut_ <<" gen particles with pdgId=" << pdgGen;
00343     return;
00344   }
00345 
00346 
00347   // fill L1 and HLT info
00348   // get objects possed by each filter
00349   edm::Handle<trigger::TriggerEventWithRefs> triggerObj;
00350   event.getByLabel(triggerobjwithrefs,triggerObj); 
00351   if(!triggerObj.isValid()) { 
00352     edm::LogWarning("EmDQM") << "triggerobjwithrefs=" << triggerobjwithrefs;
00353     return;
00354   }
00355 
00356 
00358   // Decide if this was an event of interest.               //
00359   //  Did the highest energy particles happen               //
00360   //  to have |eta| < 2.5 ?  Then continue.                 //
00362   edm::Handle< edm::View<reco::GenParticle> > genParticles;
00363   event.getByLabel("genParticles", genParticles);
00364   if(!genParticles.isValid()) { 
00365     edm::LogWarning("EmDQM") << "genParticles invalid";
00366     return;
00367   }
00368 
00369   std::vector<reco::GenParticle> allSortedGenParticles;
00370 
00371   for(edm::View<reco::GenParticle>::const_iterator currentGenParticle = genParticles->begin(); currentGenParticle != genParticles->end(); currentGenParticle++){
00372 
00373     if (  !( abs((*currentGenParticle).pdgId())==pdgGen  && (*currentGenParticle).status()==1 && (*currentGenParticle).et() > 2.0)  )  continue;
00374 
00375     reco::GenParticle tmpcand( *(currentGenParticle) );
00376     allSortedGenParticles.push_back(tmpcand);
00377   }
00378 
00379   std::sort(allSortedGenParticles.begin(), allSortedGenParticles.end(),pTGenComparator_);
00380 
00381   // Were enough high energy gen particles found?
00382   if (allSortedGenParticles.size() < gencut_) {
00383     // if no, throw event away
00384     return;
00385   }
00386 
00387   // We now have a sorted collection of all generated particles
00388   // with pdgId = pdgGen.
00389   // Loop over them to see if the top gen particles have eta within acceptance
00390   bool keepEvent = true;
00391   for (unsigned int i = 0 ; i < gencut_ ; i++ ) {
00392     bool inECALgap = fabs(allSortedGenParticles[i].eta()) > 1.4442 && fabs(allSortedGenParticles[i].eta()) < 1.556;
00393     if ( (fabs(allSortedGenParticles[i].eta()) > genEtaAcc) || inECALgap ) {
00394       //edm::LogWarning("EmDQM") << "Throwing event away. Gen particle with pdgId="<< allSortedGenParticles[i].pdgId() <<"; et="<< allSortedGenParticles[i].et() <<"; and eta="<< allSortedGenParticles[i].eta() <<" beyond acceptance.";
00395       keepEvent=false;
00396       break;
00397     }
00398   }
00399   if (!keepEvent) {
00400     return;
00401   }
00402 
00403 
00404   // It was an event worth keeping. Continue.
00405 
00407   //  Fill the bin labeled "Total"                          //
00408   //   This will be the number of events looked at.         //
00410   total->Fill(numOfHLTCollectionLabels+0.5);
00411   totalmatch->Fill(numOfHLTCollectionLabels+0.5);
00412 
00413 
00415   //               Fill generator info                      //
00417   // the gencut_ highest Et generator objects of the preselected type are our matches
00418 
00419   std::vector<reco::Particle> sortedGen;
00420   for(edm::View<reco::Candidate>::const_iterator genpart = cutCounter->begin(); genpart != cutCounter->end();genpart++){
00421     reco::Particle tmpcand(  genpart->charge(), genpart->p4(), genpart->vertex(),genpart->pdgId(),genpart->status() );
00422     sortedGen.push_back(tmpcand);
00423   }
00424   std::sort(sortedGen.begin(),sortedGen.end(),pTComparator_ );
00425 
00426   // Now the collection of gen particles is sorted by pt.
00427   // So, remove all particles from the collection so that we 
00428   // only have the top "1 thru gencut_" particles in it
00429   sortedGen.erase(sortedGen.begin()+gencut_,sortedGen.end());
00430   if (gencut_ != sortedGen.size() ){
00431     return;
00432   }
00433 
00434 
00435   for (unsigned int i = 0 ; i < gencut_ ; i++ ) {
00436     etgen ->Fill( sortedGen[i].et()  ); //validity has been implicitily checked by the cut on gencut_ above
00437     etagen->Fill( sortedGen[i].eta() );
00438     phigen->Fill( sortedGen[i].phi() );
00439   } // END of loop over Generated particles
00440   if (gencut_ >= reqNum) total->Fill(numOfHLTCollectionLabels+1.5); // this isn't really needed anymore keep for backward comp.
00441   if (gencut_ >= reqNum) totalmatch->Fill(numOfHLTCollectionLabels+1.5); // this isn't really needed anymore keep for backward comp.
00442           
00443 
00445   //            Loop over filter modules                    //
00447   for(unsigned int n=0; n < numOfHLTCollectionLabels ; n++) {
00448     // These numbers are from the Parameter Set, such as:
00449     //   theHLTOutputTypes = cms.uint32(100)
00450     switch(theHLTOutputTypes[n]) 
00451     {
00452       case trigger::TriggerL1NoIsoEG: // Non-isolated Level 1
00453         fillHistos<l1extra::L1EmParticleCollection>(triggerObj,event,n,sortedGen);break;
00454       case trigger::TriggerL1IsoEG: // Isolated Level 1
00455         fillHistos<l1extra::L1EmParticleCollection>(triggerObj,event,n,sortedGen);break;
00456       case trigger::TriggerPhoton: // Photon 
00457         fillHistos<reco::RecoEcalCandidateCollection>(triggerObj,event,n,sortedGen);break;
00458       case trigger::TriggerElectron: // Electron 
00459         fillHistos<reco::ElectronCollection>(triggerObj,event,n,sortedGen);break;
00460       case trigger::TriggerCluster: // TriggerCluster
00461         fillHistos<reco::RecoEcalCandidateCollection>(triggerObj,event,n,sortedGen);break;
00462       default: 
00463         throw(cms::Exception("Release Validation Error") << "HLT output type not implemented: theHLTOutputTypes[n]" );
00464     }
00465   } // END of loop over filter modules
00466 }
00467 
00468 
00470 // fillHistos                                                                 //
00471 //   Called by analyze method.                                                //
00473 template <class T> void EmDQM::fillHistos(edm::Handle<trigger::TriggerEventWithRefs>& triggerObj,const edm::Event& iEvent ,unsigned int n,std::vector<reco::Particle>& sortedGen)
00474 {
00475   std::vector<edm::Ref<T> > recoecalcands;
00476   if ( ( triggerObj->filterIndex(theHLTCollectionLabels[n])>=triggerObj->size() )){ // only process if available
00477     return;
00478   }
00479 
00481   //      Retrieve saved filter objects                     //
00483   triggerObj->getObjects(triggerObj->filterIndex(theHLTCollectionLabels[n]),theHLTOutputTypes[n],recoecalcands);
00484   //Danger: special case, L1 non-isolated
00485   // needs to be merged with L1 iso
00486   if (theHLTOutputTypes[n] == trigger::TriggerL1NoIsoEG){
00487     std::vector<edm::Ref<T> > isocands;
00488     triggerObj->getObjects(triggerObj->filterIndex(theHLTCollectionLabels[n]),trigger::TriggerL1IsoEG,isocands);
00489     if (isocands.size()>0) 
00490       {
00491         for (unsigned int i=0; i < isocands.size(); i++)
00492           recoecalcands.push_back(isocands[i]);
00493       }
00494   } // END of if theHLTOutputTypes == 82
00495   
00496 
00497   if (recoecalcands.size() < 1){ // stop if no object passed the previous filter
00498     return;
00499   }
00500 
00501 
00502   if (recoecalcands.size() >= reqNum ) 
00503     total->Fill(n+0.5);
00504 
00505 
00507   // check for validity                            //
00508   // prevents crash in CMSSW_3_1_0_pre6            //
00510   for (unsigned int j=0; j<recoecalcands.size(); j++){
00511     if(!( recoecalcands.at(j).isAvailable())){
00512       edm::LogError("EmDQM") << "Event content inconsistent: TriggerEventWithRefs contains invalid Refs" << std::endl << "invalid refs for: " << theHLTCollectionLabels[n].label();
00513       return;
00514     }
00515   }
00516 
00518   // Loop over the Generated Particles, and find the        //
00519   // closest HLT object match.                              //
00521   for (unsigned int i=0; i < gencut_; i++) {
00522     math::XYZVector currentGenParticleMomentum = sortedGen[i].momentum();
00523 
00524     float closestDeltaR = 0.5;
00525     int closestEcalCandIndex = -1;
00526     for (unsigned int j=0; j<recoecalcands.size(); j++) {
00527       float deltaR = DeltaR(recoecalcands[j]->momentum(),currentGenParticleMomentum);
00528 
00529       if (deltaR < closestDeltaR) {
00530         closestDeltaR = deltaR;
00531         closestEcalCandIndex = j;
00532       }
00533     }
00534 
00535     // If an HLT object was found within some delta-R
00536     // of this gen particle, store it in a histogram
00537     if ( closestEcalCandIndex >= 0 ) {
00538       histEtOfHltObjMatchToGen[n] ->Fill( recoecalcands[closestEcalCandIndex]->et()  );
00539       histEtaOfHltObjMatchToGen[n]->Fill( recoecalcands[closestEcalCandIndex]->eta() );
00540       histPhiOfHltObjMatchToGen[n]->Fill( recoecalcands[closestEcalCandIndex]->phi() );
00541       
00542       // Also store isolation info
00543       if (n+1 < numOfHLTCollectionLabels){ // can't plot beyond last
00544         if (plotiso[n+1] ){  // only plot if requested in config
00545           for (unsigned int j =  0 ; j < isoNames[n+1].size() ;j++  ){
00546             edm::Handle<edm::AssociationMap<edm::OneToValue< T , float > > > depMap; 
00547             iEvent.getByLabel(isoNames[n+1].at(j),depMap);
00548             if (depMap.isValid()){ //Map may not exist if only one candidate passes a double filter
00549               typename edm::AssociationMap<edm::OneToValue< T , float > >::const_iterator mapi = depMap->find(recoecalcands[closestEcalCandIndex]);
00550               if (mapi!=depMap->end()) {  // found candidate in isolation map! 
00551                 histEtaIsoOfHltObjMatchToGen[n+1]->Fill( recoecalcands[closestEcalCandIndex]->eta(),mapi->val);
00552                 histPhiIsoOfHltObjMatchToGen[n+1]->Fill( recoecalcands[closestEcalCandIndex]->phi(),mapi->val);
00553                 histEtIsoOfHltObjMatchToGen[n+1] ->Fill( recoecalcands[closestEcalCandIndex]->et(), mapi->val);
00554               }
00555             }
00556           }
00557         }
00558       }
00559     } // END of if closestEcalCandIndex >= 0
00560   }
00561 
00562 
00563 
00565   //  Loop over all HLT objects in this filter step, and    //
00566   //  fill histograms.                                      //
00568   //  bool foundAllMatches = false;
00569   //  unsigned int numOfHLTobjectsMatched = 0;
00570   for (unsigned int i=0; i<recoecalcands.size(); i++) {
00571     /*
00572     // See if this HLT object has a gen-level match
00573     float closestGenParticleDr = 99.0;
00574     for(unsigned int j =0; j < gencut_; j++) {
00575       math::XYZVector currentGenParticle = sortedGen[j].momentum();
00576 
00577       double currentDeltaR = DeltaR(recoecalcands[i]->momentum(),currentGenParticle);
00578       if ( currentDeltaR < closestGenParticleDr ) {
00579         closestGenParticleDr = currentDeltaR;
00580       }
00581     }
00582     // If this HLT object did not have a gen particle match, go to next HLT object
00583     if ( !(fabs(closestGenParticleDr < 0.3)) ) continue;
00584  
00585     numOfHLTobjectsMatched++;
00586     if (numOfHLTobjectsMatched >= gencut_) foundAllMatches=true;
00587     */
00588     // Fill HLT object histograms
00589     ethist[n] ->Fill(recoecalcands[i]->et() );
00590     etahist[n]->Fill(recoecalcands[i]->eta() );
00591     phihist[n]->Fill(recoecalcands[i]->phi() );
00592 
00594     //  Plot isolation variables (show the not-yet-cut        //
00595     //  isolation, i.e. associated to next filter)            //
00597     if ( n+1 < numOfHLTCollectionLabels ) { // can't plot beyond last
00598       if (plotiso[n+1]) {
00599         for (unsigned int j =  0 ; j < isoNames[n+1].size() ;j++  ){
00600           edm::Handle<edm::AssociationMap<edm::OneToValue< T , float > > > depMap; 
00601           iEvent.getByLabel(isoNames[n+1].at(j),depMap);
00602           if (depMap.isValid()){ //Map may not exist if only one candidate passes a double filter
00603             typename edm::AssociationMap<edm::OneToValue< T , float > >::const_iterator mapi = depMap->find(recoecalcands[i]);
00604             if (mapi!=depMap->end()){  // found candidate in isolation map! 
00605               etahistiso[n+1]->Fill(recoecalcands[i]->eta(),mapi->val);
00606               phihistiso[n+1]->Fill(recoecalcands[i]->phi(),mapi->val);
00607               ethistiso[n+1]->Fill(recoecalcands[i]->et(),mapi->val);
00608             }
00609           }
00610         }
00611       }
00612     } // END of if n+1 < then the number of hlt collections
00613   }
00614 
00615 
00617   //        Fill mc matched objects into histograms         //
00619   unsigned int mtachedMcParts = 0;
00620   float mindist=0.3;
00621   if(n==0) mindist=0.5; //low L1-resolution => allow wider matching 
00622   for(unsigned int i =0; i < gencut_; i++){
00623     //match generator candidate    
00624     bool matchThis= false;
00625     math::XYZVector candDir=sortedGen[i].momentum();
00626     unsigned int closest = 0;
00627     double closestDr = 1000.;
00628     for(unsigned int trigOb = 0 ; trigOb < recoecalcands.size(); trigOb++){
00629       double dr = DeltaR(recoecalcands[trigOb]->momentum(),candDir);
00630       if (dr < closestDr) {
00631         closestDr = dr;
00632         closest = trigOb;
00633       }
00634       if (closestDr > mindist) { // it's not really a "match" if it's that far away
00635         closest = -1;
00636       } else {
00637         mtachedMcParts++;
00638         matchThis = true;
00639       }
00640     }
00641     if ( !matchThis ) continue; // only plot matched candidates
00642     // fill coordinates of mc particle matching trigger object
00643     ethistmatch[n] ->Fill( sortedGen[i].et()  );
00644     etahistmatch[n]->Fill( sortedGen[i].eta() );
00645     phihistmatch[n]->Fill( sortedGen[i].phi() );
00647     //  Plot isolation variables (show the not-yet-cut        //
00648     //  isolation, i.e. associated to next filter)            //
00650     if (n+1 < numOfHLTCollectionLabels){ // can't plot beyond last
00651       if (plotiso[n+1] ){  // only plot if requested in config
00652         for (unsigned int j =  0 ; j < isoNames[n+1].size() ;j++  ){
00653           edm::Handle<edm::AssociationMap<edm::OneToValue< T , float > > > depMap; 
00654           iEvent.getByLabel(isoNames[n+1].at(j),depMap);
00655           if (depMap.isValid()){ //Map may not exist if only one candidate passes a double filter
00656             typename edm::AssociationMap<edm::OneToValue< T , float > >::const_iterator mapi = depMap->find(recoecalcands[closest]);
00657             if (mapi!=depMap->end()){  // found candidate in isolation map! 
00658               etahistisomatch[n+1]->Fill(sortedGen[i].eta(),mapi->val);
00659               phihistisomatch[n+1]->Fill(sortedGen[i].phi(),mapi->val);
00660               ethistisomatch[n+1]->Fill(sortedGen[i].et(),mapi->val);
00661             }
00662           }
00663         }
00664       }
00665     } // END of if n+1 < then the number of hlt collections
00666   }
00667   // fill total mc matched efficiency
00668   if (mtachedMcParts >= reqNum ) 
00669     totalmatch->Fill(n+0.5);
00670   
00671 
00672 }
00673 
00674 
00676 //      method called once each job just after ending the event loop          //
00678 void EmDQM::endJob(){
00679 
00680 }
00681 
00682 DEFINE_FWK_MODULE(EmDQM);