CMS 3D CMS Logo

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

#include <METAnalyzer.h>

Inheritance diagram for METAnalyzer:
thread_unsafe::DQMEDAnalyzer edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &)
 Get the analysis. More...
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 Inizialize parameters for histo binning. More...
 
void dqmBeginRun (const edm::Run &, const edm::EventSetup &)
 Initialize run-based parameters. More...
 
void endRun (const edm::Run &iRun, const edm::EventSetup &iSetup)
 Finish up a run. More...
 
void fillMESet (const edm::Event &, std::string, const reco::MET &, const reco::PFMET &, const reco::CaloMET &, std::map< std::string, MonitorElement * > &)
 
void fillMonitorElement (const edm::Event &, std::string, std::string, const reco::MET &, const reco::PFMET &, const reco::CaloMET &, std::map< std::string, MonitorElement * > &, bool)
 
void makeRatePlot (std::string, double)
 
 METAnalyzer (const edm::ParameterSet &)
 Constructor. More...
 
virtual ~METAnalyzer ()
 Destructor. More...
 
- Public Member Functions inherited from thread_unsafe::DQMEDAnalyzer
virtual void beginRun (edm::Run const &, edm::EventSetup const &) final
 
 DQMEDAnalyzer (void)
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

void bookMESet (std::string, DQMStore::IBooker &, std::map< std::string, MonitorElement * > &)
 
void bookMonitorElement (std::string, DQMStore::IBooker &, std::map< std::string, MonitorElement * > &, bool)
 

Private Attributes

std::vector< int > allTriggerDecisions_
 
std::vector< std::string > allTriggerNames_
 
edm::InputTag beamHaloSummaryTag_
 
edm::EDGetTokenT
< reco::BeamHaloSummary
beamHaloSummaryToken_
 
math::XYZPoint beamSpot_
 
bool bypassAllDCSChecks_
 
bool bypassAllPVChecks_
 
edm::EDGetTokenT
< reco::CaloJetCollection
caloJetsToken_
 
edm::EDGetTokenT
< reco::CaloMETCollection
caloMetToken_
 
edm::ParameterSet cleaningParameters_
 
JetMETDQMDCSFilterDCSFilter_
 
bool fill_met_high_level_histo
 
std::string FolderName_
 
std::vector< std::string > folderNames_
 
edm::InputTag gtTag_
 
edm::EDGetTokenT
< L1GlobalTriggerReadoutRecord
gtToken_
 
edm::InputTag hbheNoiseFilterResultTag_
 
edm::EDGetTokenT< bool > hbheNoiseFilterResultToken_
 
edm::InputTag hcalNoiseRBXCollectionTag_
 
edm::EDGetTokenT
< reco::HcalNoiseRBXCollection
HcalNoiseRBXToken_
 
MonitorElementhCaloEmEtFraction
 
MonitorElementhCaloEmEtFraction020
 
MonitorElementhCaloEmEtInEB
 
MonitorElementhCaloEmEtInEE
 
MonitorElementhCaloEmEtInHF
 
MonitorElementhCaloEmMET
 
MonitorElementhCaloEmMETPhi
 
MonitorElementhCaloEmMEx
 
MonitorElementhCaloEmMEy
 
MonitorElementhCaloEtFractionHadronic
 
MonitorElementhCaloHadEtInHB
 
MonitorElementhCaloHadEtInHE
 
MonitorElementhCaloHadEtInHF
 
MonitorElementhCaloHadEtInHO
 
MonitorElementhCaloHaMET
 
MonitorElementhCaloHaMETPhi
 
MonitorElementhCaloHaMEx
 
MonitorElementhCaloHaMEy
 
MonitorElementhCaloMaxEtInEmTowers
 
MonitorElementhCaloMaxEtInHadTowers
 
MonitorElementhCaloMETPhi020
 
double hfCalibFactor_
 
HLTConfigProvider hltConfig_
 
std::string hltPhysDec_
 
MonitorElementhMET
 
MonitorElementhMET_logx
 
MonitorElementhMETPhi
 
MonitorElementhMETRate
 
MonitorElementhMETSig
 
MonitorElementhMEx
 
MonitorElementhMExLS
 
MonitorElementhMEy
 
MonitorElementhMEyLS
 
MonitorElementhSumET
 
MonitorElementhSumET_logx
 
MonitorElementhTrigger
 
bool hTriggerLabelsIsSet_
 
edm::InputTag inputJetIDValueMap
 
bool isCaloMet_
 
bool isPFMet_
 
edm::InputTag jetCollectionLabel_
 
std::string jetCorrectionService_
 
edm::EDGetTokenT
< edm::ValueMap< reco::JetID > > 
jetID_ValueMapToken_
 
JetIDSelectionFunctor jetIDFunctorLoose
 
int LSBegin_
 
int LSEnd_
 
MonitorElementlumisecME
 
std::map< std::string,
MonitorElement * > 
map_dijet_MEs
 
MonitorElementmeChargedHadronEt
 
MonitorElementmeChargedHadronEt_profile
 
MonitorElementmeChargedHadronEtFraction
 
MonitorElementmeChargedHadronEtFraction_profile
 
MonitorElementmeElectronEt
 
MonitorElementmeElectronEt_profile
 
MonitorElementmeElectronEtFraction
 
MonitorElementmeElectronEtFraction_profile
 
MonitorElementmeHFEMEt
 
MonitorElementmeHFEMEt_profile
 
MonitorElementmeHFEMEtFraction
 
MonitorElementmeHFEMEtFraction_profile
 
MonitorElementmeHFHadronEt
 
MonitorElementmeHFHadronEt_profile
 
MonitorElementmeHFHadronEtFraction
 
MonitorElementmeHFHadronEtFraction_profile
 
MonitorElementmeMET_profile
 
MonitorElementmeMEx_profile
 
MonitorElementmeMEy_profile
 
MonitorElementmeMuonEt
 
MonitorElementmeMuonEt_profile
 
MonitorElementmeMuonEtFraction
 
MonitorElementmeMuonEtFraction_profile
 
MonitorElementmeNeutralHadronEt
 
MonitorElementmeNeutralHadronEt_profile
 
MonitorElementmeNeutralHadronEtFraction
 
MonitorElementmeNeutralHadronEtFraction_profile
 
MonitorElementmePhotonEt
 
MonitorElementmePhotonEt_profile
 
MonitorElementmePhotonEtFraction
 
MonitorElementmePhotonEtFraction_profile
 
MonitorElementmeSumET_profile
 
edm::InputTag metCollectionLabel_
 
std::string MetType_
 
std::string mOutputFile_
 
int nbinsPV_
 
double nPVMax_
 
double nPVMin_
 
int numPV_
 
bool outputMEsInRootFile
 
edm::ParameterSet parameters
 
PFJetIDSelectionFunctor pfjetIDFunctorLoose
 
edm::EDGetTokenT
< reco::PFJetCollection
pfJetsToken_
 
edm::EDGetTokenT
< reco::PFMETCollection
pfMetToken_
 
double ptThreshold_
 
bool runcosmics_
 
std::vector< int > triggerFolderDecisions_
 
std::vector
< GenericTriggerEventFlag * > 
triggerFolderEventFlag_
 
std::vector< std::vector
< std::string > > 
triggerFolderExpr_
 
std::vector< std::string > triggerFolderLabels_
 
edm::InputTag triggerResultsLabel_
 
edm::EDGetTokenT
< edm::TriggerResults
triggerResultsToken_
 
edm::VParameterSet triggerSelectedSubFolders_
 
int verbose_
 
edm::InputTag vertexTag_
 
edm::EDGetTokenT< std::vector
< reco::Vertex > > 
vertexToken_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- 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)
 

Detailed Description

Definition at line 81 of file METAnalyzer.h.

Constructor & Destructor Documentation

METAnalyzer::METAnalyzer ( const edm::ParameterSet pSet)

Constructor.

Definition at line 41 of file METAnalyzer.cc.

References PFJetIDSelectionFunctor::FIRSTDATA, edm::ParameterSet::getParameter(), PFJetIDSelectionFunctor::LOOSE, JetIDSelectionFunctor::LOOSE, Parameters::parameters, JetIDSelectionFunctor::PURE09, AlCaHLTBitMon_QueryRunRegistry::string, and triggerResultsToken_().

41  {
42  parameters = pSet;
43 
44  outputMEsInRootFile = parameters.getParameter<bool>("OutputMEsInRootFile");
46 
47  LSBegin_ = pSet.getParameter<int>("LSBegin");
48  LSEnd_ = pSet.getParameter<int>("LSEnd");
49 
51 
53  triggerResultsToken_= consumes<edm::TriggerResults>(edm::InputTag(triggerResultsLabel_));
54  jetCorrectionService_ = pSet.getParameter<std::string> ("JetCorrections");
55 
56  isCaloMet_ = (std::string("calo")==MetType_);
57  //isTCMet_ = (std::string("tc") ==MetType_);
58  isPFMet_ = (std::string("pf") ==MetType_);
59 
60  // MET information
62 
63  if(/*isTCMet_ || */isCaloMet_){
64  inputJetIDValueMap = pSet.getParameter<edm::InputTag>("InputJetIDValueMap");
65  jetID_ValueMapToken_= consumes< edm::ValueMap<reco::JetID> >(inputJetIDValueMap);
67  }
68  if(isPFMet_){
70  }
71  ptThreshold_ = parameters.getParameter<double>("ptThreshold");
72 
73 
74  if(isPFMet_){
75  pfMetToken_= consumes<reco::PFMETCollection>(edm::InputTag(metCollectionLabel_));
76  }
77  if(isCaloMet_){
78  caloMetToken_= consumes<reco::CaloMETCollection>(edm::InputTag(metCollectionLabel_));
79  }
80  //if(isTCMet_){
81  // tcMetToken_= consumes<reco::METCollection>(edm::InputTag(metCollectionLabel_));
82  //}
83 
84  fill_met_high_level_histo = parameters.getParameter<bool>("fillMetHighLevel");
85 
86  hTriggerLabelsIsSet_ = false;
87  //jet cleanup parameters
88  cleaningParameters_ = pSet.getParameter<ParameterSet>("CleaningParameters");
89 
91  //DCS
93 
94  //Vertex requirements
95  bypassAllPVChecks_ = cleaningParameters_.getParameter<bool>("bypassAllPVChecks");
96  bypassAllDCSChecks_ = cleaningParameters_.getParameter<bool>("bypassAllDCSChecks");
97  runcosmics_ = parameters.getUntrackedParameter<bool>("runcosmics");
99  vertexToken_ = consumes<std::vector<reco::Vertex> >(edm::InputTag(vertexTag_));
100 
101  //Trigger parameters
103  gtToken_= consumes<L1GlobalTriggerReadoutRecord>(edm::InputTag(gtTag_));
104 
105  //inputTrackLabel_ = parameters.getParameter<edm::InputTag>("InputTrackLabel");
106  //inputMuonLabel_ = parameters.getParameter<edm::InputTag>("InputMuonLabel");
107  //inputElectronLabel_ = parameters.getParameter<edm::InputTag>("InputElectronLabel");
108  //inputBeamSpotLabel_ = parameters.getParameter<edm::InputTag>("InputBeamSpotLabel");
109  //inputTCMETValueMap_ = parameters.getParameter<edm::InputTag>("InputTCMETValueMap");
110  //TrackToken_= consumes<edm::View <reco::Track> >(inputTrackLabel_);
111  //MuonToken_= consumes<reco::MuonCollection>(inputMuonLabel_);
112  //ElectronToken_= consumes<edm::View<reco::GsfElectron> >(inputElectronLabel_);
113  //BeamspotToken_= consumes<reco::BeamSpot>(inputBeamSpotLabel_);
114  //tcMETValueMapToken_= consumes< edm::ValueMap<reco::MuonMETCorrectionData> >(inputTCMETValueMap_);
115 
116  // Other data collections
117  jetCollectionLabel_ = parameters.getParameter<edm::InputTag>("JetCollectionLabel");
118  if (isCaloMet_) caloJetsToken_ = consumes<reco::CaloJetCollection>(jetCollectionLabel_);
119  //if (isTCMet_) jptJetsToken_ = consumes<reco::JPTJetCollection>(jetCollectionLabel_);
120  if (isPFMet_) pfJetsToken_ = consumes<reco::PFJetCollection>(jetCollectionLabel_);
121 
123  HcalNoiseRBXToken_ = consumes<reco::HcalNoiseRBXCollection>(hcalNoiseRBXCollectionTag_);
124 
125  beamHaloSummaryTag_ = parameters.getParameter<edm::InputTag>("BeamHaloSummaryLabel");
126  beamHaloSummaryToken_ = consumes<BeamHaloSummary>(beamHaloSummaryTag_);
127  hbheNoiseFilterResultTag_ = parameters.getParameter<edm::InputTag>("HBHENoiseFilterResultLabel");
129 
130  //
131  nbinsPV_ = parameters.getParameter<int>("pVBin");
132  nPVMin_ = parameters.getParameter<double>("pVMin");
133  nPVMax_ = parameters.getParameter<double>("pVMax");
134 
136  for (edm::VParameterSet::const_iterator it = triggerSelectedSubFolders_.begin(); it!= triggerSelectedSubFolders_.end(); it++) {
138  triggerFolderExpr_.push_back(it->getParameter<std::vector<std::string> >("hltPaths"));
139  triggerFolderLabels_.push_back(it->getParameter<std::string>("label"));
140  }
141 
142  cleaningParameters_ = parameters.getParameter<ParameterSet>("CleaningParameters");
143 
144  verbose_ = parameters.getParameter<int>("verbose");
145 
147 
148 }
edm::EDGetTokenT< reco::CaloMETCollection > caloMetToken_
Definition: METAnalyzer.h:154
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::InputTag gtTag_
Definition: METAnalyzer.h:141
edm::EDGetTokenT< reco::CaloJetCollection > caloJetsToken_
Definition: METAnalyzer.h:145
bool bypassAllDCSChecks_
Definition: METAnalyzer.h:240
double ptThreshold_
Definition: METAnalyzer.h:176
bool runcosmics_
Definition: METAnalyzer.h:241
bool isCaloMet_
Definition: METAnalyzer.h:403
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:33
std::string MetType_
Definition: METAnalyzer.h:130
edm::InputTag hbheNoiseFilterResultTag_
Definition: METAnalyzer.h:139
edm::EDGetTokenT< reco::PFMETCollection > pfMetToken_
Definition: METAnalyzer.h:153
bool outputMEsInRootFile
Definition: METAnalyzer.h:131
double nPVMin_
Definition: METAnalyzer.h:232
Provides a code based selection for trigger and DCS information in order to have no failing filters i...
std::string FolderName_
Definition: METAnalyzer.h:133
std::vector< GenericTriggerEventFlag * > triggerFolderEventFlag_
Definition: METAnalyzer.h:198
edm::InputTag beamHaloSummaryTag_
Definition: METAnalyzer.h:138
std::vector< std::vector< std::string > > triggerFolderExpr_
Definition: METAnalyzer.h:199
std::string mOutputFile_
Definition: METAnalyzer.h:132
edm::InputTag vertexTag_
Definition: METAnalyzer.h:140
edm::EDGetTokenT< std::vector< reco::Vertex > > vertexToken_
Definition: METAnalyzer.h:143
bool fill_met_high_level_histo
Definition: METAnalyzer.h:407
edm::InputTag metCollectionLabel_
Definition: METAnalyzer.h:135
std::string jetCorrectionService_
Definition: METAnalyzer.h:174
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
edm::ParameterSet cleaningParameters_
Definition: METAnalyzer.h:228
edm::EDGetTokenT< bool > hbheNoiseFilterResultToken_
Definition: METAnalyzer.h:149
edm::EDGetTokenT< reco::HcalNoiseRBXCollection > HcalNoiseRBXToken_
Definition: METAnalyzer.h:155
PF Jet selector for pat::Jets.
edm::EDGetTokenT< edm::TriggerResults > triggerResultsToken_
Definition: METAnalyzer.h:190
edm::EDGetTokenT< reco::PFJetCollection > pfJetsToken_
Definition: METAnalyzer.h:146
JetMETDQMDCSFilter * DCSFilter_
Definition: METAnalyzer.h:267
Jet selector for pat::Jets and for CaloJets.
PFJetIDSelectionFunctor pfjetIDFunctorLoose
Definition: METAnalyzer.h:172
double nPVMax_
Definition: METAnalyzer.h:233
edm::InputTag triggerResultsLabel_
Definition: METAnalyzer.h:189
edm::ParameterSet parameters
Definition: METAnalyzer.h:125
edm::VParameterSet triggerSelectedSubFolders_
Definition: METAnalyzer.h:197
edm::EDGetTokenT< L1GlobalTriggerReadoutRecord > gtToken_
Definition: METAnalyzer.h:144
bool bypassAllPVChecks_
Definition: METAnalyzer.h:239
edm::EDGetTokenT< edm::ValueMap< reco::JetID > > jetID_ValueMapToken_
Definition: METAnalyzer.h:169
edm::InputTag inputJetIDValueMap
Definition: METAnalyzer.h:168
edm::InputTag hcalNoiseRBXCollectionTag_
Definition: METAnalyzer.h:136
JetIDSelectionFunctor jetIDFunctorLoose
Definition: METAnalyzer.h:171
edm::EDGetTokenT< reco::BeamHaloSummary > beamHaloSummaryToken_
Definition: METAnalyzer.h:150
edm::InputTag jetCollectionLabel_
Definition: METAnalyzer.h:137
bool hTriggerLabelsIsSet_
Definition: METAnalyzer.h:277
std::vector< std::string > triggerFolderLabels_
Definition: METAnalyzer.h:200
METAnalyzer::~METAnalyzer ( )
virtual

