CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Protected Member Functions | Private Attributes
EwkElecDQM Class Reference

#include <EwkElecDQM.h>

Inheritance diagram for EwkElecDQM:
DQMOneEDAnalyzer<> edm::one::EDProducer< edm::EndRunProducer, edm::one::WatchRuns, edm::Accumulator, Args...> edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 EwkElecDQM (const edm::ParameterSet &)
 
- Public Member Functions inherited from DQMOneEDAnalyzer<>
void accumulate (edm::Event const &event, edm::EventSetup const &setup) override
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
 DQMOneEDAnalyzer ()
 
void endRun (edm::Run const &, edm::EventSetup const &) final
 
void endRunProduce (edm::Run &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::one::EDProducer< edm::EndRunProducer, edm::one::WatchRuns, edm::Accumulator, Args...>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector
< edm::ProductResolverIndex >
const & 
indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector
< edm::ProductResolverIndex >
const & 
putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex >
const & 
esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector
< ProductResolverIndexAndSkipBit >
const & 
itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const * > *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Protected Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
double calcDeltaPhi (double phi1, double phi2)
 
void dqmBeginRun (const edm::Run &, const edm::EventSetup &) override
 
void dqmEndRun (const edm::Run &, const edm::EventSetup &) override
 
- Protected Member Functions inherited from edm::ProducerBase
ProducesCollector producesCollector ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< B > consumes (edm::InputTag tag) noexcept
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes () noexcept
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag) noexcept
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
void resetItemsToGetFrom (BranchType iType)
 

Private Attributes

edm::EDGetTokenT< reco::BeamSpotbeamSpotTag_
 
MonitorElementdetainbarrel_after_
 
MonitorElementdetainbarrel_before_
 
double detainCutBarrel_
 
double detainCutEndcap_
 
MonitorElementdetainendcap_after_
 
MonitorElementdetainendcap_before_
 
MonitorElementecalisobarrel_after_
 
MonitorElementecalisobarrel_before_
 
double ecalIsoCutBarrel_
 
double ecalIsoCutEndcap_
 
MonitorElementecalisoendcap_after_
 
MonitorElementecalisoendcap_before_
 
double eJetMin_
 
edm::EDGetTokenT< edm::View
< reco::GsfElectron > > 
elecTag_
 
const std::vector< std::string > elecTrig_
 
MonitorElementeta_after_
 
MonitorElementeta_before_
 
double etaCut_
 
MonitorElementhcalisobarrel_after_
 
MonitorElementhcalisobarrel_before_
 
double hcalIsoCutBarrel_
 
double hcalIsoCutEndcap_
 
MonitorElementhcalisoendcap_after_
 
MonitorElementhcalisoendcap_before_
 
HLTPrescaleProvider hltPrescaleProvider_
 
MonitorElementinvmass_after_
 
MonitorElementinvmass_before_
 
MonitorElementinvmassPU_afterZ_
 
MonitorElementinvmassPU_before_
 
bool isValidHltConfig_
 
MonitorElementjet_et_after_
 
MonitorElementjet_et_before_
 
MonitorElementjet_eta_after_
 
MonitorElementjet_eta_before_
 
edm::InputTag jetTag_
 
edm::EDGetTokenT< edm::View
< reco::Jet > > 
jetToken_
 
MonitorElementmet_after_
 
MonitorElementmet_before_
 
bool metIncludesMuons_
 
double metMax_
 
double metMin_
 
edm::InputTag metTag_
 
edm::EDGetTokenT< edm::View
< reco::MET > > 
metToken_
 
MonitorElementmt_after_
 
MonitorElementmt_before_
 
double mtMax_
 
double mtMin_
 
unsigned int nall
 
unsigned int neid
 
MonitorElementnelectrons_after_
 
MonitorElementnelectrons_before_
 
unsigned int nGoodElectrons
 
unsigned int niso
 
int nJetMax_
 
MonitorElementnjets_after_
 
MonitorElementnjets_before_
 
MonitorElementnpvs_afterZ_
 
MonitorElementnpvs_before_
 
unsigned int nrec
 
unsigned int nsel
 
MonitorElementpt_after_
 
MonitorElementpt_before_
 
double ptCut_
 
unsigned int PUBinCount_
 
unsigned int PUMax_
 
MonitorElementsieiebarrel_after_
 
MonitorElementsieiebarrel_before_
 
double sieieCutBarrel_
 
double sieieCutEndcap_
 
MonitorElementsieieendcap_after_
 
MonitorElementsieieendcap_before_
 
MonitorElementtrig_after_
 
MonitorElementtrig_before_
 
edm::EDGetTokenT
< edm::TriggerResults
trigTag_
 
MonitorElementtrkisobarrel_after_
 
MonitorElementtrkisobarrel_before_
 
double trkIsoCutBarrel_
 
double trkIsoCutEndcap_
 
MonitorElementtrkisoendcap_after_
 
MonitorElementtrkisoendcap_before_
 
edm::EDGetTokenT< edm::View
< reco::Vertex > > 
vertexTag_
 

Additional Inherited Members

- Public Types inherited from DQMOneEDAnalyzer<>
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Attributes inherited from DQMOneEDAnalyzer<>
edm::EDPutTokenT< DQMTokenrunToken_
 

Detailed Description

DQM offline for EWK Electrons

Definition at line 27 of file EwkElecDQM.h.

Constructor & Destructor Documentation

EwkElecDQM::EwkElecDQM ( const edm::ParameterSet cfg)

Definition at line 36 of file EwkElecDQM.cc.

References isValidHltConfig_.

37  : // Input collections
39  jetTag_(cfg.getUntrackedParameter<edm::InputTag>("JetTag", edm::InputTag("sisCone5CaloJets"))),
40  // trigTag_(consumes<edm::TriggerResults>(
41  // cfg.getUntrackedParameter<edm::InputTag>(
42  // "TrigTag", edm::InputTag("TriggerResults::HLT")))),
44  cfg.getUntrackedParameter<edm::InputTag>("ElecTag", edm::InputTag("gsfElectrons")))),
45  metToken_(
48  cfg.getUntrackedParameter<edm::InputTag>("JetTag", edm::InputTag("sisCone5CaloJets")))),
50  cfg.getUntrackedParameter<edm::InputTag>("VertexTag", edm::InputTag("offlinePrimaryVertices")))),
52  consumes<reco::BeamSpot>(cfg.getUntrackedParameter<edm::InputTag>("BeamSpotTag", edm::InputTag("BeamSpot")))),
53 
54  // Main cuts
55  // muonTrig_(cfg.getUntrackedParameter<std::string> ("MuonTrig",
56  // "HLT_Mu9")),
57  // elecTrig_(cfg.getUntrackedParameter<std::vector< std::string >
58  // >("ElecTrig", "HLT_Ele10_SW_L1R")),
59  elecTrig_(cfg.getUntrackedParameter<std::vector<std::string> >("ElecTrig")),
60  // ptCut_(cfg.getUntrackedParameter<double>("PtCut", 25.)),
61  ptCut_(cfg.getUntrackedParameter<double>("PtCut", 10.)),
62  // etaCut_(cfg.getUntrackedParameter<double>("EtaCut", 2.1)),
63  etaCut_(cfg.getUntrackedParameter<double>("EtaCut", 2.4)),
64  sieieCutBarrel_(cfg.getUntrackedParameter<double>("SieieBarrel", 0.01)),
65  sieieCutEndcap_(cfg.getUntrackedParameter<double>("SieieEndcap", 0.028)),
66  detainCutBarrel_(cfg.getUntrackedParameter<double>("DetainBarrel", 0.0071)),
67  detainCutEndcap_(cfg.getUntrackedParameter<double>("DetainEndcap", 0.0066)),
68  // isRelativeIso_(cfg.getUntrackedParameter<bool>("IsRelativeIso",
69  // true)),
70  // isCombinedIso_(cfg.getUntrackedParameter<bool>("IsCombinedIso",
71  // false)),
72  // isoCut03_(cfg.getUntrackedParameter<double>("IsoCut03", 0.1)),
73  ecalIsoCutBarrel_(cfg.getUntrackedParameter<double>("EcalIsoCutBarrel", 5.7)),
74  ecalIsoCutEndcap_(cfg.getUntrackedParameter<double>("EcalIsoCutEndcap", 5.0)),
75  hcalIsoCutBarrel_(cfg.getUntrackedParameter<double>("HcalIsoCutBarrel", 8.1)),
76  hcalIsoCutEndcap_(cfg.getUntrackedParameter<double>("HcalIsoCutEndcap", 3.4)),
77  trkIsoCutBarrel_(cfg.getUntrackedParameter<double>("TrkIsoCutBarrel", 7.2)),
78  trkIsoCutEndcap_(cfg.getUntrackedParameter<double>("TrkIsoCutEndcap", 5.1)),
79  mtMin_(cfg.getUntrackedParameter<double>("MtMin", -999999)),
80  mtMax_(cfg.getUntrackedParameter<double>("MtMax", 999999.)),
81  metMin_(cfg.getUntrackedParameter<double>("MetMin", -999999.)),
82  metMax_(cfg.getUntrackedParameter<double>("MetMax", 999999.)),
83  // acopCut_(cfg.getUntrackedParameter<double>("AcopCut", 2.)),
84 
85  // Muon quality cuts
86  // dxyCut_(cfg.getUntrackedParameter<double>("DxyCut", 0.2)),
87  // normalizedChi2Cut_(cfg.getUntrackedParameter<double>("NormalizedChi2Cut",
88  // 10.)),
89  // trackerHitsCut_(cfg.getUntrackedParameter<int>("TrackerHitsCut",
90  // 11)),
91  // isAlsoTrackerMuon_(cfg.getUntrackedParameter<bool>("IsAlsoTrackerMuon",
92  // true)),
93 
94  // Z rejection
95  // ptThrForZ1_(cfg.getUntrackedParameter<double>("PtThrForZ1", 20.)),
96  // ptThrForZ2_(cfg.getUntrackedParameter<double>("PtThrForZ2", 10.)),
97 
98  // Top rejection
99  eJetMin_(cfg.getUntrackedParameter<double>("EJetMin", 999999.)),
100  nJetMax_(cfg.getUntrackedParameter<int>("NJetMax", 999999)),
101  PUMax_(cfg.getUntrackedParameter<unsigned int>("PUMax", 60)),
102  PUBinCount_(cfg.getUntrackedParameter<unsigned int>("PUBinCount", 12)),
104 // caloJetCollection_(cfg.getUntrackedParameter<edm:InputTag>("CaloJetCollection","sisCone5CaloJets"))
105 {
106  isValidHltConfig_ = false;
107 }
T getUntrackedParameter(std::string const &, T const &) const
double metMax_
Definition: EwkElecDQM.h:75
double sieieCutBarrel_
Definition: EwkElecDQM.h:57
edm::EDGetTokenT< reco::BeamSpot > beamSpotTag_
Definition: EwkElecDQM.h:49
edm::InputTag jetTag_
Definition: EwkElecDQM.h:43
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
double ecalIsoCutEndcap_
Definition: EwkElecDQM.h:67
double detainCutEndcap_
Definition: EwkElecDQM.h:60
edm::EDGetTokenT< edm::View< reco::GsfElectron > > elecTag_
Definition: EwkElecDQM.h:45
double metMin_
Definition: EwkElecDQM.h:74
double trkIsoCutEndcap_
Definition: EwkElecDQM.h:71
unsigned int PUBinCount_
Definition: EwkElecDQM.h:90
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
int nJetMax_
Definition: EwkElecDQM.h:87
double ecalIsoCutBarrel_
Definition: EwkElecDQM.h:66
double trkIsoCutBarrel_
Definition: EwkElecDQM.h:70
edm::InputTag metTag_
Definition: EwkElecDQM.h:42
edm::EDGetTokenT< edm::View< reco::Jet > > jetToken_
Definition: EwkElecDQM.h:47
unsigned int PUMax_
Definition: EwkElecDQM.h:90
double etaCut_
Definition: EwkElecDQM.h:55
double hcalIsoCutBarrel_
Definition: EwkElecDQM.h:68
bool isValidHltConfig_
Definition: EwkElecDQM.h:92
double sieieCutEndcap_
Definition: EwkElecDQM.h:58
edm::EDGetTokenT< edm::View< reco::Vertex > > vertexTag_
Definition: EwkElecDQM.h:48
edm::EDGetTokenT< edm::View< reco::MET > > metToken_
Definition: EwkElecDQM.h:46
double eJetMin_
Definition: EwkElecDQM.h:86
double mtMax_
Definition: EwkElecDQM.h:73
double detainCutBarrel_
Definition: EwkElecDQM.h:59
const std::vector< std::string > elecTrig_
Definition: EwkElecDQM.h:53
HLTPrescaleProvider hltPrescaleProvider_
Definition: EwkElecDQM.h:93
double hcalIsoCutEndcap_
Definition: EwkElecDQM.h:69
double mtMin_
Definition: EwkElecDQM.h:72
double ptCut_
Definition: EwkElecDQM.h:54

