CMS 3D CMS Logo

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

#include <L1Trigger/L1TNtuples/src/L1JetRecoTreeProducer.cc>

Inheritance diagram for L1JetRecoTreeProducer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 L1JetRecoTreeProducer (const edm::ParameterSet &)
 
 ~L1JetRecoTreeProducer () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Public Attributes

L1Analysis::L1AnalysisRecoJetDataFormatjet_data
 
L1Analysis::L1AnalysisRecoMetDataFormatmet_data
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob (void) override
 
bool caloJetID (const reco::CaloJet &jet)
 
void doCaloJetCorr (edm::Handle< reco::CaloJetCollection > caloJets, edm::Handle< reco::JetCorrector > caloJetCorr)
 
void doCaloJets (edm::Handle< reco::CaloJetCollection > caloJets)
 
void doCaloMet (edm::Handle< reco::CaloMETCollection > caloMet)
 
void doCaloMetBE (edm::Handle< reco::CaloMETCollection > caloMetBE)
 
void doPFJetCorr (edm::Handle< reco::PFJetCollection > pfJets, edm::Handle< reco::JetCorrector > pfJetCorr)
 
void doPFJets (edm::Handle< reco::PFJetCollection > pfJets)
 
void doPFMet (edm::Handle< reco::PFMETCollection > pfMet)
 
void doPFMetNoMu (edm::Handle< reco::PFMETCollection > pfMet, edm::Handle< reco::MuonCollection >)
 
void endJob () override
 
bool pfJetID (const reco::PFJet &jet)
 

Private Attributes

edm::EDGetTokenT< reco::JetCorrectorcaloJECToken_
 
bool caloJetCorrMissing_
 
bool caloJetIDMissing_
 
edm::EDGetTokenT< edm::ValueMap< reco::JetID > > caloJetIDToken_
 
bool caloJetsMissing_
 
edm::EDGetTokenT< reco::CaloJetCollectioncaloJetToken_
 
bool caloMetBEMissing_
 
edm::EDGetTokenT< reco::CaloMETCollectioncaloMetBEToken_
 
bool caloMetMissing_
 
edm::EDGetTokenT< reco::CaloMETCollectioncaloMetToken_
 
edm::Service< TFileServicefs_
 
double jetetaMax_
 
double jetptThreshold_
 
unsigned int maxCl_
 
unsigned int maxJet_
 
unsigned int maxTrk_
 
unsigned int maxVtx_
 
bool muonsMissing_
 
edm::EDGetTokenT< reco::MuonCollectionmuonToken_
 
edm::EDGetTokenT< reco::JetCorrectorpfJECToken_
 
bool pfJetCorrMissing_
 
bool pfJetsMissing_
 
edm::EDGetTokenT< reco::PFJetCollectionpfJetToken_
 
bool pfMetMissing_
 
edm::EDGetTokenT< reco::PFMETCollectionpfMetToken_
 
TTree * tree_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- 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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
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

Description: Produces tree containing reco quantities

Definition at line 60 of file L1JetRecoTreeProducer.cc.

Constructor & Destructor Documentation

L1JetRecoTreeProducer::L1JetRecoTreeProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 130 of file L1JetRecoTreeProducer.cc.

References caloJECToken_, caloJetIDToken_, caloJetToken_, caloMetBEToken_, caloMetToken_, fs_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), jet_data, jetetaMax_, jetptThreshold_, TFileService::make(), maxJet_, met_data, muonToken_, pfJECToken_, pfJetToken_, pfMetToken_, and tree_.