Destructor.

Definition at line 151 of file METAnalyzer.cc.

151  {
152  for (std::vector<GenericTriggerEventFlag *>::const_iterator it = triggerFolderEventFlag_.begin(); it!= triggerFolderEventFlag_.end(); it++) {
153  delete *it;
154  }
155  delete DCSFilter_;
156 }
std::vector< GenericTriggerEventFlag * > triggerFolderEventFlag_
Definition: METAnalyzer.h:198
JetMETDQMDCSFilter * DCSFilter_
Definition: METAnalyzer.h:267

Member Function Documentation

void METAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Get the analysis.

  • if(isTCMet_){

Implements edm::EDAnalyzer.

Definition at line 599 of file METAnalyzer.cc.

References FamosSequences_cff::caloJets, reco::JetCorrector::correction(), mvaPFMET_cff::corrector, gather_cfg::cout, reco::BeamHaloSummary::CSCTightHaloId(), cmsPerfPublish::DirName, reco::BeamHaloSummary::EcalTightHaloId(), edm::Event::getByToken(), JetCorrector::getJetCorrector(), reco::BeamHaloSummary::GlobalTightHaloId(), reco::BeamHaloSummary::HcalTightHaloId(), i, edm::HandleBase::isValid(), bTagSequences_cff::jetID, LogDebug, LogTrace, edm::EventBase::luminosityBlock(), caloMETBenchmarkGeneric_cfi::met, metname, NULL, pfJets_cff::pfJets, phi, edm::Handle< T >::product(), pileupReCalc_HLTpaths::scale, AlCaHLTBitMon_QueryRunRegistry::string, edm::triggerResults(), triggerResultsToken_(), and GoodVertex_cfg::vertexCollection.

599  {
600 
601 
602  // *** Fill lumisection ME
603  int myLuminosityBlock;
604  myLuminosityBlock = iEvent.luminosityBlock();
606  lumisecME=map_dijet_MEs["JetMET/lumisec"]; if(lumisecME && lumisecME->getRootObject()) lumisecME->Fill(myLuminosityBlock);
607  }
608 
609  if (myLuminosityBlock<LSBegin_) return;
610  if (myLuminosityBlock>LSEnd_ && LSEnd_>0) return;
611 
612  if (verbose_) std::cout << "METAnalyzer analyze" << std::endl;
613 
615 
616 
617 
618 
619  // ==========================================================
620  // Trigger information
621  //
622 // trigJetMB_=0;
623 // trigHighPtJet_=0;
624 // trigLowPtJet_=0;
625 // trigMinBias_=0;
626 // trigHighMET_=0;
627 // // _trig_LowMET=0;
628 // trigEle_=0;
629 // trigMuon_=0;
630 // trigPhysDec_=0;
631  std::vector<int> triggerFolderDecisions;
632  triggerFolderDecisions_ = std::vector<int> (triggerFolderEventFlag_.size(), 0);
633  // **** Get the TriggerResults container
635  iEvent.getByToken(triggerResultsToken_, triggerResults);
636 
637  if( triggerResults.isValid()) {
639  // Check how many HLT triggers are in triggerResults
640  int ntrigs = (*triggerResults).size();
641  if (verbose_) std::cout << "ntrigs=" << ntrigs << std::endl;
642  // If index=ntrigs, this HLT trigger doesn't exist in the HLT table for this data.
643  for (std::vector<GenericTriggerEventFlag *>::const_iterator it = triggerFolderEventFlag_.begin(); it!=triggerFolderEventFlag_.end();it++) {
644  unsigned pos = it - triggerFolderEventFlag_.begin();
645  bool fd = (*it)->accept(iEvent, iSetup);
646  triggerFolderDecisions_[pos] = fd;
647  }
648  allTriggerDecisions_.clear();
649  for (unsigned i=0;i<allTriggerNames_.size();++i) {
650  allTriggerDecisions_.push_back((*triggerResults).accept(i));
651 // std::cout<<"TR "<<(*triggerResults).size()<<" "<<(*triggerResults).accept(i)<<" "<<allTriggerNames_[i]<<std::endl;
652  }
653  }
654 
655  // ==========================================================
656  // MET information
657 
658  // **** Get the MET container
662 
663  //if(isTCMet_){
664  //iEvent.getByToken(tcMetToken_, tcmetcoll);
665  //if(!tcmetcoll.isValid()) return;
666  //}
667  if(isCaloMet_){
668  iEvent.getByToken(caloMetToken_, calometcoll);
669  if(!calometcoll.isValid()) return;
670  }
671  if(isPFMet_){
672  iEvent.getByToken(pfMetToken_, pfmetcoll);
673  if(!pfmetcoll.isValid()) return;
674  }
675 
676  const MET *met=NULL;
677  const PFMET *pfmet=NULL;
678  const CaloMET *calomet=NULL;
679  //if(isTCMet_){
680  //met=&(tcmetcoll->front());
681  //}
682  if(isPFMet_){
683  met=&(pfmetcoll->front());
684  pfmet=&(pfmetcoll->front());
685  }
686  if(isCaloMet_){
687  met=&(calometcoll->front());
688  calomet=&(calometcoll->front());
689  }
690 
691  LogTrace(metname)<<"[METAnalyzer] Call to the MET analyzer";
692 
693  // ==========================================================
694  // TCMET
695 
696  //if (/*isTCMet_ || */(isCaloMet_ && metCollectionLabel_.label() == "corMetGlobalMuons")) {
697 
698  //iEvent.getByToken(MuonToken_, muonHandle_);
699  //iEvent.getByToken(TrackToken_, trackHandle_);
700  //iEvent.getByToken(ElectronToken_, electronHandle_);
701  //iEvent.getByToken(BeamspotToken_, beamSpotHandle_);
702  //iEvent.getByToken(tcMETValueMapToken_,tcMetValueMapHandle_);
703 
704  //if(!muonHandle_.isValid()) edm::LogInfo("OutputInfo") << "falied to retrieve muon data require by MET Task";
705  //if(!trackHandle_.isValid()) edm::LogInfo("OutputInfo") << "falied to retrieve track data require by MET Task";
706  //if(!electronHandle_.isValid()) edm::LogInfo("OutputInfo") << "falied to retrieve electron data require by MET Task";
707  //if(!beamSpotHandle_.isValid()) edm::LogInfo("OutputInfo") << "falied to retrieve beam spot data require by MET Task";
708 
709  //beamSpot_ = ( beamSpotHandle_.isValid() ) ? beamSpotHandle_->position() : math::XYZPoint(0, 0, 0);
710  //}
711 
712  // ==========================================================
713  //
714 
716  iEvent.getByToken(HcalNoiseRBXToken_,HRBXCollection);
717  if (!HRBXCollection.isValid()) {
718  LogDebug("") << "METAnalyzer: Could not find HcalNoiseRBX Collection" << std::endl;
719  if (verbose_) std::cout << "METAnalyzer: Could not find HcalNoiseRBX Collection" << std::endl;
720  }
721 
722  edm::Handle<bool> HBHENoiseFilterResultHandle;
723  iEvent.getByToken(hbheNoiseFilterResultToken_, HBHENoiseFilterResultHandle);
724  bool HBHENoiseFilterResult = *HBHENoiseFilterResultHandle;
725  if (!HBHENoiseFilterResultHandle.isValid()) {
726  LogDebug("") << "METAnalyzer: Could not find HBHENoiseFilterResult" << std::endl;
727  if (verbose_) std::cout << "METAnalyzer: Could not find HBHENoiseFilterResult" << std::endl;
728  }
729 
730  // ==========================================================
731  bool bJetID = false;
732  bool bDiJetID = false;
733  // Jet ID -------------------------------------------------------
734  //
735 
739 
740  int collsize=-1;
741 
742  if (isCaloMet_){
743  iEvent.getByToken(caloJetsToken_, caloJets);
744  if (!caloJets.isValid()) {
745  LogDebug("") << "METAnalyzer: Could not find calojet product" << std::endl;
746  if (verbose_) std::cout << "METAnalyzer: Could not find calojet product" << std::endl;
747  }
748  collsize=caloJets->size();
749  }
751  //if (isTCMet_){
752  //iEvent.getByToken(jptJetsToken_, jptJets);
753  //if (!jptJets.isValid()) {
754  // LogDebug("") << "METAnalyzer: Could not find jptjet product" << std::endl;
755  // if (verbose_) std::cout << "METAnalyzer: Could not find jptjet product" << std::endl;
756  //}
757  //collsize=jptJets->size();
758  //}*/
759 
760  edm::Handle< edm::ValueMap<reco::JetID> >jetID_ValueMap_Handle;
761  if(/*isTCMet_ || */isCaloMet_){
762  if(!runcosmics_){
763  iEvent.getByToken(jetID_ValueMapToken_,jetID_ValueMap_Handle);
764  }
765  }
766 
767  if (isPFMet_){ iEvent.getByToken(pfJetsToken_, pfJets);
768  if (!pfJets.isValid()) {
769  LogDebug("") << "METAnalyzer: Could not find pfjet product" << std::endl;
770  if (verbose_) std::cout << "METAnalyzer: Could not find pfjet product" << std::endl;
771  }
772  collsize=pfJets->size();
773  }
774 
775  unsigned int ind1=-1;
776  double pt1=-1;
777  bool pass_jetID1=false;
778  unsigned int ind2=-1;
779  double pt2=-1;
780  bool pass_jetID2=false;
781 
782  //do loose jet ID-> check threshold on corrected jets
783  for (int ijet=0; ijet<collsize; ijet++) {
784  double pt_jet=-10;
785  double scale=1.;
786  bool iscleaned=false;
787  if (!jetCorrectionService_.empty()) {
789  if(isCaloMet_){
790  scale = corrector->correction((*caloJets)[ijet], iEvent, iSetup);
791  }
792  //if(isTCMet_){
793  //scale = corrector->correction((*jptJets)[ijet], iEvent, iSetup);
794  //}
795  if(isPFMet_){
796  scale = corrector->correction((*pfJets)[ijet], iEvent, iSetup);
797  }
798  }
799  if(isCaloMet_){
800  pt_jet=scale*(*caloJets)[ijet].pt();
801  if(pt_jet> ptThreshold_){
802  reco::CaloJetRef calojetref(caloJets, ijet);
803  if(!runcosmics_){
804  reco::JetID jetID = (*jetID_ValueMap_Handle)[calojetref];
805  iscleaned = jetIDFunctorLoose((*caloJets)[ijet], jetID);
806  }else{
807  iscleaned=true;
808  }
809  }
810  }
812  //if(isTCMet_){
813  //pt_jet=scale*(*jptJets)[ijet].pt();
814  //if(pt_jet> ptThreshold_){
815  // const edm::RefToBase<reco::Jet>& rawJet = (*jptJets)[ijet].getCaloJetRef();
816  // const reco::CaloJet *rawCaloJet = dynamic_cast<const reco::CaloJet*>(&*rawJet);
817  // reco::CaloJetRef const theCaloJetRef = (rawJet).castTo<reco::CaloJetRef>();
818  // if(!runcosmics_){
819  // reco::JetID jetID = (*jetID_ValueMap_Handle)[theCaloJetRef];
820  // iscleaned = jetIDFunctorLoose(*rawCaloJet, jetID);
821  // }else{
822  // iscleaned=true;
823  // }
824  //}
825  //}*/
826  if(isPFMet_){
827  pt_jet=scale*(*pfJets)[ijet].pt();
828  if(pt_jet> ptThreshold_){
829  iscleaned = pfjetIDFunctorLoose((*pfJets)[ijet]);
830  }
831  }
832  if(iscleaned){
833  bJetID=true;
834  }
835  if(pt_jet>pt1){
836  pt2=pt1;
837  ind2=ind1;
838  pass_jetID2=pass_jetID1;
839  pt1=pt_jet;
840  ind1=ijet;
841  pass_jetID1=iscleaned;
842  }else if (pt_jet>pt2){
843  pt2=pt_jet;
844  ind2=ijet;
845  pass_jetID2=iscleaned;
846  }
847  }
848  if(pass_jetID1 && pass_jetID2){
849  double dphi=-1.0;
850  if(isCaloMet_){
851  dphi=fabs((*caloJets)[ind1].phi()-(*caloJets)[ind2].phi());
852  }
854  //dphi=fabs((*jptJets)[ind1].phi()-(*jptJets)[ind2].phi());
855  //}*/
856  if(isPFMet_){
857  dphi=fabs((*pfJets)[ind1].phi()-(*pfJets)[ind2].phi());
858  }
859  if(dphi>acos(-1.)){
860  dphi=2*acos(-1.)-dphi;
861  }
862  if(dphi>2.7){
863  bDiJetID=true;
864  }
865  }
866 
867  // ==========================================================
868  // HCAL Noise filter
869 
870  bool bHBHENoiseFilter = HBHENoiseFilterResult;
871 
872  // ==========================================================
873  // Get BeamHaloSummary
874  edm::Handle<BeamHaloSummary> TheBeamHaloSummary ;
875  iEvent.getByToken(beamHaloSummaryToken_, TheBeamHaloSummary) ;
876 
877  if (!TheBeamHaloSummary.isValid()) {
878  if (verbose_) std::cout << "BeamHaloSummary doesn't exist" << std::endl;
879  }
880 
881  bool bBeamHaloID = true;
882 
883  if(!TheBeamHaloSummary.isValid()) {
884 
885  const BeamHaloSummary TheSummary = (*TheBeamHaloSummary.product() );
886 
887  if( !TheSummary.EcalTightHaloId() && !TheSummary.HcalTightHaloId() &&
888  !TheSummary.CSCTightHaloId() && !TheSummary.GlobalTightHaloId() )
889  bBeamHaloID = false;
890  }
891 
892  // ==========================================================
893  //Vertex information
894  Handle<VertexCollection> vertexHandle;
895  iEvent.getByToken(vertexToken_, vertexHandle);
896 
897  if (!vertexHandle.isValid()) {
898  LogDebug("") << "CaloMETAnalyzer: Could not find vertex collection" << std::endl;
899  if (verbose_) std::cout << "CaloMETAnalyzer: Could not find vertex collection" << std::endl;
900  }
901  numPV_ = 0;
902  if ( vertexHandle.isValid() ){
903  VertexCollection vertexCollection = *(vertexHandle.product());
904  numPV_ = vertexCollection.size();
905  }
906  bool bPrimaryVertex = (bypassAllPVChecks_ || (numPV_>0));
907  // ==========================================================
908 
910  iEvent.getByToken( gtToken_, gtReadoutRecord);
911 
912  if (!gtReadoutRecord.isValid()) {
913  LogDebug("") << "CaloMETAnalyzer: Could not find GT readout record" << std::endl;
914  if (verbose_) std::cout << "CaloMETAnalyzer: Could not find GT readout record product" << std::endl;
915  }
916  // DCS Filter
917  bool bDCSFilter = (bypassAllDCSChecks_ || DCSFilter_->filter(iEvent, iSetup));
918  // ==========================================================
919  // Reconstructed MET Information - fill MonitorElements
920 
921  for (std::vector<std::string>::const_iterator ic = folderNames_.begin();
922  ic != folderNames_.end(); ic++){
923  if ((*ic=="Uncleaned") &&(isCaloMet_ || bPrimaryVertex)) fillMESet(iEvent, DirName+"/"+*ic, *met,*pfmet,*calomet,map_dijet_MEs);
924  //take two lines out for first check
925  if ((*ic=="Cleaned") &&bDCSFilter&&bHBHENoiseFilter&&bPrimaryVertex&&bBeamHaloID&&bJetID) fillMESet(iEvent, DirName+"/"+*ic, *met,*pfmet,*calomet,map_dijet_MEs);
926  if ((*ic=="DiJet" ) &&bDCSFilter&&bHBHENoiseFilter&&bPrimaryVertex&&bBeamHaloID&&bDiJetID) fillMESet(iEvent, DirName+"/"+*ic, *met,*pfmet,*calomet,map_dijet_MEs);
927  }
928 }
#define LogDebug(id)
edm::EDGetTokenT< reco::CaloMETCollection > caloMetToken_
Definition: METAnalyzer.h:154
std::map< std::string, MonitorElement * > map_dijet_MEs
Definition: METAnalyzer.h:401
const bool EcalTightHaloId() const
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< reco::CaloJetCollection > caloJetsToken_
Definition: METAnalyzer.h:145
void fillMESet(const edm::Event &, std::string, const reco::MET &, const reco::PFMET &, const reco::CaloMET &, std::map< std::string, MonitorElement * > &)
Definition: METAnalyzer.cc:932
std::vector< int > allTriggerDecisions_
Definition: METAnalyzer.h:195
bool bypassAllDCSChecks_
Definition: METAnalyzer.h:240
const bool HcalTightHaloId() const
std::vector< int > triggerFolderDecisions_
Definition: METAnalyzer.h:201
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
double ptThreshold_
Definition: METAnalyzer.h:176
const std::string metname
bool runcosmics_
Definition: METAnalyzer.h:241
bool isCaloMet_
Definition: METAnalyzer.h:403
edm::EDGetTokenT< reco::PFMETCollection > pfMetToken_
Definition: METAnalyzer.h:153
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:59
#define NULL
Definition: scimark2.h:8
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:47
std::string FolderName_
Definition: METAnalyzer.h:133
const bool GlobalTightHaloId() const
std::vector< GenericTriggerEventFlag * > triggerFolderEventFlag_
Definition: METAnalyzer.h:198
Jet ID object.
Definition: JetID.h:16
tuple vertexCollection
void Fill(long long x)
const bool CSCTightHaloId() const
MonitorElement * lumisecME
Definition: METAnalyzer.h:281
edm::EDGetTokenT< std::vector< reco::Vertex > > vertexToken_
Definition: METAnalyzer.h:143
bool fill_met_high_level_histo
Definition: METAnalyzer.h:407
edm::InputTag metCollectionLabel_
Definition: METAnalyzer.h:135
std::string jetCorrectionService_
Definition: METAnalyzer.h:174
tuple corrector
Definition: mvaPFMET_cff.py:50
Definition: MET.h:42
edm::EDGetTokenT< bool > hbheNoiseFilterResultToken_
Definition: METAnalyzer.h:149
std::vector< std::string > folderNames_
Definition: METAnalyzer.h:269
static std::string const triggerResults
Definition: EdmProvDump.cc:41
edm::EDGetTokenT< reco::HcalNoiseRBXCollection > HcalNoiseRBXToken_
Definition: METAnalyzer.h:155
bool isValid() const
Definition: HandleBase.h:76
std::vector< std::string > allTriggerNames_
Definition: METAnalyzer.h:194
bool filter(const edm::Event &evt, const edm::EventSetup &es)
#define LogTrace(id)
edm::EDGetTokenT< edm::TriggerResults > triggerResultsToken_
Definition: METAnalyzer.h:190
edm::EDGetTokenT< reco::PFJetCollection > pfJetsToken_
Definition: METAnalyzer.h:146
JetMETDQMDCSFilter * DCSFilter_
Definition: METAnalyzer.h:267
TObject * getRootObject(void) const
T const * product() const
Definition: Handle.h:81
PFJetIDSelectionFunctor pfjetIDFunctorLoose
Definition: METAnalyzer.h:172
static const JetCorrector * getJetCorrector(const std::string &fName, const edm::EventSetup &fSetup)
retrieve corrector from the event setup. troughs exception if something is missing ...
Definition: JetCorrector.cc:50
std::string const & label() const
Definition: InputTag.h:42
edm::EDGetTokenT< L1GlobalTriggerReadoutRecord > gtToken_
Definition: METAnalyzer.h:144
bool bypassAllPVChecks_
Definition: METAnalyzer.h:239
tuple cout
Definition: gather_cfg.py:121
edm::EDGetTokenT< edm::ValueMap< reco::JetID > > jetID_ValueMapToken_
Definition: METAnalyzer.h:169
JetIDSelectionFunctor jetIDFunctorLoose
Definition: METAnalyzer.h:171
edm::EDGetTokenT< reco::BeamHaloSummary > beamHaloSummaryToken_
Definition: METAnalyzer.h:150
tuple pfJets
Definition: pfJets_cff.py:8
Definition: DDAxes.h:10
void METAnalyzer::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &   
)
overridevirtual