Member Function Documentation

void EwkElecDQM::analyze ( const edm::Event ev,
const edm::EventSetup iSet 
)
overrideprotectedvirtual

Reimplemented from DQMOneEDAnalyzer<>.

Definition at line 347 of file EwkElecDQM.cc.

References beamSpotTag_, calcDeltaPhi(), reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), detainbarrel_after_, detainbarrel_before_, detainCutBarrel_, detainCutEndcap_, detainendcap_after_, detainendcap_before_, reco::GsfElectron::dr03EcalRecHitSumEt(), reco::GsfElectron::dr03HcalTowerSumEt(), reco::GsfElectron::dr04TkSumPt(), ecalisobarrel_after_, ecalisobarrel_before_, ecalIsoCutBarrel_, ecalIsoCutEndcap_, ecalisoendcap_after_, ecalisoendcap_before_, eJetMin_, elecTag_, metsig::electron, reco::LeafCandidate::energy(), reco::LeafCandidate::et(), PVValHelper::eta, reco::LeafCandidate::eta(), eta_after_, eta_before_, etaCut_, dqm::impl::MonitorElement::Fill(), edm::Event::getByToken(), hcalisobarrel_after_, hcalisobarrel_before_, hcalIsoCutBarrel_, hcalIsoCutEndcap_, hcalisoendcap_after_, hcalisoendcap_before_, hltPrescaleProvider_, mps_fire::i, invmass_after_, invmass_before_, invmassPU_afterZ_, invmassPU_before_, PixelPluginsPhase0_cfi::isBarrel, GeomDetEnumerators::isEndcap(), reco::Vertex::isValid(), dqmiolumiharvest::j, metsig::jet, jet_et_after_, jet_et_before_, jet_eta_after_, jet_eta_before_, jetToken_, LogTrace, reco::LeafCandidate::massSqr(), objects.METAnalyzer::met, met_after_, met_before_, HLT_FULL_cff::metCollection, metMax_, metMin_, metToken_, mt_after_, mt_before_, mtMax_, mtMin_, nall, neid, nelectrons_after_, nelectrons_before_, nGoodElectrons, niso, nJetMax_, HLT_FULL_cff::njets, njets_after_, njets_before_, npvs_afterZ_, npvs_before_, nrec, nsel, reco::LeafCandidate::phi(), HLTPrescaleProvider::prescaleSet(), DiDispStaMuonMonitor_cfi::pt, reco::LeafCandidate::pt(), pt_after_, pt_before_, ptCut_, reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), sieiebarrel_after_, sieiebarrel_before_, sieieCutBarrel_, sieieCutEndcap_, sieieendcap_after_, sieieendcap_before_, reco::GsfElectron::sigmaIetaIeta(), mathSSE::sqrt(), trkisobarrel_after_, trkisobarrel_before_, trkIsoCutBarrel_, trkIsoCutEndcap_, trkisoendcap_after_, trkisoendcap_before_, GoodVertex_cfg::vertexCollection, and vertexTag_.

