CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/DQM/Physics/src/EwkElecDQM.cc

Go to the documentation of this file.
00001 #include "DQM/Physics/src/EwkElecDQM.h"
00002 
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 
00007 #include "FWCore/ServiceRegistry/interface/Service.h"
00008 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00009 
00010 #include "DQMServices/Core/interface/DQMStore.h"
00011 #include "DQMServices/Core/interface/MonitorElement.h"
00012 
00013 #include "DataFormats/TrackReco/interface/Track.h"
00014 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00015 
00016 //#include "DataFormats/MuonReco/interface/Muon.h"
00017 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h" // I guess this is the right one??
00018 // also need Fwd.h file ???
00019 #include "DataFormats/METReco/interface/MET.h"
00020 #include "DataFormats/JetReco/interface/Jet.h"
00021 #include "DataFormats/JetReco/interface/CaloJet.h"
00022 #include "DataFormats/METReco/interface/CaloMET.h"
00023 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00024 #include "DataFormats/METReco/interface/CaloMETFwd.h"
00025 
00026 #include "DataFormats/GeometryVector/interface/Phi.h"
00027 
00028 #include "FWCore/Common/interface/TriggerNames.h"
00029 #include "FWCore/Framework/interface/Event.h"
00030 #include "DataFormats/Common/interface/TriggerResults.h"
00031 
00032 #include "DataFormats/Common/interface/View.h"
00033   
00034 using namespace edm;
00035 using namespace std;
00036 using namespace reco;
00037 
00038 EwkElecDQM::EwkElecDQM( const ParameterSet & cfg ) :
00039       // Input collections
00040       trigTag_(cfg.getUntrackedParameter<edm::InputTag> ("TrigTag", edm::InputTag("TriggerResults::HLT"))),
00041       //      muonTag_(cfg.getUntrackedParameter<edm::InputTag> ("MuonTag", edm::InputTag("muons"))),
00042       elecTag_(cfg.getUntrackedParameter<edm::InputTag> ("ElecTag", edm::InputTag("gsfElectrons"))), 
00043       metTag_(cfg.getUntrackedParameter<edm::InputTag> ("METTag", edm::InputTag("met"))),
00044       //      metIncludesMuons_(cfg.getUntrackedParameter<bool> ("METIncludesMuons", false)),
00045       jetTag_(cfg.getUntrackedParameter<edm::InputTag> ("JetTag", edm::InputTag("sisCone5CaloJets"))),
00046 
00047       // Main cuts 
00048       //      muonTrig_(cfg.getUntrackedParameter<std::string> ("MuonTrig", "HLT_Mu9")),
00049       //elecTrig_(cfg.getUntrackedParameter<std::vector< std::string > >("ElecTrig", "HLT_Ele10_SW_L1R")), 
00050       elecTrig_(cfg.getUntrackedParameter<std::vector< std::string > >("ElecTrig")), 
00051       //      ptCut_(cfg.getUntrackedParameter<double>("PtCut", 25.)),
00052       ptCut_(cfg.getUntrackedParameter<double>("PtCut", 10.)),
00053       //      etaCut_(cfg.getUntrackedParameter<double>("EtaCut", 2.1)),
00054       etaCut_(cfg.getUntrackedParameter<double>("EtaCut", 2.4)),
00055       sieieCutBarrel_(cfg.getUntrackedParameter<double>("SieieBarrel", 0.01)),
00056       sieieCutEndcap_(cfg.getUntrackedParameter<double>("SieieEndcap", 0.028)),
00057       detainCutBarrel_(cfg.getUntrackedParameter<double>("DetainBarrel", 0.0071)),
00058       detainCutEndcap_(cfg.getUntrackedParameter<double>("DetainEndcap", 0.0066)),
00059       //      isRelativeIso_(cfg.getUntrackedParameter<bool>("IsRelativeIso", true)),
00060       //      isCombinedIso_(cfg.getUntrackedParameter<bool>("IsCombinedIso", false)),
00061       //      isoCut03_(cfg.getUntrackedParameter<double>("IsoCut03", 0.1)),
00062       ecalIsoCutBarrel_(cfg.getUntrackedParameter<double>("EcalIsoCutBarrel", 5.7)),
00063       ecalIsoCutEndcap_(cfg.getUntrackedParameter<double>("EcalIsoCutEndcap", 5.0)),
00064       hcalIsoCutBarrel_(cfg.getUntrackedParameter<double>("HcalIsoCutBarrel", 8.1)),
00065       hcalIsoCutEndcap_(cfg.getUntrackedParameter<double>("HcalIsoCutEndcap", 3.4)),
00066       trkIsoCutBarrel_(cfg.getUntrackedParameter<double>("TrkIsoCutBarrel", 7.2)),
00067       trkIsoCutEndcap_(cfg.getUntrackedParameter<double>("TrkIsoCutEndcap", 5.1)),
00068       mtMin_(cfg.getUntrackedParameter<double>("MtMin", -999999)),
00069       mtMax_(cfg.getUntrackedParameter<double>("MtMax", 999999.)),
00070       metMin_(cfg.getUntrackedParameter<double>("MetMin", -999999.)),
00071       metMax_(cfg.getUntrackedParameter<double>("MetMax", 999999.)),
00072       //      acopCut_(cfg.getUntrackedParameter<double>("AcopCut", 2.)),
00073 
00074       // Muon quality cuts
00075       //      dxyCut_(cfg.getUntrackedParameter<double>("DxyCut", 0.2)),
00076       //      normalizedChi2Cut_(cfg.getUntrackedParameter<double>("NormalizedChi2Cut", 10.)),
00077       //      trackerHitsCut_(cfg.getUntrackedParameter<int>("TrackerHitsCut", 11)),
00078       //      isAlsoTrackerMuon_(cfg.getUntrackedParameter<bool>("IsAlsoTrackerMuon", true)),
00079 
00080       // Z rejection
00081       //      ptThrForZ1_(cfg.getUntrackedParameter<double>("PtThrForZ1", 20.)),
00082       //      ptThrForZ2_(cfg.getUntrackedParameter<double>("PtThrForZ2", 10.)),
00083 
00084       // Top rejection
00085       eJetMin_(cfg.getUntrackedParameter<double>("EJetMin", 999999.)),
00086       nJetMax_(cfg.getUntrackedParameter<int>("NJetMax", 999999))
00087       
00088 //       caloJetCollection_(cfg.getUntrackedParameter<edm:InputTag>("CaloJetCollection","sisCone5CaloJets"))
00089 
00090 
00091 
00092 {
00093 
00094   // access to dbe
00095   theDbe = Service<DQMStore>().operator->();
00096   theDbe->setCurrentFolder("Physics/EwkElecDQM");
00097   init_histograms();
00098 
00099 }
00100 
00101 void EwkElecDQM::beginRun(const Run& r, const EventSetup&) {
00102       nall = 0;
00103       nsel = 0;
00104 
00105       nrec = 0; 
00106       neid = 0;
00107       niso = 0; 
00108 //       nhlt = 0; 
00109 //       nmet = 0;
00110 }
00111 
00112 
00113 void EwkElecDQM::beginJob() {
00114 
00115 }
00116 
00117 void EwkElecDQM::init_histograms() {
00118 
00119   char chtitle[256] = "";
00120   for (int i=0; i<2; ++i) {
00121     
00122     //             pt_before_ = theDbe->book1D("PT_BEFORECUTS","Muon transverse momentum (global muon) [GeV],100,0.,100.);
00123     //             pt_after_ = theDbe->book1D("PT_LASTCUT","Muon transverse momentum (global muon) [GeV],100,0.,100.);
00124     
00125     pt_before_ = theDbe->book1D("PT_BEFORECUTS","Electron transverse momentum [GeV]",100,0.,100.);
00126     pt_after_ = theDbe->book1D("PT_LASTCUT","Electron transverse momentum [GeV]",100,0.,100.);
00127 
00128     eta_before_ = theDbe->book1D("ETA_BEFORECUTS","Electron pseudo-rapidity",50,-2.5,2.5);
00129     eta_after_ = theDbe->book1D("ETA_LASTCUT","Electron pseudo-rapidity",50,-2.5,2.5);
00130 
00131     sieiebarrel_before_ = theDbe->book1D("SIEIEBARREL_BEFORECUTS","Electron #sigma_{i#etai#eta} (barrel)",70,0.,0.07);
00132     sieiebarrel_after_ = theDbe->book1D("SIEIEBARREL_LASTCUT","Electron #sigma_{i#etai#eta} (barrel)",70,0.,0.07);
00133 
00134     sieieendcap_before_ = theDbe->book1D("SIEIEENDCAP_BEFORECUTS","Electron #sigma_{i#etai#eta} (endcap)",70,0.,0.07);
00135     sieieendcap_after_ = theDbe->book1D("SIEIEENDCAP_LASTCUT","Electron #sigma_{i#etai#eta} (endcap)",70,0.,0.07);
00136     
00137     detainbarrel_before_ = theDbe->book1D("DETAINBARREL_BEFORECUTS","Electron #Delta#eta_{in} (barrel)",40,-0.02,0.02);
00138     detainbarrel_after_ = theDbe->book1D("DETAINBARREL_LASTCUT","Electron #Delta#eta_{in} (barrel)",40,-0.02,0.02);
00139 
00140     detainendcap_before_ = theDbe->book1D("DETAINENDCAP_BEFORECUTS","Electron #Delta#eta_{in} (endcap)",40,-0.02,0.02);
00141     detainendcap_after_ = theDbe->book1D("DETAINENDCAP_LASTCUT","Electron #Delta#eta_{in} (endcap)",40,-0.02,0.02);
00142     
00143     //             dxy_before_ = theDbe->book1D("DXY_BEFORECUTS","Muon transverse distance to beam spot [cm]",100,-0.5,0.5);
00144     //             dxy_after_ = theDbe->book1D("DXY_LASTCUT","Muon transverse distance to beam spot [cm]",100,-0.5,0.5);
00145 
00146     //             chi2_before_ = theDbe->book1D("CHI2_BEFORECUTS","Normalized Chi2, inner track fit",100,0.,100.);
00147     //             chi2_after_ = theDbe->book1D("CHI2_LASTCUT","Normalized Chi2, inner track fit",100,0.,100.);
00148 
00149     //             nhits_before_ = theDbe->book1D("NHITS_BEFORECUTS","Number of hits, inner track",40,-0.5,39.5);
00150     //             nhits_after_ = theDbe->book1D("NHITS_LASTCUT","Number of hits, inner track",40,-0.5,39.5);
00151 
00152     //             tkmu_before_ = theDbe->book1D("TKMU_BEFORECUTS","Tracker-muon flag (for global muons)",2,-0.5,1.5);
00153     //             tkmu_after_ = theDbe->book1D("TKMU_LASTCUT","Tracker-muon flag (for global muons)",2,-0.5,1.5);
00154 
00155     //             if (isRelativeIso_) {
00156     //               if (isCombinedIso_) {
00157     //                 snprintf(chtitle, 255, "Relative (combined) isolation variable");
00158     //                 iso_before_ = theDbe->book1D("ISO_BEFORECUTS","Relative (combined) isolation variable",100, 0., 1.);
00159     //                 iso_after_ = theDbe->book1D("ISO_LASTCUT","Relative (combined) isolation variable",100, 0., 1.);
00160     //               } else {
00161     //                 snprintf(chtitle, 255, "Relative (tracker) isolation variable");
00162     //                 iso_before_ = theDbe->book1D("ISO_BEFORECUTS","Relative (tracker) isolation variable",100, 0., 1.);
00163     //                 iso_after_ = theDbe->book1D("ISO_LASTCUT","Relative (tracker) isolation variable",100, 0., 1.);
00164     //               }
00165     //             } else {
00166     //               if (isCombinedIso_) {
00167     //                 iso_before_ = theDbe->book1D("ISO_BEFORECUTS","Absolute (combined) isolation variable [GeV]",100, 0., 20.);
00168     //                 iso_after_ = theDbe->book1D("ISO_LASTCUT","Absolute (combined) isolation variable [GeV]",100, 0., 20.);
00169     //               } else {
00170     //                 iso_before_ = theDbe->book1D("ISO_BEFORECUTS","Absolute (tracker) isolation variable [GeV]",100, 0., 20.);
00171     //                 iso_after_ = theDbe->book1D("ISO_LASTCUT","Absolute (tracker) isolation variable [GeV]",100, 0., 20.);
00172     //               }
00173     //             }
00174 
00175             ecalisobarrel_before_ = theDbe->book1D("ECALISOBARREL_BEFORECUTS","Absolute electron ECAL isolation variable (barrel) [GeV]",50,0.,50.);
00176             ecalisobarrel_after_ = theDbe->book1D("ECALISOBARREL_LASTCUT","Absolute electron ECAL isolation variable (barrel) [GeV]",50,0.,50.);
00177 
00178             ecalisoendcap_before_ = theDbe->book1D("ECALISOENDCAP_BEFORECUTS","Absolute electron ECAL isolation variable (endcap) [GeV]",50,0.,50.);
00179             ecalisoendcap_after_ = theDbe->book1D("ECALISOENDCAP_LASTCUT","Absolute electron ECAL isolation variable (endcap) [GeV]",50,0.,50.);
00180 
00181             hcalisobarrel_before_ = theDbe->book1D("HCALISOBARREL_BEFORECUTS","Absolute electron HCAL isolation variable (barrel) [GeV]",50,0.,50.);
00182             hcalisobarrel_after_ = theDbe->book1D("HCALISOBARREL_LASTCUT","Absolute electron HCAL isolation variable (barrel) [GeV]",50,0.,50.);
00183 
00184             hcalisoendcap_before_ = theDbe->book1D("HCALISOENDCAP_BEFORECUTS","Absolute electron HCAL isolation variable (endcap) [GeV]",50,0.,50.);
00185             hcalisoendcap_after_ = theDbe->book1D("HCALISOENDCAP_LASTCUT","Absolute electron HCAL isolation variable (endcap) [GeV]",50,0.,50.);
00186 
00187             trkisobarrel_before_ = theDbe->book1D("TRKISOBARREL_BEFORECUTS","Absolute electron track isolation variable (barrel) [GeV]",50,0.,50.);
00188             trkisobarrel_after_ = theDbe->book1D("TRKISOBARREL_LASTCUT","Absolute electron track isolation variable (barrel) [GeV]",50,0.,50.);
00189 
00190             trkisoendcap_before_ = theDbe->book1D("TRKISOENDCAP_BEFORECUTS","Absolute electron track isolation variable (endcap) [GeV]",50,0.,50.);
00191             trkisoendcap_after_ = theDbe->book1D("TRKISOENDCAP_LASTCUT","Absolute electron track isolation variable (endcap) [GeV]",50,0.,50.);
00192 
00193 //             snprintf(chtitle, 255, "Trigger response (bit %s)", muonTrig_.data());
00194 //             trig_before_ = theDbe->book1D("TRIG_BEFORECUTS",chtitle,2,-0.5,1.5);
00195 //             trig_after_ = theDbe->book1D("TRIG_LASTCUT",chtitle,2,-0.5,1.5);
00196 
00197             //snprintf(chtitle, 255, "Trigger response (bit %s)", elecTrig_.data());
00198             trig_before_ = theDbe->book1D("TRIG_BEFORECUTS","Trigger response",2,-0.5,1.5); // elecTrig_ is now a vector of strings!
00199             trig_after_ = theDbe->book1D("TRIG_LASTCUT","Trigger response",2,-0.5,1.5);
00200 
00201             invmass_before_ = theDbe->book1D("INVMASS_BEFORECUTS","Di-electron invariant mass [GeV]",100,0.,200.);
00202             invmass_after_ = theDbe->book1D("INVMASS_AFTERCUTS","Di-electron invariant mass [GeV]",100,0.,200.);
00203 
00204             nelectrons_before_ = theDbe->book1D("NELECTRONS_BEFORECUTS","Number of electrons in event",10,-0.5,9.5);
00205             nelectrons_after_ = theDbe->book1D("NELECTRONS_AFTERCUTS","Number of electrons in event",10,-0.5,9.5);
00206 
00207             snprintf(chtitle, 255, "Transverse mass (%s) [GeV]", metTag_.label().data());
00208             mt_before_ = theDbe->book1D("MT_BEFORECUTS",chtitle,150,0.,300.);
00209             mt_after_ = theDbe->book1D("MT_LASTCUT",chtitle,150,0.,300.);
00210 
00211 //             snprintf(chtitle, 255, "Transverse mass (%s) [GeV]", metTag_.label().data());
00212 //             mt_before_ = theDbe->book1D("MT_BEFORECUTS",chtitle,150,0.,300.);
00213 //             mt_after_ = theDbe->book1D("MT_LASTCUT",chtitle,150,0.,300.);
00214 
00215             snprintf(chtitle, 255, "Missing transverse energy (%s) [GeV]", metTag_.label().data());
00216             met_before_ = theDbe->book1D("MET_BEFORECUTS",chtitle,100,0.,200.);
00217             met_after_ = theDbe->book1D("MET_LASTCUT",chtitle,100,0.,200.);
00218 
00219 //             snprintf(chtitle, 255, "MU-MET (%s) acoplanarity", metTag_.label().data());
00220 //             acop_before_ = theDbe->book1D("ACOP_BEFORECUTS",chtitle,50,0.,M_PI);
00221 //             acop_after_ = theDbe->book1D("ACOP_LASTCUT",chtitle,50,0.,M_PI);
00222 
00223 //             snprintf(chtitle, 255, "Z rejection: number of muons above %.2f GeV", ptThrForZ1_);
00224 //             nz1_before_ = theDbe->book1D("NZ1_BEFORECUTS",chtitle,10,-0.5,9.5);
00225 //             nz1_after_ = theDbe->book1D("NZ1_LASTCUT",chtitle,10,-0.5,9.5);
00226 
00227 //             snprintf(chtitle, 255, "Z rejection: number of muons above %.2f GeV", ptThrForZ2_);
00228 //             nz2_before_ = theDbe->book1D("NZ2_BEFORECUTS",chtitle,10,-0.5,9.5);
00229 //             nz2_after_ = theDbe->book1D("NZ2_LASTCUT",chtitle,10,-0.5,9.5);
00230 
00231              snprintf(chtitle, 255, "Number of jets (%s) above %.2f GeV", jetTag_.label().data(), eJetMin_);
00232              njets_before_ = theDbe->book1D("NJETS_BEFORECUTS",chtitle,10,-0.5,9.5);
00233              njets_after_  = theDbe->book1D("NJETS_LASTCUT",chtitle,10,-0.5,9.5);
00234 
00235              snprintf(chtitle, 255, "Jet with highest E_{T} (%s)", jetTag_.label().data());
00236              jet_et_before_       = theDbe->book1D("JETET1_BEFORECUTS",chtitle, 20, 0., 200.0);
00237              jet_et_after_       = theDbe->book1D("JETET1_AFTERCUTS",chtitle, 20, 0., 200.0);
00238 
00239              snprintf(chtitle, 255, "Eta of Jet with highest E_{T} (%s)", jetTag_.label().data());
00240              jet_eta_before_       = theDbe->book1D("JETETA1_BEFORECUTS",chtitle, 20, -5, 5);
00241              jet_eta_after_       = theDbe->book1D("JETETA1_AFTERCUTS",chtitle, 20, -5, 5);
00242 
00243 //           snprintf(chtitle, 255, "Jet with 2nd highest E_{T} (%s)", jetTag_.label().data());
00244 //           jet2_et_before      = theDbe->book1D("JETET2_BEFORECUTS",chtitle, 20, 0., 200.0);
00245 //           jet2_et_after       = theDbe->book1D("JETET2_AFTERCUTS",chtitle, 20, 0., 200.0);
00246              
00247 
00248 
00249       }
00250 }
00251 
00252 
00253 void EwkElecDQM::endJob() {
00254 }
00255 
00256 void EwkElecDQM::endRun(const Run& r, const EventSetup&) {
00257 
00258   // overall
00259   double all = nall;
00260   double esel = nsel/all;
00261   LogVerbatim("") << "\n>>>>>> SELECTION SUMMARY BEGIN >>>>>>>>>>>>>>>";
00262   LogVerbatim("") << "Total number of events analyzed: " << nall << " [events]";
00263   LogVerbatim("") << "Total number of events selected: " << nsel << " [events]";
00264   LogVerbatim("") << "Overall efficiency:             " << "(" << setprecision(4) << esel*100. <<" +/- "<< setprecision(2) << sqrt(esel*(1-esel)/all)*100. << ")%";
00265   
00266   double erec = nrec/all;
00267   double eeid = neid/all;
00268   double eiso = niso/all;
00269 //   double ehlt = nhlt/all;
00270 //   double emet = nmet/all;
00271   
00272   // general reconstruction step??
00273   double num = nrec;
00274   double eff = erec;
00275   double err = sqrt(eff*(1-eff)/all);
00276   LogVerbatim("") << "Passing Pt/Eta/Quality cuts:    " << num << " [events], (" << setprecision(4) << eff*100. <<" +/- "<< setprecision(2) << err*100. << ")%";
00277 
00278   // electron ID step  
00279   num = neid;
00280   eff = eeid;
00281   err = sqrt(eff*(1-eff)/all);
00282   double effstep = 0.;
00283   double errstep = 0.;
00284   if (nrec>0) effstep = eeid/erec;
00285   if (nrec>0) errstep = sqrt(effstep*(1-effstep)/nrec);
00286   LogVerbatim("") << "Passing eID cuts:         " << num << " [events], (" << setprecision(4) << eff*100. <<" +/- "<< setprecision(2) << err*100. << ")%, to previous step: (" <<  setprecision(4) << effstep*100. << " +/- "<< setprecision(2) << errstep*100. <<")%";
00287   
00288   // isolation step  
00289   num = niso;
00290   eff = eiso;
00291   err = sqrt(eff*(1-eff)/all);
00292   effstep = 0.;
00293   errstep = 0.;
00294   if (neid>0) effstep = eiso/eeid;
00295   if (neid>0) errstep = sqrt(effstep*(1-effstep)/neid);
00296   LogVerbatim("") << "Passing isolation cuts:         " << num << " [events], (" << setprecision(4) << eff*100. <<" +/- "<< setprecision(2) << err*100. << ")%, to previous step: (" <<  setprecision(4) << effstep*100. << " +/- "<< setprecision(2) << errstep*100. <<")%";
00297   
00298 //   // trigger step
00299 //   num = nhlt;
00300 //   eff = ehlt;
00301 //   err = sqrt(eff*(1-eff)/all);
00302 //   effstep = 0.;
00303 //   errstep = 0.;
00304 //   if (niso>0) effstep = ehlt/eiso;
00305 //   if (niso>0) errstep = sqrt(effstep*(1-effstep)/niso);
00306 //   LogVerbatim("") << "Passing HLT criteria:           " << num << " [events], (" << setprecision(4) << eff*100. <<" +/- "<< setprecision(2) << err*100. << ")%, to previous step: (" <<  setprecision(4) << effstep*100. << " +/- "<< setprecision(2) << errstep*100. <<")%";
00307   
00308   // trigger step
00309   num = nsel;
00310   eff = esel;
00311   err = sqrt(eff*(1-eff)/all);
00312   effstep = 0.;
00313   errstep = 0.;
00314   if (niso>0) effstep = esel/eiso;
00315   if (niso>0) errstep = sqrt(effstep*(1-effstep)/niso);
00316   LogVerbatim("") << "Passing HLT criteria:           " << num << " [events], (" << setprecision(4) << eff*100. <<" +/- "<< setprecision(2) << err*100. << ")%, to previous step: (" <<  setprecision(4) << effstep*100. << " +/- "<< setprecision(2) << errstep*100. <<")%";
00317   
00318 //   // met/acoplanarity cuts 
00319 //   num = nmet;
00320 //   eff = emet;
00321 //   err = sqrt(eff*(1-eff)/all);
00322 //   effstep = 0.;
00323 //   errstep = 0.;
00324 //   if (nhlt>0) effstep = emet/ehlt;
00325 //   if (nhlt>0) errstep = sqrt(effstep*(1-effstep)/nhlt);
00326 //   LogVerbatim("") << "Passing MET/acoplanarity cuts:  " << num << " [events], (" << setprecision(4) << eff*100. <<" +/- "<< setprecision(2) << err*100. << ")%, to previous step: (" <<  setprecision(4) << effstep*100. << " +/- "<< setprecision(2) << errstep*100. <<")%";
00327   
00328 //   // Z/top selection cuts ALSO LAST STEP so "sel" for "selection"
00329 //   num = nsel;
00330 //   eff = esel;
00331 //   err = sqrt(eff*(1-eff)/all);
00332 //   effstep = 0.;
00333 //   errstep = 0.;
00334 //   if (nmet>0) effstep = esel/emet;
00335 //   if (nmet>0) errstep = sqrt(effstep*(1-effstep)/nmet);
00336 //   LogVerbatim("") << "Passing Z/top rejection cuts:   " << num << " [events], (" << setprecision(4) << eff*100. <<" +/- "<< setprecision(2) << err*100. << ")%, to previous step: (" <<  setprecision(4) << effstep*100. << " +/- "<< setprecision(2) << errstep*100. <<")%";
00337   
00338   LogVerbatim("") << ">>>>>> SELECTION SUMMARY END   >>>>>>>>>>>>>>>\n";
00339 }
00340 
00341 void EwkElecDQM::analyze (const Event & ev, const EventSetup &) {
00342       
00343       // Reset global event selection flags
00344       bool rec_sel = false;
00345       bool eid_sel = false;
00346       bool iso_sel = false;
00347 //       bool hlt_sel = false;
00348       bool all_sel = false;
00349 
00350 //       // Muon collection
00351 //       Handle<View<Muon> > muonCollection;
00352 //       if (!ev.getByLabel(muonTag_, muonCollection)) {
00353 //             LogWarning("") << ">>> Muon collection does not exist !!!";
00354 //             return;
00355 //       }
00356 //       unsigned int muonCollectionSize = muonCollection->size();
00357 
00358       // Electron collection
00359       Handle<View<GsfElectron> > electronCollection;
00360       if (!ev.getByLabel(elecTag_, electronCollection)) {
00361         //LogWarning("") << ">>> Electron collection does not exist !!!";
00362         return;
00363       }
00364       unsigned int electronCollectionSize = electronCollection->size();
00365 
00366       // Beam spot
00367       Handle<reco::BeamSpot> beamSpotHandle;
00368       if (!ev.getByLabel(InputTag("offlineBeamSpot"), beamSpotHandle)) {
00369         //LogWarning("") << ">>> No beam spot found !!!";
00370         return;
00371       }
00372         // MET
00373       double met_px = 0.;
00374       double met_py = 0.;
00375       Handle<View<MET> > metCollection;
00376       if (!ev.getByLabel(metTag_, metCollection)) {
00377         //LogWarning("") << ">>> MET collection does not exist !!!";
00378         return;
00379       }
00380       const MET& met = metCollection->at(0);
00381       met_px = met.px();
00382       met_py = met.py();
00383 //       if (!metIncludesMuons_) {
00384 //             for (unsigned int i=0; i<muonCollectionSize; i++) {
00385 //                   const Muon& mu = muonCollection->at(i);
00386 //                   if (!mu.isGlobalMuon()) continue;
00387 //                   met_px -= mu.px();
00388 //                   met_py -= mu.py();
00389 //             }
00390 //       }
00391       double met_et = sqrt(met_px*met_px+met_py*met_py);
00392       LogTrace("") << ">>> MET, MET_px, MET_py: " << met_et << ", " << met_px << ", " << met_py << " [GeV]";
00393       met_before_->Fill(met_et);
00394 
00395       // Trigger
00396       Handle<TriggerResults> triggerResults;
00397       if (!ev.getByLabel(trigTag_, triggerResults)) {
00398         //LogWarning("") << ">>> TRIGGER collection does not exist !!!";
00399         return;
00400       }
00401       const edm::TriggerNames & trigNames = ev.triggerNames(*triggerResults);
00402       bool trigger_fired = false;
00403       /*
00404       for (unsigned int i=0; i<triggerResults->size(); i++) {
00405             if (triggerResults->accept(i)) {
00406                   LogTrace("") << "Accept by: " << i << ", Trigger: " << trigNames.triggerName(i);
00407             }
00408       }
00409       */
00410 
00411       // the following gives error on CRAFT08 data where itrig1=19 (vector index out of range)
00412       /*
00413       int itrig1 = trigNames.triggerIndex(muonTrig_);
00414       if (triggerResults->accept(itrig1)) trigger_fired = true;
00415       */
00416       //suggested replacement: lm250909
00417       for (unsigned int i=0; i<triggerResults->size(); i++) 
00418         {
00419           std::string trigName = trigNames.triggerName(i);
00420           for (unsigned int j = 0; j < elecTrig_.size(); j++)
00421             {
00422               if ( trigName == elecTrig_.at(j) && triggerResults->accept(i)) 
00423                 {
00424                   trigger_fired = true;
00425                 }
00426             }
00427         }
00428 
00429 
00430       LogTrace("") << ">>> Trigger bit: " << trigger_fired << " for one of ( " ;
00431       for (unsigned int k = 0; k < elecTrig_.size(); k++)
00432         {
00433           LogTrace("") << elecTrig_.at(k) << " ";
00434         }
00435       LogTrace("") << ")";
00436       trig_before_->Fill(trigger_fired);
00437 
00438 //       // Loop to reject/control Z->mumu is done separately
00439 //       unsigned int nmuonsForZ1 = 0;
00440 //       unsigned int nmuonsForZ2 = 0;
00441 //       for (unsigned int i=0; i<muonCollectionSize; i++) {
00442 //             const Muon& mu = muonCollection->at(i);
00443 //             if (!mu.isGlobalMuon()) continue;
00444 //             double pt = mu.pt();
00445 //             if (pt>ptThrForZ1_) nmuonsForZ1++;
00446 //             if (pt>ptThrForZ2_) nmuonsForZ2++;
00447 //       }
00448 //       LogTrace("") << "> Z rejection: muons above " << ptThrForZ1_ << " [GeV]: " << nmuonsForZ1;
00449 //       LogTrace("") << "> Z rejection: muons above " << ptThrForZ2_ << " [GeV]: " << nmuonsForZ2;
00450 //       nz1_before_->Fill(nmuonsForZ1);
00451 //       nz2_before_->Fill(nmuonsForZ2);
00452       
00453       // Jet collection
00454       Handle<View<Jet> > jetCollection;
00455       if (!ev.getByLabel(jetTag_, jetCollection)) {
00456         //LogError("") << ">>> JET collection does not exist !!!";
00457         return;
00458       }      
00459       float electron_et   = -8.0;
00460       float electron_eta  = -8.0;
00461       float electron_phi  = -8.0;
00462       float electron2_et  = -9.0;
00463       float electron2_eta = -9.0;
00464       float electron2_phi = -9.0;
00465       //need to get some electron info so jets can be cleaned of them
00466       for (unsigned int i=0; i<electronCollectionSize; i++) 
00467         {
00468           const GsfElectron& elec = electronCollection->at(i);   
00469 
00470           if (i < 1)
00471             {
00472                electron_et   = elec.pt();
00473                electron_eta  = elec.eta();
00474                electron_phi  = elec.phi();
00475             }
00476           if (i == 2)
00477             {
00478                electron2_et  = elec.pt();
00479                electron2_eta = elec.eta();
00480                electron2_phi = elec.phi();
00481             }       
00482         }
00483 
00484       float jet_et    = -8.0;
00485       float jet_eta   = -8.0;
00486       int   jet_count = 0;
00487       float jet2_et   = -9.0;
00488       unsigned int jetCollectionSize = jetCollection->size();
00489       int njets = 0;
00490       for (unsigned int i=0; i<jetCollectionSize; i++) {
00491         const Jet& jet = jetCollection->at(i);
00492         
00493         float jet_current_et = jet.et();
00494 //      cout << "jet_current_et " << jet_current_et << endl;
00495         // if it overlaps with electron, it is not a jet
00496         if ( electron_et>0.0 && fabs(jet.eta()-electron_eta ) < 0.2 
00497              && calcDeltaPhi(jet.phi(), electron_phi ) < 0.2)
00498           continue;
00499         if ( electron2_et>0.0&& fabs(jet.eta()-electron2_eta) < 0.2 
00500              && calcDeltaPhi(jet.phi(), electron2_phi) < 0.2) 
00501           continue;
00502 
00503         // if it has too low Et, throw away
00504 //      if (jet_current_et < eJetMin_) continue; //Keep if only want to plot above jet cut
00505 
00506         if (jet.et()>eJetMin_) 
00507           {
00508             njets++;
00509             jet_count++;
00510           }
00511         if (jet_current_et > jet_et) 
00512           {
00513             jet2_et  = jet_et;  // 2nd highest jet get's et from current highest
00514             jet_et   = jet.et(); // current highest jet gets et from the new highest
00515             jet_eta  = jet.eta();
00516           } 
00517         else if (jet_current_et > jet2_et) 
00518           {
00519             jet2_et  = jet.et();
00520           }
00521       }
00522 
00523       //Fill After all electron cuts (or both before and after)
00524       if (jet_et>10) //don't want low energy "jets"
00525         {
00526           jet_et_before_   ->Fill(jet_et);
00527 //        jet2_et_before_  ->Fill(jet2_et);
00528           jet_eta_before_  ->Fill(jet_eta);
00529         }
00530 
00531 
00532       LogTrace("") << ">>> Total number of jets: " << jetCollectionSize;
00533       LogTrace("") << ">>> Number of jets above " << eJetMin_ << " [GeV]: " << njets;
00534       njets_before_->Fill(njets);
00535 
00536       // Start counting
00537       nall++;
00538 
00539       // Histograms per event should be done only once, so keep track of them
00540       bool hlt_hist_done = false;
00541       //bool minv_hist_done = false;
00542        bool met_hist_done = false;
00543 //       bool nz1_hist_done = false;
00544 //       bool nz2_hist_done = false;
00545       bool njets_hist_done = false;
00546 
00547       // Central selection criteria
00548       //const int NFLAGS = 13; // number of individual selection criteria
00549       const int NFLAGS = 11; // number of individual selection criteria
00550       // 0: pt cut           | rec
00551       // 1: eta cut          | rec
00552       // 2: sieie            | eid
00553       // 3: detain           | eid
00554       // 4: ecal iso         | iso
00555       // 5: hcal iso         | iso
00556       // 6: trk iso          | iso
00557       // 7: trigger fired    | hlt/all
00558       bool electron_sel[NFLAGS];
00559 
00560       // for invariant mass calculation
00561       // keep track of highest-pt electrons for initial (RECO) electrons 
00562       // and "good" electrons (passing all cuts)
00563       // first array dimension is for first or second good electron
00564       // second array dimension is for relevant quantities of good electron
00565       //    [0]: 1 for electron found or 0 for not found (if 0, no other quantities filled)
00566       //    [1]: mSqr
00567       //    [2]: E
00568       //    [3]: px
00569       //    [4]: py
00570       //    [5]: pz
00571       // inv mass = sqrt(m_1^2 + m_2^2 + 2*(E_1*E_2 - (px1*px2 + py1*py2 + pz1+pz2) ) )
00572       double electron[2][6];
00573       double goodElectron[2][6];
00574       nGoodElectrons = 0;
00575       for (unsigned int i = 0; i < 2; i++)
00576         {
00577           for (unsigned int j = 0; j < 6; j++)
00578             {
00579               electron[i][j] = 0.;
00580               goodElectron[i][j] = 0.;
00581             }
00582         }
00583 
00584       for (unsigned int i=0; i<electronCollectionSize; i++) 
00585         {
00586           for (int j=0; j<NFLAGS; ++j) 
00587             {
00588               electron_sel[j] = false;
00589             }
00590           
00591           const GsfElectron& elec = electronCollection->at(i);
00592           //if (!mu.isGlobalMuon()) continue;
00593           //if (mu.globalTrack().isNull()) continue;
00594           //if (mu.innerTrack().isNull()) continue;
00595           
00596           LogTrace("") << "> elec: processing electron number " << i << "...";
00597           //reco::TrackRef gm = mu.globalTrack();
00598           //reco::TrackRef tk = mu.innerTrack();
00599           // should have stuff for electron track?
00600           
00601           if (i < 2)
00602             {
00603               electron[i][0] = 1.;
00604               electron[i][1] = elec.massSqr();
00605               electron[i][2] = elec.energy();
00606               electron[i][3] = elec.px();
00607               electron[i][4] = elec.py();
00608               electron[i][5] = elec.pz();
00609             }
00610 
00611           // Pt,eta cuts
00612           double pt = elec.pt();
00613           double px = elec.px();
00614           double py = elec.py();
00615           double eta = elec.eta();
00616           LogTrace("") << "\t... pt, eta: " << pt << " [GeV], " << eta;;
00617           if (pt>ptCut_) electron_sel[0] = true; 
00618           if (fabs(eta)<etaCut_) electron_sel[1] = true; 
00619           
00620           bool isBarrel = false;
00621           bool isEndcap = false;
00622           if (eta < 1.4442 && eta > -1.4442)
00623             {
00624               isBarrel = true;
00625             }
00626           else if ((eta > 1.56 && eta < 2.4) || (eta < -1.56 && eta > -2.4))
00627             {
00628               isEndcap = true;
00629             }
00630           
00631           //             // d0, chi2, nhits quality cuts
00632           //             double dxy = tk->dxy(beamSpotHandle->position());
00633           //             double normalizedChi2 = gm->normalizedChi2();
00634           //             double trackerHits = tk->numberOfValidHits();
00635           //             LogTrace("") << "\t... dxy, normalizedChi2, trackerHits, isTrackerMuon?: " << dxy << " [cm], " << normalizedChi2 << ", " << trackerHits << ", " << mu.isTrackerMuon();
00636           //             if (fabs(dxy)<dxyCut_) muon_sel[2] = true; 
00637           //             if (normalizedChi2<normalizedChi2Cut_) muon_sel[3] = true; 
00638           //             if (trackerHits>=trackerHitsCut_) muon_sel[4] = true; 
00639           //             if (mu.isTrackerMuon()) muon_sel[5] = true; 
00640           
00641           pt_before_->Fill(pt);
00642           eta_before_->Fill(eta);
00643           //             dxy_before_->Fill(dxy);
00644           //             chi2_before_->Fill(normalizedChi2);
00645           //             nhits_before_->Fill(trackerHits);
00646           //             tkmu_before_->Fill(mu.isTrackerMuon());
00647           
00648           
00649           // Electron ID cuts
00650           double sieie = (double) elec.sigmaIetaIeta();
00651           double detain = (double) elec.deltaEtaSuperClusterTrackAtVtx(); // think this is detain
00652           if (sieie < sieieCutBarrel_ && isBarrel) electron_sel[2] = true; 
00653           if (sieie < sieieCutEndcap_ && isEndcap) electron_sel[2] = true; 
00654           if (detain < detainCutBarrel_ && isBarrel) electron_sel[3] = true; 
00655           if (detain < detainCutEndcap_ && isEndcap) electron_sel[3] = true; 
00656           if (isBarrel)
00657             {
00658               LogTrace("") << "\t... sieie value " << sieie << " (barrel), pass? " << electron_sel[2]; 
00659               LogTrace("") << "\t... detain value " << detain << " (barrel), pass? " << electron_sel[3]; 
00660             }
00661           else if (isEndcap)
00662             {
00663               LogTrace("") << "\t... sieie value " << sieie << " (endcap), pass? " << electron_sel[2]; 
00664               LogTrace("") << "\t... detain value " << detain << " (endcap), pass? " << electron_sel[2]; 
00665             }
00666           
00667           if (isBarrel) 
00668             {
00669               sieiebarrel_before_->Fill(sieie);
00670               detainbarrel_before_->Fill(detain);
00671             }
00672           else if (isEndcap) 
00673             {
00674               sieieendcap_before_->Fill(sieie);
00675               detainendcap_before_->Fill(detain);
00676             }
00677           
00678           
00679           
00680           // Isolation cuts
00681           //double isovar = mu.isolationR03().sumPt;
00682           double ecalisovar = elec.dr03EcalRecHitSumEt();  // picked one set! 
00683           double hcalisovar = elec.dr03HcalTowerSumEt();   // try others if 
00684           double trkisovar = elec.dr04TkSumPt();           // doesn't work 
00685           //if (isCombinedIso_) {
00686           //isovar += mu.isolationR03().emEt;
00687           //isovar += mu.isolationR03().hadEt;
00688           //}
00689           //if (isRelativeIso_) isovar /= pt;
00690           if (ecalisovar<ecalIsoCutBarrel_ && isBarrel) electron_sel[4] = true; 
00691           if (ecalisovar<ecalIsoCutEndcap_ && isEndcap) electron_sel[4] = true; 
00692           if (hcalisovar<hcalIsoCutBarrel_ && isBarrel) electron_sel[5] = true; 
00693           if (hcalisovar<hcalIsoCutEndcap_ && isEndcap) electron_sel[5] = true; 
00694           if (trkisovar<trkIsoCutBarrel_ && isBarrel) electron_sel[6] = true;   
00695           if (trkisovar<trkIsoCutEndcap_ && isEndcap) electron_sel[6] = true;   
00696           if (isBarrel)
00697             {
00698               LogTrace("") << "\t... ecal isolation value " << ecalisovar << " (barrel), pass? " << electron_sel[4]; 
00699               LogTrace("") << "\t... hcal isolation value " << hcalisovar << " (barrel), pass? " << electron_sel[5];
00700               LogTrace("") << "\t... trk isolation value " << trkisovar << " (barrel), pass? " << electron_sel[6];
00701             }
00702           else if (isEndcap)
00703             {
00704               LogTrace("") << "\t... ecal isolation value " << ecalisovar << " (endcap), pass? " << electron_sel[4]; 
00705               LogTrace("") << "\t... hcal isolation value " << hcalisovar << " (endcap), pass? " << electron_sel[5];
00706               LogTrace("") << "\t... trk isolation value " << trkisovar << " (endcap), pass? " << electron_sel[6];
00707             }
00708           
00709           //iso_before_->Fill(isovar);
00710           if (isBarrel) 
00711             {
00712               ecalisobarrel_before_->Fill(ecalisovar);
00713               hcalisobarrel_before_->Fill(hcalisovar);
00714               trkisobarrel_before_->Fill(trkisovar);
00715             }
00716           else if (isEndcap) 
00717             {
00718               ecalisoendcap_before_->Fill(ecalisovar);
00719               hcalisoendcap_before_->Fill(hcalisovar);
00720               trkisoendcap_before_->Fill(trkisovar);
00721             }
00722           
00723           
00724           // HLT 
00725           if (trigger_fired) electron_sel[7] = true; 
00726           
00727           
00728           //             // MET/MT cuts
00729           double w_et = met_et+pt;
00730           double w_px = met_px+px;
00731           double w_py = met_py+py;
00732           
00733           double massT = w_et*w_et - w_px*w_px - w_py*w_py;
00734           massT = (massT>0) ? sqrt(massT) : 0;
00735           
00736           LogTrace("") << "\t... W mass, W_et, W_px, W_py: " << massT << ", " << w_et << ", " << w_px << ", " << w_py << " [GeV]";
00737           if (massT>mtMin_ && massT<mtMax_) electron_sel[8] = true; 
00738           mt_before_->Fill(massT);
00739           if (met_et>metMin_ && met_et<metMax_) electron_sel[9] = true; 
00740           
00741           //             // Acoplanarity cuts
00742           //             Geom::Phi<double> deltaphi(mu.phi()-atan2(met_py,met_px));
00743           //             double acop = deltaphi.value();
00744           //             if (acop<0) acop = - acop;
00745           //             acop = M_PI - acop;
00746           //             LogTrace("") << "\t... acoplanarity: " << acop;
00747           //             if (acop<acopCut_) muon_sel[10] = true; 
00748           //             acop_before_->Fill(acop);
00749           
00750           //             // Remaining flags (from global event information)
00751           //             if (nmuonsForZ1<1 || nmuonsForZ2<2) muon_sel[11] = true; 
00752           if (njets<=nJetMax_) electron_sel[10] = true; 
00753           
00754           // Collect necessary flags "per electron"
00755           int flags_passed = 0;
00756           bool rec_sel_this = true;
00757           bool eid_sel_this = true;
00758           bool iso_sel_this = true;
00759           bool all_sel_this = true;
00760           for (int j=0; j<NFLAGS; ++j) 
00761             {
00762               if (electron_sel[j]) flags_passed += 1;
00763               if (j<2  && !electron_sel[j]) rec_sel_this = false;
00764               if (j<4  && !electron_sel[j]) eid_sel_this = false;
00765               if (j<7  && !electron_sel[j]) iso_sel_this = false;
00766               if (!electron_sel[j]) all_sel_this = false;
00767             }
00768           
00769           if (all_sel_this)
00770             {
00771               if (nGoodElectrons < 2)
00772                 {
00773                   goodElectron[nGoodElectrons][0] = 1.;
00774                   goodElectron[nGoodElectrons][1] = elec.massSqr();
00775                   goodElectron[nGoodElectrons][2] = elec.energy();
00776                   goodElectron[nGoodElectrons][3] = elec.px();
00777                   goodElectron[nGoodElectrons][4] = elec.py();
00778                   goodElectron[nGoodElectrons][5] = elec.pz();
00779                 }
00780               nGoodElectrons++;
00781             }
00782 
00783           //             // "rec" => pt,eta and quality cuts are satisfied
00784           //             if (rec_sel_this) rec_sel = true;
00785           //             // "iso" => "rec" AND "muon is isolated"
00786           //             if (iso_sel_this) iso_sel = true;
00787           //             // "hlt" => "iso" AND "event is triggered"
00788           //             if (hlt_sel_this) hlt_sel = true;
00789           //             // "all" => "met" AND "Z/top rejection cuts"
00790           //             if (all_sel_this) all_sel = true;
00791           
00792           // "rec" => pt,eta cuts are satisfied
00793           if (rec_sel_this) rec_sel = true;
00794           // "eid" => "rec" AND "electron passes ID"
00795           if (eid_sel_this) iso_sel = true;
00796           // "iso" => "eid" AND "electron is isolated"
00797           if (iso_sel_this) iso_sel = true;
00798           // "met" => "iso" AND "MET/MT"
00799           // "all" => "met" AND "event is triggered"
00800           if (all_sel_this) all_sel = true;
00801 
00802           // Do N-1 histograms now (and only once for global event quantities)
00803           if (flags_passed >= (NFLAGS-1)) 
00804             {
00805               if (!electron_sel[0] || flags_passed==NFLAGS) 
00806                 {
00807                   pt_after_->Fill(pt);
00808                 }
00809               if (!electron_sel[1] || flags_passed==NFLAGS) 
00810                 {
00811                   eta_after_->Fill(eta);
00812                 }
00813               if (!electron_sel[2] || flags_passed==NFLAGS) 
00814                 {
00815                   if (isBarrel)
00816                     {
00817                       sieiebarrel_after_->Fill(sieie);
00818                     }
00819                   else if (isEndcap)
00820                     {
00821                       sieieendcap_after_->Fill(sieie);
00822                     }
00823                 }
00824               if (!electron_sel[3] || flags_passed==NFLAGS) 
00825                 {
00826                   if (isBarrel)
00827                     {
00828                       detainbarrel_after_->Fill(detain);
00829                     }
00830                   else if (isEndcap)
00831                     {
00832                       detainendcap_after_->Fill(detain);
00833                     }
00834                 }
00835               if (!electron_sel[4] || flags_passed==NFLAGS) 
00836                 {
00837                   if (isBarrel)
00838                     {
00839                       ecalisobarrel_after_->Fill(ecalisovar);
00840                     }
00841                   else if (isEndcap)
00842                     {
00843                       ecalisoendcap_after_->Fill(ecalisovar);
00844                     }
00845                 }
00846               if (!electron_sel[5] || flags_passed==NFLAGS) 
00847                 {
00848                   if (isBarrel)
00849                     {
00850                       hcalisobarrel_after_->Fill(hcalisovar);
00851                     }
00852                   else if (isEndcap)
00853                     {
00854                       hcalisoendcap_after_->Fill(hcalisovar);
00855                     }
00856                 }
00857               if (!electron_sel[6] || flags_passed==NFLAGS) 
00858                 {
00859                   if (isBarrel)
00860                     {
00861                       trkisobarrel_after_->Fill(trkisovar);
00862                     }
00863                   else if (isEndcap)
00864                     {
00865                       trkisoendcap_after_->Fill(trkisovar);
00866                     }
00867                 }
00868               //                if (!electron_sel[3] || flags_passed==NFLAGS) 
00869               //                  {
00870               //                    detain_after_->Fill(detain);
00871               //                  }
00872               //                if (!electron_sel[4] || flags_passed==NFLAGS) 
00873               //                  {
00874               //                    ecaliso_after_->Fill(trackerHits);
00875               //                  }
00876               //                if (!electron_sel[5] || flags_passed==NFLAGS) 
00877               //                  {
00878               //                    tkelectr_after_->Fill(electr.isTrackerElectron());
00879               //                  }
00880               //                if (!electron_sel[6] || flags_passed==NFLAGS) 
00881               //                  {
00882               //                    iso_after_->Fill(isovar);
00883               //                  }
00884               if (!electron_sel[7] || flags_passed==NFLAGS) 
00885                 {
00886                   if (!hlt_hist_done) 
00887                     {
00888                       trig_after_->Fill(trigger_fired);
00889                     }
00890                 }
00891               hlt_hist_done = true;
00892               if (!electron_sel[8] || flags_passed==NFLAGS) 
00893                 {
00894                   mt_after_->Fill(massT);
00895                 }
00896               if (!electron_sel[9] || flags_passed==NFLAGS) 
00897                 {
00898                   if (!met_hist_done) {
00899                     met_after_->Fill(met_et);
00900                   }
00901                 }
00902               met_hist_done = true;
00903               //                   if (!muon_sel[10] || flags_passed==NFLAGS) 
00904               //                         acop_after_->Fill(acop);
00905               //                   if (!muon_sel[11] || flags_passed==NFLAGS) 
00906               //                         if (!nz1_hist_done) nz1_after_->Fill(nmuonsForZ1);
00907               //                         nz1_hist_done = true;
00908               //                   if (!muon_sel[11] || flags_passed==NFLAGS) 
00909               //                         if (!nz2_hist_done) nz2_after_->Fill(nmuonsForZ2);
00910               //                         nz2_hist_done = true;
00911               if (!electron_sel[10] || flags_passed==NFLAGS) {
00912                 if (!njets_hist_done) 
00913                   {
00914                     njets_after_->Fill(njets);
00915                     if (jet_et>10) //don't want low energy "jets"
00916                       {
00917                         jet_et_after_   ->Fill(jet_et);
00918                         jet_eta_after_  ->Fill(jet_eta);
00919                       }
00920                   }
00921               }
00922               njets_hist_done = true;
00923 
00924             } // end N-1 histos block
00925           
00926         } // end loop through electrons
00927 
00928       // inv mass = sqrt(m_1^2 + m_2^2 + 2*(E_1*E_2 - (px1*px2 + py1*py2 + pz1+pz2) ) )
00929       double invMass;
00930 
00931       nelectrons_before_->Fill(electronCollectionSize);
00932       if (electronCollectionSize > 1)
00933         {
00934           invMass = sqrt(electron[0][1] + electron[1][1] + 2*(electron[0][2]*electron[1][2] - (electron[0][3]*electron[1][3] + electron[0][4]*electron[1][4] + electron[0][5]*electron[1][5]) ) );
00935           invmass_before_->Fill(invMass);
00936         }
00937 
00938       nelectrons_after_->Fill(nGoodElectrons);
00939       if (nGoodElectrons > 1)
00940         {
00941           invMass = sqrt(goodElectron[0][1] + goodElectron[1][1] + 2*(goodElectron[0][2]*goodElectron[1][2] - (goodElectron[0][3]*goodElectron[1][3] + goodElectron[0][4]*goodElectron[1][4] + goodElectron[0][5]*goodElectron[1][5]) ) );
00942           invmass_after_->Fill(invMass);
00943         }
00944 
00945       // Collect final flags
00946       if (rec_sel) nrec++;
00947       if (eid_sel) neid++;
00948       if (iso_sel) niso++;
00949       //      if (hlt_sel) nhlt++;
00950       //      if (met_sel) nmet++;
00951 
00952       if (all_sel) {
00953             nsel++;
00954             LogTrace("") << ">>>> Event ACCEPTED";
00955       } else {
00956             LogTrace("") << ">>>> Event REJECTED";
00957       }
00958 
00959 
00960 
00961       return;
00962 
00963 }
00964 
00965 // This always returns only a positive deltaPhi
00966 double EwkElecDQM::calcDeltaPhi(double phi1, double phi2) {
00967 
00968   double deltaPhi = phi1 - phi2;
00969 
00970   if (deltaPhi < 0) deltaPhi = -deltaPhi;
00971 
00972   if (deltaPhi > 3.1415926) {
00973     deltaPhi = 2 * 3.1415926 - deltaPhi;
00974   }
00975 
00976   return deltaPhi;
00977 }
00978