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 Attributes
MVAJetIdMaker Class Reference

#include <MVAJetIdMaker.h>

Public Member Functions

 MVAJetIdMaker (const edm::ParameterSet &, edm::ConsumesCollector)
 
void SetVars (HWW &, const edm::Event &, const edm::EventSetup &)
 

Private Attributes

PileupJetIdAlgofPUJetIdAlgo
 
std::string jetCorrector_
 
edm::EDGetTokenT
< reco::PFJetCollection
PFJetCollection_
 
edm::EDGetTokenT
< reco::VertexCollection
thePVCollection_
 
edm::EDGetTokenT< double > theRhoCollection_
 

Detailed Description

Definition at line 11 of file MVAJetIdMaker.h.

Constructor & Destructor Documentation

MVAJetIdMaker::MVAJetIdMaker ( const edm::ParameterSet iConfig,
edm::ConsumesCollector  iCollector 
)

Definition at line 7 of file MVAJetIdMaker.cc.

References edm::ConsumesCollector::consumes(), fPUJetIdAlgo, edm::ParameterSet::getParameter(), jetCorrector_, PFJetCollection_, AlCaHLTBitMon_QueryRunRegistry::string, thePVCollection_, and theRhoCollection_.

7  {
8 
9  PFJetCollection_ = iCollector.consumes<reco::PFJetCollection> (iConfig.getParameter<edm::InputTag>("pfJetsInputTag"));
10  thePVCollection_ = iCollector.consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("primaryVertexInputTag"));
11  jetCorrector_ = iConfig.getParameter<std::string>("jetCorrector");
12  theRhoCollection_ = iCollector.consumes<double>(iConfig.getParameter<edm::InputTag>("rhoInputTag"));
13 
14  fPUJetIdAlgo = new PileupJetIdAlgo(iConfig.getParameter<edm::ParameterSet>("puJetIDParams"),true);
15 
16 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< reco::PFJetCollection > PFJetCollection_
Definition: MVAJetIdMaker.h:20
T getParameter(std::string const &) const
edm::EDGetTokenT< double > theRhoCollection_
Definition: MVAJetIdMaker.h:23
PileupJetIdAlgo * fPUJetIdAlgo
Definition: MVAJetIdMaker.h:24
std::string jetCorrector_
Definition: MVAJetIdMaker.h:22
edm::EDGetTokenT< reco::VertexCollection > thePVCollection_
Definition: MVAJetIdMaker.h:21
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< PFJet > PFJetCollection
collection of PFJet objects

Member Function Documentation

void MVAJetIdMaker::SetVars ( HWW hww,
const edm::Event iEvent,
const edm::EventSetup iSetup 
)

Definition at line 29 of file MVAJetIdMaker.cc.

References StoredPileupJetIdentifier::beta(), StoredPileupJetIdentifier::betaStar(), PileupJetIdAlgo::computeIdVariables(), mvaPFMET_cff::corrector, PileupJetIdentifier::d0(), StoredPileupJetIdentifier::dRMean(), StoredPileupJetIdentifier::dZ(), reco::LeafCandidate::eta(), fPUJetIdAlgo, StoredPileupJetIdentifier::frac01(), StoredPileupJetIdentifier::frac02(), StoredPileupJetIdentifier::frac03(), StoredPileupJetIdentifier::frac04(), StoredPileupJetIdentifier::frac05(), edm::Event::getByToken(), JetCorrector::getJetCorrector(), edm::EventBase::id(), reco::Vertex::isFake(), metsig::jet, reco::Jet::jetArea(), jetCorrector_, StoredPileupJetIdentifier::jetEta(), PileupJetIdentifier::jetPhi(), HWW::Load_pfjets_corr_p4(), HWW::Load_pfjets_JEC(), HWW::Load_pfjets_mvavalue(), LogDebug, PileupJetIdentifier::mva(), StoredPileupJetIdentifier::nCharged(), reco::Vertex::ndof(), StoredPileupJetIdentifier::nNeutrals(), StoredPileupJetIdentifier::nvtx(), PFJetCollection_, HWW::pfjets_JEC(), HWW::pfjets_mvavalue(), reco::Vertex::position(), reco::LeafCandidate::pt(), dt_dqm_sourceclient_common_cff::reco, reco::Jet::scaleEnergy(), thePVCollection_, and theRhoCollection_.