Inizialize parameters for histo binning.

Implements thread_unsafe::DQMEDAnalyzer.

Definition at line 159 of file METAnalyzer.cc.

References cmsPerfPublish::DirName, DQMStore::IBooker::setCurrentFolder(), and AlCaHLTBitMon_QueryRunRegistry::string.

161  {
163  ibooker.setCurrentFolder(DirName);
164 
165  if(!folderNames_.empty()){
166  folderNames_.clear();
167  }
168 
169  folderNames_.push_back("Uncleaned");
170  folderNames_.push_back("Cleaned");
171  folderNames_.push_back("DiJet");
172 
173  for (std::vector<std::string>::const_iterator ic = folderNames_.begin();
174  ic != folderNames_.end(); ic++){
175  bookMESet(DirName+"/"+*ic, ibooker,map_dijet_MEs);
176  }
177 }
std::map< std::string, MonitorElement * > map_dijet_MEs
Definition: METAnalyzer.h:401
std::string FolderName_
Definition: METAnalyzer.h:133
edm::InputTag metCollectionLabel_
Definition: METAnalyzer.h:135
void bookMESet(std::string, DQMStore::IBooker &, std::map< std::string, MonitorElement * > &)
Definition: METAnalyzer.cc:181
std::vector< std::string > folderNames_
Definition: METAnalyzer.h:269
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
std::string const & label() const
Definition: InputTag.h:42
void METAnalyzer::bookMESet ( std::string  DirName,
DQMStore::IBooker ibooker,
std::map< std::string, MonitorElement * > &  map_of_MEs 
)
private

Definition at line 181 of file METAnalyzer.cc.

References i.

182 {
183  bool bLumiSecPlot=fill_met_high_level_histo;
184  //if (DirName.find("Uncleaned")!=std::string::npos)bLumiSecPlot=true;//now defined on configlevel
185  bookMonitorElement(DirName,ibooker,map_of_MEs,bLumiSecPlot);
186 
187  if (DirName.find("Cleaned")!=std::string::npos) {
188  for (unsigned int i = 0; i<triggerFolderEventFlag_.size(); i++) {
189  if (triggerFolderEventFlag_[i]->on()) {
190  bookMonitorElement(DirName+"/"+triggerFolderLabels_[i],ibooker,map_of_MEs,false);
191  }
192  }
193  }
194 }
int i
Definition: DBlmapReader.cc:9
void bookMonitorElement(std::string, DQMStore::IBooker &, std::map< std::string, MonitorElement * > &, bool)
Definition: METAnalyzer.cc:197
std::vector< GenericTriggerEventFlag * > triggerFolderEventFlag_
Definition: METAnalyzer.h:198
bool fill_met_high_level_histo
Definition: METAnalyzer.h:407
std::vector< std::string > triggerFolderLabels_
Definition: METAnalyzer.h:200
void METAnalyzer::bookMonitorElement ( std::string  DirName,
DQMStore::IBooker ibooker,
std::map< std::string, MonitorElement * > &  map_of_MEs,
bool  bLumiSecPlot = false 
)
private

Definition at line 197 of file METAnalyzer.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::book2D(), DQMStore::IBooker::bookProfile(), gather_cfg::cout, MonitorElement::setAxisTitle(), and DQMStore::IBooker::setCurrentFolder().