347  {
348  // Reset global event selection flags
349  bool rec_sel = false;
350  bool eid_sel = false;
351  bool iso_sel = false;
352  bool all_sel = false;
353 
354  // Electron collection
355  Handle<View<GsfElectron> > electronCollection;
356  if (!ev.getByToken(elecTag_, electronCollection)) {
357  // LogWarning("") << ">>> Electron collection does not exist !!!";
358  return;
359  }
360  unsigned int electronCollectionSize = electronCollection->size();
361 
362  // Beam spot
363  Handle<reco::BeamSpot> beamSpotHandle;
364  if (!ev.getByToken(beamSpotTag_, beamSpotHandle)) {
365  // LogWarning("") << ">>> No beam spot found !!!";
366  return;
367  }
368 
369  // MET
370  double met_px = 0.;
371  double met_py = 0.;
373  if (!ev.getByToken(metToken_, metCollection)) {
374  // LogWarning("") << ">>> MET collection does not exist !!!";
375  return;
376  }
377 
378  const MET& met = metCollection->at(0);
379  met_px = met.px();
380  met_py = met.py();
381  // if (!metIncludesMuons_) {
382  // for (unsigned int i=0; i<muonCollectionSize; i++) {
383  // const Muon& mu = muonCollection->at(i);
384  // if (!mu.isGlobalMuon()) continue;
385  // met_px -= mu.px();
386  // met_py -= mu.py();
387  // }
388  // }
389  double met_et = sqrt(met_px * met_px + met_py * met_py);
390  LogTrace("") << ">>> MET, MET_px, MET_py: " << met_et << ", " << met_px << ", " << met_py << " [GeV]";
391  met_before_->Fill(met_et);
392 
393  // Vertices in the event
394  int npvCount = 0;
396  if (!ev.getByToken(vertexTag_, vertexCollection)) {
397  LogError("") << ">>> Vertex collection does not exist !!!";
398  return;
399  }
400 
401  for (unsigned int i = 0; i < vertexCollection->size(); i++) {
402  const Vertex& vertex = vertexCollection->at(i);
403  if (vertex.isValid())
404  npvCount++;
405  }
406  npvs_before_->Fill(npvCount);
407 
408  // Trigger
409  // Handle<TriggerResults> triggerResults;
410  // if (!ev.getByToken(trigTag_, triggerResults)) {
411  // LogWarning("") << ">>> TRIGGER collection does not exist !!!";
412  return;
413  // }
414  // const edm::TriggerNames& trigNames = ev.triggerNames(*triggerResults);
415  // bool trigger_fired = false;
416 
417  // HLTConfigProvider const& hltConfigProvider = hltPrescaleProvider_.hltConfigProvider();
418 
419  /* very old code
420  for (unsigned int i=0; i<triggerResults->size(); i++) {
421  if (triggerResults->accept(i)) {
422  LogTrace("") << "Accept by: " << i << ", Trigger: " <<
423  trigNames.triggerName(i);
424  }
425  }
426  */
427 
428  // the following gives error on CRAFT08 data where itrig1=19 (vector index out
429  // of range)
430  /*
431  int itrig1 = trigNames.triggerIndex(muonTrig_);
432  if (triggerResults->accept(itrig1)) trigger_fired = true;
433  */
434  // suggested replacement: lm250909
435  /* Fix buggy trigger logic
436  for (unsigned int i=0; i<triggerResults->size(); i++)
437  {
438  std::string trigName = trigNames.triggerName(i);
439  bool found=false;
440 
441 // for (unsigned int j = 0; j < elecTrig_.size(); j++)
442 // {
443 // if ( trigName == elecTrig_.at(j) && triggerResults->accept(i))
444 // {
445 // trigger_fired = true;
446 // }
447 // }
448 
449 
450  for(unsigned int index=0; index<elecTrig_.size() && found==false; index++)
451 {
452  size_t trigPath = trigName.find(elecTrig_.at(index)); // 0 if found, pos
453 if not
454  if (trigPath==0) found=true;
455  }
456  if(!found) continue;
457 
458  bool prescaled=false;
459  for (unsigned int ps= 0; ps< hltConfigProvider.prescaleSize();
460 ps++){
461  const unsigned int prescaleValue =
462 hltConfigProvider.prescaleValue(ps, trigName) ;
463  if (prescaleValue != 1) prescaled =true;
464  }
465 
466  if(triggerResults->accept(i) && !prescaled) trigger_fired=true;
467 
468  }
469  */
470 
471  // get the prescale set for this event
472  const int prescaleSet = hltPrescaleProvider_.prescaleSet(ev, iSet);
473  if (prescaleSet == -1) {
474  LogTrace("") << "Failed to determine prescaleSet\n";
475  // std::cout << "Failed to determine prescaleSet. Check cmsRun GlobalTag\n";
476  return;
477  }
478 
479  // for (unsigned int i = 0;
480  // (i < triggerResults->size()) && (trigger_fired == false); i++) {
481  // skip trigger, if it did not fire
482  // if (!triggerResults->accept(i)) continue;
483 
484  // skip trigger, if it is not on our list
485  // bool found = false;
486  // const std::string trigName = trigNames.triggerName(i);
487  // for (unsigned int index = 0; index < elecTrig_.size() && found == false;
488  // index++) {
489  // if (trigName.find(elecTrig_.at(index)) == 0) found = true;
490  // }
491  // if (!found) continue;
492 
493  // skip trigger, if it is prescaled
494  // if (hltConfigProvider.prescaleValue(prescaleSet, trigName) != 1) continue;
495 
496  // std::cout << "found unprescaled trigger that fired: " << trigName <<
497  // "\n";
498  // trigger_fired = true;
499  // }
500 
501  /* LogTrace("") << ">>> Trigger bit: " << trigger_fired << " for one of ( ";
502  for (unsigned int k = 0; k < elecTrig_.size(); k++) {
503  LogTrace("") << elecTrig_.at(k) << " ";
504  }
505  LogTrace("") << ")";
506  trig_before_->Fill(trigger_fired);
507 */
508  // Jet collection
509  Handle<View<Jet> > jetCollection;
510  if (!ev.getByToken(jetToken_, jetCollection)) {
511  // LogError("") << ">>> JET collection does not exist !!!";
512  return;
513  }
514  float electron_et = -8.0;
515  float electron_eta = -8.0;
516  float electron_phi = -8.0;
517  float electron2_et = -9.0;
518  float electron2_eta = -9.0;
519  float electron2_phi = -9.0;
520  // need to get some electron info so jets can be cleaned of them
521  for (unsigned int i = 0; i < electronCollectionSize; i++) {
522  const GsfElectron& elec = electronCollection->at(i);
523 
524  if (i < 1) {
525  electron_et = elec.pt();
526  electron_eta = elec.eta();
527  electron_phi = elec.phi();
528  }
529  if (i == 2) {
530  electron2_et = elec.pt();
531  electron2_eta = elec.eta();
532  electron2_phi = elec.phi();
533  }
534  }
535 
536  float jet_et = -8.0;
537  float jet_eta = -8.0;
538  int jet_count = 0;
539  float jet2_et = -9.0;
540  unsigned int jetCollectionSize = jetCollection->size();
541  int njets = 0;
542  for (unsigned int i = 0; i < jetCollectionSize; i++) {
543  const Jet& jet = jetCollection->at(i);
544 
545  float jet_current_et = jet.et();
546  // cout << "jet_current_et " << jet_current_et << endl;
547  // if it overlaps with electron, it is not a jet
548  if (electron_et > 0.0 && fabs(jet.eta() - electron_eta) < 0.2 && calcDeltaPhi(jet.phi(), electron_phi) < 0.2)
549  continue;
550  if (electron2_et > 0.0 && fabs(jet.eta() - electron2_eta) < 0.2 && calcDeltaPhi(jet.phi(), electron2_phi) < 0.2)
551  continue;
552 
553  // if it has too low Et, throw away
554  // if (jet_current_et < eJetMin_) continue; //Keep if only want to plot
555  // above jet cut
556 
557  if (jet.et() > eJetMin_) {
558  njets++;
559  jet_count++;
560  }
561  if (jet_current_et > jet_et) {
562  jet2_et = jet_et; // 2nd highest jet get's et from current highest
563  jet_et = jet.et(); // current highest jet gets et from the new highest
564  jet_eta = jet.eta();
565  } else if (jet_current_et > jet2_et) {
566  jet2_et = jet.et();
567  }
568  }
569 
570  // Fill After all electron cuts (or both before and after)
571  if (jet_et > 10) // don't want low energy "jets"
572  {
573  jet_et_before_->Fill(jet_et);
574  // jet2_et_before_ ->Fill(jet2_et);
575  jet_eta_before_->Fill(jet_eta);
576  }
577 
578  LogTrace("") << ">>> Total number of jets: " << jetCollectionSize;
579  LogTrace("") << ">>> Number of jets above " << eJetMin_ << " [GeV]: " << njets;
580  njets_before_->Fill(njets);
581 
582  // Start counting
583  nall++;
584 
585  // Histograms per event should be done only once, so keep track of them
586  //bool hlt_hist_done = false;
587  // bool minv_hist_done = false;
588  bool met_hist_done = false;
589  // bool nz1_hist_done = false;
590  // bool nz2_hist_done = false;
591  bool njets_hist_done = false;
592 
593  // Central selection criteria
594  // const int NFLAGS = 13; // number of individual selection criteria
595  const int NFLAGS = 10; // number of individual selection criteria
596  // 0: pt cut | rec
597  // 1: eta cut | rec
598  // 2: sieie | eid
599  // 3: detain | eid
600  // 4: ecal iso | iso
601  // 5: hcal iso | iso
602  // 6: trk iso | iso
603  // 7: trigger fired | hlt/all
604  bool electron_sel[NFLAGS];
605 
606  // for invariant mass calculation
607  // keep track of highest-pt electrons for initial (RECO) electrons
608  // and "good" electrons (passing all cuts)
609  // first array dimension is for first or second good electron
610  // second array dimension is for relevant quantities of good electron
611  // [0]: 1 for electron found or 0 for not found (if 0, no other quantities
612  // filled)
613  // [1]: mSqr
614  // [2]: E
615  // [3]: px
616  // [4]: py
617  // [5]: pz
618  // inv mass = sqrt(m_1^2 + m_2^2 + 2*(E_1*E_2 - (px1*px2 + py1*py2 + pz1+pz2)
619  // ) )
620  double electron[2][6];
621  double goodElectron[2][6];
622  nGoodElectrons = 0;
623  for (unsigned int i = 0; i < 2; i++) {
624  for (unsigned int j = 0; j < 6; j++) {
625  electron[i][j] = 0.;
626  goodElectron[i][j] = 0.;
627  }
628  }
629 
630  for (unsigned int i = 0; i < electronCollectionSize; i++) {
631  for (int j = 0; j < NFLAGS; ++j) {
632  electron_sel[j] = false;
633  }
634 
635  const GsfElectron& elec = electronCollection->at(i);
636  // if (!mu.isGlobalMuon()) continue;
637  // if (mu.globalTrack().isNull()) continue;
638  // if (mu.innerTrack().isNull()) continue;
639 
640  LogTrace("") << "> elec: processing electron number " << i << "...";
641  // reco::TrackRef gm = mu.globalTrack();
642  // reco::TrackRef tk = mu.innerTrack();
643  // should have stuff for electron track?
644 
645  if (i < 2) {
646  electron[i][0] = 1.;
647  electron[i][1] = elec.massSqr();
648  electron[i][2] = elec.energy();
649  electron[i][3] = elec.px();
650  electron[i][4] = elec.py();
651  electron[i][5] = elec.pz();
652  }
653 
654  // Pt,eta cuts
655  double pt = elec.pt();
656  double px = elec.px();
657  double py = elec.py();
658  double eta = elec.eta();
659  LogTrace("") << "\t... pt, eta: " << pt << " [GeV], " << eta;
660  ;
661  if (pt > ptCut_)
662  electron_sel[0] = true;
663  if (fabs(eta) < etaCut_)
664  electron_sel[1] = true;
665 
666  bool isBarrel = false;
667  bool isEndcap = false;
668  if (eta < 1.4442 && eta > -1.4442) {
669  isBarrel = true;
670  } else if ((eta > 1.56 && eta < 2.4) || (eta < -1.56 && eta > -2.4)) {
671  isEndcap = true;
672  }
673 
674  // // d0, chi2, nhits quality cuts
675  // double dxy = tk->dxy(beamSpotHandle->position());
676  // double normalizedChi2 = gm->normalizedChi2();
677  // double trackerHits = tk->numberOfValidHits();
678  // LogTrace("") << "\t... dxy, normalizedChi2, trackerHits,
679  // isTrackerMuon?: " << dxy << " [cm], " << normalizedChi2 << ", " <<
680  // trackerHits << ", " << mu.isTrackerMuon();
681  // if (fabs(dxy)<dxyCut_) muon_sel[2] = true;
682  // if (normalizedChi2<normalizedChi2Cut_) muon_sel[3] = true;
683  // if (trackerHits>=trackerHitsCut_) muon_sel[4] = true;
684  // if (mu.isTrackerMuon()) muon_sel[5] = true;
685 
686  pt_before_->Fill(pt);
687  eta_before_->Fill(eta);
688  // dxy_before_->Fill(dxy);
689  // chi2_before_->Fill(normalizedChi2);
690  // nhits_before_->Fill(trackerHits);
691  // tkmu_before_->Fill(mu.isTrackerMuon());
692 
693  // Electron ID cuts
694  double sieie = (double)elec.sigmaIetaIeta();
695  double detain = (double)elec.deltaEtaSuperClusterTrackAtVtx(); // think this is detain
696  if (sieie < sieieCutBarrel_ && isBarrel)
697  electron_sel[2] = true;
698  if (sieie < sieieCutEndcap_ && isEndcap)
699  electron_sel[2] = true;
700  if (detain < detainCutBarrel_ && isBarrel)
701  electron_sel[3] = true;
702  if (detain < detainCutEndcap_ && isEndcap)
703  electron_sel[3] = true;
704  if (isBarrel) {
705  LogTrace("") << "\t... sieie value " << sieie << " (barrel), pass? " << electron_sel[2];
706  LogTrace("") << "\t... detain value " << detain << " (barrel), pass? " << electron_sel[3];
707  } else if (isEndcap) {
708  LogTrace("") << "\t... sieie value " << sieie << " (endcap), pass? " << electron_sel[2];
709  LogTrace("") << "\t... detain value " << detain << " (endcap), pass? " << electron_sel[2];
710  }
711 
712  if (isBarrel) {
713  sieiebarrel_before_->Fill(sieie);
714  detainbarrel_before_->Fill(detain);
715  } else if (isEndcap) {
716  sieieendcap_before_->Fill(sieie);
717  detainendcap_before_->Fill(detain);
718  }
719 
720  // Isolation cuts
721  // double isovar = mu.isolationR03().sumPt;
722  double ecalisovar = elec.dr03EcalRecHitSumEt(); // picked one set!
723  double hcalisovar = elec.dr03HcalTowerSumEt(); // try others if
724  double trkisovar = elec.dr04TkSumPt(); // doesn't work
725  // if (isCombinedIso_) {
726  // isovar += mu.isolationR03().emEt;
727  // isovar += mu.isolationR03().hadEt;
728  //}
729  // if (isRelativeIso_) isovar /= pt;
730  if (ecalisovar < ecalIsoCutBarrel_ && isBarrel)
731  electron_sel[4] = true;
732  if (ecalisovar < ecalIsoCutEndcap_ && isEndcap)
733  electron_sel[4] = true;
734  if (hcalisovar < hcalIsoCutBarrel_ && isBarrel)
735  electron_sel[5] = true;
736  if (hcalisovar < hcalIsoCutEndcap_ && isEndcap)
737  electron_sel[5] = true;
738  if (trkisovar < trkIsoCutBarrel_ && isBarrel)
739  electron_sel[6] = true;
740  if (trkisovar < trkIsoCutEndcap_ && isEndcap)
741  electron_sel[6] = true;
742  if (isBarrel) {
743  LogTrace("") << "\t... ecal isolation value " << ecalisovar << " (barrel), pass? " << electron_sel[4];
744  LogTrace("") << "\t... hcal isolation value " << hcalisovar << " (barrel), pass? " << electron_sel[5];
745  LogTrace("") << "\t... trk isolation value " << trkisovar << " (barrel), pass? " << electron_sel[6];
746  } else if (isEndcap) {
747  LogTrace("") << "\t... ecal isolation value " << ecalisovar << " (endcap), pass? " << electron_sel[4];
748  LogTrace("") << "\t... hcal isolation value " << hcalisovar << " (endcap), pass? " << electron_sel[5];
749  LogTrace("") << "\t... trk isolation value " << trkisovar << " (endcap), pass? " << electron_sel[6];
750  }
751 
752  // iso_before_->Fill(isovar);
753  if (isBarrel) {
754  ecalisobarrel_before_->Fill(ecalisovar);
755  hcalisobarrel_before_->Fill(hcalisovar);
756  trkisobarrel_before_->Fill(trkisovar);
757  } else if (isEndcap) {
758  ecalisoendcap_before_->Fill(ecalisovar);
759  hcalisoendcap_before_->Fill(hcalisovar);
760  trkisoendcap_before_->Fill(trkisovar);
761  }
762 
763  // HLT
764  // if (trigger_fired) electron_sel[7] = true;
765 
766  // // MET/MT cuts
767  double w_et = met_et + pt;
768  double w_px = met_px + px;
769  double w_py = met_py + py;
770 
771  double massT = w_et * w_et - w_px * w_px - w_py * w_py;
772  massT = (massT > 0) ? sqrt(massT) : 0;
773 
774  LogTrace("") << "\t... W mass, W_et, W_px, W_py: " << massT << ", " << w_et << ", " << w_px << ", " << w_py
775  << " [GeV]";
776  if (massT > mtMin_ && massT < mtMax_)
777  electron_sel[7] = true;
778  mt_before_->Fill(massT);
779  if (met_et > metMin_ && met_et < metMax_)
780  electron_sel[8] = true;
781 
782  // // Acoplanarity cuts
783  // Geom::Phi<double> deltaphi(mu.phi()-atan2(met_py,met_px));
784  // double acop = deltaphi.value();
785  // if (acop<0) acop = - acop;
786  // acop = M_PI - acop;
787  // LogTrace("") << "\t... acoplanarity: " << acop;
788  // if (acop<acopCut_) muon_sel[10] = true;
789  // acop_before_->Fill(acop);
790 
791  // // Remaining flags (from global event information)
792  // if (nmuonsForZ1<1 || nmuonsForZ2<2) muon_sel[11] = true;
793  if (njets <= nJetMax_)
794  electron_sel[9] = true;
795 
796  // Collect necessary flags "per electron"
797  int flags_passed = 0;
798  bool rec_sel_this = true;
799  bool eid_sel_this = true;
800  bool iso_sel_this = true;
801  bool all_sel_this = true;
802  for (int j = 0; j < NFLAGS; ++j) {
803  if (electron_sel[j])
804  flags_passed += 1;
805  if (j < 2 && !electron_sel[j])
806  rec_sel_this = false;
807  if (j < 4 && !electron_sel[j])
808  eid_sel_this = false;
809  if (j < 7 && !electron_sel[j])
810  iso_sel_this = false;
811  if (!electron_sel[j])
812  all_sel_this = false;
813  }
814 
815  if (all_sel_this) {
816  if (nGoodElectrons < 2) {
817  goodElectron[nGoodElectrons][0] = 1.;
818  goodElectron[nGoodElectrons][1] = elec.massSqr();
819  goodElectron[nGoodElectrons][2] = elec.energy();
820  goodElectron[nGoodElectrons][3] = elec.px();
821  goodElectron[nGoodElectrons][4] = elec.py();
822  goodElectron[nGoodElectrons][5] = elec.pz();
823  }
824  nGoodElectrons++;
825  }
826 
827  // // "rec" => pt,eta and quality cuts are satisfied
828  // if (rec_sel_this) rec_sel = true;
829  // // "iso" => "rec" AND "muon is isolated"
830  // if (iso_sel_this) iso_sel = true;
831  // // "hlt" => "iso" AND "event is triggered"
832  // if (hlt_sel_this) hlt_sel = true;
833  // // "all" => "met" AND "Z/top rejection cuts"
834  // if (all_sel_this) all_sel = true;
835 
836  // "rec" => pt,eta cuts are satisfied
837  if (rec_sel_this)
838  rec_sel = true;
839  // "eid" => "rec" AND "electron passes ID"
840  if (eid_sel_this)
841  iso_sel = true;
842  // "iso" => "eid" AND "electron is isolated"
843  if (iso_sel_this)
844  iso_sel = true;
845  // "met" => "iso" AND "MET/MT"
846  // "all" => "met" AND "event is triggered"
847  if (all_sel_this)
848  all_sel = true;
849 
850  // Do N-1 histograms now (and only once for global event quantities)
851  if (flags_passed >= (NFLAGS - 1)) {
852  if (!electron_sel[0] || flags_passed == NFLAGS) {
853  pt_after_->Fill(pt);
854  }
855  if (!electron_sel[1] || flags_passed == NFLAGS) {
856  eta_after_->Fill(eta);
857  }
858  if (!electron_sel[2] || flags_passed == NFLAGS) {
859  if (isBarrel) {
860  sieiebarrel_after_->Fill(sieie);
861  } else if (isEndcap) {
862  sieieendcap_after_->Fill(sieie);
863  }
864  }
865  if (!electron_sel[3] || flags_passed == NFLAGS) {
866  if (isBarrel) {
867  detainbarrel_after_->Fill(detain);
868  } else if (isEndcap) {
869  detainendcap_after_->Fill(detain);
870  }
871  }
872  if (!electron_sel[4] || flags_passed == NFLAGS) {
873  if (isBarrel) {
874  ecalisobarrel_after_->Fill(ecalisovar);
875  } else if (isEndcap) {
876  ecalisoendcap_after_->Fill(ecalisovar);
877  }
878  }
879  if (!electron_sel[5] || flags_passed == NFLAGS) {
880  if (isBarrel) {
881  hcalisobarrel_after_->Fill(hcalisovar);
882  } else if (isEndcap) {
883  hcalisoendcap_after_->Fill(hcalisovar);
884  }
885  }
886  if (!electron_sel[6] || flags_passed == NFLAGS) {
887  if (isBarrel) {
888  trkisobarrel_after_->Fill(trkisovar);
889  } else if (isEndcap) {
890  trkisoendcap_after_->Fill(trkisovar);
891  }
892  }
893  // if (!electron_sel[3] || flags_passed==NFLAGS)
894  // {
895  // detain_after_->Fill(detain);
896  // }
897  // if (!electron_sel[4] || flags_passed==NFLAGS)
898  // {
899  // ecaliso_after_->Fill(trackerHits);
900  // }
901  // if (!electron_sel[5] || flags_passed==NFLAGS)
902  // {
903  // tkelectr_after_->Fill(electr.isTrackerElectron());
904  // }
905  // if (!electron_sel[6] || flags_passed==NFLAGS)
906  // {
907  // iso_after_->Fill(isovar);
908  // }
909  /* if (!electron_sel[7] || flags_passed == NFLAGS) {
910  if (!hlt_hist_done) {
911  trig_after_->Fill(trigger_fired);
912  }
913  }*/
914  //hlt_hist_done = true;
915  if (!electron_sel[7] || flags_passed == NFLAGS) {
916  mt_after_->Fill(massT);
917  }
918  if (!electron_sel[8] || flags_passed == NFLAGS) {
919  if (!met_hist_done) {
920  met_after_->Fill(met_et);
921  }
922  }
923  met_hist_done = true;
924  // if (!muon_sel[10] || flags_passed==NFLAGS)
925  // acop_after_->Fill(acop);
926  // if (!muon_sel[11] || flags_passed==NFLAGS)
927  // if (!nz1_hist_done)
928  // nz1_after_->Fill(nmuonsForZ1);
929  // nz1_hist_done = true;
930  // if (!muon_sel[11] || flags_passed==NFLAGS)
931  // if (!nz2_hist_done)
932  // nz2_after_->Fill(nmuonsForZ2);
933  // nz2_hist_done = true;
934  if (!electron_sel[9] || flags_passed == NFLAGS) {
935  if (!njets_hist_done) {
936  njets_after_->Fill(njets);
937  if (jet_et > 10) // don't want low energy "jets"
938  {
939  jet_et_after_->Fill(jet_et);
940  jet_eta_after_->Fill(jet_eta);
941  }
942  }
943  }
944  njets_hist_done = true;
945 
946  } // end N-1 histos block
947 
948  } // end loop through electrons
949 
950  // inv mass = sqrt(m_1^2 + m_2^2 + 2*(E_1*E_2 - (px1*px2 + py1*py2 + pz1+pz2)
951  // ) )
952  double invMass = 0;
953 
954  nelectrons_before_->Fill(electronCollectionSize);
955  if (electronCollectionSize > 1) {
956  invMass =
957  sqrt(electron[0][1] + electron[1][1] +
958  2 * (electron[0][2] * electron[1][2] - (electron[0][3] * electron[1][3] + electron[0][4] * electron[1][4] +
959  electron[0][5] * electron[1][5])));
960  invmass_before_->Fill(invMass);
961  invmassPU_before_->Fill(invMass, npvCount);
962  }
963 
965  if (nGoodElectrons > 1) {
966  invMass = sqrt(goodElectron[0][1] + goodElectron[1][1] +
967  2 * (goodElectron[0][2] * goodElectron[1][2] -
968  (goodElectron[0][3] * goodElectron[1][3] + goodElectron[0][4] * goodElectron[1][4] +
969  goodElectron[0][5] * goodElectron[1][5])));
970  invmass_after_->Fill(invMass);
971  invmassPU_afterZ_->Fill(invMass, npvCount);
972  npvs_afterZ_->Fill(npvCount);
973  }
974 
975  // Collect final flags
976  if (rec_sel)
977  nrec++;
978  if (eid_sel)
979  neid++;
980  if (iso_sel)
981  niso++;
982  // if (hlt_sel) nhlt++;
983  // if (met_sel) nmet++;
984 
985  if (all_sel) {
986  nsel++;
987  LogTrace("") << ">>>> Event ACCEPTED";
988  } else {
989  LogTrace("") << ">>>> Event REJECTED";
990  }
991 
992  return;
993 }
MonitorElement * ecalisoendcap_after_
Definition: EwkElecDQM.h:140
MonitorElement * jet_et_after_
Definition: EwkElecDQM.h:187
unsigned int nGoodElectrons
Definition: EwkElecDQM.h:104
MonitorElement * invmass_after_
Definition: EwkElecDQM.h:158
unsigned int nall
Definition: EwkElecDQM.h:95
MonitorElement * hcalisobarrel_after_
Definition: EwkElecDQM.h:143
double pz() const final
z coordinate of momentum vector
unsigned int nsel
Definition: EwkElecDQM.h:101
double metMax_
Definition: EwkElecDQM.h:75
MonitorElement * nelectrons_before_
Definition: EwkElecDQM.h:166
double pt() const final
transverse momentum
double sieieCutBarrel_
Definition: EwkElecDQM.h:57
edm::EDGetTokenT< reco::BeamSpot > beamSpotTag_
Definition: EwkElecDQM.h:49
MonitorElement * sieieendcap_before_
Definition: EwkElecDQM.h:115
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
MonitorElement * sieieendcap_after_
Definition: EwkElecDQM.h:116
Base class for all types of Jets.
Definition: Jet.h:20
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:72
MonitorElement * ecalisobarrel_before_
Definition: EwkElecDQM.h:136
MonitorElement * sieiebarrel_before_
Definition: EwkElecDQM.h:112
MonitorElement * hcalisobarrel_before_
Definition: EwkElecDQM.h:142
float dr03HcalTowerSumEt(int depth=0) const
Definition: GsfElectron.h:576
float dr04TkSumPt() const
Definition: GsfElectron.h:597
MonitorElement * detainbarrel_after_
Definition: EwkElecDQM.h:119
MonitorElement * ecalisobarrel_after_
Definition: EwkElecDQM.h:137
MonitorElement * detainendcap_before_
Definition: EwkElecDQM.h:121
Log< level::Error, false > LogError
MonitorElement * sieiebarrel_after_
Definition: EwkElecDQM.h:113
MonitorElement * detainbarrel_before_
Definition: EwkElecDQM.h:118
tuple vertexCollection
MonitorElement * mt_before_
Definition: EwkElecDQM.h:169
MonitorElement * nelectrons_after_
Definition: EwkElecDQM.h:167
#define LogTrace(id)
MonitorElement * hcalisoendcap_before_
Definition: EwkElecDQM.h:145
MonitorElement * npvs_afterZ_
Definition: EwkElecDQM.h:164
void Fill(long long x)
unsigned int niso
Definition: EwkElecDQM.h:98
double calcDeltaPhi(double phi1, double phi2)
Definition: EwkElecDQM.cc:996
MonitorElement * mt_after_
Definition: EwkElecDQM.h:170
double px() const final
x coordinate of momentum vector
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:225
double ecalIsoCutEndcap_
Definition: EwkElecDQM.h:67
float sigmaIetaIeta() const
Definition: GsfElectron.h:419
MonitorElement * jet_et_before_
Definition: EwkElecDQM.h:186
double detainCutEndcap_
Definition: EwkElecDQM.h:60
MonitorElement * njets_after_
Definition: EwkElecDQM.h:185
MonitorElement * pt_before_
Definition: EwkElecDQM.h:106
MonitorElement * pt_after_
Definition: EwkElecDQM.h:107
edm::EDGetTokenT< edm::View< reco::GsfElectron > > elecTag_
Definition: EwkElecDQM.h:45
double metMin_
Definition: EwkElecDQM.h:74
double trkIsoCutEndcap_
Definition: EwkElecDQM.h:71
double massSqr() const final
mass squared
Definition: MET.h:41
T sqrt(T t)
Definition: SSEVec.h:19
unsigned int nrec
Definition: EwkElecDQM.h:96
MonitorElement * hcalisoendcap_after_
Definition: EwkElecDQM.h:146
int nJetMax_
Definition: EwkElecDQM.h:87
double ecalIsoCutBarrel_
Definition: EwkElecDQM.h:66
double trkIsoCutBarrel_
Definition: EwkElecDQM.h:70
edm::EDGetTokenT< edm::View< reco::Jet > > jetToken_
Definition: EwkElecDQM.h:47
MonitorElement * invmass_before_
Definition: EwkElecDQM.h:157
MonitorElement * ecalisoendcap_before_
Definition: EwkElecDQM.h:139
MonitorElement * jet_eta_after_
Definition: EwkElecDQM.h:189
double py() const final
y coordinate of momentum vector
bool isEndcap(GeomDetEnumerators::SubDetector m)
MonitorElement * jet_eta_before_
Definition: EwkElecDQM.h:188
int prescaleSet(const edm::Event &iEvent, const edm::EventSetup &iSetup)
unsigned int neid
Definition: EwkElecDQM.h:97
double etaCut_
Definition: EwkElecDQM.h:55
double hcalIsoCutBarrel_
Definition: EwkElecDQM.h:68
MonitorElement * detainendcap_after_
Definition: EwkElecDQM.h:122
MonitorElement * njets_before_
Definition: EwkElecDQM.h:184
MonitorElement * invmassPU_afterZ_
Definition: EwkElecDQM.h:160
MonitorElement * trkisobarrel_before_
Definition: EwkElecDQM.h:148
double sieieCutEndcap_
Definition: EwkElecDQM.h:58
edm::EDGetTokenT< edm::View< reco::Vertex > > vertexTag_
Definition: EwkElecDQM.h:48
MonitorElement * met_before_
Definition: EwkElecDQM.h:172
MonitorElement * trkisoendcap_before_
Definition: EwkElecDQM.h:151
edm::EDGetTokenT< edm::View< reco::MET > > metToken_
Definition: EwkElecDQM.h:46
float dr03EcalRecHitSumEt() const
Definition: GsfElectron.h:559
double et() const final
transverse energy
MonitorElement * met_after_
Definition: EwkElecDQM.h:173
double eJetMin_
Definition: EwkElecDQM.h:86
MonitorElement * trkisobarrel_after_
Definition: EwkElecDQM.h:149
double mtMax_
Definition: EwkElecDQM.h:73
MonitorElement * eta_after_
Definition: EwkElecDQM.h:110
double detainCutBarrel_
Definition: EwkElecDQM.h:59
MonitorElement * npvs_before_
Definition: EwkElecDQM.h:162
MonitorElement * trkisoendcap_after_
Definition: EwkElecDQM.h:152
HLTPrescaleProvider hltPrescaleProvider_
Definition: EwkElecDQM.h:93
double phi() const final
momentum azimuthal angle
double hcalIsoCutEndcap_
Definition: EwkElecDQM.h:69
MonitorElement * eta_before_
Definition: EwkElecDQM.h:109
MonitorElement * invmassPU_before_
Definition: EwkElecDQM.h:159
double mtMin_
Definition: EwkElecDQM.h:72
double energy() const final
energy
double ptCut_
Definition: EwkElecDQM.h:54
double eta() const final
momentum pseudorapidity
void EwkElecDQM::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  ,
edm::EventSetup const &   
)
overrideprotectedvirtual

