CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
ZToMuMuGammaAnalyzer Class Reference

EgammaCoreTools. More...

#include <ZToMuMuGammaAnalyzer.h>

Inheritance diagram for ZToMuMuGammaAnalyzer:
DQMEDAnalyzer edm::one::EDProducer< edm::Accumulator, edm::EndLuminosityBlockProducer, edm::EndRunProducer, edm::one::WatchLuminosityBlocks, edm::one::WatchRuns > edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
 ZToMuMuGammaAnalyzer (const edm::ParameterSet &)
 
 ~ZToMuMuGammaAnalyzer () override
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &ev, edm::EventSetup const &es) final
 
virtual void analyze (edm::Event const &, edm::EventSetup const &)
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) override
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
 DQMEDAnalyzer (DQMEDAnalyzer const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer &&)=delete
 
void endLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &) override
 
void endLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) override
 
void endRunProduce (edm::Run &run, edm::EventSetup const &setup) override
 
 ~DQMEDAnalyzer () override=default
 
- Public Member Functions inherited from edm::one::EDProducer< edm::Accumulator, edm::EndLuminosityBlockProducer, edm::EndRunProducer, edm::one::WatchLuminosityBlocks, edm::one::WatchRuns >
 EDProducer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () 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
 
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::vector< ModuleDescription const * > &modules, 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
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

bool basicMuonSelection (const reco::Muon &m)
 
float mumuGammaInvMass (const reco::Muon &mu1, const reco::Muon &mu2, const reco::PhotonRef &pho)
 
float mumuInvMass (const reco::Muon &m1, const reco::Muon &m2)
 
bool muonSelection (const reco::Muon &m, const reco::BeamSpot &bs)
 
bool photonSelection (const reco::PhotonRef &p)
 

Private Attributes

edm::EDGetTokenT< edm::SortedCollection< EcalRecHit, edm::StrictWeakOrdering< EcalRecHit > > > barrelRecHit_token_
 
edm::EDGetTokenT< reco::BeamSpotbeamSpot_token_
 
std::stringstream currentFolder_
 
int eBin_
 
double eMax_
 
double eMin_
 
edm::EDGetTokenT< edm::SortedCollection< EcalRecHit, edm::StrictWeakOrdering< EcalRecHit > > > endcapRecHit_token_
 
int etaBin_
 
double etaMax_
 
double etaMin_
 
int etBin_
 
double etMax_
 
double etMin_
 
float farMuonEcalIso_
 
float farMuonMinPt_
 
float farMuonTrackIso_
 
std::string fName_
 
MonitorElementh1_mumuGammaInvMass_ [3]
 
MonitorElementh1_mumuInvMass_ [3]
 photon histos More...
 
MonitorElementh2_e1x5VsEt_ [3]
 
MonitorElementh2_e1x5VsEta_ [3]
 
MonitorElementh2_e2x5VsEt_ [3]
 
MonitorElementh2_e2x5VsEta_ [3]
 
MonitorElementh2_ecalSumVsEt_ [3]
 
MonitorElementh2_ecalSumVsEta_ [3]
 
MonitorElementh2_hcalSumVsEt_ [3]
 
MonitorElementh2_hcalSumVsEta_ [3]
 
MonitorElementh2_nTrackIsolHollowVsEt_ [3]
 
MonitorElementh2_nTrackIsolHollowVsEta_ [3]
 
MonitorElementh2_nTrackIsolSolidVsEt_ [3]
 
MonitorElementh2_nTrackIsolSolidVsEta_ [3]
 
MonitorElementh2_r1x5VsEt_ [3]
 
MonitorElementh2_r1x5VsEta_ [3]
 
MonitorElementh2_r2x5VsEt_ [3]
 
MonitorElementh2_r2x5VsEta_ [3]
 
MonitorElementh2_r9VsEt_ [3]
 
MonitorElementh2_r9VsEta_ [3]
 
MonitorElementh2_sigmaIetaIetaVsEta_ [3]
 
MonitorElementh2_trackPtSumHollowVsEt_ [3]
 
MonitorElementh2_trackPtSumHollowVsEta_ [3]
 
MonitorElementh2_trackPtSumSolidVsEt_ [3]
 
MonitorElementh2_trackPtSumSolidVsEta_ [3]
 
MonitorElementh_chHadIso_ [3]
 
MonitorElementh_dRPhoPFcand_ChHad_Cleaned_ [3]
 
MonitorElementh_dRPhoPFcand_ChHad_unCleaned_ [3]
 
MonitorElementh_dRPhoPFcand_NeuHad_Cleaned_ [3]
 
MonitorElementh_dRPhoPFcand_NeuHad_unCleaned_ [3]
 
MonitorElementh_dRPhoPFcand_Pho_Cleaned_ [3]
 
MonitorElementh_dRPhoPFcand_Pho_unCleaned_ [3]
 
MonitorElementh_e1x5_ [3]
 
MonitorElementh_e2x5_ [3]
 
MonitorElementh_ecalSum_ [3]
 
MonitorElementh_etOutsideMustache_ [3]
 
MonitorElementh_h1OverE_ [3]
 
MonitorElementh_h2OverE_ [3]
 
MonitorElementh_hcalSum_ [3]
 
MonitorElementh_hOverE_ [3]
 
MonitorElementh_nCluOutsideMustache_ [3]
 
MonitorElementh_newhOverE_ [3]
 
MonitorElementh_nHadIso_ [3]
 
MonitorElementh_nPho_ [3]
 
MonitorElementh_nRecoVtx_
 
MonitorElementh_nTrackIsolHollow_ [3]
 
MonitorElementh_nTrackIsolSolid_ [3]
 
MonitorElementh_pfMva_ [3]
 
MonitorElementh_phoE_ [3]
 
MonitorElementh_phoEt_ [3]
 
MonitorElementh_phoEta_ [3]
 
MonitorElementh_phoIso_ [3]
 
MonitorElementh_phoPhi_ [3]
 
MonitorElementh_phoSigmaEoverE_ [3]
 
MonitorElementh_phoSigmaIetaIeta_ [3]
 
MonitorElementh_r1x5_ [3]
 
MonitorElementh_r2x5_ [3]
 
MonitorElementh_r9_ [3]
 
MonitorElementh_scEta_ [3]
 
MonitorElementh_scPhi_ [3]
 
MonitorElementh_SumPtOverPhoPt_ChHad_Cleaned_ [3]
 
MonitorElementh_SumPtOverPhoPt_ChHad_unCleaned_ [3]
 
MonitorElementh_SumPtOverPhoPt_NeuHad_Cleaned_ [3]
 
MonitorElementh_SumPtOverPhoPt_NeuHad_unCleaned_ [3]
 
MonitorElementh_SumPtOverPhoPt_Pho_Cleaned_ [3]
 
MonitorElementh_SumPtOverPhoPt_Pho_unCleaned_ [3]
 
MonitorElementh_trackPtSumHollow_ [3]
 
MonitorElementh_trackPtSumSolid_ [3]
 
int hOverEBin_
 
double hOverEMax_
 
double hOverEMin_
 
bool makeProfiles_
 
float maxMumuGammaInvMass_
 
float maxMumuInvMass_
 
float minMumuGammaInvMass_
 
float minMumuInvMass_
 
int minPixStripHits_
 
edm::EDGetTokenT< std::vector< reco::Muon > > muon_token_
 
int muonMatches_
 
float muonMaxChi2_
 
float muonMaxDxy_
 
float muonMinPt_
 
float muonTightEta_
 
float muonTrackIso_
 
float nearMuonDr_
 
float nearMuonHcalIso_
 
int nEvt_
 
int numberBin_
 
double numberMax_
 
double numberMin_
 
edm::EDGetTokenT< reco::VertexCollectionoffline_pvToken_
 
MonitorElementp_e1x5VsEt_ [3]
 
MonitorElementp_e1x5VsEta_ [3]
 
MonitorElementp_e2x5VsEt_ [3]
 
MonitorElementp_e2x5VsEta_ [3]
 
MonitorElementp_ecalSumVsEt_ [3]
 
MonitorElementp_ecalSumVsEta_ [3]
 
MonitorElementp_hcalSumVsEt_ [3]
 
MonitorElementp_hcalSumVsEta_ [3]
 
MonitorElementp_hOverEVsEt_ [3]
 
MonitorElementp_hOverEVsEta_ [3]
 
MonitorElementp_newhOverEVsEt_ [3]
 
MonitorElementp_newhOverEVsEta_ [3]
 
MonitorElementp_nTrackIsolHollowVsEt_ [3]
 
MonitorElementp_nTrackIsolHollowVsEta_ [3]
 
MonitorElementp_nTrackIsolSolidVsEt_ [3]
 
MonitorElementp_nTrackIsolSolidVsEta_ [3]
 
MonitorElementp_phoSigmaEoverEVsNVtx_ [3]
 
MonitorElementp_r1x5VsEt_ [3]
 
MonitorElementp_r1x5VsEta_ [3]
 
MonitorElementp_r2x5VsEt_ [3]
 
MonitorElementp_r2x5VsEta_ [3]
 
MonitorElementp_r9VsEt_ [3]
 
MonitorElementp_r9VsEta_ [3]
 
MonitorElementp_sigmaIetaIetaVsEta_ [3]
 
MonitorElementp_trackPtSumHollowVsEt_ [3]
 
MonitorElementp_trackPtSumHollowVsEta_ [3]
 
MonitorElementp_trackPtSumSolidVsEt_ [3]
 
MonitorElementp_trackPtSumSolidVsEta_ [3]
 
edm::EDGetTokenT< reco::PFCandidateCollectionpfCandidates_
 
int phiBin_
 
double phiMax_
 
double phiMin_
 
edm::EDGetTokenT< std::vector< reco::Photon > > photon_token_
 
edm::EDGetTokenT< edm::ValueMap< bool > > PhotonIDLoose_token_
 
edm::EDGetTokenT< edm::ValueMap< bool > > PhotonIDTight_token_
 
edm::EDGetTokenT< edm::ValueMap< std::vector< reco::PFCandidateRef > > > photonIsoValmap_token_
 
float photonMaxEta_
 
float photonMinEt_
 
float photonTrackIso_
 
unsigned int prescaleFactor_
 
int r9Bin_
 
double r9Max_
 
double r9Min_
 
int reducedEtaBin_
 
int reducedEtBin_
 
int reducedR9Bin_
 
int reducedSumBin_
 
int sigmaIetaBin_
 
double sigmaIetaMax_
 
double sigmaIetaMin_
 
int sumBin_
 
double sumMax_
 
double sumMin_
 
edm::EDGetTokenT< trigger::TriggerEventtriggerEvent_token_
 
bool use2DHistos_
 
int validMuonHits_
 
int validPixHits_
 

Additional Inherited Members

- 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 Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
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 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)
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 

Detailed Description

EgammaCoreTools.

$Id: ZToMuMuGammaAnalyzer authors: Nancy Marinelli, U. of Notre Dame, US Nathan Kellams, U. of Notre Dame, US

$Id: ZToMuMuGammaAnalyzer authors: Nancy Marinelli, U. of Notre Dame, US Jamie Antonelli, U. of Notre Dame, US Nathan Kellams, U. of Notre Dame, US

Definition at line 99 of file ZToMuMuGammaAnalyzer.h.

Constructor & Destructor Documentation

ZToMuMuGammaAnalyzer::ZToMuMuGammaAnalyzer ( const edm::ParameterSet pset)
explicit

Definition at line 18 of file ZToMuMuGammaAnalyzer.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), and AlCaHLTBitMon_QueryRunRegistry::string.