198 {
199  if (verbose_) std::cout << "bookMonitorElement " << DirName << std::endl;
200 
201  ibooker.setCurrentFolder(DirName);
202 
203  hTrigger = ibooker.book1D("triggerResults", "triggerResults", 500, 0, 500);
204  hMEx = ibooker.book1D("MEx", "MEx", 200, -500, 500);
205  hMEy = ibooker.book1D("MEy", "MEy", 200, -500, 500);
206  hMET = ibooker.book1D("MET", "MET", 200, 0, 1000);
207  hSumET = ibooker.book1D("SumET", "SumET", 400, 0, 4000);
208  hMETSig = ibooker.book1D("METSig", "METSig", 51, 0, 51);
209  hMETPhi = ibooker.book1D("METPhi", "METPhi", 60, -3.2, 3.2);
210  hMET_logx = ibooker.book1D("MET_logx", "MET_logx", 40, -1, 7);
211  hSumET_logx = ibooker.book1D("SumET_logx", "SumET_logx", 40, -1, 7);
212 
213  hMEx ->setAxisTitle("MEx [GeV]", 1);
214  hMEy ->setAxisTitle("MEy [GeV]", 1);
215  hMET ->setAxisTitle("MET [GeV]", 1);
216  hSumET ->setAxisTitle("SumET [GeV]", 1);
217  hMETSig ->setAxisTitle("METSig", 1);
218  hMETPhi ->setAxisTitle("METPhi [rad]", 1);
219  hMET_logx ->setAxisTitle("log(MET) [GeV]", 1);
220  hSumET_logx->setAxisTitle("log(SumET) [GeV]", 1);
221 
222  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"triggerResults",hTrigger));
223  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"MEx",hMEx));
224  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"MEy",hMEy));
225  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"MET",hMET));
226  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"SumET",hSumET));
227  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"METSig",hMETSig));
228  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"METPhi",hMETPhi));
229  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"MET_logx",hMET_logx));
230  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"SumET_logx",hSumET_logx));
231 
232  // Book NPV profiles --> would some of these profiles be interesting for other MET types too
233  //----------------------------------------------------------------------------
234  meMEx_profile = ibooker.bookProfile("MEx_profile", "met.px()", nbinsPV_, nPVMin_, nPVMax_, 200, -500, 500);
235  meMEy_profile = ibooker.bookProfile("MEy_profile", "met.py()", nbinsPV_, nPVMin_, nPVMax_, 200, -500, 500);
236  meMET_profile = ibooker.bookProfile("MET_profile", "met.pt()", nbinsPV_, nPVMin_, nPVMax_, 200, 0, 1000);
237  meSumET_profile = ibooker.bookProfile("SumET_profile", "met.sumEt()", nbinsPV_, nPVMin_, nPVMax_, 400, 0, 4000);
238  // Set NPV profiles x-axis title
239  //----------------------------------------------------------------------------
240  meMEx_profile ->setAxisTitle("nvtx", 1);
241  meMEy_profile ->setAxisTitle("nvtx", 1);
242  meMET_profile ->setAxisTitle("nvtx", 1);
243  meSumET_profile->setAxisTitle("nvtx", 1);
244 
245  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"MEx_profile",meMEx_profile));
246  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"MEy_profile",meMEy_profile));
247  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"MET_profile",meMET_profile));
248  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"SumET_profile",meSumET_profile));
249 
250 
251  if(isCaloMet_){
252  hCaloMaxEtInEmTowers = ibooker.book1D("CaloMaxEtInEmTowers", "CaloMaxEtInEmTowers" ,100,0,2000);
253  hCaloMaxEtInEmTowers->setAxisTitle("Et(Max) in EM Tower [GeV]",1);
254  hCaloMaxEtInHadTowers = ibooker.book1D("CaloMaxEtInHadTowers", "CaloMaxEtInHadTowers" ,100,0,2000);
255  hCaloMaxEtInHadTowers->setAxisTitle("Et(Max) in Had Tower [GeV]",1);
256 
257  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloMaxEtInEmTowers",hCaloMaxEtInEmTowers));
258  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloMaxEtInHadTowers",hCaloMaxEtInHadTowers));
259 
260  hCaloHadEtInHB = ibooker.book1D("CaloHadEtInHB","CaloHadEtInHB",100,0,2000);
261  hCaloHadEtInHB->setAxisTitle("Had Et [GeV]",1);
262  hCaloHadEtInHO = ibooker.book1D("CaloHadEtInHO","CaloHadEtInHO",25,0,500);
263  hCaloHadEtInHO->setAxisTitle("Had Et [GeV]",1);
264  hCaloHadEtInHE = ibooker.book1D("CaloHadEtInHE","CaloHadEtInHE",100,0,2000);
265  hCaloHadEtInHE->setAxisTitle("Had Et [GeV]",1);
266  hCaloHadEtInHF = ibooker.book1D("CaloHadEtInHF","CaloHadEtInHF",50,0,1000);
267  hCaloHadEtInHF->setAxisTitle("Had Et [GeV]",1);
268  hCaloEmEtInHF = ibooker.book1D("CaloEmEtInHF" ,"CaloEmEtInHF" ,25,0,500);
269  hCaloEmEtInHF->setAxisTitle("EM Et [GeV]",1);
270  hCaloEmEtInEE = ibooker.book1D("CaloEmEtInEE" ,"CaloEmEtInEE" ,50,0,1000);
271  hCaloEmEtInEE->setAxisTitle("EM Et [GeV]",1);
272  hCaloEmEtInEB = ibooker.book1D("CaloEmEtInEB" ,"CaloEmEtInEB" ,100,0,2000);
273  hCaloEmEtInEB->setAxisTitle("EM Et [GeV]",1);
274 
275  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloHadEtInHO",hCaloHadEtInHO));
276  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloHadEtInHF",hCaloHadEtInHF));
277  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloHadEtInHE",hCaloHadEtInHE));
278  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloHadEtInHB",hCaloHadEtInHB));
279  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloEmEtInHF",hCaloEmEtInHF));
280  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloEmEtInEE",hCaloEmEtInEE));
281  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloEmEtInEB",hCaloEmEtInEB));
282 
283  hCaloMETPhi020 = ibooker.book1D("CaloMETPhi020", "CaloMETPhi020", 60, -3.2, 3.2);
284  hCaloMETPhi020 ->setAxisTitle("METPhi [rad] (MET>20 GeV)", 1);
285 
286  //hCaloMaxEtInEmTowers = ibooker.book1D("CaloMaxEtInEmTowers", "CaloMaxEtInEmTowers" ,100,0,2000);
287  //hCaloMaxEtInEmTowers->setAxisTitle("Et(Max) in EM Tower [GeV]",1);
288  //hCaloMaxEtInHadTowers = ibooker.book1D("CaloMaxEtInHadTowers", "CaloMaxEtInHadTowers" ,100,0,2000);
289  //hCaloMaxEtInHadTowers->setAxisTitle("Et(Max) in Had Tower [GeV]",1);
290  hCaloEtFractionHadronic = ibooker.book1D("CaloEtFractionHadronic","CaloEtFractionHadronic",100,0,1);
291  hCaloEtFractionHadronic->setAxisTitle("Hadronic Et Fraction",1);
292  hCaloEmEtFraction = ibooker.book1D("CaloEmEtFraction", "CaloEmEtFraction" ,100,0,1);
293  hCaloEmEtFraction->setAxisTitle("EM Et Fraction",1);
294 
295  //hCaloEmEtFraction002 = ibooker.book1D("CaloEmEtFraction002", "CaloEmEtFraction002" ,100,0,1);
296  //hCaloEmEtFraction002->setAxisTitle("EM Et Fraction (MET>2 GeV)",1);
297  //hCaloEmEtFraction010 = ibooker.book1D("CaloEmEtFraction010", "CaloEmEtFraction010" ,100,0,1);
298  //hCaloEmEtFraction010->setAxisTitle("EM Et Fraction (MET>10 GeV)",1);
299  hCaloEmEtFraction020 = ibooker.book1D("CaloEmEtFraction020", "CaloEmEtFraction020" ,100,0,1);
300  hCaloEmEtFraction020->setAxisTitle("EM Et Fraction (MET>20 GeV)",1);
301 
302  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloMETPhi020",hCaloMETPhi020));
303  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloEtFractionHadronic",hCaloEtFractionHadronic));
304  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloEmEtFraction", hCaloEmEtFraction));
305  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloEmEtFraction020",hCaloEmEtFraction020));
306 
307 
308 
309  //if (metCollectionLabel_.label() == "corMetGlobalMuons" ) {
310  //hCalomuPt = ibooker.book1D("CalomuonPt", "CalomuonPt", 50, 0, 500);
311  //hCalomuEta = ibooker.book1D("CalomuonEta", "CalomuonEta", 60, -3.0, 3.0);
312  //hCalomuNhits = ibooker.book1D("CalomuonNhits", "CalomuonNhits", 50, 0, 50);
313  //hCalomuChi2 = ibooker.book1D("CalomuonNormalizedChi2", "CalomuonNormalizedChi2", 20, 0, 20);
314  //hCalomuD0 = ibooker.book1D("CalomuonD0", "CalomuonD0", 50, -1, 1);
315  //hCaloMExCorrection = ibooker.book1D("CaloMExCorrection", "CaloMExCorrection", 100, -500.0,500.0);
316  //hCaloMEyCorrection = ibooker.book1D("CaloMEyCorrection", "CaloMEyCorrection", 100, -500.0,500.0);
317  //hCaloMuonCorrectionFlag = ibooker.book1D("CaloCorrectionFlag","CaloCorrectionFlag", 5, -0.5, 4.5);
318 
319  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CalomuonPt",hCalomuPt));
320  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CalomuonEta",hCalomuEta));
321  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CalomuonNhit",hCalomuNhits));
322  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CalomuonNormalizedChi2",hCalomuChi2));
323  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloMExCorrection",hCaloMExCorrection));
324  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloMEyCorrection",hCaloMEyCorrection));
325  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CaloCorrectionFlag",hCaloMuonCorrectionFlag));
326  //}
327 
328  }
329 
330  if(isPFMet_){
331  mePhotonEtFraction = ibooker.book1D("PfPhotonEtFraction", "pfmet.photonEtFraction()", 50, 0, 1);
332  mePhotonEt = ibooker.book1D("PfPhotonEt", "pfmet.photonEt()", 100, 0, 1000);
333  meNeutralHadronEtFraction = ibooker.book1D("PfNeutralHadronEtFraction", "pfmet.neutralHadronEtFraction()", 50, 0, 1);
334  meNeutralHadronEt = ibooker.book1D("PfNeutralHadronEt", "pfmet.neutralHadronEt()", 100, 0, 1000);
335  meElectronEtFraction = ibooker.book1D("PfElectronEtFraction", "pfmet.electronEtFraction()", 50, 0, 1);
336  meElectronEt = ibooker.book1D("PfElectronEt", "pfmet.electronEt()", 100, 0, 1000);
337  meChargedHadronEtFraction = ibooker.book1D("PfChargedHadronEtFraction", "pfmet.chargedHadronEtFraction()", 50, 0, 1);
338  meChargedHadronEt = ibooker.book1D("PfChargedHadronEt", "pfmet.chargedHadronEt()", 100, 0, 1000);
339  meMuonEtFraction = ibooker.book1D("PfMuonEtFraction", "pfmet.muonEtFraction()", 50, 0, 1);
340  meMuonEt = ibooker.book1D("PfMuonEt", "pfmet.muonEt()", 100, 0, 1000);
341  meHFHadronEtFraction = ibooker.book1D("PfHFHadronEtFraction", "pfmet.HFHadronEtFraction()", 50, 0, 1);
342  meHFHadronEt = ibooker.book1D("PfHFHadronEt", "pfmet.HFHadronEt()", 100, 0, 1000);
343  meHFEMEtFraction = ibooker.book1D("PfHFEMEtFraction", "pfmet.HFEMEtFraction()", 50, 0, 1);
344  meHFEMEt = ibooker.book1D("PfHFEMEt", "pfmet.HFEMEt()", 100, 0, 1000);
345 
346  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfPhotonEtFraction" ,mePhotonEtFraction));
347  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfPhotonEt" ,mePhotonEt));
348  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfNeutralHadronEtFraction",meNeutralHadronEtFraction));
349  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfNeutralHadronEt" ,meNeutralHadronEt));
350  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfElectronEtFraction" ,meElectronEtFraction));
351  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfElectronEt" ,meElectronEt));
352  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfChargedHadronEtFraction",meChargedHadronEtFraction));
353  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfChargedHadronEt" ,meChargedHadronEt));
354  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfMuonEtFraction" ,meMuonEtFraction));
355  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfMuonEt" ,meMuonEt));
356  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfHFHadronEtFraction" ,meHFHadronEtFraction));
357  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfHFHadronEt" ,meHFHadronEt));
358  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfHFEMEtFraction" ,meHFEMEtFraction));
359  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfHFEMEt" ,meHFEMEt));
360 
361  mePhotonEtFraction_profile = ibooker.bookProfile("PfPhotonEtFraction_profile", "pfmet.photonEtFraction()", nbinsPV_, nPVMin_, nPVMax_, 50, 0, 1);
362  mePhotonEt_profile = ibooker.bookProfile("PfPhotonEt_profile", "pfmet.photonEt()", nbinsPV_, nPVMin_, nPVMax_, 100, 0, 1000);
363  meNeutralHadronEtFraction_profile = ibooker.bookProfile("PfNeutralHadronEtFraction_profile", "pfmet.neutralHadronEtFraction()", nbinsPV_, nPVMin_, nPVMax_, 50, 0, 1);
364  meNeutralHadronEt_profile = ibooker.bookProfile("PfNeutralHadronEt_profile", "pfmet.neutralHadronEt()", nbinsPV_, nPVMin_, nPVMax_, 100, 0, 1000);
365  meElectronEtFraction_profile = ibooker.bookProfile("PfElectronEtFraction_profile", "pfmet.electronEtFraction()", nbinsPV_, nPVMin_, nPVMax_, 50, 0, 1);
366  meElectronEt_profile = ibooker.bookProfile("PfElectronEt_profile", "pfmet.electronEt()", nbinsPV_, nPVMin_, nPVMax_, 100, 0, 1000);
367  meChargedHadronEtFraction_profile = ibooker.bookProfile("PfChargedHadronEtFraction_profile", "pfmet.chargedHadronEtFraction()", nbinsPV_, nPVMin_, nPVMax_, 50, 0, 1);
368  meChargedHadronEt_profile = ibooker.bookProfile("PfChargedHadronEt_profile", "pfmet.chargedHadronEt()", nbinsPV_, nPVMin_, nPVMax_, 100, 0, 1000);
369  meMuonEtFraction_profile = ibooker.bookProfile("PfMuonEtFraction_profile", "pfmet.muonEtFraction()", nbinsPV_, nPVMin_, nPVMax_, 50, 0, 1);
370  meMuonEt_profile = ibooker.bookProfile("PfMuonEt_profile", "pfmet.muonEt()", nbinsPV_, nPVMin_, nPVMax_, 100, 0, 1000);
371  meHFHadronEtFraction_profile = ibooker.bookProfile("PfHFHadronEtFraction_profile", "pfmet.HFHadronEtFraction()", nbinsPV_, nPVMin_, nPVMax_, 50, 0, 1);
372  meHFHadronEt_profile = ibooker.bookProfile("PfHFHadronEt_profile", "pfmet.HFHadronEt()", nbinsPV_, nPVMin_, nPVMax_, 100, 0, 1000);
373  meHFEMEtFraction_profile = ibooker.bookProfile("PfHFEMEtFraction_profile", "pfmet.HFEMEtFraction()", nbinsPV_, nPVMin_, nPVMax_, 50, 0, 1);
374  meHFEMEt_profile = ibooker.bookProfile("PfHFEMEt_profile", "pfmet.HFEMEt()", nbinsPV_, nPVMin_, nPVMax_, 100, 0, 1000);
375 
377  mePhotonEt_profile ->setAxisTitle("nvtx", 1);
381  meElectronEt_profile ->setAxisTitle("nvtx", 1);
385  meMuonEt_profile ->setAxisTitle("nvtx", 1);
387  meHFHadronEt_profile ->setAxisTitle("nvtx", 1);
389  meHFEMEt_profile ->setAxisTitle("nvtx", 1);
390 
391  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfPhotonEtFraction_profile" ,mePhotonEtFraction_profile));
392  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfPhotonEt_profile" ,mePhotonEt_profile));
393  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfNeutralHadronEtFraction_profile" ,meNeutralHadronEtFraction_profile));
394  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfNeutralHadronEt_profile" ,meNeutralHadronEt_profile));
395  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfElectronEtFraction_profile" ,meElectronEtFraction_profile));
396  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfElectronEt_profile" ,meElectronEt_profile));
397  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfChargedHadronEtFraction_profile" ,meChargedHadronEtFraction_profile));
398  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfChargedHadronEt_profile" ,meChargedHadronEt_profile));
399  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfMuonEtFraction_profile" ,meMuonEtFraction_profile));
400  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfMuonEt_profile" ,meMuonEt_profile));
401  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfHFHadronEtFraction_profile" ,meHFHadronEtFraction_profile));
402  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfHFHadronEt_profile" ,meHFHadronEt_profile));
403  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfHFEMEtFraction_profile" ,meHFEMEtFraction_profile));
404  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"PfHFEMEt_profile" ,meHFEMEt_profile));
405  }
406 
407  if (isCaloMet_){
408  if (fill_met_high_level_histo){//now configurable in python file
409  hMExLS = ibooker.book2D("MExLS","MEx_LS",200,-200,200,250,0.,2500.);
410  hMExLS->setAxisTitle("MEx [GeV]",1);
411  hMExLS->setAxisTitle("Lumi Section",2);
412  hMExLS->getTH2F()->SetOption("colz");
413  hMEyLS = ibooker.book2D("MEyLS","MEy_LS",200,-200,200,250,0.,2500.);
414  hMEyLS->setAxisTitle("MEy [GeV]",1);
415  hMEyLS->setAxisTitle("Lumi Section",2);
416  hMEyLS->getTH2F()->SetOption("colz");
417  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"MExLS",hMExLS));
418  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"MEyLS",hMEyLS));
419  }
420  }
421 
422  //if (isTCMet_) {
423  //htrkPt = ibooker.book1D("trackPt", "trackPt", 50, 0, 500);
424  //htrkEta = ibooker.book1D("trackEta", "trackEta", 60, -3.0, 3.0);
425  //htrkNhits = ibooker.book1D("trackNhits", "trackNhits", 50, 0, 50);
426  //htrkChi2 = ibooker.book1D("trackNormalizedChi2", "trackNormalizedChi2", 20, 0, 20);
427  //htrkD0 = ibooker.book1D("trackD0", "trackd0", 50, -1, 1);
428  //helePt = ibooker.book1D("electronPt", "electronPt", 50, 0, 500);
429  //heleEta = ibooker.book1D("electronEta", "electronEta", 60, -3.0, 3.0);
430  //heleHoE = ibooker.book1D("electronHoverE", "electronHoverE", 25, 0, 0.5);
431  //hmuPt = ibooker.book1D("muonPt", "muonPt", 50, 0, 500);
432  //hmuEta = ibooker.book1D("muonEta", "muonEta", 60, -3.0, 3.0);
433  //hmuNhits = ibooker.book1D("muonNhits", "muonNhits", 50, 0, 50);
434  //hmuChi2 = ibooker.book1D("muonNormalizedChi2", "muonNormalizedChi2", 20, 0, 20);
435  //hmuD0 = ibooker.book1D("muonD0", "muonD0", 50, -1, 1);
436 
437  //hMExCorrection = ibooker.book1D("MExCorrection", "MExCorrection", 100, -500.0,500.0);
438  //hMEyCorrection = ibooker.book1D("MEyCorrection", "MEyCorrection", 100, -500.0,500.0);
439  //hMuonCorrectionFlag = ibooker.book1D("CorrectionFlag","CorrectionFlag", 5, -0.5, 4.5);
440 
441  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"trackPt" ,htrkPt));
442  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"trackEta" ,htrkEta));
443  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"trackNhits",htrkNhits));
444  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"trackNormalizedChi2" ,htrkChi2));
445  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"trackD0" ,htrkD0));
446  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"electronPt" ,helePt));
447  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"electronEta" ,heleEta));
448  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"electronHoverE",heleHoE));
449  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"muonPt" ,hmuPt));
450  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"muonEta" ,hmuEta));
451  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"muonNhits",hmuNhits));
452  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"muonNormalizedChi2" ,hmuChi2));
453  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"muonD0" ,hmuD0));
454  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"MExCorrection" ,hMExCorrection));
455  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"MEyCorrection" ,hMEyCorrection));
456  //map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"CorrectionFlag" ,hMuonCorrectionFlag));
457  //}
458 
459  hMETRate = ibooker.book1D("METRate", "METRate", 200, 0, 1000);
460  map_of_MEs.insert(std::pair<std::string,MonitorElement*>(DirName+"/"+"METRate",hMETRate));
461 
462 
463  ibooker.setCurrentFolder("JetMET");
464  lumisecME = ibooker.book1D("lumisec", "lumisec", 2500, 0., 2500.);
465  map_of_MEs.insert(std::pair<std::string,MonitorElement*>("JetMET/lumisec",lumisecME));
466 
467 }
MonitorElement * meElectronEt_profile
Definition: METAnalyzer.h:391
MonitorElement * hMExLS
Definition: METAnalyzer.h:292
MonitorElement * hSumET
Definition: METAnalyzer.h:290
MonitorElement * meNeutralHadronEtFraction
Definition: METAnalyzer.h:365
MonitorElement * meNeutralHadronEtFraction_profile
Definition: METAnalyzer.h:388
MonitorElement * meElectronEtFraction_profile
Definition: METAnalyzer.h:390
MonitorElement * hMEy
Definition: METAnalyzer.h:285
MonitorElement * hCaloMETPhi020
Definition: METAnalyzer.h:299
MonitorElement * mePhotonEt
Definition: METAnalyzer.h:364
MonitorElement * meChargedHadronEtFraction_profile
Definition: METAnalyzer.h:392
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:157
MonitorElement * hCaloHadEtInHF
Definition: METAnalyzer.h:313
bool isCaloMet_
Definition: METAnalyzer.h:403
MonitorElement * hMEx
Definition: METAnalyzer.h:284
MonitorElement * hCaloHadEtInHE
Definition: METAnalyzer.h:312
MonitorElement * hCaloEtFractionHadronic
Definition: METAnalyzer.h:303
MonitorElement * hMEyLS
Definition: METAnalyzer.h:293
double nPVMin_
Definition: METAnalyzer.h:232
MonitorElement * meHFHadronEtFraction
Definition: METAnalyzer.h:373
MonitorElement * hMET_logx
Definition: METAnalyzer.h:295
MonitorElement * hCaloEmEtInHF
Definition: METAnalyzer.h:314
MonitorElement * meHFEMEtFraction_profile
Definition: METAnalyzer.h:398
MonitorElement * hCaloHadEtInHO
Definition: METAnalyzer.h:311
MonitorElement * meChargedHadronEtFraction
Definition: METAnalyzer.h:369
MonitorElement * meMEx_profile
Definition: METAnalyzer.h:381
MonitorElement * hMET
Definition: METAnalyzer.h:288
MonitorElement * hTrigger
Definition: METAnalyzer.h:282
MonitorElement * lumisecME
Definition: METAnalyzer.h:281
MonitorElement * hCaloMaxEtInHadTowers
Definition: METAnalyzer.h:302
bool fill_met_high_level_histo
Definition: METAnalyzer.h:407
MonitorElement * hMETPhi
Definition: METAnalyzer.h:289
MonitorElement * hCaloEmEtInEE
Definition: METAnalyzer.h:315
MonitorElement * meMuonEt
Definition: METAnalyzer.h:372
MonitorElement * meHFEMEt_profile
Definition: METAnalyzer.h:399
MonitorElement * hSumET_logx
Definition: METAnalyzer.h:296
MonitorElement * meMuonEt_profile
Definition: METAnalyzer.h:395
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * hCaloEmEtFraction
Definition: METAnalyzer.h:304
MonitorElement * meChargedHadronEt
Definition: METAnalyzer.h:370
MonitorElement * hMETSig
Definition: METAnalyzer.h:287
MonitorElement * meChargedHadronEt_profile
Definition: METAnalyzer.h:393
MonitorElement * meMuonEtFraction
Definition: METAnalyzer.h:371
MonitorElement * mePhotonEtFraction
Definition: METAnalyzer.h:363
MonitorElement * hCaloMaxEtInEmTowers
Definition: METAnalyzer.h:301
MonitorElement * mePhotonEt_profile
Definition: METAnalyzer.h:387
MonitorElement * meElectronEt
Definition: METAnalyzer.h:368
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
MonitorElement * meMET_profile
Definition: METAnalyzer.h:383
MonitorElement * meHFHadronEt
Definition: METAnalyzer.h:374
MonitorElement * hCaloEmEtInEB
Definition: METAnalyzer.h:316
MonitorElement * hCaloHadEtInHB
Definition: METAnalyzer.h:310
double nPVMax_
Definition: METAnalyzer.h:233
MonitorElement * meHFHadronEtFraction_profile
Definition: METAnalyzer.h:396
MonitorElement * meHFEMEtFraction
Definition: METAnalyzer.h:375
MonitorElement * hCaloEmEtFraction020
Definition: METAnalyzer.h:308
MonitorElement * meHFEMEt
Definition: METAnalyzer.h:376
MonitorElement * meMEy_profile
Definition: METAnalyzer.h:382
tuple cout
Definition: gather_cfg.py:121
MonitorElement * meElectronEtFraction
Definition: METAnalyzer.h:367
MonitorElement * meSumET_profile
Definition: METAnalyzer.h:384
TH2F * getTH2F(void) const
MonitorElement * meHFHadronEt_profile
Definition: METAnalyzer.h:397
MonitorElement * hMETRate
Definition: METAnalyzer.h:226
MonitorElement * meNeutralHadronEt_profile
Definition: METAnalyzer.h:389
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
MonitorElement * mePhotonEtFraction_profile
Definition: METAnalyzer.h:386
MonitorElement * meMuonEtFraction_profile
Definition: METAnalyzer.h:394
MonitorElement * meNeutralHadronEt
Definition: METAnalyzer.h:366
void METAnalyzer::dqmBeginRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
virtual