1.5); // elecTrig_ is now a vector of strings!

Implements DQMOneEDAnalyzer<>.

Definition at line 128 of file EwkElecDQM.cc.

References dqm::implementation::IBooker::book1D(), dqm::implementation::IBooker::book2D(), detainbarrel_after_, detainbarrel_before_, detainendcap_after_, detainendcap_before_, ecalisobarrel_after_, ecalisobarrel_before_, ecalisoendcap_after_, ecalisoendcap_before_, eJetMin_, eta_after_, eta_before_, hcalisobarrel_after_, hcalisobarrel_before_, hcalisoendcap_after_, hcalisoendcap_before_, invmass_after_, invmass_before_, invmassPU_afterZ_, invmassPU_before_, jet_et_after_, jet_et_before_, jet_eta_after_, jet_eta_before_, jetTag_, edm::InputTag::label(), met_after_, met_before_, metTag_, mt_after_, mt_before_, nelectrons_after_, nelectrons_before_, njets_after_, njets_before_, npvs_afterZ_, npvs_before_, pt_after_, pt_before_, PUBinCount_, PUMax_, dqm::implementation::NavigatorBase::setCurrentFolder(), sieiebarrel_after_, sieiebarrel_before_, sieieendcap_after_, sieieendcap_before_, trkisobarrel_after_, trkisobarrel_before_, trkisoendcap_after_, and trkisoendcap_before_.