19 {
20  fName_ = pset.getParameter<std::string>("analyzerName");
21  prescaleFactor_ = pset.getUntrackedParameter<int>("prescaleFactor",1);
22  use2DHistos_ = pset.getParameter<bool>("use2DHistos");
23  makeProfiles_ = pset.getParameter<bool>("makeProfiles");
24 
25  triggerEvent_token_ = consumes<trigger::TriggerEvent>(pset.getParameter<edm::InputTag>("triggerEvent"));
26  offline_pvToken_ = consumes<reco::VertexCollection>(pset.getUntrackedParameter<edm::InputTag>("offlinePV", edm::InputTag("offlinePrimaryVertices")));
27  photon_token_ = consumes<vector<reco::Photon> >(pset.getParameter<edm::InputTag>("phoProducer"));
28  muon_token_ = consumes<vector<reco::Muon> >(pset.getParameter<edm::InputTag>("muonProducer"));
29  pfCandidates_ = consumes<reco::PFCandidateCollection>(pset.getParameter<edm::InputTag>("pfCandidates"));
30  photonIsoValmap_token_ = consumes<edm::ValueMap<std::vector<reco::PFCandidateRef> > >(pset.getParameter<edm::InputTag>("particleBasedIso"));
31  barrelRecHit_token_ = consumes<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > >(pset.getParameter<edm::InputTag>("barrelRecHitProducer"));
32  endcapRecHit_token_ = consumes<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > >(pset.getParameter<edm::InputTag>("endcapRecHitProducer"));
33  beamSpot_token_ = consumes<reco::BeamSpot>(pset.getParameter<edm::InputTag>("beamSpot"));
34 
35  nEvt_=0;
36 
37  // Muon selection
38  muonMinPt_ = pset.getParameter<double>("muonMinPt");
39  minPixStripHits_ = pset.getParameter<int>("minPixStripHits");
40  muonMaxChi2_ = pset.getParameter<double>("muonMaxChi2");
41  muonMaxDxy_ = pset.getParameter<double>("muonMaxDxy");
42  muonMatches_ = pset.getParameter<int>("muonMatches");
43  validPixHits_ = pset.getParameter<int>("validPixHits");
44  validMuonHits_ = pset.getParameter<int>("validMuonHits");
45  muonTrackIso_ = pset.getParameter<double>("muonTrackIso");
46  muonTightEta_ = pset.getParameter<double>("muonTightEta");
47  // Dimuon selection
48  minMumuInvMass_ = pset.getParameter<double>("minMumuInvMass");
49  maxMumuInvMass_ = pset.getParameter<double>("maxMumuInvMass");
50  // Photon selection
51  photonMinEt_ = pset.getParameter<double>("photonMinEt");
52  photonMaxEta_ = pset.getParameter<double>("photonMaxEta");
53  photonTrackIso_ = pset.getParameter<double>("photonTrackIso");
54  // mumuGamma selection
55  nearMuonDr_ = pset.getParameter<double>("nearMuonDr");
56  nearMuonHcalIso_ = pset.getParameter<double>("nearMuonHcalIso");
57  farMuonEcalIso_ = pset.getParameter<double>("farMuonEcalIso");
58  farMuonTrackIso_ = pset.getParameter<double>("farMuonTrackIso");
59  farMuonMinPt_ = pset.getParameter<double>("farMuonMinPt");
60  minMumuGammaInvMass_ = pset.getParameter<double>("minMumuGammaInvMass");
61  maxMumuGammaInvMass_ = pset.getParameter<double>("maxMumuGammaInvMass");
62 
63  // Histogram parameters
64  eMin_ = pset.getParameter<double>("eMin");
65  eMax_ = pset.getParameter<double>("eMax");
66  eBin_ = pset.getParameter<int>("eBin");
67 
68  etMin_ = pset.getParameter<double>("etMin");
69  etMax_ = pset.getParameter<double>("etMax");
70  etBin_ = pset.getParameter<int>("etBin");
71 
72  sumMin_ = pset.getParameter<double>("sumMin");
73  sumMax_ = pset.getParameter<double>("sumMax");
74  sumBin_ = pset.getParameter<int>("sumBin");
75 
76  etaMin_ = pset.getParameter<double>("etaMin");
77  etaMax_ = pset.getParameter<double>("etaMax");
78  etaBin_ = pset.getParameter<int>("etaBin");
79 
80  phiMin_ = pset.getParameter<double>("phiMin");
81  phiMax_ = pset.getParameter<double>("phiMax");
82  phiBin_ = pset.getParameter<int>("phiBin");
83 
84  r9Min_ = pset.getParameter<double>("r9Min");
85  r9Max_ = pset.getParameter<double>("r9Max");
86  r9Bin_ = pset.getParameter<int>("r9Bin");
87 
88  hOverEMin_ = pset.getParameter<double>("hOverEMin");
89  hOverEMax_ = pset.getParameter<double>("hOverEMax");
90  hOverEBin_ = pset.getParameter<int>("hOverEBin");
91 
92  numberMin_ = pset.getParameter<double>("numberMin");
93  numberMax_ = pset.getParameter<double>("numberMax");
94  numberBin_ = pset.getParameter<int>("numberBin");
95 
96  sigmaIetaMin_ = pset.getParameter<double>("sigmaIetaMin");
97  sigmaIetaMax_ = pset.getParameter<double>("sigmaIetaMax");
98  sigmaIetaBin_ = pset.getParameter<int>("sigmaIetaBin");
99 
100  reducedEtBin_ = etBin_/4;
103  reducedR9Bin_ = r9Bin_/4;
104 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::SortedCollection< EcalRecHit, edm::StrictWeakOrdering< EcalRecHit > > > barrelRecHit_token_
edm::EDGetTokenT< reco::PFCandidateCollection > pfCandidates_
edm::EDGetTokenT< std::vector< reco::Photon > > photon_token_
edm::EDGetTokenT< edm::SortedCollection< EcalRecHit, edm::StrictWeakOrdering< EcalRecHit > > > endcapRecHit_token_
edm::EDGetTokenT< trigger::TriggerEvent > triggerEvent_token_
edm::EDGetTokenT< reco::VertexCollection > offline_pvToken_
edm::EDGetTokenT< edm::ValueMap< std::vector< reco::PFCandidateRef > > > photonIsoValmap_token_
edm::EDGetTokenT< reco::BeamSpot > beamSpot_token_
edm::EDGetTokenT< std::vector< reco::Muon > > muon_token_
ZToMuMuGammaAnalyzer::~ZToMuMuGammaAnalyzer ( )
override

Definition at line 106 of file ZToMuMuGammaAnalyzer.cc.

107 {
108 }

Member Function Documentation

void ZToMuMuGammaAnalyzer::analyze ( const edm::Event e,
const edm::EventSetup esup 
)
override

uncleaned

Definition at line 421 of file ZToMuMuGammaAnalyzer.cc.

References begin, deltaR(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, reco::PFCandidate::e, reco::MuonIsolation::emEt, edm::ValueMap< T >::end(), HcalObjRepresent::Fill(), trigger::TriggerEvent::filterKeys(), trigger::TriggerEvent::filterTag(), reco::PFCandidate::gamma, edm::Event::getByToken(), reco::PFCandidate::h, reco::PFCandidate::h0, reco::MuonIsolation::hadEt, mps_fire::i, edm::EventBase::id(), reco::Muon::isolationR03(), edm::HandleBase::isValid(), edm::InputTag::label(), diffTwoXMLs::label, HiRecoMuon_cff::muonCollection, softMuonTagInfos_cfi::muonSelection, correctedPhotons_cfi::photonCollection, edm::Handle< T >::product(), reco::LeafCandidate::pt(), trigger::TriggerEvent::sizeFilters(), createPayload::skip, reco::MuonIsolation::sumPt, and parallelization::uint().

422 {
423  using namespace edm;
424 
425  if (nEvt_% prescaleFactor_ ) return;
426  nEvt_++;
427  LogInfo("ZToMuMuGammaAnalyzer") << "ZToMuMuGammaAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n";
428 
429  // Get the trigger results
430  bool validTriggerEvent=true;
431  edm::Handle<trigger::TriggerEvent> triggerEventHandle;
432  trigger::TriggerEvent triggerEvent;
433  e.getByToken(triggerEvent_token_,triggerEventHandle);
434  if(!triggerEventHandle.isValid()) {
435  edm::LogInfo("PhotonAnalyzer") << "Error! Can't get the product: triggerEvent_token_" << endl;
436  validTriggerEvent=false;
437  }
438  if(validTriggerEvent) triggerEvent = *(triggerEventHandle.product());
439 
440  // Get the reconstructed photons
441  // bool validPhotons=true;
442  Handle<reco::PhotonCollection> photonHandle;
444  e.getByToken(photon_token_ , photonHandle);
445  if ( !photonHandle.isValid()) {
446  edm::LogInfo("ZToMuMuGammaAnalyzer") << "Error! Can't get the product: photon_token_" << endl;
447  //validPhotons=false;
448  }
449  // if(validPhotons) photonCollection = *(photonHandle.product());
450 
451  // Get the PF refined cluster collection
452  Handle<reco::PFCandidateCollection> pfCandidateHandle;
453  e.getByToken(pfCandidates_,pfCandidateHandle);
454  if (!pfCandidateHandle.isValid()) {
455  edm::LogError("PhotonValidator") << "Error! Can't get the product pfCandidates "<< std::endl ;
456  }
457 
458  edm::Handle<edm::ValueMap<std::vector<reco::PFCandidateRef> > > phoToParticleBasedIsoMapHandle;
459  edm::ValueMap<std::vector<reco::PFCandidateRef> > phoToParticleBasedIsoMap;
460  if ( fName_ == "zmumugammaGedValidation") {
461  e.getByToken(photonIsoValmap_token_,phoToParticleBasedIsoMapHandle);
462  // e.getByLabel("particleBasedIsolation",valueMapPhoPFCandIso_,phoToParticleBasedIsoMapHandle);
463  if ( ! phoToParticleBasedIsoMapHandle.isValid()) {
464  edm::LogInfo("PhotonValidator") << "Error! Can't get the product: valueMap photons to particle based iso " << std::endl;
465 
466  }
467  phoToParticleBasedIsoMap = *(phoToParticleBasedIsoMapHandle.product());
468  }
469 
470  // Get the reconstructed muons
471  bool validMuons=true;
472  Handle<reco::MuonCollection> muonHandle;
474  e.getByToken(muon_token_, muonHandle);
475  if ( !muonHandle.isValid()) {
476  edm::LogInfo("ZToMuMuGammaAnalyzer") << "Error! Can't get the product: muon_token_" << endl;
477  validMuons=false;
478  }
479  if(validMuons) muonCollection = *(muonHandle.product());
480 
481  // Get the beam spot
483  e.getByToken(beamSpot_token_, bsHandle);
484  if (!bsHandle.isValid()) {
485  edm::LogError("TrackerOnlyConversionProducer") << "Error! Can't get the product primary Vertex Collection "<< "\n";
486  return;
487  }
488  const reco::BeamSpot &thebs = *bsHandle.product();
489 
490  //Prepare list of photon-related HLT filter names
491  vector<int> Keys;
492  for(uint filterIndex=0;filterIndex<triggerEvent.sizeFilters();++filterIndex){ //loop over all trigger filters in event (i.e. filters passed)
493  string label = triggerEvent.filterTag(filterIndex).label();
494  if(label.find( "Photon" ) != string::npos ) { //get photon-related filters
495  for(uint filterKeyIndex=0;filterKeyIndex<triggerEvent.filterKeys(filterIndex).size();++filterKeyIndex){ //loop over keys to objects passing this filter
496  Keys.push_back(triggerEvent.filterKeys(filterIndex)[filterKeyIndex]); //add keys to a vector for later reference
497  }
498  }
499  }
500 
501  // sort Keys vector in ascending order
502  // and erases duplicate entries from the vector
503  sort(Keys.begin(),Keys.end());
504  for ( uint i=0; i<Keys.size(); ) {
505  if (i!=(Keys.size()-1)) {
506  if (Keys[i]==Keys[i+1]) {
507  Keys.erase(Keys.begin()+i+1);
508  } else {
509  ++i;
510  }
511  } else {
512  ++i;
513  }
514  }
515 
517  e.getByToken(offline_pvToken_, vtxH);
518  h_nRecoVtx_ ->Fill (float(vtxH->size()));
519 
520  //photon counters
521  int nPho = 0;
522  int nPhoBarrel = 0;
523  int nPhoEndcap = 0;
524 
526  if ( muonCollection.size() < 2 ) return;
527 
528  for( reco::MuonCollection::const_iterator iMu = muonCollection.begin(); iMu != muonCollection.end(); iMu++) {
529  if ( !basicMuonSelection (*iMu) ) continue;
530 
531  for( reco::MuonCollection::const_iterator iMu2 = iMu+1; iMu2 != muonCollection.end(); iMu2++) {
532  if ( !basicMuonSelection (*iMu2) ) continue;
533  if ( iMu->charge()*iMu2->charge() > 0) continue;
534 
535  if ( !muonSelection(*iMu,thebs) && !muonSelection(*iMu2,thebs) ) continue;
536 
537  float mumuMass = mumuInvMass(*iMu,*iMu2) ;
538  if ( mumuMass < minMumuInvMass_ || mumuMass > maxMumuInvMass_ ) continue;
539 
540  h1_mumuInvMass_[0] -> Fill (mumuMass);
541 
542  if ( photonHandle->empty() ) continue;
543 
544  reco::Muon nearMuon;
545  reco::Muon farMuon;
546  for(unsigned int iPho=0; iPho < photonHandle->size(); iPho++) {
547  reco::PhotonRef aPho(reco::PhotonRef(photonHandle, iPho));
548  //
549  double dr1 = deltaR((*iMu).eta(), aPho->eta(), (*iMu).phi(), aPho->phi());
550  double dr2 = deltaR((*iMu2).eta(), aPho->eta(), (*iMu2).phi(), aPho->phi());
551  double drNear = dr1;
552  if (dr1 < dr2) {
553  nearMuon =*iMu ; farMuon = *iMu2; drNear = dr1;
554  } else {
555  nearMuon = *iMu2; farMuon = *iMu; drNear = dr2;
556  }
557  //
558  if ( nearMuon.isolationR03().hadEt > nearMuonHcalIso_ ) continue;
559  if ( farMuon.isolationR03().sumPt > farMuonTrackIso_ ) continue;
560  if ( farMuon.isolationR03().emEt > farMuonEcalIso_ ) continue;
561  if ( farMuon.pt() < farMuonMinPt_ ) continue;
562  if ( drNear > nearMuonDr_) continue;
563  //
564  if ( !photonSelection (aPho) ) continue;
565  float mumuGammaMass = mumuGammaInvMass(*iMu,*iMu2,aPho) ;
566  if ( mumuGammaMass < minMumuGammaInvMass_ || mumuGammaMass > maxMumuGammaInvMass_ ) continue;
567  //
568  //counter: number of photons
569  int iDet=0;
570  if ( aPho->isEB() || aPho->isEE() ) {
571  nPho++;
572  }
573  if ( aPho->isEB() ) {
574  iDet=1;
575  nPhoBarrel++;
576  }
577  if ( aPho->isEE() ) {
578  iDet=2;
579  nPhoEndcap++;
580  }
581 
582  //PHOTON RELATED HISTOGRAMS
583  h1_mumuGammaInvMass_[0] ->Fill (mumuGammaMass);
584  h1_mumuGammaInvMass_[iDet] ->Fill (mumuGammaMass);
585  //ENERGY
586  h_phoE_[0] ->Fill (aPho->energy());
587  h_phoSigmaEoverE_[0] ->Fill( aPho->getCorrectedEnergyError(aPho->getCandidateP4type())/aPho->energy() );
588  h_phoEt_[0] ->Fill (aPho->et());
589  h_phoE_[iDet] ->Fill (aPho->energy());
590  h_phoSigmaEoverE_[iDet] ->Fill( aPho->getCorrectedEnergyError(aPho->getCandidateP4type())/aPho->energy() );
591  p_phoSigmaEoverEVsNVtx_[iDet] ->Fill( float(vtxH->size()), aPho->getCorrectedEnergyError(aPho->getCandidateP4type())/aPho->energy() );
592  h_phoEt_[iDet] ->Fill (aPho->et());
593  //GEOMETRICAL
594  h_phoEta_[0] ->Fill (aPho->eta());
595  h_phoPhi_[0] ->Fill (aPho->phi());
596  h_scEta_[0] ->Fill (aPho->superCluster()->eta());
597  h_scPhi_[0] ->Fill (aPho->superCluster()->phi());
598  h_phoEta_[iDet] ->Fill (aPho->eta());
599  h_phoPhi_[iDet] ->Fill (aPho->phi());
600  h_scEta_[iDet] ->Fill (aPho->superCluster()->eta());
601  h_scPhi_[iDet] ->Fill (aPho->superCluster()->phi());
602  //SHOWER SHAPE
603  h_r9_[0] ->Fill (aPho->r9());
604  h_e1x5_[0]->Fill(aPho->e1x5());
605  h_e2x5_[0]->Fill(aPho->e2x5());
606  h_r1x5_[0]->Fill(aPho->r1x5());
607  h_r2x5_[0]->Fill(aPho->r2x5());
608  h_phoSigmaIetaIeta_[0] ->Fill(aPho->sigmaIetaIeta());
609  //
610  h_r9_[iDet] ->Fill (aPho->r9());
611  h_e1x5_[iDet] ->Fill(aPho->e1x5());
612  h_e2x5_[iDet] ->Fill(aPho->e2x5());
613  h_r1x5_[iDet] ->Fill( aPho->r1x5());
614  h_r2x5_[iDet] ->Fill(aPho->r2x5());
615  h_phoSigmaIetaIeta_[iDet] ->Fill(aPho->sigmaIetaIeta());
616  //TRACK ISOLATION
617  h_nTrackIsolSolid_[0] ->Fill(aPho->nTrkSolidConeDR04());
618  h_nTrackIsolHollow_[0] ->Fill(aPho->nTrkHollowConeDR04());
619  h_trackPtSumSolid_[0] ->Fill(aPho->trkSumPtSolidConeDR04());
620  h_trackPtSumHollow_[0] ->Fill(aPho->trkSumPtSolidConeDR04());
621  h_nTrackIsolSolid_[iDet] ->Fill(aPho->nTrkSolidConeDR04());
622  h_nTrackIsolHollow_[iDet] ->Fill(aPho->nTrkHollowConeDR04());
623  h_trackPtSumSolid_[iDet] ->Fill(aPho->trkSumPtSolidConeDR04());
624  h_trackPtSumHollow_[iDet] ->Fill(aPho->trkSumPtSolidConeDR04());
625  //CALORIMETER ISOLATION
626  h_ecalSum_[0] ->Fill(aPho->ecalRecHitSumEtConeDR04());
627  h_hcalSum_[0] ->Fill(aPho->hcalTowerSumEtConeDR04());
628  h_hOverE_[0] ->Fill(aPho->hadTowOverEm());
629  h_h1OverE_[0] ->Fill(aPho->hadTowDepth1OverEm());
630  h_h2OverE_[0] ->Fill(aPho->hadTowDepth2OverEm());
631  h_newhOverE_[0]->Fill( aPho->hadTowOverEm());
632  h_ecalSum_[iDet] ->Fill(aPho->ecalRecHitSumEtConeDR04());
633  h_hcalSum_[iDet] ->Fill(aPho->hcalTowerSumEtConeDR04());
634  h_hOverE_[iDet] ->Fill(aPho->hadTowOverEm());
635  h_h1OverE_[iDet] ->Fill(aPho->hadTowDepth1OverEm());
636  h_h2OverE_[iDet] ->Fill(aPho->hadTowDepth2OverEm());
637  h_newhOverE_[iDet]->Fill( aPho->hadTowOverEm());
638  // Isolation from particle flow
639  h_chHadIso_[0]-> Fill (aPho->chargedHadronIso());
640  h_nHadIso_[0]-> Fill (aPho->neutralHadronIso());
641  h_phoIso_[0]-> Fill (aPho->photonIso());
642  h_nCluOutsideMustache_[0]->Fill(float(aPho->nClusterOutsideMustache()));
643  h_etOutsideMustache_[0]->Fill(aPho->etOutsideMustache());
644  h_pfMva_[0]->Fill(aPho->pfMVA());
645  h_chHadIso_[iDet]-> Fill (aPho->chargedHadronIso());
646  h_nHadIso_[iDet]-> Fill (aPho->neutralHadronIso());
647  h_phoIso_[iDet]-> Fill (aPho->photonIso());
648  h_nCluOutsideMustache_[iDet]->Fill(float(aPho->nClusterOutsideMustache()));
649  h_etOutsideMustache_[iDet]->Fill(aPho->etOutsideMustache());
650  h_pfMva_[iDet]->Fill(aPho->pfMVA());
651 
653  if ( fName_ == "zmumugammaGedValidation") {
654 
655  float SumPtIsoValCh = 0.;
656  float SumPtIsoValNh = 0.;
657  float SumPtIsoValPh = 0.;
658 
659  float SumPtIsoValCleanCh = 0.;
660  float SumPtIsoValCleanNh = 0.;
661  float SumPtIsoValCleanPh = 0.;
662 
663  for(unsigned int lCand=0; lCand < pfCandidateHandle->size(); lCand++) {
664  reco::PFCandidateRef pfCandRef(reco::PFCandidateRef(pfCandidateHandle,lCand));
665  float dR= deltaR(aPho->eta(), aPho->phi(),pfCandRef->eta(), pfCandRef->phi());
666  if ( dR<0.4) {
668  reco::PFCandidate::ParticleType type = pfCandRef->particleId();
669  if ( type == reco::PFCandidate::e ) continue;
670  if ( type == reco::PFCandidate::gamma && pfCandRef->mva_nothing_gamma() > 0.) continue;
671 
672  if( type == reco::PFCandidate::h ) {
673  SumPtIsoValCh += pfCandRef->pt();
676  }
677  if( type == reco::PFCandidate::h0 ) {
678  SumPtIsoValNh += pfCandRef->pt();
681  }
682  if( type == reco::PFCandidate::gamma ) {
683  SumPtIsoValPh += pfCandRef->pt();
686  }
688  bool skip=false;
689  for( std::vector<reco::PFCandidateRef>::const_iterator i = phoToParticleBasedIsoMap[aPho].begin(); i != phoToParticleBasedIsoMap[aPho].end(); ++i ) {
690  // std::cout << " PhotonValidator PfCand pt " << pfCandRef->pt() << " id " <<pfCandRef->particleId() << " and in the map " << (*i)->pt() << " type " << (*i)->particleId() << std::endl;
691  if ( (*i) == pfCandRef ) {
692  skip=true;
693  }
694  } // loop over the PFCandidates flagged as overlapping with the photon
695 
696  if ( skip ) continue;
697  if( type == reco::PFCandidate::h ) {
698  SumPtIsoValCleanCh += pfCandRef->pt();
701  }
702  if( type == reco::PFCandidate::h0 ) {
703  SumPtIsoValCleanNh += pfCandRef->pt();
706  }
707  if( type == reco::PFCandidate::gamma ) {
708  SumPtIsoValCleanPh += pfCandRef->pt();
710  h_dRPhoPFcand_Pho_Cleaned_[iDet]->Fill(dR);
711  }
712  } // dr=0.4
713  } // loop over all PF Candidates
714 
715  h_SumPtOverPhoPt_ChHad_Cleaned_[0]->Fill(SumPtIsoValCleanCh/aPho->pt());
716  h_SumPtOverPhoPt_NeuHad_Cleaned_[0]->Fill(SumPtIsoValCleanNh/aPho->pt());
717  h_SumPtOverPhoPt_Pho_Cleaned_[0]->Fill(SumPtIsoValCleanPh/aPho->pt());
718  h_SumPtOverPhoPt_ChHad_unCleaned_[0]->Fill(SumPtIsoValCh/aPho->pt());
719  h_SumPtOverPhoPt_NeuHad_unCleaned_[0]->Fill(SumPtIsoValNh/aPho->pt());
720  h_SumPtOverPhoPt_Pho_unCleaned_[0]->Fill(SumPtIsoValPh/aPho->pt());
721  //
722  h_SumPtOverPhoPt_ChHad_Cleaned_[iDet]->Fill(SumPtIsoValCleanCh/aPho->pt());
723  h_SumPtOverPhoPt_NeuHad_Cleaned_[iDet]->Fill(SumPtIsoValCleanNh/aPho->pt());
724  h_SumPtOverPhoPt_Pho_Cleaned_[iDet]->Fill(SumPtIsoValCleanPh/aPho->pt());
725  h_SumPtOverPhoPt_ChHad_unCleaned_[iDet]->Fill(SumPtIsoValCh/aPho->pt());
726  h_SumPtOverPhoPt_NeuHad_unCleaned_[iDet]->Fill(SumPtIsoValNh/aPho->pt());
727  h_SumPtOverPhoPt_Pho_unCleaned_[iDet]->Fill(SumPtIsoValPh/aPho->pt());
728  } // only for zmumugammaGedValidation
729 
730  if ( makeProfiles_ ) {
731  p_r9VsEt_[iDet] ->Fill (aPho->et(),aPho->r9());
732  p_r9VsEta_[0]->Fill (aPho->eta(),aPho->r9());
733  p_e1x5VsEt_[iDet] ->Fill(aPho->et(), aPho->e1x5());
734  p_e1x5VsEta_[0]->Fill(aPho->eta(),aPho->e1x5());
735  p_e2x5VsEt_[iDet] ->Fill(aPho->et(), aPho->e2x5());
736  p_e2x5VsEta_[0]->Fill(aPho->eta(),aPho->e2x5());
737  p_r1x5VsEt_[iDet] ->Fill(aPho->et(), aPho->r1x5());
738  p_r1x5VsEta_[0]->Fill(aPho->eta(),aPho->r1x5());
739  p_r2x5VsEt_[iDet] ->Fill(aPho->et(), aPho->r2x5());
740  p_r2x5VsEta_[0]->Fill(aPho->eta(),aPho->r2x5());
741  //
742  p_sigmaIetaIetaVsEta_[0] ->Fill(aPho->eta(),aPho->sigmaIetaIeta());
743  p_nTrackIsolSolidVsEt_[iDet] ->Fill(aPho->et(), aPho->nTrkSolidConeDR04());
744  p_nTrackIsolSolidVsEta_[0] ->Fill(aPho->eta(),aPho->nTrkSolidConeDR04());
745  p_nTrackIsolHollowVsEt_[iDet] ->Fill(aPho->et(), aPho->nTrkHollowConeDR04());
746  p_nTrackIsolHollowVsEta_[0] ->Fill(aPho->eta(),aPho->nTrkHollowConeDR04());
747  p_trackPtSumSolidVsEt_[iDet] ->Fill(aPho->et(), aPho->trkSumPtSolidConeDR04());
748  p_trackPtSumSolidVsEta_[0] ->Fill(aPho->eta(),aPho->trkSumPtSolidConeDR04());
749  p_trackPtSumHollowVsEt_[iDet] ->Fill(aPho->et(), aPho->trkSumPtHollowConeDR04());
750  p_trackPtSumHollowVsEta_[0] ->Fill(aPho->eta(),aPho->trkSumPtHollowConeDR04());
751  //
752  p_ecalSumVsEt_[iDet] ->Fill(aPho->et(), aPho->ecalRecHitSumEtConeDR04());
753  p_ecalSumVsEta_[0] ->Fill(aPho->eta(),aPho->ecalRecHitSumEtConeDR04());
754  p_hcalSumVsEt_[iDet] ->Fill(aPho->et(), aPho->hcalTowerSumEtConeDR04());
755  p_hcalSumVsEta_[0] ->Fill(aPho->eta(),aPho->hcalTowerSumEtConeDR04());
756  p_hOverEVsEt_[iDet] ->Fill(aPho->et(), aPho->hadTowOverEm());
757  p_hOverEVsEta_[0] ->Fill(aPho->eta(),aPho->hadTowOverEm());
758 
759  }
760 
762  if(use2DHistos_){
763  //SHOWER SHAPE
764  h2_r9VsEt_[iDet] ->Fill (aPho->et(),aPho->r9());
765  h2_r9VsEta_[0]->Fill (aPho->eta(),aPho->r9());
766  h2_e1x5VsEt_[iDet] ->Fill(aPho->et(), aPho->e1x5());
767  h2_e1x5VsEta_[0]->Fill(aPho->eta(),aPho->e1x5());
768  h2_e2x5VsEta_[0]->Fill(aPho->eta(),aPho->e2x5());
769  h2_e2x5VsEt_[iDet] ->Fill(aPho->et(), aPho->e2x5());
770  h2_r1x5VsEta_[0]->Fill(aPho->eta(),aPho->r1x5());
771  h2_r1x5VsEt_[iDet] ->Fill(aPho->et(), aPho->r1x5());
772  h2_r2x5VsEt_[iDet] ->Fill(aPho->et(), aPho->r2x5());
773  h2_r2x5VsEta_[0]->Fill(aPho->eta(),aPho->r2x5());
774  h2_sigmaIetaIetaVsEta_[0] ->Fill(aPho->eta(),aPho->sigmaIetaIeta());
775  //TRACK ISOLATION
776  h2_nTrackIsolSolidVsEt_[0] ->Fill(aPho->et(), aPho->nTrkSolidConeDR04());
777  h2_nTrackIsolSolidVsEta_[0] ->Fill(aPho->eta(),aPho->nTrkSolidConeDR04());
778  h2_nTrackIsolHollowVsEt_[0] ->Fill(aPho->et(), aPho->nTrkHollowConeDR04());
779  h2_nTrackIsolHollowVsEta_[0] ->Fill(aPho->eta(),aPho->nTrkHollowConeDR04());
780  h2_trackPtSumSolidVsEt_[0] ->Fill(aPho->et(), aPho->trkSumPtSolidConeDR04());
781  h2_trackPtSumSolidVsEta_[0] ->Fill(aPho->eta(),aPho->trkSumPtSolidConeDR04());
782  h2_trackPtSumHollowVsEt_[0] ->Fill(aPho->et(), aPho->trkSumPtHollowConeDR04());
783  h2_trackPtSumHollowVsEta_[0] ->Fill(aPho->eta(),aPho->trkSumPtHollowConeDR04());
784  //CALORIMETER ISOLATION
785  h2_ecalSumVsEt_[iDet] ->Fill(aPho->et(), aPho->ecalRecHitSumEtConeDR04());
786  h2_ecalSumVsEta_[0] ->Fill(aPho->eta(),aPho->ecalRecHitSumEtConeDR04());
787  h2_hcalSumVsEt_[iDet] ->Fill(aPho->et(), aPho->hcalTowerSumEtConeDR04());
788  h2_hcalSumVsEta_[0] ->Fill(aPho->eta(),aPho->hcalTowerSumEtConeDR04());
789  }
790  } //end photon loop
791 
792  h_nPho_[0] ->Fill (float(nPho));
793  h_nPho_[1] ->Fill (float(nPhoBarrel));
794  h_nPho_[2] ->Fill (float(nPhoEndcap));
795  } //end inner muon loop
796  } //end outer muon loop
797 }//End of Analyze method
MonitorElement * h_dRPhoPFcand_Pho_unCleaned_[3]
type
Definition: HCALResponse.h:21
MonitorElement * p_hOverEVsEta_[3]
float hadEt
hcal sum-Et
Definition: MuonIsolation.h:9
MonitorElement * h_pfMva_[3]
MonitorElement * h_SumPtOverPhoPt_ChHad_unCleaned_[3]
MonitorElement * h_dRPhoPFcand_Pho_Cleaned_[3]
MonitorElement * h2_nTrackIsolSolidVsEt_[3]
ParticleType
particle types
Definition: PFCandidate.h:45
MonitorElement * p_sigmaIetaIetaVsEta_[3]
MonitorElement * p_r1x5VsEta_[3]
const_iterator end() const
Definition: ValueMap.h:209
float sumPt
sum-pt of tracks
Definition: MuonIsolation.h:7
bool muonSelection(const reco::Muon &m, const reco::BeamSpot &bs)
MonitorElement * h_phoSigmaIetaIeta_[3]
The single EDProduct to be saved for each event (AOD case)
Definition: TriggerEvent.h:25
trigger::size_type sizeFilters() const
Definition: TriggerEvent.h:135
MonitorElement * h2_r2x5VsEta_[3]
MonitorElement * h_newhOverE_[3]
edm::EDGetTokenT< reco::PFCandidateCollection > pfCandidates_
MonitorElement * h2_trackPtSumHollowVsEt_[3]
MonitorElement * h_SumPtOverPhoPt_NeuHad_unCleaned_[3]
edm::EDGetTokenT< std::vector< reco::Photon > > photon_token_
MonitorElement * p_nTrackIsolHollowVsEta_[3]
MonitorElement * h2_r9VsEt_[3]
MonitorElement * h_nCluOutsideMustache_[3]
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
MonitorElement * h2_ecalSumVsEt_[3]
MonitorElement * h_nTrackIsolSolid_[3]
MonitorElement * h_SumPtOverPhoPt_Pho_Cleaned_[3]
MonitorElement * h_trackPtSumHollow_[3]
MonitorElement * h2_nTrackIsolSolidVsEta_[3]
MonitorElement * p_ecalSumVsEt_[3]
MonitorElement * h_phoEta_[3]
MonitorElement * h_nPho_[3]
MonitorElement * h_chHadIso_[3]
const Keys & filterKeys(trigger::size_type index) const
Definition: TriggerEvent.h:111
MonitorElement * h_phoSigmaEoverE_[3]
MonitorElement * h_dRPhoPFcand_NeuHad_Cleaned_[3]
MonitorElement * h_r1x5_[3]
MonitorElement * p_phoSigmaEoverEVsNVtx_[3]
MonitorElement * h2_sigmaIetaIetaVsEta_[3]
MonitorElement * h_dRPhoPFcand_NeuHad_unCleaned_[3]
double pt() const final
transverse momentum
MonitorElement * h2_hcalSumVsEt_[3]
MonitorElement * h_hcalSum_[3]
MonitorElement * h_phoPhi_[3]
MonitorElement * p_trackPtSumHollowVsEt_[3]
MonitorElement * h2_ecalSumVsEta_[3]
edm::EDGetTokenT< trigger::TriggerEvent > triggerEvent_token_
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
MonitorElement * h_hOverE_[3]
void Fill(long long x)
MonitorElement * h_trackPtSumSolid_[3]
MonitorElement * p_hcalSumVsEt_[3]
MonitorElement * p_e1x5VsEt_[3]
MonitorElement * h2_r2x5VsEt_[3]
MonitorElement * h_h1OverE_[3]
MonitorElement * h_SumPtOverPhoPt_ChHad_Cleaned_[3]
MonitorElement * h2_e1x5VsEta_[3]
MonitorElement * h1_mumuGammaInvMass_[3]
MonitorElement * h2_trackPtSumSolidVsEt_[3]
MonitorElement * h_phoE_[3]
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
MonitorElement * p_r2x5VsEta_[3]
MonitorElement * h_nRecoVtx_
float emEt
ecal sum-Et
Definition: MuonIsolation.h:8
MonitorElement * p_e2x5VsEt_[3]
MonitorElement * p_hcalSumVsEta_[3]
float mumuInvMass(const reco::Muon &m1, const reco::Muon &m2)
float mumuGammaInvMass(const reco::Muon &mu1, const reco::Muon &mu2, const reco::PhotonRef &pho)
MonitorElement * h_phoEt_[3]
MonitorElement * p_ecalSumVsEta_[3]
bool isValid() const
Definition: HandleBase.h:74
MonitorElement * p_r9VsEta_[3]
MonitorElement * h_nHadIso_[3]
edm::EDGetTokenT< reco::VertexCollection > offline_pvToken_
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
MonitorElement * p_trackPtSumSolidVsEta_[3]
MonitorElement * h1_mumuInvMass_[3]
photon histos
MonitorElement * h2_trackPtSumSolidVsEta_[3]
const edm::InputTag filterTag(trigger::size_type index) const
Definition: TriggerEvent.h:103
MonitorElement * h2_hcalSumVsEta_[3]
MonitorElement * h_dRPhoPFcand_ChHad_unCleaned_[3]
bool photonSelection(const reco::PhotonRef &p)
MonitorElement * h_scPhi_[3]
MonitorElement * h2_e1x5VsEt_[3]
T const * product() const
Definition: Handle.h:81
MonitorElement * h_scEta_[3]
edm::EDGetTokenT< edm::ValueMap< std::vector< reco::PFCandidateRef > > > photonIsoValmap_token_
std::vector< size_type > Keys
std::vector< Photon > PhotonCollection
collectin of Photon objects
Definition: PhotonFwd.h:9
MonitorElement * h_nTrackIsolHollow_[3]
MonitorElement * p_r2x5VsEt_[3]
def uint(string)
MonitorElement * h_SumPtOverPhoPt_Pho_unCleaned_[3]
MonitorElement * h_SumPtOverPhoPt_NeuHad_Cleaned_[3]
MonitorElement * h_etOutsideMustache_[3]
MonitorElement * h_e1x5_[3]
MonitorElement * h2_nTrackIsolHollowVsEt_[3]
MonitorElement * h2_nTrackIsolHollowVsEta_[3]
edm::EDGetTokenT< reco::BeamSpot > beamSpot_token_
MonitorElement * h2_e2x5VsEta_[3]
std::string const & label() const
Definition: InputTag.h:36
MonitorElement * h2_trackPtSumHollowVsEta_[3]
edm::EventID id() const
Definition: EventBase.h:60
MonitorElement * h_r9_[3]
#define begin
Definition: vmac.h:32
HLT enums.
MonitorElement * p_r1x5VsEt_[3]
MonitorElement * p_nTrackIsolHollowVsEt_[3]
MonitorElement * p_r9VsEt_[3]
bool basicMuonSelection(const reco::Muon &m)
MonitorElement * p_trackPtSumSolidVsEt_[3]
MonitorElement * p_e2x5VsEta_[3]
MonitorElement * h_e2x5_[3]
MonitorElement * h2_r1x5VsEta_[3]
MonitorElement * p_trackPtSumHollowVsEta_[3]
MonitorElement * h_dRPhoPFcand_ChHad_Cleaned_[3]
MonitorElement * p_hOverEVsEt_[3]
MonitorElement * h2_r1x5VsEt_[3]
MonitorElement * p_nTrackIsolSolidVsEt_[3]
MonitorElement * p_e1x5VsEta_[3]
MonitorElement * p_nTrackIsolSolidVsEta_[3]
MonitorElement * h2_r9VsEta_[3]
const MuonIsolation & isolationR03() const
Definition: Muon.h:162
edm::EDGetTokenT< std::vector< reco::Muon > > muon_token_
MonitorElement * h_ecalSum_[3]
MonitorElement * h_phoIso_[3]
MonitorElement * h_h2OverE_[3]
MonitorElement * h_r2x5_[3]
MonitorElement * h2_e2x5VsEt_[3]
bool ZToMuMuGammaAnalyzer::basicMuonSelection ( const reco::Muon m)
private

Definition at line 799 of file ZToMuMuGammaAnalyzer.cc.

References reco::LeafCandidate::eta(), reco::Muon::globalTrack(), reco::Muon::innerTrack(), reco::Muon::isGlobalMuon(), edm::Ref< C, T, F >::isNonnull(), reco::LeafCandidate::pt(), and mps_fire::result.

799  {
800  bool result=true;
801  if (!mu.innerTrack().isNonnull()) result=false;
802  if (!mu.globalTrack().isNonnull()) result=false;
803  if ( !mu.isGlobalMuon() ) result=false;
804  if ( mu.pt() < muonMinPt_ ) result=false;
805  if ( fabs(mu.eta())>2.4 ) result=false;
806 
807  int pixHits=0;
808  int tkHits=0;
809  if ( mu.innerTrack().isNonnull() ) {
810  pixHits=mu.innerTrack()->hitPattern().numberOfValidPixelHits();
811  tkHits=mu.innerTrack()->hitPattern().numberOfValidStripHits();
812  }
813 
814  if ( pixHits+tkHits < minPixStripHits_ ) result=false;
815 
816  return result;
817 }
const int mu
Definition: Constants.h:22
void ZToMuMuGammaAnalyzer::bookHistograms ( DQMStore::IBooker iBooker,
edm::Run const &  ,
edm::EventSetup const &   
)
overridevirtual

Implements DQMEDAnalyzer.

Definition at line 110 of file ZToMuMuGammaAnalyzer.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::book2D(), DQMStore::IBooker::bookProfile(), and DQMStore::IBooker::setCurrentFolder().

113 {
115  iBooker.setCurrentFolder("Egamma/"+fName_+"/ZToMuMuGamma");
116 
117  h1_mumuInvMass_[0] = iBooker.book1D("mumuInvMass","Two muon invariant mass: M (GeV)",etBin_,etMin_,etMax_);
118  h1_mumuGammaInvMass_[0] = iBooker.book1D("mumuGammaInvMass","Two-muon plus gamma invariant mass: M (GeV)",etBin_,etMin_,etMax_);
119  h1_mumuGammaInvMass_[1] = iBooker.book1D("mumuGammaInvMassBarrel","Two-muon plus gamma invariant mass: M (GeV)",etBin_,etMin_,etMax_);
120  h1_mumuGammaInvMass_[2] = iBooker.book1D("mumuGammaInvMassEndcap","Two-muon plus gamma invariant mass: M (GeV)",etBin_,etMin_,etMax_);
121 
124  h_nRecoVtx_ = iBooker.book1D("nOfflineVtx","# of Offline Vertices",80, -0.5, 79.5);
125 
126  //ENERGY
127  h_phoE_[0] = iBooker.book1D("phoE","Energy;E (GeV)",eBin_,eMin_,eMax_);
128  h_phoSigmaEoverE_[0] = iBooker.book1D("phoSigmaEoverE","All Ecal: #sigma_{E}/E;#sigma_{E}/E",eBin_,eMin_,eMax_);
129  h_phoEt_[0] = iBooker.book1D("phoEt","E_{T};E_{T} (GeV)", etBin_,etMin_,etMax_);
130 
131  //NUMBER OF PHOTONS
132  h_nPho_[0] = iBooker.book1D("nPho", "Number of Photons per Event;# #gamma", numberBin_,numberMin_,numberMax_);
133 
134  //GEOMETRICAL
135  h_phoEta_[0] = iBooker.book1D("phoEta", "#eta;#eta",etaBin_,etaMin_,etaMax_);
136  h_phoPhi_[0] = iBooker.book1D("phoPhi", "#phi;#phi",phiBin_,phiMin_,phiMax_);
137 
138  h_scEta_[0] = iBooker.book1D("scEta", "SuperCluster #eta;#eta",etaBin_,etaMin_,etaMax_);
139  h_scPhi_[0] = iBooker.book1D("scPhi", "SuperCluster #phi;#phi",phiBin_,phiMin_,phiMax_);
140 
141  //SHOWER SHAPE
142  h_r9_[0] = iBooker.book1D("r9","R9;R9",r9Bin_,r9Min_, r9Max_);
143  h_e1x5_[0] = iBooker.book1D("e1x5","E1x5;E1X5 (GeV)",reducedEtBin_,etMin_,etMax_);
144  h_e2x5_[0] = iBooker.book1D("e2x5","E2x5;E2X5 (GeV)",reducedEtBin_,etMin_,etMax_);
145  h_r1x5_[0] = iBooker.book1D("r1x5","r1x5;r1X5 (GeV)",reducedEtBin_,etMin_,etMax_);
146  h_r2x5_[0] = iBooker.book1D("r2x5","r2x5;r2X5 (GeV)",reducedEtBin_,etMin_,etMax_);
147  h_phoSigmaIetaIeta_[0] = iBooker.book1D("phoSigmaIetaIeta","#sigma_{i#etai#eta};#sigma_{i#etai#eta}",sigmaIetaBin_,sigmaIetaMin_,sigmaIetaMax_);
148  //TRACK ISOLATION
149  h_nTrackIsolSolid_[0] = iBooker.book1D("nIsoTracksSolid","Number Of Tracks in the Solid Iso Cone;# tracks",numberBin_,numberMin_,numberMax_);
150  h_nTrackIsolHollow_[0] = iBooker.book1D("nIsoTracksHollow","Number Of Tracks in the Hollow Iso Cone;# tracks",numberBin_,numberMin_,numberMax_);
151  h_trackPtSumSolid_[0] = iBooker.book1D("isoPtSumSolid","Track P_{T} Sum in the Solid Iso Cone;P_{T} (GeV)",sumBin_,sumMin_,sumMax_);
152  h_trackPtSumHollow_[0] = iBooker.book1D("isoPtSumHollow","Track P_{T} Sum in the Hollow Iso Cone;P_{T} (GeV)",sumBin_,sumMin_,sumMax_);
153  //CALORIMETER ISOLATION VARIABLES
154  h_ecalSum_[0] = iBooker.book1D("ecalSum","Ecal Sum in the Iso Cone;E (GeV)",sumBin_,sumMin_,sumMax_);
155  h_hcalSum_[0] = iBooker.book1D("hcalSum","Hcal Sum in the Iso Cone;E (GeV)",sumBin_,sumMin_,sumMax_);
156  h_hOverE_[0] = iBooker.book1D("hOverE","H/E;H/E",hOverEBin_,hOverEMin_,hOverEMax_);
157  h_h1OverE_[0] = iBooker.book1D("h1OverE","H/E for Depth 1;H/E",hOverEBin_,hOverEMin_,hOverEMax_);
158  h_h2OverE_[0] = iBooker.book1D("h2OverE","H/E for Depth 2;H/E",hOverEBin_,hOverEMin_,hOverEMax_);
159  string histname = "newhOverE";
160  h_newhOverE_[0] = iBooker.book1D(histname+"All", "new H/E: All Ecal",100,0., 0.1) ;
161  // Information from Particle Flow
162  histname = "chargedHadIso";
163  h_chHadIso_[0]= iBooker.book1D(histname+"All", "PF chargedHadIso: All Ecal",etBin_,etMin_,20.);
164  histname = "neutralHadIso";
165  h_nHadIso_[0]= iBooker.book1D(histname+"All", "PF neutralHadIso: All Ecal",etBin_,etMin_,20.);
166  histname = "photonIso";
167  h_phoIso_[0]= iBooker.book1D(histname+"All", "PF photonIso: All Ecal",etBin_,etMin_,20.);
168  histname = "nCluOutMustache";
169  h_nCluOutsideMustache_[0]= iBooker.book1D(histname+"All", "PF number of clusters outside Mustache: All Ecal",50,0.,50.);
170  histname = "etOutMustache";
171  h_etOutsideMustache_[0]= iBooker.book1D(histname+"All", "PF et outside Mustache: All Ecal",etBin_,etMin_,20.);
172  histname = "pfMVA";
173  h_pfMva_[0]= iBooker.book1D(histname+"All", "PF MVA output: All Ecal",50,-1.,2.);
175  histname = "SumPtOverPhoPt_ChHad_Cleaned";
176  h_SumPtOverPhoPt_ChHad_Cleaned_[0]= iBooker.book1D(histname+"All", "Pf Cand Sum Pt Over photon pt Charged Hadrons: All Ecal",etBin_,etMin_,2.);
177  histname = "SumPtOverPhoPt_NeuHad_Cleaned";
178  h_SumPtOverPhoPt_NeuHad_Cleaned_[0]= iBooker.book1D(histname+"All", "Pf Cand Sum Pt Over photon pt Neutral Hadrons: All Ecal",etBin_,etMin_,2.);
179  histname = "SumPtOverPhoPt_Pho_Cleaned";
180  h_SumPtOverPhoPt_Pho_Cleaned_[0]= iBooker.book1D(histname+"All", "Pf Cand Sum Pt Over photon pt Photons Hadrons: All Ecal",etBin_,etMin_,2.);
181  histname = "dRPhoPFcand_ChHad_Cleaned";
182  h_dRPhoPFcand_ChHad_Cleaned_[0]= iBooker.book1D(histname+"All", "dR(pho,cand) Charged Hadrons : All Ecal",etBin_,etMin_,0.7);
183  histname = "dRPhoPFcand_NeuHad_Cleaned";
184  h_dRPhoPFcand_NeuHad_Cleaned_[0]= iBooker.book1D(histname+"All", "dR(pho,cand) Neutral Hadrons : All Ecal",etBin_,etMin_,0.7);
185  histname = "dRPhoPFcand_Pho_Cleaned";
186  h_dRPhoPFcand_Pho_Cleaned_[0]= iBooker.book1D(histname+"All", "dR(pho,cand) Photons : All Ecal",etBin_,etMin_,0.7);
187  //
188  histname = "SumPtOverPhoPt_ChHad_unCleaned";
189  h_SumPtOverPhoPt_ChHad_unCleaned_[0]= iBooker.book1D(histname+"All", "Pf Cand Sum Pt Over photon pt Charged Hadrons : All Ecal",etBin_,etMin_,2.);
190  histname = "SumPtOverPhoPt_NeuHad_unCleaned";
191  h_SumPtOverPhoPt_NeuHad_unCleaned_[0]= iBooker.book1D(histname+"All", "Pf Cand Sum Pt Over photon pt Neutral Hadrons : All Ecal",etBin_,etMin_,2.);
192  histname = "SumPtOverPhoPt_Pho_unCleaned";
193  h_SumPtOverPhoPt_Pho_unCleaned_[0]= iBooker.book1D(histname+"All", "Pf Cand Sum Pt Over photon pt Photons: All Ecal",etBin_,etMin_,2.);
194  histname = "dRPhoPFcand_ChHad_unCleaned";
195  h_dRPhoPFcand_ChHad_unCleaned_[0]= iBooker.book1D(histname+"All", "dR(pho,cand) Charged Hadrons : All Ecal",etBin_,etMin_,0.7);
196  histname = "dRPhoPFcand_NeuHad_unCleaned";
197  h_dRPhoPFcand_NeuHad_unCleaned_[0]= iBooker.book1D(histname+"All", "dR(pho,cand) Neutral Hadrons : All Ecal",etBin_,etMin_,0.7);
198  histname = "dRPhoPFcand_Pho_unCleaned";
199  h_dRPhoPFcand_Pho_unCleaned_[0]= iBooker.book1D(histname+"All", "dR(pho,cand) Photons: All Ecal",etBin_,etMin_,0.7);
200 
201  // NUMBER OF PHOTONS
202  h_nPho_[1] = iBooker.book1D("nPhoBarrel","Number of Photons per Event;# #gamma", numberBin_,numberMin_,numberMax_);
203  h_nPho_[2] = iBooker.book1D("nPhoEndcap","Number of Photons per Event;# #gamma", numberBin_,numberMin_,numberMax_);
204  //EB ENERGY
205  h_phoE_[1] = iBooker.book1D("phoEBarrel","Energy for Barrel;E (GeV)",eBin_,eMin_,eMax_);
206  h_phoSigmaEoverE_[1] = iBooker.book1D("phoSigmaEoverEBarrel","Barrel: #sigma_E/E;#sigma_{E}/E",eBin_,eMin_,eMax_);
207  h_phoEt_[1] = iBooker.book1D("phoEtBarrel","E_{T};E_{T} (GeV)", etBin_,etMin_,etMax_);
208  //EE ENERGY
209  h_phoEt_[2] = iBooker.book1D("phoEtEndcap","E_{T};E_{T} (GeV)", etBin_,etMin_,etMax_);
210  h_phoE_[2] = iBooker.book1D("phoEEndcap","Energy for Endcap;E (GeV)",eBin_,eMin_,eMax_);
211  h_phoSigmaEoverE_[2] = iBooker.book1D("phoSigmaEoverEEndcap","Endcap: #sigma_{E}/E;#sigma_{E}/E",eBin_,eMin_,eMax_);
212  //EB GEOMETRICAL
213  h_phoEta_[1] = iBooker.book1D("phoEtaBarrel","#eta;#eta",etaBin_,etaMin_,etaMax_);
214  h_phoPhi_[1] = iBooker.book1D("phoPhiBarrel","#phi;#phi",phiBin_,phiMin_,phiMax_);
215  h_scEta_[1] = iBooker.book1D("scEtaBarrel","SuperCluster #eta;#eta",etaBin_,etaMin_,etaMax_);
216  h_scPhi_[1] = iBooker.book1D("scPhiBarrel","SuperCluster #phi;#phi",phiBin_,phiMin_,phiMax_);
217  //EE GEOMETRICAL
218  h_phoEta_[2] = iBooker.book1D("phoEtaEndcap","#eta;#eta",etaBin_,etaMin_,etaMax_);
219  h_phoPhi_[2] = iBooker.book1D("phoPhiEndcap","#phi;#phi",phiBin_,phiMin_,phiMax_);
220  h_scEta_[2] = iBooker.book1D("scEtaEndcap","SuperCluster #eta;#eta",etaBin_,etaMin_,etaMax_);
221  h_scPhi_[2] = iBooker.book1D("scPhiEndcap","SuperCluster #phi;#phi",phiBin_,phiMin_,phiMax_);
222  //SHOWER SHAPES
223  h_r9_[1] = iBooker.book1D("r9Barrel","R9;R9",r9Bin_,r9Min_, r9Max_);
224  h_r9_[2] = iBooker.book1D("r9Endcap","R9;R9",r9Bin_,r9Min_, r9Max_);
225  h_e1x5_[1] = iBooker.book1D("e1x5Barrel","E1x5;E1X5 (GeV)",reducedEtBin_,etMin_,etMax_);
226  h_e1x5_[2] = iBooker.book1D("e1x5Endcap","E1x5;E1X5 (GeV)",reducedEtBin_,etMin_,etMax_);
227  h_e2x5_[1] = iBooker.book1D("e2x5Barrel","E2x5;E2X5 (GeV)",reducedEtBin_,etMin_,etMax_);
228  h_e2x5_[2] = iBooker.book1D("e2x5Endcap","E2x5;E2X5 (GeV)",reducedEtBin_,etMin_,etMax_);
229  h_r1x5_[1] = iBooker.book1D("r1x5Barrel","r1x5;r1X5 (GeV)",reducedEtBin_,etMin_,etMax_);
230  h_r1x5_[2] = iBooker.book1D("r1x5Endcap","r1x5;r1X5 (GeV)",reducedEtBin_,etMin_,etMax_);
231  h_r2x5_[1] = iBooker.book1D("r2x5Barrel","r2x5;r2X5 (GeV)",reducedEtBin_,etMin_,etMax_);
232  h_r2x5_[2] = iBooker.book1D("r2x5Endcap","r2x5;r2X5 (GeV)",reducedEtBin_,etMin_,etMax_);
233  h_phoSigmaIetaIeta_[1] = iBooker.book1D("phoSigmaIetaIetaBarrel","#sigma_{i#etai#eta};#sigma_{i#etai#eta}",sigmaIetaBin_,sigmaIetaMin_,sigmaIetaMax_);
234  h_phoSigmaIetaIeta_[2] = iBooker.book1D("phoSigmaIetaIetaEndcap","#sigma_{i#etai#eta};#sigma_{i#etai#eta}",sigmaIetaBin_,sigmaIetaMin_,sigmaIetaMax_);
235  // TRACK ISOLATION
236  h_nTrackIsolSolid_[1] = iBooker.book1D("nIsoTracksSolidBarrel","Number Of Tracks in the Solid Iso Cone;# tracks",numberBin_,numberMin_,numberMax_);
237  h_nTrackIsolSolid_[2] = iBooker.book1D("nIsoTracksSolidEndcap","Number Of Tracks in the Solid Iso Cone;# tracks",numberBin_,numberMin_,numberMax_);
238  h_nTrackIsolHollow_[1] = iBooker.book1D("nIsoTracksHollowBarrel","Number Of Tracks in the Hollow Iso Cone;# tracks",numberBin_,numberMin_,numberMax_);
239  h_nTrackIsolHollow_[2] = iBooker.book1D("nIsoTracksHollowEndcap","Number Of Tracks in the Hollow Iso Cone;# tracks",numberBin_,numberMin_,numberMax_);
240  h_trackPtSumSolid_[1] = iBooker.book1D("isoPtSumSolidBarrel","Track P_{T} Sum in the Solid Iso Cone;P_{T} (GeV)",sumBin_,sumMin_,sumMax_);
241  h_trackPtSumSolid_[2] = iBooker.book1D("isoPtSumSolidEndcap","Track P_{T} Sum in the Solid Iso Cone;P_{T} (GeV)",sumBin_,sumMin_,sumMax_);
242  h_trackPtSumHollow_[1] = iBooker.book1D("isoPtSumHollowBarrel","Track P_{T} Sum in the Hollow Iso Cone;P_{T} (GeV)",sumBin_,sumMin_,sumMax_);
243  h_trackPtSumHollow_[2] = iBooker.book1D("isoPtSumHollowEndcap","Track P_{T} Sum in the Hollow Iso Cone;P_{T} (GeV)",sumBin_,sumMin_,sumMax_);
244  // CALORIMETER ISOLATION VARIABLES
245  h_ecalSum_[1] = iBooker.book1D("ecalSumBarrel","Ecal Sum in the Iso Cone;E (GeV)",sumBin_,sumMin_,sumMax_);
246  h_ecalSum_[2] = iBooker.book1D("ecalSumEndcap","Ecal Sum in the Iso Cone;E (GeV)",sumBin_,sumMin_,sumMax_);
247  h_hcalSum_[1] = iBooker.book1D("hcalSumBarrel","Hcal Sum in the Iso Cone;E (GeV)",sumBin_,sumMin_,sumMax_);
248  h_hcalSum_[2] = iBooker.book1D("hcalSumEndcap","Hcal Sum in the Iso Cone;E (GeV)",sumBin_,sumMin_,sumMax_);
249  //H/E
250  // EB
251  h_hOverE_[1] = iBooker.book1D("hOverEBarrel","H/E;H/E",hOverEBin_,hOverEMin_,hOverEMax_);
252  h_h1OverE_[1] = iBooker.book1D("h1OverEBarrel","H/E for Depth 1;H/E",hOverEBin_,hOverEMin_,hOverEMax_);
253  h_h2OverE_[1] = iBooker.book1D("h2OverEBarrel","H/E for Depth 2;H/E",hOverEBin_,hOverEMin_,hOverEMax_);
254  histname = "newhOverE";
255  h_newhOverE_[1] = iBooker.book1D(histname+"Barrel", "new H/E: Barrel",100,0., 0.1) ;
256  //EE
257  h_hOverE_[2] = iBooker.book1D("hOverEEndcap","H/E;H/E",hOverEBin_,hOverEMin_,hOverEMax_);
258  h_h1OverE_[2] = iBooker.book1D("h1OverEEndcap","H/E for Depth 1;H/E",hOverEBin_,hOverEMin_,hOverEMax_);
259  h_h2OverE_[2] = iBooker.book1D("h2OverEEndcap","H/E for Depth 2;H/E",hOverEBin_,hOverEMin_,hOverEMax_);
260  histname = "newhOverE";
261  h_newhOverE_[2] = iBooker.book1D(histname+"Endcap", "new H/E: Endcap",100,0., 0.1) ;
262  // Information from Particle Flow
263  histname = "chargedHadIso";
264  h_chHadIso_[1]= iBooker.book1D(histname+"Barrel", "PF chargedHadIso: Barrel",etBin_,etMin_,20.);
265  h_chHadIso_[2]= iBooker.book1D(histname+"Endcap", "PF chargedHadIso: Endcap",etBin_,etMin_,20.);
266  histname = "neutralHadIso";
267  h_nHadIso_[1]= iBooker.book1D(histname+"Barrel", "PF neutralHadIso: Barrel",etBin_,etMin_,20.);
268  h_nHadIso_[2]= iBooker.book1D(histname+"Endcap", "PF neutralHadIso: Endcap",etBin_,etMin_,20.);
269  histname = "photonIso";
270  h_phoIso_[1]= iBooker.book1D(histname+"Barrel", "PF photonIso: Barrel",etBin_,etMin_,20.);
271  h_phoIso_[2]= iBooker.book1D(histname+"Endcap", "PF photonIso: Endcap",etBin_,etMin_,20.);
272  histname = "nCluOutMustache";
273  h_nCluOutsideMustache_[1]= iBooker.book1D(histname+"Barrel", "PF number of clusters outside Mustache: Barrel",50,0.,50.);
274  h_nCluOutsideMustache_[2]= iBooker.book1D(histname+"Endcap", "PF number of clusters outside Mustache: Endcap",50,0.,50.);
275  histname = "etOutMustache";
276  h_etOutsideMustache_[1]= iBooker.book1D(histname+"Barrel", "PF et outside Mustache: Barrel",etBin_,etMin_,20.);
277  h_etOutsideMustache_[2]= iBooker.book1D(histname+"Endcap", "PF et outside Mustache: Endcap",etBin_,etMin_,20.);
278  histname = "pfMVA";
279  h_pfMva_[1]= iBooker.book1D(histname+"Barrel", "PF MVA output: Barrel",50,-1.,2.);
280  h_pfMva_[2]= iBooker.book1D(histname+"Endcap", "PF MVA output: Endcap",50,-1,2.);
282  histname = "SumPtOverPhoPt_ChHad_Cleaned";
283  h_SumPtOverPhoPt_ChHad_Cleaned_[1]= iBooker.book1D(histname+"Barrel","PF Cand Sum Pt Over photon pt Charged Hadrons: Barrel",etBin_,etMin_,2.);
284  h_SumPtOverPhoPt_ChHad_Cleaned_[2]= iBooker.book1D(histname+"Endcap","PF Cand Sum Pt Over photon pt Charged Hadrons: Endcap",etBin_,etMin_,2.);
285  histname = "SumPtOverPhoPt_NeuHad_Cleaned";
286  h_SumPtOverPhoPt_NeuHad_Cleaned_[1]= iBooker.book1D(histname+"Barrel","PF Cand Sum Pt Over photon pt Neutral Hadrons: Barrel",etBin_,etMin_,2.);
287  h_SumPtOverPhoPt_NeuHad_Cleaned_[2]= iBooker.book1D(histname+"Endcap","PF Cand Sum Pt Over photon pt Neutral Hadrons: Endcap",etBin_,etMin_,2.);
288  histname = "SumPtOverPhoPt_Pho_Cleaned";
289  h_SumPtOverPhoPt_Pho_Cleaned_[1]= iBooker.book1D(histname+"Barrel","PF Cand Sum Pt Over photon pt Photons Hadrons: Barrel",etBin_,etMin_,2.);
290  h_SumPtOverPhoPt_Pho_Cleaned_[2]= iBooker.book1D(histname+"Endcap","PF Cand Sum Pt Over photon pt Photons Hadrons: Endcap",etBin_,etMin_,2.);
291  histname = "dRPhoPFcand_ChHad_Cleaned";
292  h_dRPhoPFcand_ChHad_Cleaned_[1]= iBooker.book1D(histname+"Barrel","dR(pho,cand) Charged Hadrons : Barrel",etBin_,etMin_,0.7);
293  h_dRPhoPFcand_ChHad_Cleaned_[2]= iBooker.book1D(histname+"Endcap","dR(pho,cand) Charged Hadrons : Endcap",etBin_,etMin_,0.7);
294  histname = "dRPhoPFcand_NeuHad_Cleaned";
295  h_dRPhoPFcand_NeuHad_Cleaned_[1]= iBooker.book1D(histname+"Barrel","dR(pho,cand) Neutral Hadrons : Barrel",etBin_,etMin_,0.7);
296  h_dRPhoPFcand_NeuHad_Cleaned_[2]= iBooker.book1D(histname+"Endcap","dR(pho,cand) Neutral Hadrons : Endcap",etBin_,etMin_,0.7);
297  histname = "dRPhoPFcand_Pho_Cleaned";
298  h_dRPhoPFcand_Pho_Cleaned_[1]= iBooker.book1D(histname+"Barrel","dR(pho,cand) Photons : Barrel",etBin_,etMin_,0.7);
299  h_dRPhoPFcand_Pho_Cleaned_[2]= iBooker.book1D(histname+"Endcap","dR(pho,cand) Photons : Endcap",etBin_,etMin_,0.7);
300  //
301  histname = "SumPtOverPhoPt_ChHad_unCleaned";
302  h_SumPtOverPhoPt_ChHad_unCleaned_[1]= iBooker.book1D(histname+"Barrel","PF Cand Sum Pt Over photon pt Charged Hadrons: Barrel",etBin_,etMin_,2.);
303  h_SumPtOverPhoPt_ChHad_unCleaned_[2]= iBooker.book1D(histname+"Endcap","PF Cand Sum Pt Over photon pt Charged Hadrons: Endcap",etBin_,etMin_,2.);
304  histname = "SumPtOverPhoPt_NeuHad_unCleaned";
305  h_SumPtOverPhoPt_NeuHad_unCleaned_[1]= iBooker.book1D(histname+"Barrel","PF Cand Sum Pt Over photon pt Neutral Hadrons: Barrel",etBin_,etMin_,2.);
306  h_SumPtOverPhoPt_NeuHad_unCleaned_[2]= iBooker.book1D(histname+"Endcap","PF Cand Sum Pt Over photon pt Neutral Hadrons: Endcap",etBin_,etMin_,2.);
307  histname = "SumPtOverPhoPt_Pho_unCleaned";
308  h_SumPtOverPhoPt_Pho_unCleaned_[1]= iBooker.book1D(histname+"Barrel","PF Cand Sum Pt Over photon pt Photons: Barrel",etBin_,etMin_,2.);
309  h_SumPtOverPhoPt_Pho_unCleaned_[2]= iBooker.book1D(histname+"Endcap","PF Cand Sum Pt Over photon pt Photons: Endcap",etBin_,etMin_,2.);
310  histname = "dRPhoPFcand_ChHad_unCleaned";
311  h_dRPhoPFcand_ChHad_unCleaned_[1]= iBooker.book1D(histname+"Barrel","dR(pho,cand) Charged Hadrons : Barrel",etBin_,etMin_,0.7);
312  h_dRPhoPFcand_ChHad_unCleaned_[2]= iBooker.book1D(histname+"Endcap","dR(pho,cand) Charged Hadrons : Endcap",etBin_,etMin_,0.7);
313  histname = "dRPhoPFcand_NeuHad_unCleaned";
314  h_dRPhoPFcand_NeuHad_unCleaned_[1]= iBooker.book1D(histname+"Barrel","dR(pho,cand) Neutral Hadrons : Barrel",etBin_,etMin_,0.7);
315  h_dRPhoPFcand_NeuHad_unCleaned_[2]= iBooker.book1D(histname+"Endcap","dR(pho,cand) Neutral Hadrons : Endcap",etBin_,etMin_,0.7);
316  histname = "dRPhoPFcand_Pho_unCleaned";
317  h_dRPhoPFcand_Pho_unCleaned_[1]= iBooker.book1D(histname+"Barrel","dR(pho,cand) Photons: Barrel",etBin_,etMin_,0.7);
318  h_dRPhoPFcand_Pho_unCleaned_[2]= iBooker.book1D(histname+"Endcap","dR(pho,cand) Photons: Endcap",etBin_,etMin_,0.7);
319 
321  if(makeProfiles_){
322  // p_r9VsEt_[0] = iBooker.bookProfile("r9VsEt","Avg R9 vs E_{T};E_{T} (GeV);R9",etBin_,etMin_,etMax_,r9Bin_,r9Min_,r9Max_);
323  p_r9VsEt_[1] = iBooker.bookProfile("r9VsEtBarrel","Avg R9 vs E_{T};E_{T} (GeV);R9",etBin_,etMin_,etMax_,r9Bin_,r9Min_,r9Max_);
324  p_r9VsEt_[2] = iBooker.bookProfile("r9VsEtEndcap","Avg R9 vs E_{T};E_{T} (GeV);R9",etBin_,etMin_,etMax_,r9Bin_,r9Min_,r9Max_);
325  p_r9VsEta_[0] = iBooker.bookProfile("r9VsEta","Avg R9 vs #eta;#eta;R9",etaBin_,etaMin_,etaMax_,r9Bin_,r9Min_,r9Max_);
326  //
327  p_sigmaIetaIetaVsEta_[0] = iBooker.bookProfile("sigmaIetaIetaVsEta","Avg #sigma_{i#etai#eta} vs #eta;#eta;#sigma_{i#etai#eta}",etaBin_,etaMin_,etaMax_,sigmaIetaBin_,sigmaIetaMin_,sigmaIetaMax_);
328  p_e1x5VsEt_[1] = iBooker.bookProfile("e1x5VsEtBarrel","Avg E1x5 vs E_{T};E_{T} (GeV);E1X5 (GeV)",etBin_,etMin_,etMax_,etBin_,etMin_,etMax_);
329  p_e1x5VsEt_[2] = iBooker.bookProfile("e1x5VsEtEndcap","Avg E1x5 vs E_{T};E_{T} (GeV);E1X5 (GeV)",etBin_,etMin_,etMax_,etBin_,etMin_,etMax_);
330  p_e1x5VsEta_[0] = iBooker.bookProfile("e1x5VsEta","Avg E1x5 vs #eta;#eta;E1X5 (GeV)",etaBin_,etaMin_,etaMax_,etBin_,etMin_,etMax_);
331  p_e2x5VsEt_[1] = iBooker.bookProfile("e2x5VsEtBarrel","Avg E2x5 vs E_{T};E_{T} (GeV);E2X5 (GeV)",etBin_,etMin_,etMax_,etBin_,etMin_,etMax_);
332  p_e2x5VsEt_[2] = iBooker.bookProfile("e2x5VsEtEndcap","Avg E2x5 vs E_{T};E_{T} (GeV);E2X5 (GeV)",etBin_,etMin_,etMax_,etBin_,etMin_,etMax_);
333  p_e2x5VsEta_[0] = iBooker.bookProfile("e2x5VsEta","Avg E2x5 vs #eta;#eta;E2X5 (GeV)",etaBin_,etaMin_,etaMax_,etBin_,etMin_,etMax_);
334  p_r1x5VsEt_[1] = iBooker.bookProfile("r1x5VsEtBarrel","Avg R1x5 vs E_{T};E_{T} (GeV);R1X5",etBin_,etMin_,etMax_,r9Bin_,r9Min_,r9Max_);
335  p_r1x5VsEt_[2] = iBooker.bookProfile("r1x5VsEtEndcap","Avg R1x5 vs E_{T};E_{T} (GeV);R1X5",etBin_,etMin_,etMax_,r9Bin_,r9Min_,r9Max_);
336  p_r1x5VsEta_[0] = iBooker.bookProfile("r1x5VsEta","Avg R1x5 vs #eta;#eta;R1X5",etaBin_,etaMin_,etaMax_,r9Bin_,r9Min_,r9Max_);
337  p_r2x5VsEt_[1] = iBooker.bookProfile("r2x5VsEtBarrel","Avg R2x5 vs E_{T};E_{T} (GeV);R2X5",etBin_,etMin_,etMax_,r9Bin_,r9Min_,r9Max_);
338  p_r2x5VsEt_[2] = iBooker.bookProfile("r2x5VsEtEndcap","Avg R2x5 vs E_{T};E_{T} (GeV);R2X5",etBin_,etMin_,etMax_,r9Bin_,r9Min_,r9Max_);
339  p_r2x5VsEta_[0] = iBooker.bookProfile("r2x5VsEta","Avg R2x5 vs #eta;#eta;R2X5",etaBin_,etaMin_,etaMax_,r9Bin_,r9Min_,r9Max_);
340  p_nTrackIsolSolidVsEt_[1] = iBooker.bookProfile("nIsoTracksSolidVsEtBarrel","Avg Number Of Tracks in the Solid Iso Cone vs E_{T};E_{T};# tracks",etBin_,etMin_,etMax_,numberBin_,numberMin_,numberMax_);
341  p_nTrackIsolSolidVsEt_[2] = iBooker.bookProfile("nIsoTracksSolidVsEtEndcap","Avg Number Of Tracks in the Solid Iso Cone vs E_{T};E_{T};# tracks",etBin_,etMin_,etMax_,numberBin_,numberMin_,numberMax_);
342  p_nTrackIsolSolidVsEta_[0] = iBooker.bookProfile("nIsoTracksSolidVsEta","Avg Number Of Tracks in the Solid Iso Cone vs #eta;#eta;# tracks",etaBin_,etaMin_, etaMax_,numberBin_,numberMin_,numberMax_);
343  p_nTrackIsolHollowVsEt_[1] = iBooker.bookProfile("nIsoTracksHollowVsEtBarrel","Avg Number Of Tracks in the Hollow Iso Cone vs E_{T};E_{T};# tracks",etBin_,etMin_,etMax_,numberBin_,numberMin_,numberMax_);
344  p_nTrackIsolHollowVsEt_[2] = iBooker.bookProfile("nIsoTracksHollowVsEtEndcap","Avg Number Of Tracks in the Hollow Iso Cone vs E_{T};E_{T};# tracks",etBin_,etMin_,etMax_,numberBin_,numberMin_,numberMax_);
345  p_nTrackIsolHollowVsEta_[0] = iBooker.bookProfile("nIsoTracksHollowVsEta","Avg Number Of Tracks in the Hollow Iso Cone vs #eta;#eta;# tracks",etaBin_,etaMin_, etaMax_,numberBin_,numberMin_,numberMax_);
346  p_trackPtSumSolidVsEt_[1] = iBooker.bookProfile("isoPtSumSolidVsEtBarrel","Avg Track P_{T} Sum in the Solid Iso Cone vs E_{T};E_{T} (GeV);P_{T} (GeV)",etBin_,etMin_,etMax_,sumBin_,sumMin_,sumMax_);
347  p_trackPtSumSolidVsEt_[2] = iBooker.bookProfile("isoPtSumSolidVsEtEndcap","Avg Track P_{T} Sum in the Solid Iso Cone vs E_{T};E_{T} (GeV);P_{T} (GeV)",etBin_,etMin_,etMax_,sumBin_,sumMin_,sumMax_);
348  p_trackPtSumSolidVsEta_[0] = iBooker.bookProfile("isoPtSumSolidVsEta","Avg Track P_{T} Sum in the Solid Iso Cone vs #eta;#eta;P_{T} (GeV)",etaBin_,etaMin_, etaMax_,sumBin_,sumMin_,sumMax_);
349  p_trackPtSumHollowVsEt_[1] = iBooker.bookProfile("isoPtSumHollowVsEtBarrel","Avg Track P_{T} Sum in the Hollow Iso Cone vs E_{T};E_{T} (GeV);P_{T} (GeV)",etBin_,etMin_,etMax_,sumBin_,sumMin_,sumMax_);
350  p_trackPtSumHollowVsEt_[2] = iBooker.bookProfile("isoPtSumHollowVsEtEndcap","Avg Track P_{T} Sum in the Hollow Iso Cone vs E_{T};E_{T} (GeV);P_{T} (GeV)",etBin_,etMin_,etMax_,sumBin_,sumMin_,sumMax_);
351  p_trackPtSumHollowVsEta_[0] = iBooker.bookProfile("isoPtSumHollowVsEta","Avg Track P_{T} Sum in the Hollow Iso Cone vs #eta;#eta;P_{T} (GeV)",etaBin_,etaMin_, etaMax_,sumBin_,sumMin_,sumMax_);
352  p_ecalSumVsEt_[1] = iBooker.bookProfile("ecalSumVsEtBarrel","Avg Ecal Sum in the Iso Cone vs E_{T};E_{T} (GeV);E (GeV)",etBin_,etMin_, etMax_,sumBin_,sumMin_,sumMax_);
353  p_ecalSumVsEt_[2] = iBooker.bookProfile("ecalSumVsEtEndcap","Avg Ecal Sum in the Iso Cone vs E_{T};E_{T} (GeV);E (GeV)",etBin_,etMin_, etMax_,sumBin_,sumMin_,sumMax_);
354  p_ecalSumVsEta_[0] = iBooker.bookProfile("ecalSumVsEta","Avg Ecal Sum in the Iso Cone vs #eta;#eta;E (GeV)",etaBin_,etaMin_, etaMax_,sumBin_,sumMin_,sumMax_);
355  p_hcalSumVsEt_[1] = iBooker.bookProfile("hcalSumVsEtBarrel","Avg Hcal Sum in the Iso Cone vs E_{T};E_{T} (GeV);E (GeV)",etBin_,etMin_, etMax_,sumBin_,sumMin_,sumMax_);
356  p_hcalSumVsEt_[2] = iBooker.bookProfile("hcalSumVsEtEndcap","Avg Hcal Sum in the Iso Cone vs E_{T};E_{T} (GeV);E (GeV)",etBin_,etMin_, etMax_,sumBin_,sumMin_,sumMax_);
357  p_hcalSumVsEta_[0] = iBooker.bookProfile("hcalSumVsEta","Avg Hcal Sum in the Iso Cone vs #eta;#eta;E (GeV)",etaBin_,etaMin_, etaMax_,sumBin_,sumMin_,sumMax_);
358  p_hOverEVsEt_[1] = iBooker.bookProfile("p_hOverEVsEtBarrel","Avg H/E vs Et;E_{T} (GeV);H/E",etBin_,etMin_,etMax_,hOverEBin_,hOverEMin_,hOverEMax_);
359  p_hOverEVsEt_[2] = iBooker.bookProfile("p_hOverEVsEtEndcap","Avg H/E vs Et;E_{T} (GeV);H/E",etBin_,etMin_,etMax_,hOverEBin_,hOverEMin_,hOverEMax_);
360  p_hOverEVsEta_[0] = iBooker.bookProfile("p_hOverEVsEta","Avg H/E vs #eta;#eta;H/E",etaBin_,etaMin_,etaMax_,hOverEBin_,hOverEMin_,hOverEMax_);
361 
362  // sigmaE/E
363  histname = "sigmaEoverEVsNVtx";
364  p_phoSigmaEoverEVsNVtx_[1] = iBooker.bookProfile(histname+"Barrel","Photons #sigma_{E}/E vs N_{vtx}: Barrel; N_{vtx}; #sigma_{E}/E ",80, -0.5, 79.5, 100,0., 0.08, "");
365  p_phoSigmaEoverEVsNVtx_[2] = iBooker.bookProfile(histname+"Endcap","Photons #sigma_{E}/E vs N_{vtx}: Endcap; N_{vtx}; #sigma_{E}/E",80, -0.5, 79.5, 100,0., 0.08, "");
366  }
367 
369  if(use2DHistos_){
370  //SHOWER SHAPE
371  //r9
372  h2_r9VsEt_[0] = iBooker.book2D("r9VsEt2D","R9 vs E_{T};E_{T} (GeV);R9",reducedEtBin_,etMin_,etMax_,reducedR9Bin_,r9Min_,r9Max_);
373  h2_r9VsEt_[1] = iBooker.book2D("r9VsEt2DBarrel","R9 vs E_{T};E_{T} (GeV);R9",reducedEtBin_,etMin_,etMax_,reducedR9Bin_,r9Min_,r9Max_);
374  h2_r9VsEt_[2] = iBooker.book2D("r9VsEt2DEndcap","R9 vs E_{T};E_{T} (GeV);R9",reducedEtBin_,etMin_,etMax_,reducedR9Bin_,r9Min_,r9Max_);
375  h2_r9VsEta_[0] = iBooker.book2D("r9VsEta2D","R9 vs #eta;#eta;R9",reducedEtaBin_,etaMin_,etaMax_,reducedR9Bin_,r9Min_,r9Max_);
376  //sigmaIetaIeta
377  h2_sigmaIetaIetaVsEta_[0] = iBooker.book2D("sigmaIetaIetaVsEta2D","#sigma_{i#etai#eta} vs #eta;#eta;#sigma_{i#etai#eta}",reducedEtaBin_,etaMin_,etaMax_,sigmaIetaBin_,sigmaIetaMin_,sigmaIetaMax_);
378  //e1x5
379  h2_e1x5VsEt_[0] = iBooker.book2D("e1x5VsEt2D","E1x5 vs E_{T};E_{T} (GeV);E1X5 (GeV)",reducedEtBin_,etMin_,etMax_,reducedEtBin_,etMin_,etMax_);
380  h2_e1x5VsEt_[1] = iBooker.book2D("e1x5VsEt2DBarrel","E1x5 vs E_{T};E_{T} (GeV);E1X5 (GeV)",reducedEtBin_,etMin_,etMax_,reducedEtBin_,etMin_,etMax_);
381  h2_e1x5VsEt_[2] = iBooker.book2D("e1x5VsEt2DEndcap","E1x5 vs E_{T};E_{T} (GeV);E1X5 (GeV)",reducedEtBin_,etMin_,etMax_,reducedEtBin_,etMin_,etMax_);
382  h2_e1x5VsEta_[0] = iBooker.book2D("e1x5VsEta2D","E1x5 vs #eta;#eta;E1X5 (GeV)",reducedEtaBin_,etaMin_,etaMax_,reducedEtBin_,etMin_,etMax_);
383  //e2x5
384  h2_e2x5VsEt_[0] = iBooker.book2D("e2x5VsEt2D","E2x5 vs E_{T};E_{T} (GeV);E2X5 (GeV)",reducedEtBin_,etMin_,etMax_,reducedEtBin_,etMin_,etMax_);
385  h2_e2x5VsEt_[1] = iBooker.book2D("e2x5VsEt2DBarrel","E2x5 vs E_{T};E_{T} (GeV);E2X5 (GeV)",reducedEtBin_,etMin_,etMax_,reducedEtBin_,etMin_,etMax_);
386  h2_e2x5VsEt_[2] = iBooker.book2D("e2x5VsEt2DEndcap","E2x5 vs E_{T};E_{T} (GeV);E2X5 (GeV)",reducedEtBin_,etMin_,etMax_,reducedEtBin_,etMin_,etMax_);
387  h2_e2x5VsEta_[0] = iBooker.book2D("e2x5VsEta2D","E2x5 vs #eta;#eta;E2X5 (GeV)",reducedEtaBin_,etaMin_,etaMax_,reducedEtBin_,etMin_,etMax_);
388  //r1x5
389  h2_r1x5VsEt_[0] = iBooker.book2D("r1x5VsEt2D","R1x5 vs E_{T};E_{T} (GeV);R1X5",reducedEtBin_,etMin_,etMax_,reducedR9Bin_,r9Min_,r9Max_);
390  h2_r1x5VsEt_[1] = iBooker.book2D("r1x5VsEt2DBarrel","R1x5 vs E_{T};E_{T} (GeV);R1X5",reducedEtBin_,etMin_,etMax_,reducedR9Bin_,r9Min_,r9Max_);
391  h2_r1x5VsEt_[2] = iBooker.book2D("r1x5VsEt2DEndcap","R1x5 vs E_{T};E_{T} (GeV);R1X5",reducedEtBin_,etMin_,etMax_,reducedR9Bin_,r9Min_,r9Max_);
392  h2_r1x5VsEta_[0] = iBooker.book2D("r1x5VsEta2D","R1x5 vs #eta;#eta;R1X5",reducedEtaBin_,etaMin_,etaMax_,reducedR9Bin_,r9Min_,r9Max_);
393  //r2x5
394  h2_r2x5VsEt_[0] = iBooker.book2D("r2x5VsEt2D","R2x5 vs E_{T};E_{T} (GeV);R2X5",reducedEtBin_,etMin_,etMax_,reducedR9Bin_,r9Min_,r9Max_);
395  h2_r2x5VsEt_[1] = iBooker.book2D("r2x5VsEt2DBarrel","R2x5 vs E_{T};E_{T} (GeV);R2X5",reducedEtBin_,etMin_,etMax_,reducedR9Bin_,r9Min_,r9Max_);
396  h2_r2x5VsEt_[2] = iBooker.book2D("r2x5VsEt2DEndcap","R2x5 vs E_{T};E_{T} (GeV);R2X5",reducedEtBin_,etMin_,etMax_,reducedR9Bin_,r9Min_,r9Max_);
397  h2_r2x5VsEta_[0] = iBooker.book2D("r2x5VsEta2D","R2x5 vs #eta;#eta;R2X5",reducedEtaBin_,etaMin_,etaMax_,reducedR9Bin_,r9Min_,r9Max_);
398  //TRACK ISOLATION
399  //nTrackIsolSolid
400  h2_nTrackIsolSolidVsEt_[0] = iBooker.book2D("nIsoTracksSolidVsEt2D","Number Of Tracks in the Solid Iso Cone vs E_{T};E_{T};# tracks",reducedEtBin_,etMin_, etMax_,numberBin_,numberMin_,numberMax_);
401  h2_nTrackIsolSolidVsEta_[0] = iBooker.book2D("nIsoTracksSolidVsEta2D","Number Of Tracks in the Solid Iso Cone vs #eta;#eta;# tracks",reducedEtaBin_,etaMin_, etaMax_,numberBin_,numberMin_,numberMax_);
402  //nTrackIsolHollow
403  h2_nTrackIsolHollowVsEt_[0] = iBooker.book2D("nIsoTracksHollowVsEt2D","Number Of Tracks in the Hollow Iso Cone vs E_{T};E_{T};# tracks",reducedEtBin_,etMin_, etMax_,numberBin_,numberMin_,numberMax_);
404  h2_nTrackIsolHollowVsEta_[0] = iBooker.book2D("nIsoTracksHollowVsEta2D","Number Of Tracks in the Hollow Iso Cone vs #eta;#eta;# tracks",reducedEtaBin_,etaMin_, etaMax_,numberBin_,numberMin_,numberMax_);
405  //trackPtSumSolid
406  h2_trackPtSumSolidVsEt_[0] = iBooker.book2D("isoPtSumSolidVsEt2D","Track P_{T} Sum in the Solid Iso Cone;E_{T} (GeV);P_{T} (GeV)",reducedEtBin_,etMin_, etMax_,reducedSumBin_,sumMin_,sumMax_);
407  h2_trackPtSumSolidVsEta_[0] = iBooker.book2D("isoPtSumSolidVsEta2D","Track P_{T} Sum in the Solid Iso Cone;#eta;P_{T} (GeV)",reducedEtaBin_,etaMin_, etaMax_,reducedSumBin_,sumMin_,sumMax_);
408  //trackPtSumHollow
409  h2_trackPtSumHollowVsEt_[0] = iBooker.book2D("isoPtSumHollowVsEt2D","Track P_{T} Sum in the Hollow Iso Cone;E_{T} (GeV);P_{T} (GeV)",reducedEtBin_,etMin_, etMax_,reducedSumBin_,sumMin_,sumMax_);
410  h2_trackPtSumHollowVsEta_[0] = iBooker.book2D("isoPtSumHollowVsEta2D","Track P_{T} Sum in the Hollow Iso Cone;#eta;P_{T} (GeV)",reducedEtaBin_,etaMin_, etaMax_,reducedSumBin_,sumMin_,sumMax_);
411  //CALORIMETER ISOLATION VARIABLES
412  //ecal sum
413  h2_ecalSumVsEt_[0] = iBooker.book2D("ecalSumVsEt2D","Ecal Sum in the Iso Cone;E_{T} (GeV);E (GeV)",reducedEtBin_,etMin_, etMax_,reducedSumBin_,sumMin_,sumMax_);
414  h2_ecalSumVsEta_[0] = iBooker.book2D("ecalSumVsEta2D","Ecal Sum in the Iso Cone;#eta;E (GeV)",reducedEtaBin_,etaMin_, etaMax_,reducedSumBin_,sumMin_,sumMax_);
415  //hcal sum
416  h2_hcalSumVsEt_[0] = iBooker.book2D("hcalSumVsEt2D","Hcal Sum in the Iso Cone;E_{T} (GeV);E (GeV)",reducedEtBin_,etMin_, etMax_,reducedSumBin_,sumMin_,sumMax_);
417  h2_hcalSumVsEta_[0] = iBooker.book2D("hcalSumVsEta2D","Hcal Sum in the Iso Cone;#eta;E (GeV)",reducedEtaBin_,etaMin_, etaMax_,reducedSumBin_,sumMin_,sumMax_);
418  }
419 }
MonitorElement * h_dRPhoPFcand_Pho_unCleaned_[3]
MonitorElement * p_hOverEVsEta_[3]
MonitorElement * h_pfMva_[3]
MonitorElement * h_SumPtOverPhoPt_ChHad_unCleaned_[3]
MonitorElement * h_dRPhoPFcand_Pho_Cleaned_[3]
MonitorElement * h2_nTrackIsolSolidVsEt_[3]
MonitorElement * p_sigmaIetaIetaVsEta_[3]
MonitorElement * p_r1x5VsEta_[3]
MonitorElement * h_phoSigmaIetaIeta_[3]
MonitorElement * h2_r2x5VsEta_[3]
MonitorElement * h_newhOverE_[3]
MonitorElement * h2_trackPtSumHollowVsEt_[3]
MonitorElement * h_SumPtOverPhoPt_NeuHad_unCleaned_[3]
MonitorElement * p_nTrackIsolHollowVsEta_[3]
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:160
MonitorElement * h2_r9VsEt_[3]
MonitorElement * h_nCluOutsideMustache_[3]
MonitorElement * h2_ecalSumVsEt_[3]
MonitorElement * h_nTrackIsolSolid_[3]
MonitorElement * h_SumPtOverPhoPt_Pho_Cleaned_[3]
MonitorElement * h_trackPtSumHollow_[3]
MonitorElement * h2_nTrackIsolSolidVsEta_[3]
MonitorElement * p_ecalSumVsEt_[3]
MonitorElement * h_phoEta_[3]
MonitorElement * h_nPho_[3]
MonitorElement * h_chHadIso_[3]
MonitorElement * h_phoSigmaEoverE_[3]
MonitorElement * h_dRPhoPFcand_NeuHad_Cleaned_[3]
MonitorElement * h_r1x5_[3]
MonitorElement * p_phoSigmaEoverEVsNVtx_[3]
MonitorElement * h2_sigmaIetaIetaVsEta_[3]
MonitorElement * h_dRPhoPFcand_NeuHad_unCleaned_[3]
MonitorElement * h2_hcalSumVsEt_[3]
MonitorElement * h_hcalSum_[3]
MonitorElement * h_phoPhi_[3]
MonitorElement * p_trackPtSumHollowVsEt_[3]
MonitorElement * h2_ecalSumVsEta_[3]
MonitorElement * h_hOverE_[3]
MonitorElement * h_trackPtSumSolid_[3]
MonitorElement * p_hcalSumVsEt_[3]
MonitorElement * p_e1x5VsEt_[3]
MonitorElement * h2_r2x5VsEt_[3]
MonitorElement * h_h1OverE_[3]
MonitorElement * h_SumPtOverPhoPt_ChHad_Cleaned_[3]
MonitorElement * h2_e1x5VsEta_[3]
MonitorElement * h1_mumuGammaInvMass_[3]
MonitorElement * h2_trackPtSumSolidVsEt_[3]
MonitorElement * h_phoE_[3]
MonitorElement * p_r2x5VsEta_[3]
MonitorElement * h_nRecoVtx_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
MonitorElement * p_e2x5VsEt_[3]
MonitorElement * p_hcalSumVsEta_[3]
MonitorElement * h_phoEt_[3]
MonitorElement * p_ecalSumVsEta_[3]
MonitorElement * p_r9VsEta_[3]
MonitorElement * h_nHadIso_[3]
MonitorElement * p_trackPtSumSolidVsEta_[3]
MonitorElement * h1_mumuInvMass_[3]
photon histos
MonitorElement * h2_trackPtSumSolidVsEta_[3]
MonitorElement * h2_hcalSumVsEta_[3]
MonitorElement * h_dRPhoPFcand_ChHad_unCleaned_[3]
MonitorElement * h_scPhi_[3]
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * h2_e1x5VsEt_[3]
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
MonitorElement * h_scEta_[3]
MonitorElement * h_nTrackIsolHollow_[3]
MonitorElement * p_r2x5VsEt_[3]
MonitorElement * h_SumPtOverPhoPt_Pho_unCleaned_[3]
MonitorElement * h_SumPtOverPhoPt_NeuHad_Cleaned_[3]
MonitorElement * h_etOutsideMustache_[3]
MonitorElement * h_e1x5_[3]
MonitorElement * h2_nTrackIsolHollowVsEt_[3]
MonitorElement * h2_nTrackIsolHollowVsEta_[3]
MonitorElement * h2_e2x5VsEta_[3]
MonitorElement * h2_trackPtSumHollowVsEta_[3]
MonitorElement * h_r9_[3]
MonitorElement * p_r1x5VsEt_[3]
MonitorElement * p_nTrackIsolHollowVsEt_[3]
MonitorElement * p_r9VsEt_[3]
MonitorElement * p_trackPtSumSolidVsEt_[3]
MonitorElement * p_e2x5VsEta_[3]
MonitorElement * h_e2x5_[3]
MonitorElement * h2_r1x5VsEta_[3]
MonitorElement * p_trackPtSumHollowVsEta_[3]
MonitorElement * h_dRPhoPFcand_ChHad_Cleaned_[3]
MonitorElement * p_hOverEVsEt_[3]
MonitorElement * h2_r1x5VsEt_[3]
MonitorElement * p_nTrackIsolSolidVsEt_[3]
MonitorElement * p_e1x5VsEta_[3]
MonitorElement * p_nTrackIsolSolidVsEta_[3]
MonitorElement * h2_r9VsEta_[3]
MonitorElement * h_ecalSum_[3]
MonitorElement * h_phoIso_[3]
MonitorElement * h_h2OverE_[3]
MonitorElement * h_r2x5_[3]
MonitorElement * h2_e2x5VsEt_[3]
float ZToMuMuGammaAnalyzer::mumuGammaInvMass ( const reco::Muon mu1,
const reco::Muon mu2,
const reco::PhotonRef pho 
)
private

Definition at line 869 of file ZToMuMuGammaAnalyzer.cc.

References reco::LeafCandidate::p4(), and mathSSE::sqrt().

869  {
870  math::XYZTLorentzVector p12 = mu1.p4()+mu2.p4()+pho->p4() ;
871  float Mass2 = p12.Dot(p12) ;
872  float invMass = sqrt(Mass2) ;
873  return invMass ;
874 }
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
T sqrt(T t)
Definition: SSEVec.h:18
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
float ZToMuMuGammaAnalyzer::mumuInvMass ( const reco::Muon m1,
const reco::Muon m2 
)
private

Definition at line 862 of file ZToMuMuGammaAnalyzer.cc.

References reco::LeafCandidate::p4(), and mathSSE::sqrt().

862  {
863  math::XYZTLorentzVector p12 = mu1.p4()+mu2.p4() ;
864  float mumuMass2 = p12.Dot(p12) ;
865  float invMass = sqrt(mumuMass2) ;
866  return invMass ;
867 }
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
T sqrt(T t)
Definition: SSEVec.h:18
bool ZToMuMuGammaAnalyzer::muonSelection ( const reco::Muon m,
const reco::BeamSpot bs 
)
private

Definition at line 819 of file ZToMuMuGammaAnalyzer.cc.

References reco::LeafCandidate::eta(), reco::Muon::globalTrack(), reco::Muon::isolationR03(), reco::Muon::isTrackerMuon(), reco::Muon::numberOfMatches(), mps_fire::result, reco::MuonIsolation::sumPt, and reco::Muon::track().

819  {
820  bool result=true;
821  if ( mu.globalTrack()->normalizedChi2() > muonMaxChi2_ ) result=false;
822  if ( fabs( mu.globalTrack()->dxy(beamSpot)) > muonMaxDxy_ ) result=false;
823  if ( mu.numberOfMatches() < muonMatches_ ) result=false;
824 
825  if ( mu.track()-> hitPattern().numberOfValidPixelHits() < validPixHits_ ) result=false;
826  if ( mu.globalTrack()->hitPattern().numberOfValidMuonHits() < validMuonHits_ ) result=false;
827  if ( !mu.isTrackerMuon() ) result=false;
828  // track isolation
829  if ( mu.isolationR03().sumPt > muonTrackIso_ ) result=false;
830  if ( fabs(mu.eta())> muonTightEta_ ) result=false;
831 
832  return result;
833 }
const int mu
Definition: Constants.h:22
bool ZToMuMuGammaAnalyzer::photonSelection ( const reco::PhotonRef p)
private

remove after moriond if (EtCorrEcalIso>4.0) result=false;

remove after moriond if (EtCorrEcalIso>50.0) result=false;

Definition at line 835 of file ZToMuMuGammaAnalyzer.cc.

References funct::false, and mps_fire::result.

835  {
836  bool result=true;
837  if ( pho->pt() < photonMinEt_ ) result=false;
838  if ( fabs(pho->eta())> photonMaxEta_ ) result=false;
839  if ( pho->isEBEEGap() ) result=false;
840 
841  double EtCorrHcalIso = pho->hcalTowerSumEtConeDR03() - 0.005*pho->pt();
842  double EtCorrTrkIso = pho->trkSumPtHollowConeDR03() - 0.002*pho->pt();
843 
844  if (pho->r9() <=0.9) {
845  if (pho->isEB() && (pho->hadTowOverEm()>0.075 || pho->sigmaIetaIeta() > 0.014)) result=false;
846  if (pho->isEE() && (pho->hadTowOverEm()>0.075 || pho->sigmaIetaIeta() > 0.034)) result=false;
848  if (EtCorrHcalIso>4.0) result=false;
849  if (EtCorrTrkIso>4.0) result=false ;
850  if ( pho->chargedHadronIso() > 4 ) result=false;
851  } else {
852  if (pho->isEB() && (pho->hadTowOverEm()>0.082 || pho->sigmaIetaIeta() > 0.014)) result=false;
853  if (pho->isEE() && (pho->hadTowOverEm()>0.075 || pho->sigmaIetaIeta() > 0.034)) result=false;
855  if (EtCorrHcalIso>50.0) result=false;
856  if (EtCorrTrkIso>50.0) result=false;
857  if ( pho->chargedHadronIso() > 4 ) result=false;
858  }
859  return result;
860 }

Member Data Documentation

edm::EDGetTokenT<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > ZToMuMuGammaAnalyzer::barrelRecHit_token_
private

Definition at line 112 of file ZToMuMuGammaAnalyzer.h.

edm::EDGetTokenT<reco::BeamSpot> ZToMuMuGammaAnalyzer::beamSpot_token_
private

Definition at line 115 of file ZToMuMuGammaAnalyzer.h.

std::stringstream ZToMuMuGammaAnalyzer::currentFolder_
private

Definition at line 124 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::eBin_
private

Definition at line 157 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::eMax_
private

Definition at line 156 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::eMin_
private

Definition at line 155 of file ZToMuMuGammaAnalyzer.h.

edm::EDGetTokenT<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > ZToMuMuGammaAnalyzer::endcapRecHit_token_
private

Definition at line 113 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::etaBin_
private

Definition at line 169 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::etaMax_
private

Definition at line 168 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::etaMin_
private

Definition at line 167 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::etBin_
private

Definition at line 161 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::etMax_
private

Definition at line 160 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::etMin_
private

Definition at line 159 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::farMuonEcalIso_
private

Definition at line 148 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::farMuonMinPt_
private

Definition at line 150 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::farMuonTrackIso_
private

Definition at line 149 of file ZToMuMuGammaAnalyzer.h.

std::string ZToMuMuGammaAnalyzer::fName_
private

Definition at line 120 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h1_mumuGammaInvMass_[3]
private

Definition at line 205 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h1_mumuInvMass_[3]
private

photon histos

Definition at line 204 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_e1x5VsEt_[3]
private

Definition at line 228 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_e1x5VsEta_[3]
private

Definition at line 226 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_e2x5VsEt_[3]
private

Definition at line 234 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_e2x5VsEta_[3]
private

Definition at line 232 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_ecalSumVsEt_[3]
private

Definition at line 278 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_ecalSumVsEta_[3]
private

Definition at line 280 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_hcalSumVsEt_[3]
private

Definition at line 284 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_hcalSumVsEta_[3]
private

Definition at line 286 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_nTrackIsolHollowVsEt_[3]
private

Definition at line 260 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_nTrackIsolHollowVsEta_[3]
private

Definition at line 262 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_nTrackIsolSolidVsEt_[3]
private

Definition at line 254 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_nTrackIsolSolidVsEta_[3]
private

Definition at line 256 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_r1x5VsEt_[3]
private

Definition at line 240 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_r1x5VsEta_[3]
private

Definition at line 238 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_r2x5VsEt_[3]
private

Definition at line 246 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_r2x5VsEta_[3]
private

Definition at line 244 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_r9VsEt_[3]
private

Definition at line 220 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_r9VsEta_[3]
private

Definition at line 222 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_sigmaIetaIetaVsEta_[3]
private

Definition at line 250 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_trackPtSumHollowVsEt_[3]
private

Definition at line 272 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_trackPtSumHollowVsEta_[3]
private

Definition at line 274 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_trackPtSumSolidVsEt_[3]
private

Definition at line 266 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h2_trackPtSumSolidVsEta_[3]
private

Definition at line 268 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_chHadIso_[3]
private

Definition at line 300 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_dRPhoPFcand_ChHad_Cleaned_[3]
private

Definition at line 308 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_dRPhoPFcand_ChHad_unCleaned_[3]
private

Definition at line 311 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_dRPhoPFcand_NeuHad_Cleaned_[3]
private

Definition at line 309 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_dRPhoPFcand_NeuHad_unCleaned_[3]
private

Definition at line 312 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_dRPhoPFcand_Pho_Cleaned_[3]
private

Definition at line 310 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_dRPhoPFcand_Pho_unCleaned_[3]
private

Definition at line 313 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_e1x5_[3]
private

Definition at line 225 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_e2x5_[3]
private

Definition at line 231 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_ecalSum_[3]
private

Definition at line 277 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_etOutsideMustache_[3]
private

Definition at line 305 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_h1OverE_[3]
private

Definition at line 292 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_h2OverE_[3]
private

Definition at line 293 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_hcalSum_[3]
private

Definition at line 283 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_hOverE_[3]
private

Definition at line 289 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_nCluOutsideMustache_[3]
private

Definition at line 304 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_newhOverE_[3]
private

Definition at line 295 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_nHadIso_[3]
private

Definition at line 301 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_nPho_[3]
private

Definition at line 212 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_nRecoVtx_
private

Definition at line 202 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_nTrackIsolHollow_[3]
private

Definition at line 259 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_nTrackIsolSolid_[3]
private

Definition at line 253 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_pfMva_[3]
private

Definition at line 306 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_phoE_[3]
private

Definition at line 207 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_phoEt_[3]
private

Definition at line 210 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_phoEta_[3]
private

Definition at line 214 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_phoIso_[3]
private

Definition at line 302 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_phoPhi_[3]
private

Definition at line 215 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_phoSigmaEoverE_[3]
private

Definition at line 208 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_phoSigmaIetaIeta_[3]
private

Definition at line 249 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_r1x5_[3]
private

Definition at line 237 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_r2x5_[3]
private

Definition at line 243 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_r9_[3]
private

Definition at line 219 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_scEta_[3]
private

Definition at line 216 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_scPhi_[3]
private

Definition at line 217 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_SumPtOverPhoPt_ChHad_Cleaned_[3]
private

Definition at line 314 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_SumPtOverPhoPt_ChHad_unCleaned_[3]
private

Definition at line 317 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_SumPtOverPhoPt_NeuHad_Cleaned_[3]
private

Definition at line 315 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_SumPtOverPhoPt_NeuHad_unCleaned_[3]
private

Definition at line 318 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_SumPtOverPhoPt_Pho_Cleaned_[3]
private

Definition at line 316 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_SumPtOverPhoPt_Pho_unCleaned_[3]
private

Definition at line 319 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_trackPtSumHollow_[3]
private

Definition at line 271 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::h_trackPtSumSolid_[3]
private

Definition at line 265 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::hOverEBin_
private

Definition at line 181 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::hOverEMax_
private

Definition at line 180 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::hOverEMin_
private

Definition at line 179 of file ZToMuMuGammaAnalyzer.h.

bool ZToMuMuGammaAnalyzer::makeProfiles_
private

Definition at line 122 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::maxMumuGammaInvMass_
private

Definition at line 152 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::maxMumuInvMass_
private

Definition at line 139 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::minMumuGammaInvMass_
private

Definition at line 151 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::minMumuInvMass_
private

Definition at line 138 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::minPixStripHits_
private

Definition at line 129 of file ZToMuMuGammaAnalyzer.h.

edm::EDGetTokenT<std::vector<reco::Muon> > ZToMuMuGammaAnalyzer::muon_token_
private

Definition at line 109 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::muonMatches_
private

Definition at line 132 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::muonMaxChi2_
private

Definition at line 130 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::muonMaxDxy_
private

Definition at line 131 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::muonMinPt_
private

Definition at line 128 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::muonTightEta_
private

Definition at line 136 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::muonTrackIso_
private

Definition at line 135 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::nearMuonDr_
private

Definition at line 146 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::nearMuonHcalIso_
private

Definition at line 147 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::nEvt_
private

Definition at line 125 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::numberBin_
private

Definition at line 185 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::numberMax_
private

Definition at line 184 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::numberMin_
private

Definition at line 183 of file ZToMuMuGammaAnalyzer.h.

edm::EDGetTokenT<reco::VertexCollection> ZToMuMuGammaAnalyzer::offline_pvToken_
private

Definition at line 117 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_e1x5VsEt_[3]
private

Definition at line 229 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_e1x5VsEta_[3]
private

Definition at line 227 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_e2x5VsEt_[3]
private

Definition at line 235 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_e2x5VsEta_[3]
private

Definition at line 233 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_ecalSumVsEt_[3]
private

Definition at line 279 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_ecalSumVsEta_[3]
private

Definition at line 281 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_hcalSumVsEt_[3]
private

Definition at line 285 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_hcalSumVsEta_[3]
private

Definition at line 287 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_hOverEVsEt_[3]
private

Definition at line 290 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_hOverEVsEta_[3]
private

Definition at line 291 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_newhOverEVsEt_[3]
private

Definition at line 297 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_newhOverEVsEta_[3]
private

Definition at line 296 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_nTrackIsolHollowVsEt_[3]
private

Definition at line 261 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_nTrackIsolHollowVsEta_[3]
private

Definition at line 263 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_nTrackIsolSolidVsEt_[3]
private

Definition at line 255 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_nTrackIsolSolidVsEta_[3]
private

Definition at line 257 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_phoSigmaEoverEVsNVtx_[3]
private

Definition at line 209 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_r1x5VsEt_[3]
private

Definition at line 241 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_r1x5VsEta_[3]
private

Definition at line 239 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_r2x5VsEt_[3]
private

Definition at line 247 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_r2x5VsEta_[3]
private

Definition at line 245 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_r9VsEt_[3]
private

Definition at line 221 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_r9VsEta_[3]
private

Definition at line 223 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_sigmaIetaIetaVsEta_[3]
private

Definition at line 251 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_trackPtSumHollowVsEt_[3]
private

Definition at line 273 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_trackPtSumHollowVsEta_[3]
private

Definition at line 275 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_trackPtSumSolidVsEt_[3]
private

Definition at line 267 of file ZToMuMuGammaAnalyzer.h.

MonitorElement* ZToMuMuGammaAnalyzer::p_trackPtSumSolidVsEta_[3]
private

Definition at line 269 of file ZToMuMuGammaAnalyzer.h.

edm::EDGetTokenT<reco::PFCandidateCollection> ZToMuMuGammaAnalyzer::pfCandidates_
private

Definition at line 116 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::phiBin_
private

Definition at line 173 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::phiMax_
private

Definition at line 172 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::phiMin_
private

Definition at line 171 of file ZToMuMuGammaAnalyzer.h.

edm::EDGetTokenT<std::vector<reco::Photon> > ZToMuMuGammaAnalyzer::photon_token_
private

Definition at line 108 of file ZToMuMuGammaAnalyzer.h.

edm::EDGetTokenT<edm::ValueMap<bool> > ZToMuMuGammaAnalyzer::PhotonIDLoose_token_
private

Definition at line 110 of file ZToMuMuGammaAnalyzer.h.

edm::EDGetTokenT<edm::ValueMap<bool> > ZToMuMuGammaAnalyzer::PhotonIDTight_token_
private

Definition at line 111 of file ZToMuMuGammaAnalyzer.h.

edm::EDGetTokenT<edm::ValueMap<std::vector<reco::PFCandidateRef> > > ZToMuMuGammaAnalyzer::photonIsoValmap_token_
private

Definition at line 118 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::photonMaxEta_
private

Definition at line 142 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::photonMinEt_
private

Definition at line 141 of file ZToMuMuGammaAnalyzer.h.

float ZToMuMuGammaAnalyzer::photonTrackIso_
private

Definition at line 143 of file ZToMuMuGammaAnalyzer.h.

unsigned int ZToMuMuGammaAnalyzer::prescaleFactor_
private

Definition at line 123 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::r9Bin_
private

Definition at line 177 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::r9Max_
private

Definition at line 176 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::r9Min_
private

Definition at line 175 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::reducedEtaBin_
private

Definition at line 192 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::reducedEtBin_
private

Definition at line 191 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::reducedR9Bin_
private

Definition at line 194 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::reducedSumBin_
private

Definition at line 193 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::sigmaIetaBin_
private

Definition at line 189 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::sigmaIetaMax_
private

Definition at line 188 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::sigmaIetaMin_
private

Definition at line 187 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::sumBin_
private

Definition at line 165 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::sumMax_
private

Definition at line 164 of file ZToMuMuGammaAnalyzer.h.

double ZToMuMuGammaAnalyzer::sumMin_
private

Definition at line 163 of file ZToMuMuGammaAnalyzer.h.

edm::EDGetTokenT<trigger::TriggerEvent> ZToMuMuGammaAnalyzer::triggerEvent_token_
private

Definition at line 114 of file ZToMuMuGammaAnalyzer.h.

bool ZToMuMuGammaAnalyzer::use2DHistos_
private

Definition at line 121 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::validMuonHits_
private

Definition at line 134 of file ZToMuMuGammaAnalyzer.h.

int ZToMuMuGammaAnalyzer::validPixHits_
private

Definition at line 133 of file ZToMuMuGammaAnalyzer.h.