Initialize run-based parameters.

Reimplemented from thread_unsafe::DQMEDAnalyzer.

Definition at line 470 of file METAnalyzer.cc.

References gather_cfg::cout, i, Parameters::parameters, and AlCaHLTBitMon_QueryRunRegistry::string.

471 {
472 
473 // std::cout << "Run " << iRun.run() << " hltconfig.init "
474 // << hltConfig_.init(iRun,iSetup,triggerResultsLabel_.process(),changed_) << " length: "<<hltConfig_.triggerNames().size()<<" changed "<<changed_<<std::endl;
475  bool changed(true);
476  if (hltConfig_.init(iRun,iSetup,triggerResultsLabel_.process(),changed)) {
477  if (changed) {
478  hltConfig_.dump("ProcessName");
479  hltConfig_.dump("GlobalTag");
480  hltConfig_.dump("TableName");
481 // hltConfig_.dump("Streams");
482 // hltConfig_.dump("Datasets");
483 // hltConfig_.dump("PrescaleTable");
484 // hltConfig_.dump("ProcessPSet");
485  }
486  } else {
487  if (verbose_) std::cout << "HLTEventAnalyzerAOD::analyze:"
488  << " config extraction failure with process name "
489  << triggerResultsLabel_.process() << std::endl;
490  }
491 
492  allTriggerNames_.clear();
493  for (unsigned i = 0; i<hltConfig_.size();i++) {
495  }
496 // std::cout<<"Length: "<<allTriggerNames_.size()<<std::endl;
497 
499  for ( std::vector<GenericTriggerEventFlag *>::const_iterator it = triggerFolderEventFlag_.begin(); it!= triggerFolderEventFlag_.end(); it++) {
500  int pos = it - triggerFolderEventFlag_.begin();
501  if ((*it)->on()) {
502  (*it)->initRun( iRun, iSetup );
503  if (triggerSelectedSubFolders_[pos].exists(std::string("hltDBKey"))) {
504 // std::cout<<"Looking for hltDBKey for"<<triggerFolderLabels_[pos]<<std::endl;
505  if ((*it)->expressionsFromDB((*it)->hltDBKey(), iSetup)[0] != "CONFIG_ERROR")
506  triggerFolderExpr_[pos] = (*it)->expressionsFromDB((*it)->hltDBKey(), iSetup);
507  }
508 // for (unsigned j = 0; j<triggerFolderExpr_[pos].size(); j++) std::cout<<"pos "<<pos<<" "<<triggerFolderLabels_[pos]<<" triggerFolderExpr_"<<triggerFolderExpr_[pos][j]<<std::endl;
509  }
510  }
511 }
unsigned int size() const
number of trigger paths in trigger table
void dump(const std::string &what) const
Dumping config info to cout.
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
HLTConfigProvider hltConfig_
Definition: METAnalyzer.h:188
const std::string & triggerName(unsigned int triggerIndex) const
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:33
std::vector< GenericTriggerEventFlag * > triggerFolderEventFlag_
Definition: METAnalyzer.h:198
std::vector< std::vector< std::string > > triggerFolderExpr_
Definition: METAnalyzer.h:199
std::vector< std::string > allTriggerNames_
Definition: METAnalyzer.h:194
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
edm::InputTag triggerResultsLabel_
Definition: METAnalyzer.h:189
edm::ParameterSet parameters
Definition: METAnalyzer.h:125
std::string const & process() const
Definition: InputTag.h:46
edm::VParameterSet triggerSelectedSubFolders_
Definition: METAnalyzer.h:197
tuple cout
Definition: gather_cfg.py:121
void METAnalyzer::endRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
virtual

Finish up a run.

Reimplemented from edm::EDAnalyzer.

Definition at line 514 of file METAnalyzer.cc.

References cmsPerfPublish::DirName, TrackerOfflineValidation_Dqm_cff::dirName, MonitorElement::getRootObject(), MonitorElement::getTH1F(), i, and AlCaHLTBitMon_QueryRunRegistry::string.

515 {
516 
517  //
518  //--- Check the time length of the Run from the lumi section plots
519 
520 
521  TH1F* tlumisec;
522 
523  MonitorElement *meLumiSec = map_dijet_MEs["aaa"];
524  meLumiSec = map_dijet_MEs["JetMET/lumisec"];
525 
526  int totlsec=0;
527  int totlssecsum=0;
528  double totltime=0.;
529  if (meLumiSec && meLumiSec->getRootObject() ) {
530  tlumisec = meLumiSec->getTH1F();
531  //check overflow bin (if we have more than 2500 LS in a run)
532  //lumisec is filled every time the analyze section is processed
533  //we know an LS is present only once in a run: normalize how many events we had on average
534  //if lumi fluctuates strongly might be unreliable for overflow bin though
535  for (int i=0; i< (tlumisec->GetNbinsX()); i++){
536  if (tlumisec->GetBinContent(i)!=0){
537  totlsec+=1;
538  totlssecsum+=tlumisec->GetBinContent(i);
539  }
540  }
541  int num_per_ls=(double)totlssecsum/(double)totlsec;
542  totlsec=totlsec+tlumisec->GetBinContent(tlumisec->GetNbinsX()+1)/(double)num_per_ls;
543  totltime = double(totlsec*90); // one lumi sec ~ 90 (sec)
544  }
545 
546  if (totltime==0.) totltime=1.;
547 
549  //dbe_->setCurrentFolder(dirName);
550 
551 
552 
553  //below is the original METAnalyzer formulation
554 
555  for (std::vector<std::string>::const_iterator ic = folderNames_.begin(); ic != folderNames_.end(); ic++) {
557  DirName = dirName+*ic;
558  makeRatePlot(DirName,totltime);
559  for ( std::vector<GenericTriggerEventFlag *>::const_iterator it = triggerFolderEventFlag_.begin(); it!= triggerFolderEventFlag_.end(); it++) {
560  int pos = it - triggerFolderEventFlag_.begin();
561  if ((*it)->on()) {
562  makeRatePlot(DirName+"/"+triggerFolderLabels_[pos],totltime);
563  }
564  }
565  }
566 
567 }
std::map< std::string, MonitorElement * > map_dijet_MEs
Definition: METAnalyzer.h:401
int i
Definition: DBlmapReader.cc:9
std::string FolderName_
Definition: METAnalyzer.h:133
std::vector< GenericTriggerEventFlag * > triggerFolderEventFlag_
Definition: METAnalyzer.h:198
void makeRatePlot(std::string, double)
Definition: METAnalyzer.cc:571
edm::InputTag metCollectionLabel_
Definition: METAnalyzer.h:135
std::vector< std::string > folderNames_
Definition: METAnalyzer.h:269
TObject * getRootObject(void) const
TH1F * getTH1F(void) const
std::string const & label() const
Definition: InputTag.h:42
std::vector< std::string > triggerFolderLabels_
Definition: METAnalyzer.h:200
void METAnalyzer::fillMESet ( const edm::Event iEvent,
std::string  DirName,
const reco::MET met,
const reco::PFMET pfmet,
const reco::CaloMET calomet,
std::map< std::string, MonitorElement * > &  map_of_MEs 
)

Definition at line 932 of file METAnalyzer.cc.

References i, and AlCaHLTBitMon_QueryRunRegistry::string.