128  {
129  ibooker.setCurrentFolder("Physics/EwkElecDQM");
130 
131  char chtitle[256] = "";
132 
133  pt_before_ = ibooker.book1D("PT_BEFORECUTS", "Electron transverse momentum [GeV]", 100, 0., 100.);
134  pt_after_ = ibooker.book1D("PT_LASTCUT", "Electron transverse momentum [GeV]", 100, 0., 100.);
135 
136  eta_before_ = ibooker.book1D("ETA_BEFORECUTS", "Electron pseudo-rapidity", 50, -2.5, 2.5);
137  eta_after_ = ibooker.book1D("ETA_LASTCUT", "Electron pseudo-rapidity", 50, -2.5, 2.5);
138 
139  sieiebarrel_before_ = ibooker.book1D("SIEIEBARREL_BEFORECUTS", "Electron #sigma_{i#etai#eta} (barrel)", 70, 0., 0.07);
140  sieiebarrel_after_ = ibooker.book1D("SIEIEBARREL_LASTCUT", "Electron #sigma_{i#etai#eta} (barrel)", 70, 0., 0.07);
141 
142  sieieendcap_before_ = ibooker.book1D("SIEIEENDCAP_BEFORECUTS", "Electron #sigma_{i#etai#eta} (endcap)", 70, 0., 0.07);
143  sieieendcap_after_ = ibooker.book1D("SIEIEENDCAP_LASTCUT", "Electron #sigma_{i#etai#eta} (endcap)", 70, 0., 0.07);
144 
146  ibooker.book1D("DETAINBARREL_BEFORECUTS", "Electron #Delta#eta_{in} (barrel)", 40, -0.02, 0.02);
147  detainbarrel_after_ = ibooker.book1D("DETAINBARREL_LASTCUT", "Electron #Delta#eta_{in} (barrel)", 40, -0.02, 0.02);
148 
150  ibooker.book1D("DETAINENDCAP_BEFORECUTS", "Electron #Delta#eta_{in} (endcap)", 40, -0.02, 0.02);
151  detainendcap_after_ = ibooker.book1D("DETAINENDCAP_LASTCUT", "Electron #Delta#eta_{in} (endcap)", 40, -0.02, 0.02);
152 
153  ecalisobarrel_before_ = ibooker.book1D(
154  "ECALISOBARREL_BEFORECUTS", "Absolute electron ECAL isolation variable (barrel) [GeV]", 50, 0., 50.);
156  ibooker.book1D("ECALISOBARREL_LASTCUT", "Absolute electron ECAL isolation variable (barrel) [GeV]", 50, 0., 50.);
157 
158  ecalisoendcap_before_ = ibooker.book1D(
159  "ECALISOENDCAP_BEFORECUTS", "Absolute electron ECAL isolation variable (endcap) [GeV]", 50, 0., 50.);
161  ibooker.book1D("ECALISOENDCAP_LASTCUT", "Absolute electron ECAL isolation variable (endcap) [GeV]", 50, 0., 50.);
162 
163  hcalisobarrel_before_ = ibooker.book1D(
164  "HCALISOBARREL_BEFORECUTS", "Absolute electron HCAL isolation variable (barrel) [GeV]", 50, 0., 50.);
166  ibooker.book1D("HCALISOBARREL_LASTCUT", "Absolute electron HCAL isolation variable (barrel) [GeV]", 50, 0., 50.);
167 
168  hcalisoendcap_before_ = ibooker.book1D(
169  "HCALISOENDCAP_BEFORECUTS", "Absolute electron HCAL isolation variable (endcap) [GeV]", 50, 0., 50.);
171  ibooker.book1D("HCALISOENDCAP_LASTCUT", "Absolute electron HCAL isolation variable (endcap) [GeV]", 50, 0., 50.);
172 
173  trkisobarrel_before_ = ibooker.book1D(
174  "TRKISOBARREL_BEFORECUTS", "Absolute electron track isolation variable (barrel) [GeV]", 50, 0., 50.);
176  ibooker.book1D("TRKISOBARREL_LASTCUT", "Absolute electron track isolation variable (barrel) [GeV]", 50, 0., 50.);
177 
178  trkisoendcap_before_ = ibooker.book1D(
179  "TRKISOENDCAP_BEFORECUTS", "Absolute electron track isolation variable (endcap) [GeV]", 50, 0., 50.);
181  ibooker.book1D("TRKISOENDCAP_LASTCUT", "Absolute electron track isolation variable (endcap) [GeV]", 50, 0., 50.);
182 
183  // trig_before_ = ibooker.book1D("TRIG_BEFORECUTS", "Trigger response", 2, -0.5,
185  // trig_after_ = ibooker.book1D("TRIG_LASTCUT", "Trigger response", 2, -0.5, 1.5);
186 
187  invmass_before_ = ibooker.book1D("INVMASS_BEFORECUTS", "Di-electron invariant mass [GeV]", 100, 0., 200.);
188  invmass_after_ = ibooker.book1D("INVMASS_AFTERCUTS", "Di-electron invariant mass [GeV]", 100, 0., 200.);
189 
190  invmassPU_before_ = ibooker.book2D("INVMASS_PU_BEFORECUTS",
191  "Di-electron invariant mass [GeV] vs PU; mass [GeV]; PU count",
192  100,
193  0.,
194  200.,
195  PUBinCount_,
196  -0.5,
197  PUMax_ + 0.5);
198  invmassPU_afterZ_ = ibooker.book2D("INVMASS_PU_AFTERZCUTS",
199  "Di-electron invariant mass [GeV] vs PU; mass [GeV]; PU count",
200  100,
201  0.,
202  200.,
203  PUBinCount_,
204  -0.5,
205  PUMax_ + 0.5);
206 
207  npvs_before_ =
208  ibooker.book1D("NPVs_BEFORECUTS", "Number of Valid Primary Vertices; nGoodPVs", PUMax_ + 1, -0.5, PUMax_ + 0.5);
209 
210  npvs_afterZ_ =
211  ibooker.book1D("NPVs_AFTERZCUTS", "Number of Valid Primary Vertices; nGoodPVs", PUMax_ + 1, -0.5, PUMax_ + 0.5);
212 
213  nelectrons_before_ = ibooker.book1D("NELECTRONS_BEFORECUTS", "Number of electrons in event", 10, -0.5, 9.5);
214  nelectrons_after_ = ibooker.book1D("NELECTRONS_AFTERCUTS", "Number of electrons in event", 10, -0.5, 9.5);
215 
216  snprintf(chtitle, 255, "Transverse mass (%s) [GeV]", metTag_.label().data());
217  mt_before_ = ibooker.book1D("MT_BEFORECUTS", chtitle, 150, 0., 300.);
218  mt_after_ = ibooker.book1D("MT_LASTCUT", chtitle, 150, 0., 300.);
219 
220  snprintf(chtitle, 255, "Missing transverse energy (%s) [GeV]", metTag_.label().data());
221  met_before_ = ibooker.book1D("MET_BEFORECUTS", chtitle, 100, 0., 200.);
222  met_after_ = ibooker.book1D("MET_LASTCUT", chtitle, 100, 0., 200.);
223 
224  snprintf(chtitle, 255, "Number of jets (%s) above %.2f GeV", jetTag_.label().data(), eJetMin_);
225  njets_before_ = ibooker.book1D("NJETS_BEFORECUTS", chtitle, 10, -0.5, 9.5);
226  njets_after_ = ibooker.book1D("NJETS_LASTCUT", chtitle, 10, -0.5, 9.5);
227 
228  snprintf(chtitle, 255, "Jet with highest E_{T} (%s)", jetTag_.label().data());
229  jet_et_before_ = ibooker.book1D("JETET1_BEFORECUTS", chtitle, 20, 0., 200.0);
230  jet_et_after_ = ibooker.book1D("JETET1_AFTERCUTS", chtitle, 20, 0., 200.0);
231 
232  snprintf(chtitle, 255, "Eta of Jet with highest E_{T} (%s)", jetTag_.label().data());
233  jet_eta_before_ = ibooker.book1D("JETETA1_BEFORECUTS", chtitle, 20, -5, 5);
234  jet_eta_after_ = ibooker.book1D("JETETA1_AFTERCUTS", chtitle, 20, -5, 5);
235 }
MonitorElement * ecalisoendcap_after_
Definition: EwkElecDQM.h:140
MonitorElement * jet_et_after_
Definition: EwkElecDQM.h:187
MonitorElement * invmass_after_
Definition: EwkElecDQM.h:158
MonitorElement * hcalisobarrel_after_
Definition: EwkElecDQM.h:143
MonitorElement * nelectrons_before_
Definition: EwkElecDQM.h:166
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
MonitorElement * sieieendcap_before_
Definition: EwkElecDQM.h:115
MonitorElement * sieieendcap_after_
Definition: EwkElecDQM.h:116
MonitorElement * ecalisobarrel_before_
Definition: EwkElecDQM.h:136
edm::InputTag jetTag_
Definition: EwkElecDQM.h:43
MonitorElement * sieiebarrel_before_
Definition: EwkElecDQM.h:112
MonitorElement * hcalisobarrel_before_
Definition: EwkElecDQM.h:142
MonitorElement * detainbarrel_after_
Definition: EwkElecDQM.h:119
MonitorElement * ecalisobarrel_after_
Definition: EwkElecDQM.h:137
MonitorElement * detainendcap_before_
Definition: EwkElecDQM.h:121
MonitorElement * sieiebarrel_after_
Definition: EwkElecDQM.h:113
MonitorElement * detainbarrel_before_
Definition: EwkElecDQM.h:118
MonitorElement * mt_before_
Definition: EwkElecDQM.h:169
MonitorElement * nelectrons_after_
Definition: EwkElecDQM.h:167
MonitorElement * hcalisoendcap_before_
Definition: EwkElecDQM.h:145
MonitorElement * npvs_afterZ_
Definition: EwkElecDQM.h:164
MonitorElement * mt_after_
Definition: EwkElecDQM.h:170
MonitorElement * jet_et_before_
Definition: EwkElecDQM.h:186
MonitorElement * njets_after_
Definition: EwkElecDQM.h:185
MonitorElement * pt_before_
Definition: EwkElecDQM.h:106
MonitorElement * pt_after_
Definition: EwkElecDQM.h:107
unsigned int PUBinCount_
Definition: EwkElecDQM.h:90
MonitorElement * hcalisoendcap_after_
Definition: EwkElecDQM.h:146
edm::InputTag metTag_
Definition: EwkElecDQM.h:42
MonitorElement * invmass_before_
Definition: EwkElecDQM.h:157
MonitorElement * ecalisoendcap_before_
Definition: EwkElecDQM.h:139
unsigned int PUMax_
Definition: EwkElecDQM.h:90
MonitorElement * jet_eta_after_
Definition: EwkElecDQM.h:189
MonitorElement * jet_eta_before_
Definition: EwkElecDQM.h:188
MonitorElement * detainendcap_after_
Definition: EwkElecDQM.h:122
MonitorElement * njets_before_
Definition: EwkElecDQM.h:184
MonitorElement * invmassPU_afterZ_
Definition: EwkElecDQM.h:160
MonitorElement * trkisobarrel_before_
Definition: EwkElecDQM.h:148
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
MonitorElement * met_before_
Definition: EwkElecDQM.h:172
MonitorElement * trkisoendcap_before_
Definition: EwkElecDQM.h:151
std::string const & label() const
Definition: InputTag.h:36
MonitorElement * met_after_
Definition: EwkElecDQM.h:173
double eJetMin_
Definition: EwkElecDQM.h:86
MonitorElement * trkisobarrel_after_
Definition: EwkElecDQM.h:149
MonitorElement * eta_after_
Definition: EwkElecDQM.h:110
MonitorElement * npvs_before_
Definition: EwkElecDQM.h:162
MonitorElement * trkisoendcap_after_
Definition: EwkElecDQM.h:152
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * eta_before_
Definition: EwkElecDQM.h:109
MonitorElement * invmassPU_before_
Definition: EwkElecDQM.h:159
double EwkElecDQM::calcDeltaPhi ( double  phi1,
double  phi2 
)
protected

