CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MVAJetIdMaker.cc
Go to the documentation of this file.
4 
6 
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 
13  fPUJetIdAlgo = new PileupJetIdAlgo(iConfig);
14 
15 }
16 
17 bool passPFLooseId(const reco::PFJet *iJet) {
18  if(iJet->energy()== 0) return false;
19  if(iJet->neutralHadronEnergy()/iJet->energy() > 0.99) return false;
20  if(iJet->neutralEmEnergy()/iJet->energy() > 0.99) return false;
21  if(iJet->nConstituents() < 2) return false;
22  if(iJet->chargedHadronEnergy()/iJet->energy() <= 0 && fabs(iJet->eta()) < 2.4 ) return false;
23  if(iJet->chargedEmEnergy()/iJet->energy() > 0.99 && fabs(iJet->eta()) < 2.4 ) return false;
24  if(iJet->chargedMultiplicity() < 1 && fabs(iJet->eta()) < 2.4 ) return false;
25  return true;
26 }
27 
28 void MVAJetIdMaker::SetVars(HWW& hww, const edm::Event& iEvent, const edm::EventSetup& iSetup){
29 
30  using namespace std;
31  using namespace edm;
32  using namespace reco;
33 
34  hww.Load_pfjets_corr_p4();
36  hww.Load_pfjets_JEC();
37 
38  bool validToken;
39 
40  //Uncorrected Jets
41  Handle<PFJetCollection> lHUCJets;
42  validToken = iEvent.getByToken(PFJetCollection_, lHUCJets);
43  if(!validToken) return;
44  PFJetCollection lUCJets = *lHUCJets;
45 
46  // vertices
48  validToken = iEvent.getByToken(thePVCollection_, lHVertices);
49  if(!validToken) return;
50  VertexCollection lVertices = *lHVertices;
51 
52  const JetCorrector* corrector=0;
54  std::vector<reco::PFJet> lCJets;
55  for(reco::PFJetCollection::const_iterator jet=lUCJets.begin(); jet!=lUCJets.end(); ++jet){
56 
57  reco::PFJet tempJet = *jet;
58  tempJet.scaleEnergy(corrector ? corrector->correction(*jet, iEvent, iSetup) : 1.);
59  lCJets.push_back(tempJet);
60 
61  }
62 
63 
64  // select good vertices
65  // make new collection to put into computeIdVariables(...)
66  VertexCollection lGoodVertices;
67  for(int ivtx = 0; ivtx < (int)lVertices.size(); ivtx++) {
68 
69  const Vertex *vtx = &(lVertices.at(ivtx));
70  if( vtx->isFake() ) continue;
71  if( vtx->ndof()<=4 ) continue;
72  if( vtx->position().Rho()>2.0 ) continue;
73  if( fabs(vtx->position().Z())>24.0 ) continue;
74  lGoodVertices.push_back(*vtx);
75  }
76 
77  // loop over jets
78  for(int i0 = 0; i0 < (int) lUCJets.size(); i0++) { // uncorrected jets collection
79  const PFJet *pUCJet = &(lUCJets.at(i0));
80  for(int i1 = 0; i1 < (int) lCJets.size(); i1++) { // corrected jets collection
81  const PFJet *pCJet = &(lCJets.at(i1));
82  if( pUCJet->jetArea() != pCJet->jetArea() ) continue;
83  if( fabs(pUCJet->eta() - pCJet->eta()) > 0.001 ) continue;
84  if( pUCJet->pt() < 0.0 ) continue;
85  double lJec = pCJet ->pt()/pUCJet->pt();
86 
87  // calculate mva value only when there are good vertices
88  // otherwise store -999
89  if( lGoodVertices.size()>0 ) {
90  PileupJetIdentifier lPUJetId = fPUJetIdAlgo->computeIdVariables(pCJet,lJec,&lGoodVertices[0],lGoodVertices,true);
91  hww.pfjets_mvavalue() .push_back( lPUJetId.mva() );
92  hww.pfjets_JEC() .push_back( lJec );
93 
94  // print out MVA inputs
95  if(0) {
96 
97  LogDebug("MVAJetIdMaker")
98  << "Debug Jet MVA : "
99  << iEvent.id() << " : "
100  << lPUJetId.nvtx() << " "
101  << pCJet->pt() << " "
102  << lPUJetId.jetEta() << " "
103  << lPUJetId.jetPhi() << " "
104  << lPUJetId.d0() << " "
105  << lPUJetId.dZ() << " "
106  << lPUJetId.beta() << " "
107  << lPUJetId.betaStar() << " "
108  << lPUJetId.nCharged() << " "
109  << lPUJetId.nNeutrals() << " "
110  << lPUJetId.dRMean() << " "
111  << lPUJetId.frac01() << " "
112  << lPUJetId.frac02() << " "
113  << lPUJetId.frac03() << " "
114  << lPUJetId.frac04() << " "
115  << lPUJetId.frac05()
116  << " === : === "
117  << lPUJetId.mva();
118  }
119  }
120  else
121 
122  hww.pfjets_mvavalue() .push_back( -999. );
123 
124  break;
125 
126  } // lCJets
127  } // lUCJets
128 }
#define LogDebug(id)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< reco::PFJetCollection > PFJetCollection_
Definition: MVAJetIdMaker.h:20
T getParameter(std::string const &) const
virtual void scaleEnergy(double fScale)
scale energy of the jet
float chargedEmEnergy() const
chargedEmEnergy
Definition: PFJet.h:142
virtual float pt() const
transverse momentum
const float & frac05() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
PileupJetIdAlgo * fPUJetIdAlgo
Definition: MVAJetIdMaker.h:23
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
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:155
const Point & position() const
position
Definition: Vertex.h:106
Jets made from PFObjects.
Definition: PFJet.h:21
float neutralEmEnergy() const
neutralEmEnergy
Definition: PFJet.h:150
virtual double energy() const
energy
int iEvent
Definition: GenABIO.cc:230
tuple corrector
Definition: mvaPFMET_cff.py:50
void Load_pfjets_JEC()
Definition: HWW.cc:1373
virtual float eta() const
momentum pseudorapidity
const float & beta() const
PileupJetIdentifier computeIdVariables(const reco::Jet *jet, float jec, const reco::Vertex *, const reco::VertexCollection &, bool calculateMva=false)
double ndof() const
Definition: Vertex.h:102
Definition: HWW.h:12
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
std::vector< float > & pfjets_mvavalue()
Definition: HWW.cc:785
void SetVars(HWW &, const edm::Event &, const edm::EventSetup &)
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:56
std::vector< PFJet > PFJetCollection
collection of PFJet objects
virtual float jetArea() const
get jet area
Definition: Jet.h:105
virtual int nConstituents() const
of constituents
Definition: Jet.h:65
std::vector< float > & pfjets_JEC()
Definition: HWW.cc:781
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: PFJet.h:102
MVAJetIdMaker(const edm::ParameterSet &, edm::ConsumesCollector)
Definition: MVAJetIdMaker.cc:7
const float & betaStar() const
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:22
const float & dRMean() const
const float & d0() const
const float & nCharged() const
bool passPFLooseId(const reco::PFJet *iJet)
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:98
const float & frac02() const
math::PtEtaPhiELorentzVectorF LorentzVector
const float & jetEta() const