934 {
935 
936  //dbe_->setCurrentFolder(DirName);
937 
938  bool bLumiSecPlot=fill_met_high_level_histo;
939  //if (DirName.find("Uncleaned")) bLumiSecPlot=true; //now done on configlevel
940  fillMonitorElement(iEvent, DirName, std::string(""), met, pfmet, calomet, map_of_MEs,bLumiSecPlot);
941  if (DirName.find("Cleaned")) {
942  for (unsigned i = 0; i<triggerFolderLabels_.size(); i++) {
943  if (triggerFolderDecisions_[i]) fillMonitorElement(iEvent, DirName, triggerFolderLabels_[i], met, pfmet, calomet, map_of_MEs, false);
944  }
945  }
946 
947  if (DirName.find("DiJet")) {
948  for (unsigned i = 0; i<triggerFolderLabels_.size(); i++) {
949  if (triggerFolderDecisions_[i]) fillMonitorElement(iEvent, DirName, triggerFolderLabels_[i], met, pfmet, calomet, map_of_MEs, false);
950  }
951  }
952 
953 
954 
955 
956 
957 
958 
959 // if (trigJetMB_)
960 // fillMonitorElement(iEvent,DirName,"",met,pfmet,calomet, bLumiSecPlot);
961 // if (trigHighPtJet_)
962 // fillMonitorElement(iEvent,DirName,"HighPtJet",met,pfmet,calomet,false);
963 // if (trigLowPtJet_)
964 // fillMonitorElement(iEvent,DirName,"LowPtJet",met,pfmet,calomet,false);
965 // if (trigMinBias_)
966 // fillMonitorElement(iEvent,DirName,"MinBias",met,pfmet,calomet,false);
967 // if (trigHighMET_)
968 // fillMonitorElement(iEvent,DirName,"HighMET",met,pfmet,calomet,false);
969 // // if (_trig_LowMET)
970 // // fillMonitorElement(iEvent,DirName,"LowMET",met,pfmet,calomet,false);
971 // if (trigEle_)
972 // fillMonitorElement(iEvent,DirName,"Ele",met,pfmet,calomet,false);
973 // if (trigMuon_)
974 // fillMonitorElement(iEvent,DirName,"Muon",met,pfmet,calomet,false);
975 }
int i
Definition: DBlmapReader.cc:9
std::vector< int > triggerFolderDecisions_
Definition: METAnalyzer.h:201
bool fill_met_high_level_histo
Definition: METAnalyzer.h:407
void fillMonitorElement(const edm::Event &, std::string, std::string, const reco::MET &, const reco::PFMET &, const reco::CaloMET &, std::map< std::string, MonitorElement * > &, bool)
Definition: METAnalyzer.cc:978
std::vector< std::string > triggerFolderLabels_
Definition: METAnalyzer.h:200
void METAnalyzer::fillMonitorElement ( const edm::Event iEvent,
std::string  DirName,
std::string  subFolderName,
const reco::MET met,
const reco::PFMET pfmet,
const reco::CaloMET calomet,
std::map< std::string, MonitorElement * > &  map_of_MEs,
bool  bLumiSecPlot 
)
  • if (isTCMet_) {

Definition at line 978 of file METAnalyzer.cc.

References reco::PFMET::chargedHadronEt(), reco::PFMET::chargedHadronEtFraction(), reco::PFMET::electronEt(), reco::PFMET::electronEtFraction(), reco::CaloMET::emEtFraction(), reco::CaloMET::emEtInEB(), reco::CaloMET::emEtInEE(), reco::CaloMET::emEtInHF(), reco::CaloMET::etFractionHadronic(), HcalObjRepresent::Fill(), reco::CaloMET::hadEtInHB(), reco::CaloMET::hadEtInHE(), reco::CaloMET::hadEtInHF(), reco::CaloMET::hadEtInHO(), reco::PFMET::HFEMEt(), reco::PFMET::HFEMEtFraction(), reco::PFMET::HFHadronEt(), reco::PFMET::HFHadronEtFraction(), i, edm::EventBase::luminosityBlock(), reco::CaloMET::maxEtInEmTowers(), reco::CaloMET::maxEtInHadTowers(), reco::MET::mEtSig(), reco::PFMET::muonEt(), reco::PFMET::muonEtFraction(), reco::PFMET::neutralHadronEt(), reco::PFMET::neutralHadronEtFraction(), reco::LeafCandidate::phi(), reco::PFMET::photonEt(), reco::PFMET::photonEtFraction(), reco::LeafCandidate::pt(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), and reco::MET::sumEt().

981 {
982 
983 // if (subFolderName=="HighPtJet") {
984 // if (!selectHighPtJetEvent(iEvent)) return;
985 // }
986 // else if (subFolderName=="LowPtJet") {
987 // if (!selectLowPtJetEvent(iEvent)) return;
988 // }
989 // else if (subFolderName=="HighMET") {
990 // if (met.pt()<highMETThreshold_) return;
991 // }
992 // // else if (subFolderName=="LowMET") {
993 // // if (met.pt()<_lowMETThreshold) return;
994 // // }
995 // else if (subFolderName=="Ele") {
996 // if (!selectWElectronEvent(iEvent)) return;
997 // }
998 // else if (subFolderName=="Muon") {
999 // if (!selectWMuonEvent(iEvent)) return;
1000 // }
1001 
1002 // Reconstructed MET Information
1003  double SumET = met.sumEt();
1004  double METSig = met.mEtSig();
1005  //double Ez = met.e_longitudinal();
1006  double MET = met.pt();
1007  double MEx = met.px();
1008  double MEy = met.py();
1009  double METPhi = met.phi();
1010  //
1011  int myLuminosityBlock;
1012  myLuminosityBlock = iEvent.luminosityBlock();
1013  //
1014 
1015  if (subFolderName!="") DirName = DirName +"/"+subFolderName;
1016 
1017 
1018  hTrigger = map_of_MEs[DirName+"/triggerResults"];
1019  // std::cout<<"Hello"<<c++<<":"<<hTrigger <<std::endl;//":"<< hTrigger->getRootObject()<<std::endl;
1020  if (hTrigger && hTrigger->getRootObject()) {
1021  // std::cout<<"Hello"<<c++<<std::endl;
1022  for (unsigned int i = 0; i<allTriggerDecisions_.size();i++){
1023  // std::cout<<"Hello"<<c++<<":"<<i<<":"<< allTriggerDecisions_[i]<<":"<<allTriggerDecisions_[i]<<std::endl;
1024  if(i<(unsigned int)hTrigger->getNbinsX()){
1026  if (!hTriggerLabelsIsSet_) {
1027  hTrigger->setBinLabel(i+1, allTriggerNames_[i]);//Can't be done in beginJob (no trigger list). Can't be done in beginRun (would have to anticipate folder structure).FIXME doesn't work
1028  }
1029  }
1030  }
1031  if (!hTriggerLabelsIsSet_) for (int i = allTriggerDecisions_.size(); i<hTrigger->getNbinsX();i++){
1032  hTrigger->setBinLabel(i+1, "");//Can't be done in beginJob (no trigger list). Can't be done in beginRun (would have to anticipate folder structure).
1033  }
1034  hTriggerLabelsIsSet_ = true;
1035  // std::cout<<"Filling decision "<<allTriggerNames_[i]<<" "<<allTriggerDecisions_[i]<<std::endl;
1036  }
1037 
1038 
1039  hMEx = map_of_MEs[DirName+"/"+"MEx"]; if (hMEx && hMEx->getRootObject()) hMEx ->Fill(MEx);
1040  hMEy = map_of_MEs[DirName+"/"+"MEy"]; if (hMEy && hMEy->getRootObject()) hMEy ->Fill(MEy);
1041  hMET = map_of_MEs[DirName+"/"+"MET"]; if (hMET && hMET->getRootObject()) hMET ->Fill(MET);
1042  hMETPhi = map_of_MEs[DirName+"/"+"METPhi"]; if (hMETPhi && hMETPhi->getRootObject()) hMETPhi ->Fill(METPhi);
1043  hSumET = map_of_MEs[DirName+"/"+"SumET"]; if (hSumET && hSumET->getRootObject()) hSumET ->Fill(SumET);
1044  hMETSig = map_of_MEs[DirName+"/"+"METSig"]; if (hMETSig && hMETSig->getRootObject()) hMETSig ->Fill(METSig);
1045  hMET_logx = map_of_MEs[DirName+"/"+"MET_logx"]; if (hMET_logx && hMET_logx->getRootObject()) hMET_logx->Fill(log10(MET));
1046  hSumET_logx = map_of_MEs[DirName+"/"+"SumET_logx"]; if (hSumET_logx && hSumET_logx->getRootObject()) hSumET_logx->Fill(log10(SumET));
1047 
1048  // Fill NPV profiles
1049  //--------------------------------------------------------------------------
1050  meMEx_profile = map_of_MEs[DirName + "/MEx_profile"];
1051  meMEy_profile = map_of_MEs[DirName + "/MEy_profile"];
1052  meMET_profile = map_of_MEs[DirName + "/MET_profile"];
1053  meSumET_profile = map_of_MEs[DirName + "/SumET_profile"];
1054 
1055  if (meMEx_profile && meMEx_profile ->getRootObject()) meMEx_profile ->Fill(numPV_, MEx);
1056  if (meMEy_profile && meMEy_profile ->getRootObject()) meMEy_profile ->Fill(numPV_, MEy);
1057  if (meMET_profile && meMET_profile ->getRootObject()) meMET_profile ->Fill(numPV_, MET);
1059 
1060  if(isCaloMet_){
1061  //const reco::CaloMETCollection *calometcol = calometcoll.product();
1062  //const reco::CaloMET *calomet;
1063  //calomet = &(calometcol->front());
1064 
1065  double caloEtFractionHadronic = calomet.etFractionHadronic();
1066  double caloEmEtFraction = calomet.emEtFraction();
1067 
1068  double caloMaxEtInEMTowers = calomet.maxEtInEmTowers();
1069  double caloMaxEtInHadTowers = calomet.maxEtInHadTowers();
1070 
1071  double caloHadEtInHB = calomet.hadEtInHB();
1072  double caloHadEtInHO = calomet.hadEtInHO();
1073  double caloHadEtInHE = calomet.hadEtInHE();
1074  double caloHadEtInHF = calomet.hadEtInHF();
1075  double caloEmEtInEB = calomet.emEtInEB();
1076  double caloEmEtInEE = calomet.emEtInEE();
1077  double caloEmEtInHF = calomet.emEtInHF();
1078 
1079  hCaloMaxEtInEmTowers = map_of_MEs[DirName+"/"+"CaloMaxEtInEmTowers"]; if (hCaloMaxEtInEmTowers && hCaloMaxEtInEmTowers->getRootObject()) hCaloMaxEtInEmTowers->Fill(caloMaxEtInEMTowers);
1080  hCaloMaxEtInHadTowers = map_of_MEs[DirName+"/"+"CaloMaxEtInHadTowers"]; if (hCaloMaxEtInHadTowers && hCaloMaxEtInHadTowers->getRootObject()) hCaloMaxEtInHadTowers->Fill(caloMaxEtInHadTowers);
1081 
1082  hCaloHadEtInHB = map_of_MEs[DirName+"/"+"CaloHadEtInHB"]; if (hCaloHadEtInHB && hCaloHadEtInHB->getRootObject()) hCaloHadEtInHB->Fill(caloHadEtInHB);
1083  hCaloHadEtInHO = map_of_MEs[DirName+"/"+"CaloHadEtInHO"]; if (hCaloHadEtInHO && hCaloHadEtInHO->getRootObject()) hCaloHadEtInHO->Fill(caloHadEtInHO);
1084  hCaloHadEtInHE = map_of_MEs[DirName+"/"+"CaloHadEtInHE"]; if (hCaloHadEtInHE && hCaloHadEtInHE->getRootObject()) hCaloHadEtInHE->Fill(caloHadEtInHE);
1085  hCaloHadEtInHF = map_of_MEs[DirName+"/"+"CaloHadEtInHF"]; if (hCaloHadEtInHF && hCaloHadEtInHF->getRootObject()) hCaloHadEtInHF->Fill(caloHadEtInHF);
1086  hCaloEmEtInEB = map_of_MEs[DirName+"/"+"CaloEmEtInEB"]; if (hCaloEmEtInEB && hCaloEmEtInEB->getRootObject()) hCaloEmEtInEB->Fill(caloEmEtInEB);
1087  hCaloEmEtInEE = map_of_MEs[DirName+"/"+"CaloEmEtInEE"]; if (hCaloEmEtInEE && hCaloEmEtInEE->getRootObject()) hCaloEmEtInEE->Fill(caloEmEtInEE);
1088  hCaloEmEtInHF = map_of_MEs[DirName+"/"+"CaloEmEtInHF"]; if (hCaloEmEtInHF && hCaloEmEtInHF->getRootObject()) hCaloEmEtInHF->Fill(caloEmEtInHF);
1089 
1090  hCaloMETPhi020 = map_of_MEs[DirName+"/"+"CaloMETPhi020"]; if (MET> 20. && hCaloMETPhi020 && hCaloMETPhi020->getRootObject()) { hCaloMETPhi020->Fill(METPhi);}
1091 
1092 
1093  hCaloEtFractionHadronic = map_of_MEs[DirName+"/"+"CaloEtFractionHadronic"]; if (hCaloEtFractionHadronic && hCaloEtFractionHadronic->getRootObject()) hCaloEtFractionHadronic->Fill(caloEtFractionHadronic);
1094  hCaloEmEtFraction = map_of_MEs[DirName+"/"+"CaloEmEtFraction"]; if (hCaloEmEtFraction && hCaloEmEtFraction->getRootObject()) hCaloEmEtFraction->Fill(caloEmEtFraction);
1095  hCaloEmEtFraction020 = map_of_MEs[DirName+"/"+"CaloEmEtFraction020"]; if (MET> 20. && hCaloEmEtFraction020 && hCaloEmEtFraction020->getRootObject()) hCaloEmEtFraction020->Fill(caloEmEtFraction);
1096  //if (metCollectionLabel_.label() == "corMetGlobalMuons" ) {
1097 
1098  //for( reco::MuonCollection::const_iterator muonit = muonHandle_->begin(); muonit != muonHandle_->end(); muonit++ ) {
1099  // const reco::TrackRef siTrack = muonit->innerTrack();
1100  // hCalomuPt = map_of_MEs[DirName+"/"+"CalomuonPt"];
1101  // if (hCalomuPt && hCalomuPt->getRootObject()) hCalomuPt->Fill( muonit->p4().pt() );
1102  // hCalomuEta = map_of_MEs[DirName+"/"+"CalomuonEta"]; if (hCalomuEta && hCalomuEta->getRootObject()) hCalomuEta->Fill( muonit->p4().eta() );
1103  // hCalomuNhits = map_of_MEs[DirName+"/"+"CalomuonNhits"]; if (hCalomuNhits && hCalomuNhits->getRootObject()) hCalomuNhits->Fill( siTrack.isNonnull() ? siTrack->numberOfValidHits() : -999 );
1104  // hCalomuChi2 = map_of_MEs[DirName+"/"+"CalomuonNormalizedChi2"]; if (hCalomuChi2 && hCalomuChi2->getRootObject()) hCalomuChi2->Fill( siTrack.isNonnull() ? siTrack->chi2()/siTrack->ndof() : -999 );
1105  // double d0 = siTrack.isNonnull() ? -1 * siTrack->dxy( beamSpot_) : -999;
1106  // hCalomuD0 = map_of_MEs[DirName+"/"+"CalomuonD0"]; if (hCalomuD0 && hCalomuD0->getRootObject()) hCalomuD0->Fill( d0 );
1107  //}
1108 
1109  //const unsigned int nMuons = muonHandle_->size();
1110  //for( unsigned int mus = 0; mus < nMuons; mus++ ) {
1111  // reco::MuonRef muref( muonHandle_, mus);
1112  // reco::MuonMETCorrectionData muCorrData = (*tcMetValueMapHandle_)[muref];
1113  // hCaloMExCorrection = map_of_MEs[DirName+"/"+"CaloMExCorrection"]; if (hCaloMExCorrection && hCaloMExCorrection->getRootObject()) hCaloMExCorrection-> Fill(muCorrData.corrY());
1114  // hCaloMEyCorrection = map_of_MEs[DirName+"/"+"CaloMEyCorrection"]; if (hCaloMEyCorrection && hCaloMEyCorrection->getRootObject()) hCaloMEyCorrection-> Fill(muCorrData.corrX());
1115  // hCaloMuonCorrectionFlag = map_of_MEs[DirName+"/"+"CaloMuonCorrectionFlag"]; if (hCaloMuonCorrectionFlag && hCaloMuonCorrectionFlag->getRootObject()) hCaloMuonCorrectionFlag-> Fill(muCorrData.type());
1116  //}
1117  //}
1118  }
1119 
1120  if(isPFMet_){
1121  // **** Get the MET container
1122  //const PFMETCollection *pfmetcol = pfmetcoll.product();
1123  //const PFMET *pfmet;
1124  //pfmet = &(pfmetcol->front());
1125 
1126  // PFMET getters
1127  //----------------------------------------------------------------------------
1128  double pfPhotonEtFraction = pfmet.photonEtFraction();
1129  double pfPhotonEt = pfmet.photonEt();
1130  double pfNeutralHadronEtFraction = pfmet.neutralHadronEtFraction();
1131  double pfNeutralHadronEt = pfmet.neutralHadronEt();
1132  double pfElectronEtFraction = pfmet.electronEtFraction();
1133  double pfElectronEt = pfmet.electronEt();
1134  double pfChargedHadronEtFraction = pfmet.chargedHadronEtFraction();
1135  double pfChargedHadronEt = pfmet.chargedHadronEt();
1136  double pfMuonEtFraction = pfmet.muonEtFraction();
1137  double pfMuonEt = pfmet.muonEt();
1138  double pfHFHadronEtFraction = pfmet.HFHadronEtFraction();
1139  double pfHFHadronEt = pfmet.HFHadronEt();
1140  double pfHFEMEtFraction = pfmet.HFEMEtFraction();
1141  double pfHFEMEt = pfmet.HFEMEt();
1142 
1143  mePhotonEtFraction = map_of_MEs[DirName + "/PfPhotonEtFraction"];
1144  mePhotonEt = map_of_MEs[DirName + "/PfPhotonEt"];
1145  meNeutralHadronEtFraction = map_of_MEs[DirName + "/PfNeutralHadronEtFraction"];
1146  meNeutralHadronEt = map_of_MEs[DirName + "/PfNeutralHadronEt"];
1147  meElectronEtFraction = map_of_MEs[DirName + "/PfElectronEtFraction"];
1148  meElectronEt = map_of_MEs[DirName + "/PfElectronEt"];
1149  meChargedHadronEtFraction = map_of_MEs[DirName + "/PfChargedHadronEtFraction"];
1150  meChargedHadronEt = map_of_MEs[DirName + "/PfChargedHadronEt"];
1151  meMuonEtFraction = map_of_MEs[DirName + "/PfMuonEtFraction"];
1152  meMuonEt = map_of_MEs[DirName + "/PfMuonEt"];
1153  meHFHadronEtFraction = map_of_MEs[DirName + "/PfHFHadronEtFraction"];
1154  meHFHadronEt = map_of_MEs[DirName + "/PfHFHadronEt"];
1155  meHFEMEtFraction = map_of_MEs[DirName + "/PfHFEMEtFraction"];
1156  meHFEMEt = map_of_MEs[DirName + "/PfHFEMEt"];
1157 
1158  if (mePhotonEtFraction && mePhotonEtFraction ->getRootObject()) mePhotonEtFraction ->Fill(pfPhotonEtFraction);
1159  if (mePhotonEt && mePhotonEt ->getRootObject()) mePhotonEt ->Fill(pfPhotonEt);
1161  if (meNeutralHadronEt && meNeutralHadronEt ->getRootObject()) meNeutralHadronEt ->Fill(pfNeutralHadronEt);
1162  if (meElectronEtFraction && meElectronEtFraction ->getRootObject()) meElectronEtFraction ->Fill(pfElectronEtFraction);
1163  if (meElectronEt && meElectronEt ->getRootObject()) meElectronEt ->Fill(pfElectronEt);
1165  if (meChargedHadronEt && meChargedHadronEt ->getRootObject()) meChargedHadronEt ->Fill(pfChargedHadronEt);
1166  if (meMuonEtFraction && meMuonEtFraction ->getRootObject()) meMuonEtFraction ->Fill(pfMuonEtFraction);
1167  if (meMuonEt && meMuonEt ->getRootObject()) meMuonEt ->Fill(pfMuonEt);
1168  if (meHFHadronEtFraction && meHFHadronEtFraction ->getRootObject()) meHFHadronEtFraction ->Fill(pfHFHadronEtFraction);
1169  if (meHFHadronEt && meHFHadronEt ->getRootObject()) meHFHadronEt ->Fill(pfHFHadronEt);
1170  if (meHFEMEtFraction && meHFEMEtFraction ->getRootObject()) meHFEMEtFraction ->Fill(pfHFEMEtFraction);
1171  if (meHFEMEt && meHFEMEt ->getRootObject()) meHFEMEt ->Fill(pfHFEMEt);
1172 
1173  //NPV profiles
1174 
1175  mePhotonEtFraction_profile = map_of_MEs[DirName + "/PfPhotonEtFraction_profile"];
1176  mePhotonEt_profile = map_of_MEs[DirName + "/PfPhotonEt_profile"];
1177  meNeutralHadronEtFraction_profile = map_of_MEs[DirName + "/PfNeutralHadronEtFraction_profile"];
1178  meNeutralHadronEt_profile = map_of_MEs[DirName + "/PfNeutralHadronEt_profile"];
1179  meElectronEtFraction_profile = map_of_MEs[DirName + "/PfElectronEtFraction_profile"];
1180  meElectronEt_profile = map_of_MEs[DirName + "/PfElectronEt_profile"];
1181  meChargedHadronEtFraction_profile = map_of_MEs[DirName + "/PfChargedHadronEtFraction_profile"];
1182  meChargedHadronEt_profile = map_of_MEs[DirName + "/PfChargedHadronEt_profile"];
1183  meMuonEtFraction_profile = map_of_MEs[DirName + "/PfMuonEtFraction_profile"];
1184  meMuonEt_profile = map_of_MEs[DirName + "/PfMuonEt_profile"];
1185  meHFHadronEtFraction_profile = map_of_MEs[DirName + "/PfHFHadronEtFraction_profile"];
1186  meHFHadronEt_profile = map_of_MEs[DirName + "/PfHFHadronEt_profile"];
1187  meHFEMEtFraction_profile = map_of_MEs[DirName + "/PfHFEMEtFraction_profile"];
1188  meHFEMEt_profile = map_of_MEs[DirName + "/PfHFEMEt_profile"];
1189 
1190  if (mePhotonEtFraction_profile && mePhotonEtFraction_profile ->getRootObject()) mePhotonEtFraction_profile ->Fill(numPV_, pfPhotonEtFraction);
1191  if (mePhotonEt_profile && mePhotonEt_profile ->getRootObject()) mePhotonEt_profile ->Fill(numPV_, pfPhotonEt);
1193  if (meNeutralHadronEt_profile && meNeutralHadronEt_profile ->getRootObject()) meNeutralHadronEt_profile ->Fill(numPV_, pfNeutralHadronEt);
1194  if (meElectronEtFraction_profile && meElectronEtFraction_profile ->getRootObject()) meElectronEtFraction_profile ->Fill(numPV_, pfElectronEtFraction);
1195  if (meElectronEt_profile && meElectronEt_profile ->getRootObject()) meElectronEt_profile ->Fill(numPV_, pfElectronEt);
1197  if (meChargedHadronEt_profile && meChargedHadronEt_profile ->getRootObject()) meChargedHadronEt_profile ->Fill(numPV_, pfChargedHadronEt);
1198  if (meMuonEtFraction_profile && meMuonEtFraction_profile ->getRootObject()) meMuonEtFraction_profile ->Fill(numPV_, pfMuonEtFraction);
1199  if (meMuonEt_profile && meMuonEt_profile ->getRootObject()) meMuonEt_profile ->Fill(numPV_, pfMuonEt);
1200  if (meHFHadronEtFraction_profile && meHFHadronEtFraction_profile ->getRootObject()) meHFHadronEtFraction_profile ->Fill(numPV_, pfHFHadronEtFraction);
1201  if (meHFHadronEt_profile && meHFHadronEt_profile ->getRootObject()) meHFHadronEt_profile ->Fill(numPV_, pfHFHadronEt);
1202  if (meHFEMEtFraction_profile && meHFEMEtFraction_profile ->getRootObject()) meHFEMEtFraction_profile ->Fill(numPV_, pfHFEMEtFraction);
1203  if (meHFEMEt_profile && meHFEMEt_profile ->getRootObject()) meHFEMEt_profile ->Fill(numPV_, pfHFEMEt);
1204  }
1205 
1206  if (isCaloMet_){
1207  //if (bLumiSecPlot){//get from config level
1209  hMExLS = map_of_MEs[DirName+"/"+"MExLS"]; if (hMExLS && hMExLS->getRootObject()) hMExLS->Fill(MEx,myLuminosityBlock);
1210  hMEyLS = map_of_MEs[DirName+"/"+"MEyLS"]; if (hMEyLS && hMEyLS->getRootObject()) hMEyLS->Fill(MEy,myLuminosityBlock);
1211  }
1212  }
1213 
1216 
1217  //if(trackHandle_.isValid()) {
1218  // for( edm::View<reco::Track>::const_iterator trkit = trackHandle_->begin(); trkit != trackHandle_->end(); trkit++ ) {
1219  // htrkPt = map_of_MEs[DirName+"/"+"trackPt"]; if (htrkPt && htrkPt->getRootObject()) htrkPt->Fill( trkit->pt() );
1220  // htrkEta = map_of_MEs[DirName+"/"+"trackEta"]; if (htrkEta && htrkEta->getRootObject()) htrkEta->Fill( trkit->eta() );
1221  // htrkNhits = map_of_MEs[DirName+"/"+"trackNhits"]; if (htrkNhits && htrkNhits->getRootObject()) htrkNhits->Fill( trkit->numberOfValidHits() );
1222  // htrkChi2 = map_of_MEs[DirName+"/"+"trackNormalizedChi2"];
1223  // if (htrkChi2 && htrkChi2->getRootObject()) htrkChi2->Fill( trkit->chi2() / trkit->ndof() );
1224  // double d0 = -1 * trkit->dxy( beamSpot_ );
1225  // htrkD0 = map_of_MEs[DirName+"/"+"trackD0"]; if (htrkD0 && htrkD0->getRootObject()) htrkD0->Fill( d0 );
1226  // }
1227  //}else{if (verbose_) std::cout<<"tracks not valid"<<std::endl;}
1228 
1229  //if(electronHandle_.isValid()) {
1230  // for( edm::View<reco::GsfElectron>::const_iterator eleit = electronHandle_->begin(); eleit != electronHandle_->end(); eleit++ ) {
1231  // helePt = map_of_MEs[DirName+"/"+"electronPt"]; if (helePt && helePt->getRootObject()) helePt->Fill( eleit->p4().pt() );
1232  // heleEta = map_of_MEs[DirName+"/"+"electronEta"]; if (heleEta && heleEta->getRootObject()) heleEta->Fill( eleit->p4().eta() );
1233  // heleHoE = map_of_MEs[DirName+"/"+"electronHoverE"]; if (heleHoE && heleHoE->getRootObject()) heleHoE->Fill( eleit->hadronicOverEm() );
1234  // }
1235  //}else{
1236  // if (verbose_) std::cout<<"electrons not valid"<<std::endl;
1237  //}
1238 
1239  //if(muonHandle_.isValid()) {
1240  // for( reco::MuonCollection::const_iterator muonit = muonHandle_->begin(); muonit != muonHandle_->end(); muonit++ ) {
1241  // const reco::TrackRef siTrack = muonit->innerTrack();
1242  // hmuPt = map_of_MEs[DirName+"/"+"muonPt"]; if (hmuPt && hmuPt->getRootObject()) hmuPt ->Fill( muonit->p4().pt() );
1243  // hmuEta = map_of_MEs[DirName+"/"+"muonEta"]; if (hmuEta && hmuEta->getRootObject()) hmuEta ->Fill( muonit->p4().eta() );
1244  // hmuNhits = map_of_MEs[DirName+"/"+"muonNhits"]; if (hmuNhits && hmuNhits->getRootObject()) hmuNhits->Fill( siTrack.isNonnull() ? siTrack->numberOfValidHits() : -999 );
1245  // hmuChi2 = map_of_MEs[DirName+"/"+"muonNormalizedChi2"]; if (hmuChi2 && hmuChi2->getRootObject()) hmuChi2 ->Fill( siTrack.isNonnull() ? siTrack->chi2()/siTrack->ndof() : -999 );
1246  // double d0 = siTrack.isNonnull() ? -1 * siTrack->dxy( beamSpot_) : -999;
1247  // hmuD0 = map_of_MEs[DirName+"/"+"muonD0"]; if (hmuD0 && hmuD0->getRootObject()) hmuD0->Fill( d0 );
1248  // }
1249  // const unsigned int nMuons = muonHandle_->size();
1250  // for( unsigned int mus = 0; mus < nMuons; mus++ ) {
1251  // reco::MuonRef muref( muonHandle_, mus);
1252  // reco::MuonMETCorrectionData muCorrData = (*tcMetValueMapHandle_)[muref];
1253  // hMExCorrection = map_of_MEs[DirName+"/"+"MExCorrection"]; if (hMExCorrection && hMExCorrection->getRootObject()) hMExCorrection-> Fill(muCorrData.corrY());
1254  // hMEyCorrection = map_of_MEs[DirName+"/"+"MEyCorrection"]; if (hMEyCorrection && hMEyCorrection->getRootObject()) hMEyCorrection-> Fill(muCorrData.corrX());
1255  //hMuonCorrectionFlag = map_of_MEs[DirName+"/"+"CorrectionFlag"]; if (hMuonCorrectionFlag && hMuonCorrectionFlag->getRootObject()) hMuonCorrectionFlag-> Fill(muCorrData.type());
1256  // }
1257  // }else{
1258  //if (verbose_) std::cout<<"muons not valid"<<std::endl;
1259  // }
1260  // }*/
1261 }
MonitorElement * meElectronEt_profile
Definition: METAnalyzer.h:391
MonitorElement * hMExLS
Definition: METAnalyzer.h:292
int i
Definition: DBlmapReader.cc:9
double hadEtInHE() const
Definition: CaloMET.h:51
MonitorElement * hSumET
Definition: METAnalyzer.h:290
MonitorElement * meNeutralHadronEtFraction
Definition: METAnalyzer.h:365
MonitorElement * meNeutralHadronEtFraction_profile
Definition: METAnalyzer.h:388
MonitorElement * meElectronEtFraction_profile
Definition: METAnalyzer.h:390
std::vector< int > allTriggerDecisions_
Definition: METAnalyzer.h:195
double HFEMEtFraction() const
Definition: PFMET.h:47
double hadEtInHF() const
Definition: CaloMET.h:53
MonitorElement * hMEy
Definition: METAnalyzer.h:285
MonitorElement * hCaloMETPhi020
Definition: METAnalyzer.h:299
virtual float pt() const
transverse momentum
MonitorElement * mePhotonEt
Definition: METAnalyzer.h:364
MonitorElement * meChargedHadronEtFraction_profile
Definition: METAnalyzer.h:392
MonitorElement * hCaloHadEtInHF
Definition: METAnalyzer.h:313
bool isCaloMet_
Definition: METAnalyzer.h:403
virtual float phi() const
momentum azimuthal angle
double maxEtInHadTowers() const
Definition: CaloMET.h:40
double neutralHadronEtFraction() const
Definition: PFMET.h:32
MonitorElement * hMEx
Definition: METAnalyzer.h:284
double muonEt() const
Definition: PFMET.h:42
MonitorElement * hCaloHadEtInHE
Definition: METAnalyzer.h:312
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:59
MonitorElement * hCaloEtFractionHadronic
Definition: METAnalyzer.h:303
MonitorElement * hMEyLS
Definition: METAnalyzer.h:293
MonitorElement * meHFHadronEtFraction
Definition: METAnalyzer.h:373
MonitorElement * hMET_logx
Definition: METAnalyzer.h:295
MonitorElement * hCaloEmEtInHF
Definition: METAnalyzer.h:314
double maxEtInEmTowers() const
Definition: CaloMET.h:38
MonitorElement * meHFEMEtFraction_profile
Definition: METAnalyzer.h:398
MonitorElement * hCaloHadEtInHO
Definition: METAnalyzer.h:311
MonitorElement * meChargedHadronEtFraction
Definition: METAnalyzer.h:369
MonitorElement * meMEx_profile
Definition: METAnalyzer.h:381
MonitorElement * hMET
Definition: METAnalyzer.h:288
MonitorElement * hTrigger
Definition: METAnalyzer.h:282
void Fill(long long x)
MonitorElement * hCaloMaxEtInHadTowers
Definition: METAnalyzer.h:302
double mEtSig() const
Definition: MET.h:58
double HFHadronEtFraction() const
Definition: PFMET.h:44
double sumEt() const
Definition: MET.h:56
double muonEtFraction() const
Definition: PFMET.h:41
bool fill_met_high_level_histo
Definition: METAnalyzer.h:407
MonitorElement * hMETPhi
Definition: METAnalyzer.h:289
MonitorElement * hCaloEmEtInEE
Definition: METAnalyzer.h:315
MonitorElement * meMuonEt
Definition: METAnalyzer.h:372
double HFHadronEt() const
Definition: PFMET.h:45
MonitorElement * meHFEMEt_profile
Definition: METAnalyzer.h:399
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
Definition: MET.h:42
MonitorElement * hSumET_logx
Definition: METAnalyzer.h:296
double photonEtFraction() const
Definition: PFMET.h:29
MonitorElement * meMuonEt_profile
Definition: METAnalyzer.h:395
MonitorElement * hCaloEmEtFraction
Definition: METAnalyzer.h:304
double emEtInEB() const
Definition: CaloMET.h:55
MonitorElement * meChargedHadronEt
Definition: METAnalyzer.h:370
MonitorElement * hMETSig
Definition: METAnalyzer.h:287
MonitorElement * meChargedHadronEt_profile
Definition: METAnalyzer.h:393
MonitorElement * meMuonEtFraction
Definition: METAnalyzer.h:371
std::vector< std::string > allTriggerNames_
Definition: METAnalyzer.h:194
double HFEMEt() const
Definition: PFMET.h:48
MonitorElement * mePhotonEtFraction
Definition: METAnalyzer.h:363
double hadEtInHO() const
Definition: CaloMET.h:49
MonitorElement * hCaloMaxEtInEmTowers
Definition: METAnalyzer.h:301
MonitorElement * mePhotonEt_profile
Definition: METAnalyzer.h:387
double electronEt() const
Definition: PFMET.h:36
TObject * getRootObject(void) const
MonitorElement * meElectronEt
Definition: METAnalyzer.h:368
double etFractionHadronic() const
Definition: CaloMET.h:42
virtual double px() const
x coordinate of momentum vector
MonitorElement * meMET_profile
Definition: METAnalyzer.h:383
double photonEt() const
Definition: PFMET.h:30
MonitorElement * meHFHadronEt
Definition: METAnalyzer.h:374
MonitorElement * hCaloEmEtInEB
Definition: METAnalyzer.h:316
MonitorElement * hCaloHadEtInHB
Definition: METAnalyzer.h:310
MonitorElement * meHFHadronEtFraction_profile
Definition: METAnalyzer.h:396
MonitorElement * meHFEMEtFraction
Definition: METAnalyzer.h:375
double emEtInEE() const
Definition: CaloMET.h:57
double electronEtFraction() const
Definition: PFMET.h:35
double chargedHadronEtFraction() const
Definition: PFMET.h:38
MonitorElement * hCaloEmEtFraction020
Definition: METAnalyzer.h:308
int getNbinsX(void) const
get # of bins in X-axis
MonitorElement * meHFEMEt
Definition: METAnalyzer.h:376
double emEtInHF() const
Definition: CaloMET.h:59
MonitorElement * meMEy_profile
Definition: METAnalyzer.h:382
double neutralHadronEt() const
Definition: PFMET.h:33
MonitorElement * meElectronEtFraction
Definition: METAnalyzer.h:367
MonitorElement * meSumET_profile
Definition: METAnalyzer.h:384
bool hTriggerLabelsIsSet_
Definition: METAnalyzer.h:277
MonitorElement * meHFHadronEt_profile
Definition: METAnalyzer.h:397
MonitorElement * meNeutralHadronEt_profile
Definition: METAnalyzer.h:389
MonitorElement * mePhotonEtFraction_profile
Definition: METAnalyzer.h:386
MonitorElement * meMuonEtFraction_profile
Definition: METAnalyzer.h:394
double hadEtInHB() const
Definition: CaloMET.h:47
virtual double py() const
y coordinate of momentum vector
double emEtFraction() const
Definition: CaloMET.h:45
double chargedHadronEt() const
Definition: PFMET.h:39
MonitorElement * meNeutralHadronEt
Definition: METAnalyzer.h:366
void METAnalyzer::makeRatePlot ( std::string  DirName,
double  totltime 
)