Definition at line 996 of file EwkElecDQM.cc.

References srCondWrite_cfg::deltaPhi.

Referenced by analyze().

996  {
997  double deltaPhi = phi1 - phi2;
998 
999  if (deltaPhi < 0)
1000  deltaPhi = -deltaPhi;
1001 
1002  if (deltaPhi > 3.1415926) {
1003  deltaPhi = 2 * 3.1415926 - deltaPhi;
1004  }
1005 
1006  return deltaPhi;
1007 }
void EwkElecDQM::dqmBeginRun ( const edm::Run iRun,
const edm::EventSetup iSet 
)
overrideprotectedvirtual

Reimplemented from DQMOneEDAnalyzer<>.

Definition at line 109 of file EwkElecDQM.cc.

References hltPrescaleProvider_, HLTPrescaleProvider::init(), isValidHltConfig_, LogTrace, nall, neid, niso, nrec, and nsel.

109  {
110  nall = 0;
111  nsel = 0;
112 
113  nrec = 0;
114  neid = 0;
115  niso = 0;
116  // nhlt = 0;
117  // nmet = 0;
118 
119  // passed as parameter to HLTConfigProvider::init(), not yet used
120  bool isConfigChanged = false;
121  // isValidHltConfig_ could be used to short-circuit analyze() in case of
122  // problems
123  isValidHltConfig_ = hltPrescaleProvider_.init(iRun, iSet, "HLT", isConfigChanged);
124 
125  LogTrace("") << "isValidHltConfig_=" << isValidHltConfig_ << "\n";
126 }
unsigned int nall
Definition: EwkElecDQM.h:95
unsigned int nsel
Definition: EwkElecDQM.h:101
#define LogTrace(id)
unsigned int niso
Definition: EwkElecDQM.h:98
unsigned int nrec
Definition: EwkElecDQM.h:96
unsigned int neid
Definition: EwkElecDQM.h:97
bool isValidHltConfig_
Definition: EwkElecDQM.h:92
HLTPrescaleProvider hltPrescaleProvider_
Definition: EwkElecDQM.h:93
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
void EwkElecDQM::dqmEndRun ( const edm::Run r,
const edm::EventSetup  
)
overrideprotectedvirtual