130  :
131  pfJetsMissing_(false),
132  pfJetCorrMissing_(false),
133  caloJetCorrMissing_(false),
134  caloJetsMissing_(false),
135  caloJetIDMissing_(false),
136  pfMetMissing_(false),
137  caloMetMissing_(false),
138  caloMetBEMissing_(false),
139  muonsMissing_(false)
140 {
141 
142  caloJetToken_ = consumes<reco::CaloJetCollection>(iConfig.getUntrackedParameter("caloJetToken",edm::InputTag("ak4CaloJets")));
143  pfJetToken_ = consumes<reco::PFJetCollection>(iConfig.getUntrackedParameter("pfJetToken",edm::InputTag("ak4PFJetsCHS")));
144  caloJetIDToken_ = consumes<edm::ValueMap<reco::JetID> >(iConfig.getUntrackedParameter("caloJetIDToken",edm::InputTag("ak4JetID")));
145  pfJECToken_ = consumes<reco::JetCorrector>(iConfig.getUntrackedParameter<edm::InputTag>("pfJECToken",edm::InputTag("ak4PFCHSL1FastL2L3ResidualCorrector")));
146  caloJECToken_ = consumes<reco::JetCorrector>(iConfig.getUntrackedParameter<edm::InputTag>("caloJECToken",edm::InputTag("ak4CaloL1FastL2L3ResidualCorrector")));
147 
148  pfMetToken_ = consumes<reco::PFMETCollection>(iConfig.getUntrackedParameter("pfMetToken",edm::InputTag("pfMetT1")));
149  caloMetToken_ = consumes<reco::CaloMETCollection>(iConfig.getUntrackedParameter("caloMetToken",edm::InputTag("caloMet")));
150  caloMetBEToken_ = consumes<reco::CaloMETCollection>(iConfig.getUntrackedParameter("caloMetBEToken",edm::InputTag("caloMetBE")));
151 
152  muonToken_ = consumes<reco::MuonCollection>(iConfig.getUntrackedParameter("muonToken",edm::InputTag("muons")));
153 
154  jetptThreshold_ = iConfig.getParameter<double> ("jetptThreshold");
155  jetetaMax_ = iConfig.getParameter<double> ("jetetaMax");
156  maxJet_ = iConfig.getParameter<unsigned int>("maxJet");
157 
160 
161  // set up output
162  tree_=fs_->make<TTree>("JetRecoTree", "JetRecoTree");
163  tree_->Branch("Jet", "L1Analysis::L1AnalysisRecoJetDataFormat", &jet_data, 32000, 3);
164  tree_->Branch("Sums", "L1Analysis::L1AnalysisRecoMetDataFormat", &met_data, 32000, 3);
165 
166 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
L1Analysis::L1AnalysisRecoMetDataFormat * met_data
edm::Service< TFileService > fs_
edm::EDGetTokenT< reco::CaloMETCollection > caloMetToken_
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
edm::EDGetTokenT< reco::JetCorrector > pfJECToken_
L1Analysis::L1AnalysisRecoJetDataFormat * jet_data
edm::EDGetTokenT< edm::ValueMap< reco::JetID > > caloJetIDToken_
edm::EDGetTokenT< reco::PFJetCollection > pfJetToken_
edm::EDGetTokenT< reco::PFMETCollection > pfMetToken_
edm::EDGetTokenT< reco::MuonCollection > muonToken_
edm::EDGetTokenT< reco::JetCorrector > caloJECToken_
edm::EDGetTokenT< reco::CaloJetCollection > caloJetToken_
edm::EDGetTokenT< reco::CaloMETCollection > caloMetBEToken_
L1JetRecoTreeProducer::~L1JetRecoTreeProducer ( )
override

Definition at line 169 of file L1JetRecoTreeProducer.cc.

170 {
171 
172  // do anything here that needs to be done at desctruction time
173  // (e.g. close files, deallocate resources etc.)
174 
175 }

Member Function Documentation

void L1JetRecoTreeProducer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 183 of file L1JetRecoTreeProducer.cc.

References caloJECToken_, caloJetCorrMissing_, caloJetIDMissing_, caloJetIDToken_, isolatedTracks_cfi::caloJets, caloJetsMissing_, caloJetToken_, CaloMET_cfi::caloMet, CaloMET_cfi::caloMetBE, caloMetBEMissing_, caloMetBEToken_, caloMetMissing_, caloMetToken_, doCaloJetCorr(), doCaloJets(), doCaloMet(), doCaloMetBE(), doPFJetCorr(), doPFJets(), doPFMet(), doPFMetNoMu(), edm::Event::getByToken(), edm::HandleBase::isValid(), jet_data, met_data, extraflags_cff::muons, muonsMissing_, muonToken_, L1Analysis::L1AnalysisRecoJetDataFormat::nCaloJets, L1Analysis::L1AnalysisRecoJetDataFormat::nJets, pfJECToken_, pfJetCorrMissing_, pfJetBenchmark_cfi::pfJets, pfJetsMissing_, pfJetToken_, RecoPFMET_cff::pfMet, pfMetMissing_, pfMetToken_, L1Analysis::L1AnalysisRecoJetDataFormat::Reset(), L1Analysis::L1AnalysisRecoMetDataFormat::Reset(), and tree_.

184 {
185 
186  jet_data->Reset();
187  met_data->Reset();
188 
189  // get jets
191  iEvent.getByToken(pfJetToken_, pfJets);
192 
193  // get calo jets
195  iEvent.getByToken(caloJetToken_, caloJets);
196 
197  //get sums
199  iEvent.getByToken(pfMetToken_, pfMet);
200 
201  // get jet ID
203  iEvent.getByToken(caloJetIDToken_,jetsID);
204 
206  iEvent.getByToken(pfJECToken_, pfJetCorr);
207 
209  iEvent.getByToken(caloJECToken_, caloJetCorr);
210 
212  iEvent.getByToken(caloMetToken_, caloMet);
213 
215  iEvent.getByToken(caloMetBEToken_, caloMetBE);
216 
217  // get muons
219  iEvent.getByToken(muonToken_, muons);
220 
221 
222  if (pfJets.isValid()) {
223 
224  jet_data->nJets=0;
225 
226  doPFJets(pfJets);
227 
228  }
229  else {
230  if (!pfJetsMissing_) {edm::LogWarning("MissingProduct") << "PFJets not found. Branch will not be filled" << std::endl;}
231  pfJetsMissing_ = true;
232  }
233 
234  if (pfJetCorr.isValid()) {
235 
236  doPFJetCorr(pfJets,pfJetCorr);
237 
238  }
239  else {
240  if (!pfJetCorrMissing_) {edm::LogWarning("MissingProduct") << "PF Jet Corrector not found. Branch will not be filled" << std::endl;}
241  pfJetCorrMissing_ = true;
242  }
243 
244  if (caloJets.isValid()) {
245 
246  jet_data->nCaloJets=0;
247 
248  doCaloJets(caloJets);
249 
250  }
251  else {
252  if (!caloJetsMissing_) {edm::LogWarning("MissingProduct") << "Calo Jets not found. Branch will not be filled" << std::endl;}
253  caloJetsMissing_ = true;
254  }
255 
256  if (caloJetCorr.isValid()) {
257 
258  doCaloJetCorr(caloJets,caloJetCorr);
259 
260  }
261  else {
262  if (!caloJetCorrMissing_) {edm::LogWarning("MissingProduct") << "Calo Jet Corrector not found. Branch will not be filled" << std::endl;}
263  caloJetCorrMissing_ = true;
264  }
265 
266  if (!jetsID.isValid()){
267 
268  if (!caloJetIDMissing_) {edm::LogWarning("MissingProduct") << "Calo Jet ID not found. Branch will not be filled" << std::endl;}
269  caloJetIDMissing_ = true;
270  }
271 
272  if (pfMet.isValid()) {
273 
274  doPFMet(pfMet);
275 
276  if (muons.isValid()) {
277 
278  doPFMetNoMu(pfMet,muons);
279 
280  }
281  else {
282  if (!muonsMissing_) {edm::LogWarning("MissingProduct") << "Muons not found. PFMetNoMu branch will not be filled" << std::endl;}
283  muonsMissing_ = true;
284  }
285  }
286  else{
287  if (!pfMetMissing_) {edm::LogWarning("MissingProduct") << "PFMet not found. Branch will not be filled" << std::endl;}
288  pfMetMissing_ = true;
289  }
290 
291  if (caloMet.isValid()) {
292 
293  doCaloMet(caloMet);
294 
295  }
296  else {
297  if (!caloMetMissing_) {edm::LogWarning("MissingProduct") << "CaloMet not found. Branch will not be filled" << std::endl;}
298  caloMetMissing_ = true;
299  }
300 
301  if (caloMetBE.isValid()) {
302 
303  doCaloMetBE(caloMetBE);
304 
305  }
306  else {
307  if (!caloMetBEMissing_) {edm::LogWarning("MissingProduct") << "CaloMetBE not found. Branch will not be filled" << std::endl;}
308  caloMetBEMissing_ = true;
309  }
310 
311  tree_->Fill();
312 
313 }
void doPFMet(edm::Handle< reco::PFMETCollection > pfMet)
L1Analysis::L1AnalysisRecoMetDataFormat * met_data
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void doPFJetCorr(edm::Handle< reco::PFJetCollection > pfJets, edm::Handle< reco::JetCorrector > pfJetCorr)
edm::EDGetTokenT< reco::CaloMETCollection > caloMetToken_
void doPFJets(edm::Handle< reco::PFJetCollection > pfJets)
caloMetBE
____________________________________________________________________________||
Definition: CaloMET_cfi.py:19
edm::EDGetTokenT< reco::JetCorrector > pfJECToken_
void doCaloMetBE(edm::Handle< reco::CaloMETCollection > caloMetBE)
L1Analysis::L1AnalysisRecoJetDataFormat * jet_data
edm::EDGetTokenT< edm::ValueMap< reco::JetID > > caloJetIDToken_
edm::EDGetTokenT< reco::PFJetCollection > pfJetToken_
void doCaloJetCorr(edm::Handle< reco::CaloJetCollection > caloJets, edm::Handle< reco::JetCorrector > caloJetCorr)
bool isValid() const
Definition: HandleBase.h:74
caloMet
____________________________________________________________________________||
Definition: CaloMET_cfi.py:4
void doPFMetNoMu(edm::Handle< reco::PFMETCollection > pfMet, edm::Handle< reco::MuonCollection >)
void doCaloMet(edm::Handle< reco::CaloMETCollection > caloMet)
edm::EDGetTokenT< reco::PFMETCollection > pfMetToken_
void doCaloJets(edm::Handle< reco::CaloJetCollection > caloJets)
edm::EDGetTokenT< reco::MuonCollection > muonToken_
edm::EDGetTokenT< reco::JetCorrector > caloJECToken_
edm::EDGetTokenT< reco::CaloJetCollection > caloJetToken_
edm::EDGetTokenT< reco::CaloMETCollection > caloMetBEToken_
void L1JetRecoTreeProducer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 590 of file L1JetRecoTreeProducer.cc.

591 {
592 }
bool L1JetRecoTreeProducer::caloJetID ( const reco::CaloJet jet)
private

Definition at line 578 of file L1JetRecoTreeProducer.cc.

References tmp.

Referenced by doCaloJetCorr(), and doCaloJets().

578  {
579 
580  bool tmp = true;
581 
582 
583  return tmp;
584 
585 }
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
void L1JetRecoTreeProducer::doCaloJetCorr ( edm::Handle< reco::CaloJetCollection caloJets,
edm::Handle< reco::JetCorrector caloJetCorr 
)
private

Definition at line 443 of file L1JetRecoTreeProducer.cc.

References L1Analysis::L1AnalysisRecoJetDataFormat::caloCorrFactor, L1Analysis::L1AnalysisRecoJetDataFormat::caloEtCorr, L1Analysis::L1AnalysisRecoMetDataFormat::caloHt, caloJetID(), caloJetIDMissing_, reco::JetCorrector::correction(), jet_data, jetetaMax_, jetptThreshold_, maxJet_, met_data, and edm::Handle< T >::product().

Referenced by analyze().

443  {
444 
445 
446  float caloCorrFactor = 1.;
447  unsigned int nCaloJets = 0;
448 
449  met_data->caloHt = 0;
450 
451  for( auto it=caloJets->begin();
452  it!=caloJets->end() && nCaloJets < maxJet_;
453  ++it) {
454 
455  if(!caloJetIDMissing_)
456  if(!caloJetID(*it)) continue;
457 
458  caloCorrFactor = caloJetCorr.product()->correction(*it);
459 
460  jet_data->caloEtCorr.push_back(it->et()*caloCorrFactor);
461  jet_data->caloCorrFactor.push_back(caloCorrFactor);
462 
463  nCaloJets++;
464 
465  if (it->pt()*caloCorrFactor > jetptThreshold_ && fabs(it->eta())<jetetaMax_) {
466  met_data->caloHt += it->pt()*caloCorrFactor;
467  }
468 
469  }
470 }
L1Analysis::L1AnalysisRecoMetDataFormat * met_data
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:49
L1Analysis::L1AnalysisRecoJetDataFormat * jet_data
T const * product() const
Definition: Handle.h:74
bool caloJetID(const reco::CaloJet &jet)
void L1JetRecoTreeProducer::doCaloJets ( edm::Handle< reco::CaloJetCollection caloJets)
private

Definition at line 317 of file L1JetRecoTreeProducer.cc.

References L1Analysis::L1AnalysisRecoJetDataFormat::caloE, L1Analysis::L1AnalysisRecoJetDataFormat::caloEt, L1Analysis::L1AnalysisRecoJetDataFormat::caloEta, caloJetID(), caloJetIDMissing_, L1Analysis::L1AnalysisRecoJetDataFormat::caloPhi, L1Analysis::L1AnalysisRecoJetDataFormat::eEmEB, L1Analysis::L1AnalysisRecoJetDataFormat::eEmEE, L1Analysis::L1AnalysisRecoJetDataFormat::eEMF, L1Analysis::L1AnalysisRecoJetDataFormat::eEmHF, L1Analysis::L1AnalysisRecoJetDataFormat::eHadHB, L1Analysis::L1AnalysisRecoJetDataFormat::eHadHE, L1Analysis::L1AnalysisRecoJetDataFormat::eHadHF, L1Analysis::L1AnalysisRecoJetDataFormat::eHadHO, L1Analysis::L1AnalysisRecoJetDataFormat::eMaxEcalTow, L1Analysis::L1AnalysisRecoJetDataFormat::eMaxHcalTow, jet_data, maxJet_, L1Analysis::L1AnalysisRecoJetDataFormat::n60, L1Analysis::L1AnalysisRecoJetDataFormat::nCaloJets, and L1Analysis::L1AnalysisRecoJetDataFormat::towerArea.

Referenced by analyze().

317  {
318 
319  for( auto it=caloJets->begin();
320  it!=caloJets->end() && jet_data->nCaloJets < maxJet_;
321  ++it) {
322 
323  if(!caloJetIDMissing_)
324  if(!caloJetID(*it)) continue;
325 
326  jet_data->caloEt.push_back(it->et());
327  jet_data->caloEta.push_back(it->eta());
328  jet_data->caloPhi.push_back(it->phi());
329  jet_data->caloE.push_back(it->energy());
330 
331  jet_data->eEMF.push_back(it->emEnergyFraction());
332  jet_data->eEmEB.push_back(it->emEnergyInEB());
333  jet_data->eEmEE.push_back(it->emEnergyInEE());
334  jet_data->eEmHF.push_back(it->emEnergyInHF());
335  jet_data->eHadHB.push_back(it->hadEnergyInHB());
336  jet_data->eHadHE.push_back(it->hadEnergyInHE());
337  jet_data->eHadHO.push_back(it->hadEnergyInHO());
338  jet_data->eHadHF.push_back(it->hadEnergyInHF());
339  jet_data->eMaxEcalTow.push_back(it->maxEInEmTowers());
340  jet_data->eMaxHcalTow.push_back(it->maxEInHadTowers());
341  jet_data->towerArea.push_back(it->towersArea());
342  jet_data->n60.push_back(it->n60());
343 
344  jet_data->nCaloJets++;
345 
346  }
347 
348 
349 }
L1Analysis::L1AnalysisRecoJetDataFormat * jet_data
bool caloJetID(const reco::CaloJet &jet)
void L1JetRecoTreeProducer::doCaloMet ( edm::Handle< reco::CaloMETCollection caloMet)
private

Definition at line 523 of file L1JetRecoTreeProducer.cc.

References L1Analysis::L1AnalysisRecoMetDataFormat::caloMet, L1Analysis::L1AnalysisRecoMetDataFormat::caloMetPhi, L1Analysis::L1AnalysisRecoMetDataFormat::caloSumEt, reco::LeafCandidate::et(), met_data, reco::LeafCandidate::phi(), edm::Handle< T >::product(), and reco::MET::sumEt().

Referenced by analyze().

523  {
524 
525  const reco::CaloMETCollection *metCol = caloMet.product();
526  const reco::CaloMET theMet = metCol->front();
527 
528  met_data->caloMet = theMet.et();
529  met_data->caloMetPhi = theMet.phi();
530  met_data->caloSumEt = theMet.sumEt();
531 
532 }
L1Analysis::L1AnalysisRecoMetDataFormat * met_data
double sumEt() const
Definition: MET.h:56
double et() const final
transverse energy
std::vector< reco::CaloMET > CaloMETCollection
collection of CaloMET objects
T const * product() const
Definition: Handle.h:74
double phi() const final
momentum azimuthal angle
void L1JetRecoTreeProducer::doCaloMetBE ( edm::Handle< reco::CaloMETCollection caloMetBE)
private

Definition at line 535 of file L1JetRecoTreeProducer.cc.

References L1Analysis::L1AnalysisRecoMetDataFormat::caloMetBE, L1Analysis::L1AnalysisRecoMetDataFormat::caloMetPhiBE, L1Analysis::L1AnalysisRecoMetDataFormat::caloSumEtBE, reco::LeafCandidate::et(), met_data, reco::LeafCandidate::phi(), edm::Handle< T >::product(), and reco::MET::sumEt().

Referenced by analyze().

535  {
536 
537  const reco::CaloMETCollection *metCol = caloMetBE.product();
538  const reco::CaloMET theMet = metCol->front();
539 
540  met_data->caloMetBE = theMet.et();
541  met_data->caloMetPhiBE = theMet.phi();
542  met_data->caloSumEtBE = theMet.sumEt();
543 
544 }
L1Analysis::L1AnalysisRecoMetDataFormat * met_data
double sumEt() const
Definition: MET.h:56
double et() const final
transverse energy
std::vector< reco::CaloMET > CaloMETCollection
collection of CaloMET objects
T const * product() const
Definition: Handle.h:74
double phi() const final
momentum azimuthal angle
void L1JetRecoTreeProducer::doPFJetCorr ( edm::Handle< reco::PFJetCollection pfJets,
edm::Handle< reco::JetCorrector pfJetCorr 
)
private

Definition at line 397 of file L1JetRecoTreeProducer.cc.

References reco::JetCorrector::correction(), L1Analysis::L1AnalysisRecoJetDataFormat::corrFactor, L1Analysis::L1AnalysisRecoJetDataFormat::etCorr, L1Analysis::L1AnalysisRecoMetDataFormat::Ht, jet_data, jetetaMax_, jetptThreshold_, maxJet_, met_data, L1Analysis::L1AnalysisRecoMetDataFormat::mHt, L1Analysis::L1AnalysisRecoMetDataFormat::mHtPhi, pfJetID(), and edm::Handle< T >::product().

Referenced by analyze().

397  {
398 
399 
400  float corrFactor = 1.;
401  unsigned int nJets = 0;
402 
403  float mHx = 0;
404  float mHy = 0;
405 
406  met_data->Ht = 0;
407  met_data->mHt = -999.;
408  met_data->mHtPhi = -999.;
409 
410 
411  for( auto it=pfJets->begin();
412  it!=pfJets->end() && nJets < maxJet_;
413  ++it) {
414 
415  if(!pfJetID(*it)) continue;
416 
417  corrFactor = pfJetCorr.product()->correction(*it);
418 
419  jet_data->etCorr.push_back(it->et()*corrFactor);
420  jet_data->corrFactor.push_back(corrFactor);
421 
422  nJets++;
423 
424  if (it->pt()*corrFactor > jetptThreshold_ && fabs(it->eta())<jetetaMax_) {
425  mHx += -1.*it->px()*corrFactor;
426  mHy += -1.*it->py()*corrFactor;
427  met_data->Ht += it->pt()*corrFactor;
428  }
429 
430  }
431 
432  TVector2 tv2 = TVector2(mHx,mHy);
433  met_data->mHt = tv2.Mod();
434  met_data->mHtPhi = tv2.Phi();
435 
436 
437 }
bool pfJetID(const reco::PFJet &jet)
L1Analysis::L1AnalysisRecoMetDataFormat * met_data
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:49
L1Analysis::L1AnalysisRecoJetDataFormat * jet_data
T const * product() const
Definition: Handle.h:74
void L1JetRecoTreeProducer::doPFJets ( edm::Handle< reco::PFJetCollection pfJets)
private

Definition at line 353 of file L1JetRecoTreeProducer.cc.

References L1Analysis::L1AnalysisRecoJetDataFormat::cemef, L1Analysis::L1AnalysisRecoJetDataFormat::chef, L1Analysis::L1AnalysisRecoJetDataFormat::chMult, L1Analysis::L1AnalysisRecoJetDataFormat::cmef, L1Analysis::L1AnalysisRecoJetDataFormat::cMult, L1Analysis::L1AnalysisRecoJetDataFormat::e, L1Analysis::L1AnalysisRecoJetDataFormat::eef, L1Analysis::L1AnalysisRecoJetDataFormat::elMult, L1Analysis::L1AnalysisRecoJetDataFormat::et, L1Analysis::L1AnalysisRecoJetDataFormat::eta, L1Analysis::L1AnalysisRecoJetDataFormat::hfemef, L1Analysis::L1AnalysisRecoJetDataFormat::hfemMult, L1Analysis::L1AnalysisRecoJetDataFormat::hfhef, L1Analysis::L1AnalysisRecoJetDataFormat::hfhMult, jet_data, maxJet_, L1Analysis::L1AnalysisRecoJetDataFormat::mef, L1Analysis::L1AnalysisRecoJetDataFormat::muMult, L1Analysis::L1AnalysisRecoJetDataFormat::nemef, L1Analysis::L1AnalysisRecoJetDataFormat::nhef, L1Analysis::L1AnalysisRecoJetDataFormat::nhMult, L1Analysis::L1AnalysisRecoJetDataFormat::nJets, L1Analysis::L1AnalysisRecoJetDataFormat::nMult, L1Analysis::L1AnalysisRecoJetDataFormat::pef, pfJetID(), L1Analysis::L1AnalysisRecoJetDataFormat::phi, and L1Analysis::L1AnalysisRecoJetDataFormat::phMult.

Referenced by analyze().

353  {
354 
355 
356  for( auto it=pfJets->begin();
357  it!=pfJets->end() && jet_data->nJets < maxJet_;
358  ++it) {
359 
360  if(!pfJetID(*it)) continue;
361 
362 
363  jet_data->et.push_back(it->et());
364  jet_data->eta.push_back(it->eta());
365  jet_data->phi.push_back(it->phi());
366  jet_data->e.push_back(it->energy());
367 
368  jet_data->chef.push_back(it->chargedHadronEnergyFraction());
369  jet_data->nhef.push_back(it->neutralHadronEnergyFraction());
370  jet_data->pef.push_back(it->photonEnergyFraction());
371  jet_data->eef.push_back(it->electronEnergyFraction());
372  jet_data->mef.push_back(it->muonEnergyFraction());
373  jet_data->hfhef.push_back(it->HFHadronEnergyFraction());
374  jet_data->hfemef.push_back(it->HFEMEnergyFraction());
375  jet_data->chMult.push_back(it->chargedHadronMultiplicity());
376  jet_data->nhMult.push_back(it->neutralHadronMultiplicity());
377  jet_data->phMult.push_back(it->photonMultiplicity());
378  jet_data->elMult.push_back(it->electronMultiplicity());
379  jet_data->muMult.push_back(it->muonMultiplicity());
380  jet_data->hfhMult.push_back(it->HFHadronMultiplicity());
381  jet_data->hfemMult.push_back(it->HFEMMultiplicity());
382 
383  jet_data->cemef.push_back(it->chargedEmEnergyFraction());
384  jet_data->cmef.push_back(it->chargedMuEnergyFraction());
385  jet_data->nemef.push_back(it->neutralEmEnergyFraction());
386  jet_data->cMult.push_back(it->chargedMultiplicity());
387  jet_data->nMult.push_back(it->neutralMultiplicity());
388 
389  jet_data->nJets++;
390 
391  }
392 
393 }
bool pfJetID(const reco::PFJet &jet)
L1Analysis::L1AnalysisRecoJetDataFormat * jet_data
void L1JetRecoTreeProducer::doPFMet ( edm::Handle< reco::PFMETCollection pfMet)
private

Definition at line 473 of file L1JetRecoTreeProducer.cc.

References reco::LeafCandidate::et(), L1Analysis::L1AnalysisRecoMetDataFormat::met, met_data, L1Analysis::L1AnalysisRecoMetDataFormat::metPhi, L1Analysis::L1AnalysisRecoMetDataFormat::metPx, L1Analysis::L1AnalysisRecoMetDataFormat::metPy, reco::LeafCandidate::phi(), edm::Handle< T >::product(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::MET::sumEt(), and L1Analysis::L1AnalysisRecoMetDataFormat::sumEt.

Referenced by analyze().

473  {
474 
475  const reco::PFMETCollection *metCol = pfMet.product();
476  const reco::PFMET theMet = metCol->front();
477 
478  met_data->met = theMet.et();
479  met_data->metPhi = theMet.phi();
480  met_data->sumEt = theMet.sumEt();
481  met_data->metPx = theMet.px();
482  met_data->metPy = theMet.py();
483 
484 }
L1Analysis::L1AnalysisRecoMetDataFormat * met_data
double px() const final
x coordinate of momentum vector
double sumEt() const
Definition: MET.h:56
double et() const final
transverse energy
T const * product() const
Definition: Handle.h:74
double py() const final
y coordinate of momentum vector
std::vector< reco::PFMET > PFMETCollection
collection of PFMET objects
double phi() const final
momentum azimuthal angle
void L1JetRecoTreeProducer::doPFMetNoMu ( edm::Handle< reco::PFMETCollection pfMet,
edm::Handle< reco::MuonCollection muons 
)
private

Definition at line 487 of file L1JetRecoTreeProducer.cc.

References reco::LeafCandidate::et(), met_data, L1Analysis::L1AnalysisRecoMetDataFormat::pfMetNoMu, L1Analysis::L1AnalysisRecoMetDataFormat::pfMetNoMuPhi, L1Analysis::L1AnalysisRecoMetDataFormat::pfMetNoMuPx, L1Analysis::L1AnalysisRecoMetDataFormat::pfMetNoMuPy, reco::LeafCandidate::phi(), edm::Handle< T >::product(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), and reco::LeafCandidate::setP4().

Referenced by analyze().

487  {
488 
489  const reco::PFMETCollection *metCol = pfMet.product();
490  const reco::PFMET theMet = metCol->front();
491  reco::PFMET thePFMetNoMu = metCol->front();
492 
493  double pfMetNoMuPx = theMet.px();
494  double pfMetNoMuPy = theMet.py();
495 
496  double muPx(0.), muPy(0.);
497 
498  for( auto it=muons->begin();
499  it!=muons->end(); ++it) {
500  if(it->isPFMuon()){
501  muPx += it->px();
502  muPy += it->py();
503  }
504  }
505 
506  pfMetNoMuPx += muPx;
507  pfMetNoMuPy += muPy;
508 
509  math::XYZTLorentzVector pfMetNoMuP4(pfMetNoMuPx,pfMetNoMuPy,0,hypot(pfMetNoMuPx,pfMetNoMuPy));
510 
511 
512  thePFMetNoMu.setP4(pfMetNoMuP4);
513 
514  met_data->pfMetNoMu = thePFMetNoMu.et();
515  met_data->pfMetNoMuPhi = thePFMetNoMu.phi();
516  met_data->pfMetNoMuPx = thePFMetNoMu.px();
517  met_data->pfMetNoMuPy = thePFMetNoMu.py();
518 
519 }
L1Analysis::L1AnalysisRecoMetDataFormat * met_data
double px() const final
x coordinate of momentum vector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double et() const final
transverse energy
T const * product() const
Definition: Handle.h:74
double py() const final
y coordinate of momentum vector
std::vector< reco::PFMET > PFMETCollection
collection of PFMET objects
double phi() const final
momentum azimuthal angle
void setP4(const LorentzVector &p4) final
set 4-momentum
void L1JetRecoTreeProducer::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 596 of file L1JetRecoTreeProducer.cc.

References DEFINE_FWK_MODULE.

596  {
597 }
bool L1JetRecoTreeProducer::pfJetID ( const reco::PFJet jet)
private

Definition at line 547 of file L1JetRecoTreeProducer.cc.

References reco::PFJet::chargedEmEnergyFraction(), reco::PFJet::chargedHadronEnergyFraction(), reco::PFJet::chargedMultiplicity(), reco::LeafCandidate::eta(), reco::PFJet::muonEnergyFraction(), reco::PFJet::neutralEmEnergyFraction(), reco::PFJet::neutralHadronEnergyFraction(), reco::PFJet::neutralMultiplicity(), and tmp.

Referenced by doPFJetCorr(), and doPFJets().

547  {
548 
549  bool tmp = true;
550  if (fabs(jet.eta()) < 2.7) {
551  tmp &= jet.neutralHadronEnergyFraction() < 0.9 ;
552  tmp &= jet.neutralEmEnergyFraction() < 0.9 ;
553  tmp &= (jet.chargedMultiplicity() + jet.neutralMultiplicity()) > 1 ;
554  tmp &= jet.muonEnergyFraction() < 0.8 ;
555  tmp &= jet.chargedHadronEnergyFraction() > 0.0 ;
556  tmp &= jet.chargedMultiplicity() > 0 ;
557  tmp &= jet.chargedEmEnergyFraction() < 0.9 ;
558  }
559  if (fabs(jet.eta()) > 2.7 && fabs(jet.eta()) < 3.0){
560  tmp &= jet.neutralEmEnergyFraction() > 0.01 ;
561  tmp &= jet.neutralHadronEnergyFraction() < 0.98 ;
562  tmp &= jet.neutralMultiplicity() > 2 ;
563  }
564  if (fabs(jet.eta()) > 3.0) {
565  tmp &= jet.neutralEmEnergyFraction() < 0.9 ;
566  tmp &= jet.neutralMultiplicity() > 10 ;
567  }
568 
569  // our custom selection
570  //tmp &= jet.muonMultiplicity() == 0;
571  //tmp &= jet.electronMultiplicity() == 0;
572 
573  return tmp;
574 
575 }
double eta() const final
momentum pseudorapidity
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction
Definition: PFJet.h:100
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:155
float neutralEmEnergyFraction() const
neutralEmEnergyFraction
Definition: PFJet.h:152
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction
Definition: PFJet.h:104
int neutralMultiplicity() const
neutralMultiplicity
Definition: PFJet.h:157
float chargedEmEnergyFraction() const
chargedEmEnergyFraction
Definition: PFJet.h:144
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
float muonEnergyFraction() const
muonEnergyFraction
Definition: PFJet.h:116

Member Data Documentation

edm::EDGetTokenT<reco::JetCorrector> L1JetRecoTreeProducer::caloJECToken_
private

Definition at line 101 of file L1JetRecoTreeProducer.cc.

Referenced by analyze(), and L1JetRecoTreeProducer().

bool L1JetRecoTreeProducer::caloJetCorrMissing_
private

Definition at line 119 of file L1JetRecoTreeProducer.cc.

Referenced by analyze().

bool L1JetRecoTreeProducer::caloJetIDMissing_
private

Definition at line 121 of file L1JetRecoTreeProducer.cc.

Referenced by analyze(), doCaloJetCorr(), and doCaloJets().

edm::EDGetTokenT<edm::ValueMap<reco::JetID> > L1JetRecoTreeProducer::caloJetIDToken_
private

Definition at line 99 of file L1JetRecoTreeProducer.cc.

Referenced by analyze(), and L1JetRecoTreeProducer().

bool L1JetRecoTreeProducer::caloJetsMissing_
private

Definition at line 120 of file L1JetRecoTreeProducer.cc.

Referenced by analyze().

edm::EDGetTokenT<reco::CaloJetCollection> L1JetRecoTreeProducer::caloJetToken_
private

Definition at line 98 of file L1JetRecoTreeProducer.cc.

Referenced by analyze(), and L1JetRecoTreeProducer().

bool L1JetRecoTreeProducer::caloMetBEMissing_
private

Definition at line 124 of file L1JetRecoTreeProducer.cc.

Referenced by analyze().

edm::EDGetTokenT<reco::CaloMETCollection> L1JetRecoTreeProducer::caloMetBEToken_
private

Definition at line 105 of file L1JetRecoTreeProducer.cc.

Referenced by analyze(), and L1JetRecoTreeProducer().

bool L1JetRecoTreeProducer::caloMetMissing_
private

Definition at line 123 of file L1JetRecoTreeProducer.cc.

Referenced by analyze().

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

Definition at line 104 of file L1JetRecoTreeProducer.cc.

Referenced by analyze(), and L1JetRecoTreeProducer().

edm::Service<TFileService> L1JetRecoTreeProducer::fs_
private

Definition at line 91 of file L1JetRecoTreeProducer.cc.

Referenced by L1JetRecoTreeProducer().

L1Analysis::L1AnalysisRecoJetDataFormat* L1JetRecoTreeProducer::jet_data
double L1JetRecoTreeProducer::jetetaMax_
private

Definition at line 112 of file L1JetRecoTreeProducer.cc.

Referenced by doCaloJetCorr(), doPFJetCorr(), and L1JetRecoTreeProducer().

double L1JetRecoTreeProducer::jetptThreshold_
private

Definition at line 111 of file L1JetRecoTreeProducer.cc.

Referenced by doCaloJetCorr(), doPFJetCorr(), and L1JetRecoTreeProducer().

unsigned int L1JetRecoTreeProducer::maxCl_
private

Definition at line 113 of file L1JetRecoTreeProducer.cc.

unsigned int L1JetRecoTreeProducer::maxJet_
private
unsigned int L1JetRecoTreeProducer::maxTrk_
private

Definition at line 116 of file L1JetRecoTreeProducer.cc.

unsigned int L1JetRecoTreeProducer::maxVtx_
private

Definition at line 115 of file L1JetRecoTreeProducer.cc.

L1Analysis::L1AnalysisRecoMetDataFormat* L1JetRecoTreeProducer::met_data
bool L1JetRecoTreeProducer::muonsMissing_
private

Definition at line 126 of file L1JetRecoTreeProducer.cc.

Referenced by analyze().

edm::EDGetTokenT<reco::MuonCollection> L1JetRecoTreeProducer::muonToken_
private

Definition at line 107 of file L1JetRecoTreeProducer.cc.

Referenced by analyze(), and L1JetRecoTreeProducer().

edm::EDGetTokenT<reco::JetCorrector> L1JetRecoTreeProducer::pfJECToken_
private

Definition at line 100 of file L1JetRecoTreeProducer.cc.

Referenced by analyze(), and L1JetRecoTreeProducer().

bool L1JetRecoTreeProducer::pfJetCorrMissing_
private

Definition at line 118 of file L1JetRecoTreeProducer.cc.

Referenced by analyze().

bool L1JetRecoTreeProducer::pfJetsMissing_
private

Definition at line 110 of file L1JetRecoTreeProducer.cc.

Referenced by analyze().

edm::EDGetTokenT<reco::PFJetCollection> L1JetRecoTreeProducer::pfJetToken_
private

Definition at line 97 of file L1JetRecoTreeProducer.cc.

Referenced by analyze(), and L1JetRecoTreeProducer().

bool L1JetRecoTreeProducer::pfMetMissing_
private

Definition at line 122 of file L1JetRecoTreeProducer.cc.

Referenced by analyze().

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

Definition at line 103 of file L1JetRecoTreeProducer.cc.

Referenced by analyze(), and L1JetRecoTreeProducer().

TTree* L1JetRecoTreeProducer::tree_
private

Definition at line 94 of file L1JetRecoTreeProducer.cc.

Referenced by analyze(), and L1JetRecoTreeProducer().