Definition at line 571 of file METAnalyzer.cc.

References MonitorElement::getRootObject(), MonitorElement::getTH1F(), i, and MonitorElement::setBinContent().

572 {
573 
574  //dbe_->setCurrentFolder(DirName);
575  MonitorElement *meMET = map_dijet_MEs[DirName+"/"+"MET"];
576  MonitorElement *mMETRate = map_dijet_MEs[DirName+"/"+"METRate"];
577 
578  TH1F* tMET;
579  TH1F* tMETRate;
580 
581  if ( meMET && mMETRate){
582  if ( meMET->getRootObject() && mMETRate->getRootObject()) {
583  tMET = meMET->getTH1F();
584 
585  // Integral plot & convert number of events to rate (hz)
586  tMETRate = (TH1F*) tMET->Clone("METRateHist");
587  for (int i = tMETRate->GetNbinsX()-1; i>=0; i--){
588  mMETRate->setBinContent(i+1,tMETRate->GetBinContent(i+2)+tMET->GetBinContent(i+1));
589  }
590  for (int i = 0; i<tMETRate->GetNbinsX(); i++){
591  mMETRate->setBinContent(i+1,tMETRate->GetBinContent(i+1)/double(totltime));
592  }
593  }
594  }
595 
596 }
std::map< std::string, MonitorElement * > map_dijet_MEs
Definition: METAnalyzer.h:401
int i
Definition: DBlmapReader.cc:9
void setBinContent(int binx, double content)
set content of bin (1-D)
TObject * getRootObject(void) const
TH1F * getTH1F(void) const