Reimplemented from DQMOneEDAnalyzer<>.

Definition at line 237 of file EwkElecDQM.cc.

References python.cmstools::all(), submitPVValidationJobs::err, nall, neid, niso, nrec, nsel, pileupDistInMC::num, and mathSSE::sqrt().

237  {
238  // overall
239  double all = nall;
240  double esel = nsel / all;
241  LogVerbatim("") << "\n>>>>>> SELECTION SUMMARY BEGIN >>>>>>>>>>>>>>>";
242  LogVerbatim("") << "Total number of events analyzed: " << nall << " [events]";
243  LogVerbatim("") << "Total number of events selected: " << nsel << " [events]";
244  LogVerbatim("") << "Overall efficiency: "
245  << "(" << setprecision(4) << esel * 100. << " +/- " << setprecision(2)
246  << sqrt(esel * (1 - esel) / all) * 100. << ")%";
247 
248  double erec = nrec / all;
249  double eeid = neid / all;
250  double eiso = niso / all;
251  // double ehlt = nhlt/all;
252  // double emet = nmet/all;
253 
254  // general reconstruction step??
255  double num = nrec;
256  double eff = erec;
257  double err = sqrt(eff * (1 - eff) / all);
258  LogVerbatim("") << "Passing Pt/Eta/Quality cuts: " << num << " [events], (" << setprecision(4) << eff * 100.
259  << " +/- " << setprecision(2) << err * 100. << ")%";
260 
261  // electron ID step
262  num = neid;
263  eff = eeid;
264  err = sqrt(eff * (1 - eff) / all);
265  double effstep = 0.;
266  double errstep = 0.;
267  if (nrec > 0)
268  effstep = eeid / erec;
269  if (nrec > 0)
270  errstep = sqrt(effstep * (1 - effstep) / nrec);
271  LogVerbatim("") << "Passing eID cuts: " << num << " [events], (" << setprecision(4) << eff * 100. << " +/- "
272  << setprecision(2) << err * 100. << ")%, to previous step: (" << setprecision(4) << effstep * 100.
273  << " +/- " << setprecision(2) << errstep * 100. << ")%";
274 
275  // isolation step
276  num = niso;
277  eff = eiso;
278  err = sqrt(eff * (1 - eff) / all);
279  effstep = 0.;
280  errstep = 0.;
281  if (neid > 0)
282  effstep = eiso / eeid;
283  if (neid > 0)
284  errstep = sqrt(effstep * (1 - effstep) / neid);
285  LogVerbatim("") << "Passing isolation cuts: " << num << " [events], (" << setprecision(4) << eff * 100.
286  << " +/- " << setprecision(2) << err * 100. << ")%, to previous step: (" << setprecision(4)
287  << effstep * 100. << " +/- " << setprecision(2) << errstep * 100. << ")%";
288 
289  // // trigger step
290  // num = nhlt;
291  // eff = ehlt;
292  // err = sqrt(eff*(1-eff)/all);
293  // effstep = 0.;
294  // errstep = 0.;
295  // if (niso>0) effstep = ehlt/eiso;
296  // if (niso>0) errstep = sqrt(effstep*(1-effstep)/niso);
297  // LogVerbatim("") << "Passing HLT criteria: " << num << "
298  // [events], (" << setprecision(4) << eff*100. <<" +/- "<< setprecision(2) <<
299  // err*100. << ")%, to previous step: (" << setprecision(4) << effstep*100.
300  // << " +/- "<< setprecision(2) << errstep*100. <<")%";
301 
302  // trigger step
303  num = nsel;
304  eff = esel;
305  err = sqrt(eff * (1 - eff) / all);
306  effstep = 0.;
307  errstep = 0.;
308  if (niso > 0)
309  effstep = esel / eiso;
310  if (niso > 0)
311  errstep = sqrt(effstep * (1 - effstep) / niso);
312  LogVerbatim("") << "Passing HLT criteria: " << num << " [events], (" << setprecision(4) << eff * 100.
313  << " +/- " << setprecision(2) << err * 100. << ")%, to previous step: (" << setprecision(4)
314  << effstep * 100. << " +/- " << setprecision(2) << errstep * 100. << ")%";
315 
316  // // met/acoplanarity cuts
317  // num = nmet;
318  // eff = emet;
319  // err = sqrt(eff*(1-eff)/all);
320  // effstep = 0.;
321  // errstep = 0.;
322  // if (nhlt>0) effstep = emet/ehlt;
323  // if (nhlt>0) errstep = sqrt(effstep*(1-effstep)/nhlt);
324  // LogVerbatim("") << "Passing MET/acoplanarity cuts: " << num << "
325  // [events], (" << setprecision(4) << eff*100. <<" +/- "<< setprecision(2) <<
326  // err*100. << ")%, to previous step: (" << setprecision(4) << effstep*100.
327  // << " +/- "<< setprecision(2) << errstep*100. <<")%";
328 
329  // // Z/top selection cuts ALSO LAST STEP so "sel" for "selection"
330  // num = nsel;
331  // eff = esel;
332  // err = sqrt(eff*(1-eff)/all);
333  // effstep = 0.;
334  // errstep = 0.;
335  // if (nmet>0) effstep = esel/emet;
336  // if (nmet>0) errstep = sqrt(effstep*(1-effstep)/nmet);
337  // LogVerbatim("") << "Passing Z/top rejection cuts: " << num << "
338  // [events], (" << setprecision(4) << eff*100. <<" +/- "<< setprecision(2) <<
339  // err*100. << ")%, to previous step: (" << setprecision(4) << effstep*100.
340  // << " +/- "<< setprecision(2) << errstep*100. <<")%";
341 
342  LogVerbatim("") << ">>>>>> SELECTION SUMMARY END >>>>>>>>>>>>>>>\n";
343 }
Log< level::Info, true > LogVerbatim
unsigned int nall
Definition: EwkElecDQM.h:95
unsigned int nsel
Definition: EwkElecDQM.h:101
unsigned int niso
Definition: EwkElecDQM.h:98
T sqrt(T t)
Definition: SSEVec.h:19
unsigned int nrec
Definition: EwkElecDQM.h:96
def all
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
unsigned int neid
Definition: EwkElecDQM.h:97