Referenced by HWWAnalyzer::analyze().

29  {
30 
31  using namespace std;
32  using namespace edm;
33  using namespace reco;
34 
35  hww.Load_pfjets_corr_p4();
37  hww.Load_pfjets_JEC();
38 
39  bool validToken;
40 
41  //Uncorrected Jets
42  Handle<PFJetCollection> lHUCJets;
43  validToken = iEvent.getByToken(PFJetCollection_, lHUCJets);
44  if(!validToken) return;
45  PFJetCollection lUCJets = *lHUCJets;
46 
47  // vertices
49  validToken = iEvent.getByToken(thePVCollection_, lHVertices);
50  if(!validToken) return;
51  VertexCollection lVertices = *lHVertices;
52 
53  const JetCorrector* corrector=0;
55 
56  Handle<double> lHrho;
57  validToken = iEvent.getByToken(theRhoCollection_, lHrho);
58  if(!validToken) return;
59  double lrho = *lHrho;
60 
61  std::vector<reco::PFJet> lCJets;
62  for(reco::PFJetCollection::const_iterator jet=lUCJets.begin(); jet!=lUCJets.end(); ++jet){
63 
64  reco::PFJet tempJet = *jet;
65  tempJet.scaleEnergy(corrector ? corrector->correction(*jet, iEvent, iSetup) : 1.);
66  lCJets.push_back(tempJet);
67 
68  }
69 
70 
71  // select good vertices
72  // make new collection to put into computeIdVariables(...)
73  VertexCollection lGoodVertices;
74  for(int ivtx = 0; ivtx < (int)lVertices.size(); ivtx++) {
75 
76  const Vertex *vtx = &(lVertices.at(ivtx));
77  if( vtx->isFake() ) continue;
78  if( vtx->ndof()<=4 ) continue;
79  if( vtx->position().Rho()>2.0 ) continue;
80  if( fabs(vtx->position().Z())>24.0 ) continue;
81  lGoodVertices.push_back(*vtx);
82  }
83 
84  // loop over jets
85  for(int i0 = 0; i0 < (int) lUCJets.size(); i0++) { // uncorrected jets collection
86  const PFJet *pUCJet = &(lUCJets.at(i0));
87  for(int i1 = 0; i1 < (int) lCJets.size(); i1++) { // corrected jets collection
88  const PFJet *pCJet = &(lCJets.at(i1));
89  if( pUCJet->jetArea() != pCJet->jetArea() ) continue;
90  if( fabs(pUCJet->eta() - pCJet->eta()) > 0.001 ) continue;
91  if( pUCJet->pt() < 0.0 ) continue;
92  double lJec = pCJet ->pt()/pUCJet->pt();
93 
94  // calculate mva value only when there are good vertices
95  // otherwise store -999
96  if( lGoodVertices.size()>0 ) {
97  PileupJetIdentifier lPUJetId = fPUJetIdAlgo->computeIdVariables(pCJet,lJec,&lGoodVertices[0],lGoodVertices,lrho);
98  hww.pfjets_mvavalue() .push_back( lPUJetId.mva() );
99  hww.pfjets_JEC() .push_back( lJec );
100 
101  // print out MVA inputs
102  if(0) {
103 
104  LogDebug("MVAJetIdMaker")
105  << "Debug Jet MVA : "
106  << iEvent.id() << " : "
107  << lPUJetId.nvtx() << " "
108  << pCJet->pt() << " "
109  << lPUJetId.jetEta() << " "
110  << lPUJetId.jetPhi() << " "
111  << lPUJetId.d0() << " "
112  << lPUJetId.dZ() << " "
113  << lPUJetId.beta() << " "
114  << lPUJetId.betaStar() << " "
115  << lPUJetId.nCharged() << " "
116  << lPUJetId.nNeutrals() << " "
117  << lPUJetId.dRMean() << " "
118  << lPUJetId.frac01() << " "
119  << lPUJetId.frac02() << " "
120  << lPUJetId.frac03() << " "
121  << lPUJetId.frac04() << " "
122  << lPUJetId.frac05()
123  << " === : === "
124  << lPUJetId.mva();
125  }
126  }
127  else
128 
129  hww.pfjets_mvavalue() .push_back( -999. );
130 
131  break;
132 
133  } // lCJets
134  } // lUCJets
135 }
#define LogDebug(id)
edm::EDGetTokenT< reco::PFJetCollection > PFJetCollection_
Definition: MVAJetIdMaker.h:20
edm::EDGetTokenT< double > theRhoCollection_
Definition: MVAJetIdMaker.h:23
virtual void scaleEnergy(double fScale)
scale energy of the jet
const float & frac05() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
PileupJetIdAlgo * fPUJetIdAlgo
Definition: MVAJetIdMaker.h:24
const float & mva() const
std::string jetCorrector_
Definition: MVAJetIdMaker.h:22
edm::EDGetTokenT< reco::VertexCollection > thePVCollection_
Definition: MVAJetIdMaker.h:21
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const float & frac03() const
void Load_pfjets_mvavalue()
Definition: HWW.cc:1376
void Load_pfjets_corr_p4()
Definition: HWW.cc:1367
const Point & position() const
position
Definition: Vertex.h:106
Jets made from PFObjects.
Definition: PFJet.h:21
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
tuple corrector
Definition: mvaPFMET_cff.py:88
void Load_pfjets_JEC()
Definition: HWW.cc:1373
const float & beta() const
double ndof() const
Definition: Vertex.h:102
const float & jetPhi() const
const float & nNeutrals() const
const float & dZ() const
const float & nvtx() const
const float & frac04() const
const float & frac01() const
bool isFake() const
Definition: Vertex.h:64
PileupJetIdentifier computeIdVariables(const reco::Jet *jet, float jec, const reco::Vertex *, const reco::VertexCollection &, double rho)
std::vector< float > & pfjets_mvavalue()
Definition: HWW.cc:785
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
edm::EventID id() const
Definition: EventBase.h:60
std::vector< PFJet > PFJetCollection
collection of PFJet objects
virtual float jetArea() const
get jet area
Definition: Jet.h:105
std::vector< float > & pfjets_JEC()
Definition: HWW.cc:781
const float & betaStar() const
const float & dRMean() const
const float & d0() const
const float & nCharged() const
const float & frac02() const
const float & jetEta() const

Member Data Documentation

PileupJetIdAlgo* MVAJetIdMaker::fPUJetIdAlgo
private

Definition at line 24 of file MVAJetIdMaker.h.

Referenced by MVAJetIdMaker(), and SetVars().

std::string MVAJetIdMaker::jetCorrector_
private

Definition at line 22 of file MVAJetIdMaker.h.

Referenced by MVAJetIdMaker(), and SetVars().

edm::EDGetTokenT<reco::PFJetCollection> MVAJetIdMaker::PFJetCollection_
private

Definition at line 20 of file MVAJetIdMaker.h.

Referenced by MVAJetIdMaker(), and SetVars().

edm::EDGetTokenT<reco::VertexCollection> MVAJetIdMaker::thePVCollection_
private

Definition at line 21 of file MVAJetIdMaker.h.

Referenced by MVAJetIdMaker(), and SetVars().

edm::EDGetTokenT<double> MVAJetIdMaker::theRhoCollection_
private

Definition at line 23 of file MVAJetIdMaker.h.

Referenced by MVAJetIdMaker(), and SetVars().