Member Data Documentation

std::vector< int > METAnalyzer::allTriggerDecisions_
private

Definition at line 195 of file METAnalyzer.h.

std::vector<std::string > METAnalyzer::allTriggerNames_
private

Definition at line 194 of file METAnalyzer.h.

edm::InputTag METAnalyzer::beamHaloSummaryTag_
private

Definition at line 138 of file METAnalyzer.h.

edm::EDGetTokenT<reco::BeamHaloSummary> METAnalyzer::beamHaloSummaryToken_
private

Definition at line 150 of file METAnalyzer.h.

math::XYZPoint METAnalyzer::beamSpot_
private

Definition at line 271 of file METAnalyzer.h.

bool METAnalyzer::bypassAllDCSChecks_
private

Definition at line 240 of file METAnalyzer.h.

bool METAnalyzer::bypassAllPVChecks_
private

Definition at line 239 of file METAnalyzer.h.

edm::EDGetTokenT<reco::CaloJetCollection> METAnalyzer::caloJetsToken_
private

Definition at line 145 of file METAnalyzer.h.

edm::EDGetTokenT<reco::CaloMETCollection> METAnalyzer::caloMetToken_
private

Definition at line 154 of file METAnalyzer.h.

edm::ParameterSet METAnalyzer::cleaningParameters_
private

Definition at line 228 of file METAnalyzer.h.

JetMETDQMDCSFilter* METAnalyzer::DCSFilter_
private

Definition at line 267 of file METAnalyzer.h.

bool METAnalyzer::fill_met_high_level_histo
private

Definition at line 407 of file METAnalyzer.h.

std::string METAnalyzer::FolderName_
private

Definition at line 133 of file METAnalyzer.h.

std::vector<std::string> METAnalyzer::folderNames_
private

Definition at line 269 of file METAnalyzer.h.

edm::InputTag METAnalyzer::gtTag_
private

Definition at line 141 of file METAnalyzer.h.

edm::EDGetTokenT<L1GlobalTriggerReadoutRecord> METAnalyzer::gtToken_
private

Definition at line 144 of file METAnalyzer.h.

edm::InputTag METAnalyzer::hbheNoiseFilterResultTag_
private

Definition at line 139 of file METAnalyzer.h.

edm::EDGetTokenT<bool> METAnalyzer::hbheNoiseFilterResultToken_
private

Definition at line 149 of file METAnalyzer.h.

edm::InputTag METAnalyzer::hcalNoiseRBXCollectionTag_
private

Definition at line 136 of file METAnalyzer.h.

edm::EDGetTokenT<reco::HcalNoiseRBXCollection> METAnalyzer::HcalNoiseRBXToken_
private

Definition at line 155 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloEmEtFraction
private

Definition at line 304 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloEmEtFraction020
private

Definition at line 308 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloEmEtInEB
private

Definition at line 316 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloEmEtInEE
private

Definition at line 315 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloEmEtInHF
private

Definition at line 314 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloEmMET
private

Definition at line 321 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloEmMETPhi
private

Definition at line 322 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloEmMEx
private

Definition at line 318 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloEmMEy
private

Definition at line 319 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloEtFractionHadronic
private

Definition at line 303 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloHadEtInHB
private

Definition at line 310 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloHadEtInHE
private

Definition at line 312 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloHadEtInHF
private

Definition at line 313 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloHadEtInHO
private

Definition at line 311 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloHaMET
private

Definition at line 328 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloHaMETPhi
private

Definition at line 329 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloHaMEx
private

Definition at line 325 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloHaMEy
private

Definition at line 326 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloMaxEtInEmTowers
private

Definition at line 301 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloMaxEtInHadTowers
private

Definition at line 302 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hCaloMETPhi020
private

Definition at line 299 of file METAnalyzer.h.

double METAnalyzer::hfCalibFactor_
private

Definition at line 264 of file METAnalyzer.h.

HLTConfigProvider METAnalyzer::hltConfig_
private

Definition at line 188 of file METAnalyzer.h.

std::string METAnalyzer::hltPhysDec_
private

Definition at line 229 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hMET
private

Definition at line 288 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hMET_logx
private

Definition at line 295 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hMETPhi
private

Definition at line 289 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hMETRate
private

Definition at line 226 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hMETSig
private

Definition at line 287 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hMEx
private

Definition at line 284 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hMExLS
private

Definition at line 292 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hMEy
private

Definition at line 285 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hMEyLS
private

Definition at line 293 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hSumET
private

Definition at line 290 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hSumET_logx
private

Definition at line 296 of file METAnalyzer.h.

MonitorElement* METAnalyzer::hTrigger
private

Definition at line 282 of file METAnalyzer.h.

bool METAnalyzer::hTriggerLabelsIsSet_
private

Definition at line 277 of file METAnalyzer.h.

edm::InputTag METAnalyzer::inputJetIDValueMap
private

Definition at line 168 of file METAnalyzer.h.

bool METAnalyzer::isCaloMet_
private

Definition at line 403 of file METAnalyzer.h.

bool METAnalyzer::isPFMet_
private

Definition at line 405 of file METAnalyzer.h.

edm::InputTag METAnalyzer::jetCollectionLabel_
private

Definition at line 137 of file METAnalyzer.h.

std::string METAnalyzer::jetCorrectionService_
private

Definition at line 174 of file METAnalyzer.h.

edm::EDGetTokenT<edm::ValueMap <reco::JetID> > METAnalyzer::jetID_ValueMapToken_
private

Definition at line 169 of file METAnalyzer.h.

JetIDSelectionFunctor METAnalyzer::jetIDFunctorLoose
private

Definition at line 171 of file METAnalyzer.h.

int METAnalyzer::LSBegin_
private

Definition at line 236 of file METAnalyzer.h.

int METAnalyzer::LSEnd_
private

Definition at line 237 of file METAnalyzer.h.

MonitorElement* METAnalyzer::lumisecME
private

Definition at line 281 of file METAnalyzer.h.

std::map< std::string,MonitorElement* > METAnalyzer::map_dijet_MEs
private

Definition at line 401 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meChargedHadronEt
private

Definition at line 370 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meChargedHadronEt_profile
private

Definition at line 393 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meChargedHadronEtFraction
private

Definition at line 369 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meChargedHadronEtFraction_profile
private

Definition at line 392 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meElectronEt
private

Definition at line 368 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meElectronEt_profile
private

Definition at line 391 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meElectronEtFraction
private

Definition at line 367 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meElectronEtFraction_profile
private

Definition at line 390 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meHFEMEt
private

Definition at line 376 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meHFEMEt_profile
private

Definition at line 399 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meHFEMEtFraction
private

Definition at line 375 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meHFEMEtFraction_profile
private

Definition at line 398 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meHFHadronEt
private

Definition at line 374 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meHFHadronEt_profile
private

Definition at line 397 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meHFHadronEtFraction
private

Definition at line 373 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meHFHadronEtFraction_profile
private

Definition at line 396 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meMET_profile
private

Definition at line 383 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meMEx_profile
private

Definition at line 381 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meMEy_profile
private

Definition at line 382 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meMuonEt
private

Definition at line 372 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meMuonEt_profile
private

Definition at line 395 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meMuonEtFraction
private

Definition at line 371 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meMuonEtFraction_profile
private

Definition at line 394 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meNeutralHadronEt
private

Definition at line 366 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meNeutralHadronEt_profile
private

Definition at line 389 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meNeutralHadronEtFraction
private

Definition at line 365 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meNeutralHadronEtFraction_profile
private

Definition at line 388 of file METAnalyzer.h.

MonitorElement* METAnalyzer::mePhotonEt
private

Definition at line 364 of file METAnalyzer.h.

MonitorElement* METAnalyzer::mePhotonEt_profile
private

Definition at line 387 of file METAnalyzer.h.

MonitorElement* METAnalyzer::mePhotonEtFraction
private

Definition at line 363 of file METAnalyzer.h.

MonitorElement* METAnalyzer::mePhotonEtFraction_profile
private

Definition at line 386 of file METAnalyzer.h.

MonitorElement* METAnalyzer::meSumET_profile
private

Definition at line 384 of file METAnalyzer.h.

edm::InputTag METAnalyzer::metCollectionLabel_
private

Definition at line 135 of file METAnalyzer.h.

std::string METAnalyzer::MetType_
private

Definition at line 130 of file METAnalyzer.h.

std::string METAnalyzer::mOutputFile_
private

Definition at line 132 of file METAnalyzer.h.

int METAnalyzer::nbinsPV_
private

Definition at line 231 of file METAnalyzer.h.

double METAnalyzer::nPVMax_
private

Definition at line 233 of file METAnalyzer.h.

double METAnalyzer::nPVMin_
private

Definition at line 232 of file METAnalyzer.h.

int METAnalyzer::numPV_
private

Definition at line 258 of file METAnalyzer.h.

bool METAnalyzer::outputMEsInRootFile
private

Definition at line 131 of file METAnalyzer.h.

edm::ParameterSet METAnalyzer::parameters
private
PFJetIDSelectionFunctor METAnalyzer::pfjetIDFunctorLoose
private

Definition at line 172 of file METAnalyzer.h.

edm::EDGetTokenT<reco::PFJetCollection> METAnalyzer::pfJetsToken_
private

Definition at line 146 of file METAnalyzer.h.

edm::EDGetTokenT<reco::PFMETCollection> METAnalyzer::pfMetToken_
private

Definition at line 153 of file METAnalyzer.h.

double METAnalyzer::ptThreshold_
private

Definition at line 176 of file METAnalyzer.h.

bool METAnalyzer::runcosmics_
private

Definition at line 241 of file METAnalyzer.h.

std::vector<int> METAnalyzer::triggerFolderDecisions_
private

Definition at line 201 of file METAnalyzer.h.

std::vector<GenericTriggerEventFlag *> METAnalyzer::triggerFolderEventFlag_
private

Definition at line 198 of file METAnalyzer.h.

std::vector<std::vector<std::string> > METAnalyzer::triggerFolderExpr_
private

Definition at line 199 of file METAnalyzer.h.

std::vector<std::string > METAnalyzer::triggerFolderLabels_
private

Definition at line 200 of file METAnalyzer.h.

edm::InputTag METAnalyzer::triggerResultsLabel_
private

Definition at line 189 of file METAnalyzer.h.

edm::EDGetTokenT<edm::TriggerResults> METAnalyzer::triggerResultsToken_
private

Definition at line 190 of file METAnalyzer.h.

edm::VParameterSet METAnalyzer::triggerSelectedSubFolders_
private

Definition at line 197 of file METAnalyzer.h.

int METAnalyzer::verbose_
private

Definition at line 127 of file METAnalyzer.h.

edm::InputTag METAnalyzer::vertexTag_
private

Definition at line 140 of file METAnalyzer.h.

edm::EDGetTokenT<std::vector<reco::Vertex> > METAnalyzer::vertexToken_
private

Definition at line 143 of file METAnalyzer.h.