Member Data Documentation

edm::EDGetTokenT<reco::BeamSpot> EwkElecDQM::beamSpotTag_
private

Definition at line 49 of file EwkElecDQM.h.

Referenced by analyze().

MonitorElement* EwkElecDQM::detainbarrel_after_
private

Definition at line 119 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::detainbarrel_before_
private

Definition at line 118 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

double EwkElecDQM::detainCutBarrel_
private

Definition at line 59 of file EwkElecDQM.h.

Referenced by analyze().

double EwkElecDQM::detainCutEndcap_
private

Definition at line 60 of file EwkElecDQM.h.

Referenced by analyze().

MonitorElement* EwkElecDQM::detainendcap_after_
private

Definition at line 122 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::detainendcap_before_
private

Definition at line 121 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::ecalisobarrel_after_
private

Definition at line 137 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::ecalisobarrel_before_
private

Definition at line 136 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

double EwkElecDQM::ecalIsoCutBarrel_
private

Definition at line 66 of file EwkElecDQM.h.

Referenced by analyze().

double EwkElecDQM::ecalIsoCutEndcap_
private

Definition at line 67 of file EwkElecDQM.h.

Referenced by analyze().

MonitorElement* EwkElecDQM::ecalisoendcap_after_
private

Definition at line 140 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::ecalisoendcap_before_
private

Definition at line 139 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

double EwkElecDQM::eJetMin_
private

Definition at line 86 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

edm::EDGetTokenT<edm::View<reco::GsfElectron> > EwkElecDQM::elecTag_
private

Definition at line 45 of file EwkElecDQM.h.

Referenced by analyze().

const std::vector<std::string> EwkElecDQM::elecTrig_
private

Definition at line 53 of file EwkElecDQM.h.

MonitorElement* EwkElecDQM::eta_after_
private

Definition at line 110 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::eta_before_
private

Definition at line 109 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

double EwkElecDQM::etaCut_
private

Definition at line 55 of file EwkElecDQM.h.

Referenced by analyze().

MonitorElement* EwkElecDQM::hcalisobarrel_after_
private

Definition at line 143 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::hcalisobarrel_before_
private

Definition at line 142 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

double EwkElecDQM::hcalIsoCutBarrel_
private

Definition at line 68 of file EwkElecDQM.h.

Referenced by analyze().

double EwkElecDQM::hcalIsoCutEndcap_
private

Definition at line 69 of file EwkElecDQM.h.

Referenced by analyze().

MonitorElement* EwkElecDQM::hcalisoendcap_after_
private

Definition at line 146 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::hcalisoendcap_before_
private

Definition at line 145 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

HLTPrescaleProvider EwkElecDQM::hltPrescaleProvider_
private

Definition at line 93 of file EwkElecDQM.h.

Referenced by analyze(), and dqmBeginRun().

MonitorElement* EwkElecDQM::invmass_after_
private

Definition at line 158 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::invmass_before_
private

Definition at line 157 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::invmassPU_afterZ_
private

Definition at line 160 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::invmassPU_before_
private

Definition at line 159 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

bool EwkElecDQM::isValidHltConfig_
private

Definition at line 92 of file EwkElecDQM.h.

Referenced by dqmBeginRun(), and EwkElecDQM().

MonitorElement* EwkElecDQM::jet_et_after_
private

Definition at line 187 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::jet_et_before_
private

Definition at line 186 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::jet_eta_after_
private

Definition at line 189 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::jet_eta_before_
private

Definition at line 188 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

edm::InputTag EwkElecDQM::jetTag_
private

Definition at line 43 of file EwkElecDQM.h.

Referenced by bookHistograms().

edm::EDGetTokenT<edm::View<reco::Jet> > EwkElecDQM::jetToken_
private

Definition at line 47 of file EwkElecDQM.h.

Referenced by analyze().

MonitorElement* EwkElecDQM::met_after_
private

Definition at line 173 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::met_before_
private

Definition at line 172 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

bool EwkElecDQM::metIncludesMuons_
private

Definition at line 50 of file EwkElecDQM.h.

double EwkElecDQM::metMax_
private

Definition at line 75 of file EwkElecDQM.h.

Referenced by analyze().

double EwkElecDQM::metMin_
private

Definition at line 74 of file EwkElecDQM.h.

Referenced by analyze().

edm::InputTag EwkElecDQM::metTag_
private

Definition at line 42 of file EwkElecDQM.h.

Referenced by bookHistograms().

edm::EDGetTokenT<edm::View<reco::MET> > EwkElecDQM::metToken_
private

Definition at line 46 of file EwkElecDQM.h.

Referenced by analyze().

MonitorElement* EwkElecDQM::mt_after_
private

Definition at line 170 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::mt_before_
private

Definition at line 169 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

double EwkElecDQM::mtMax_
private

Definition at line 73 of file EwkElecDQM.h.

Referenced by analyze().

double EwkElecDQM::mtMin_
private

Definition at line 72 of file EwkElecDQM.h.

Referenced by analyze().

unsigned int EwkElecDQM::nall
private

Definition at line 95 of file EwkElecDQM.h.

Referenced by analyze(), dqmBeginRun(), and dqmEndRun().

unsigned int EwkElecDQM::neid
private

Definition at line 97 of file EwkElecDQM.h.

Referenced by analyze(), dqmBeginRun(), and dqmEndRun().

MonitorElement* EwkElecDQM::nelectrons_after_
private

Definition at line 167 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::nelectrons_before_
private

Definition at line 166 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

unsigned int EwkElecDQM::nGoodElectrons
private

Definition at line 104 of file EwkElecDQM.h.

Referenced by analyze().

unsigned int EwkElecDQM::niso
private

Definition at line 98 of file EwkElecDQM.h.

Referenced by analyze(), dqmBeginRun(), and dqmEndRun().

int EwkElecDQM::nJetMax_
private

Definition at line 87 of file EwkElecDQM.h.

Referenced by analyze().

MonitorElement* EwkElecDQM::njets_after_
private

Definition at line 185 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::njets_before_
private

Definition at line 184 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::npvs_afterZ_
private

Definition at line 164 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::npvs_before_
private

Definition at line 162 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

unsigned int EwkElecDQM::nrec
private

Definition at line 96 of file EwkElecDQM.h.

Referenced by analyze(), dqmBeginRun(), and dqmEndRun().

unsigned int EwkElecDQM::nsel
private

Definition at line 101 of file EwkElecDQM.h.

Referenced by analyze(), dqmBeginRun(), and dqmEndRun().

MonitorElement* EwkElecDQM::pt_after_
private

Definition at line 107 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::pt_before_
private

Definition at line 106 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

double EwkElecDQM::ptCut_
private
unsigned int EwkElecDQM::PUBinCount_
private

Definition at line 90 of file EwkElecDQM.h.

Referenced by bookHistograms().

unsigned int EwkElecDQM::PUMax_
private

Definition at line 90 of file EwkElecDQM.h.

Referenced by bookHistograms().

MonitorElement* EwkElecDQM::sieiebarrel_after_
private

Definition at line 113 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::sieiebarrel_before_
private

Definition at line 112 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

double EwkElecDQM::sieieCutBarrel_
private

Definition at line 57 of file EwkElecDQM.h.

Referenced by analyze().

double EwkElecDQM::sieieCutEndcap_
private

Definition at line 58 of file EwkElecDQM.h.

Referenced by analyze().

MonitorElement* EwkElecDQM::sieieendcap_after_
private

Definition at line 116 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::sieieendcap_before_
private

Definition at line 115 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::trig_after_
private

Definition at line 155 of file EwkElecDQM.h.

MonitorElement* EwkElecDQM::trig_before_
private

Definition at line 154 of file EwkElecDQM.h.

edm::EDGetTokenT<edm::TriggerResults> EwkElecDQM::trigTag_
private

Definition at line 44 of file EwkElecDQM.h.

MonitorElement* EwkElecDQM::trkisobarrel_after_
private

Definition at line 149 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::trkisobarrel_before_
private

Definition at line 148 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

double EwkElecDQM::trkIsoCutBarrel_
private

Definition at line 70 of file EwkElecDQM.h.

Referenced by analyze().

double EwkElecDQM::trkIsoCutEndcap_
private

Definition at line 71 of file EwkElecDQM.h.

Referenced by analyze().

MonitorElement* EwkElecDQM::trkisoendcap_after_
private

Definition at line 152 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* EwkElecDQM::trkisoendcap_before_
private

Definition at line 151 of file EwkElecDQM.h.

Referenced by analyze(), and bookHistograms().

edm::EDGetTokenT<edm::View<reco::Vertex> > EwkElecDQM::vertexTag_
private

Definition at line 48 of file EwkElecDQM.h.

Referenced by analyze().