CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
LeptonJetVarProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PhysicsTools/NanoAOD
4 // Class: LeptonJetVarProducer
5 //
13 //
14 // Original Author: Marco Peruzzi
15 // Created: Tue, 05 Sep 2017 12:24:38 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 
22 // user include files
25 
28 
31 
36 
39 
41 
42 //
43 // class declaration
44 //
45 
46 template <typename T>
48 public:
50  : srcJet_(consumes<edm::View<pat::Jet>>(iConfig.getParameter<edm::InputTag>("srcJet"))),
51  srcLep_(consumes<edm::View<T>>(iConfig.getParameter<edm::InputTag>("srcLep"))),
52  srcVtx_(consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("srcVtx"))) {
53  produces<edm::ValueMap<float>>("ptRatio");
54  produces<edm::ValueMap<float>>("ptRel");
55  produces<edm::ValueMap<float>>("jetNDauChargedMVASel");
56  produces<edm::ValueMap<reco::CandidatePtr>>("jetForLepJetVar");
57  }
58  ~LeptonJetVarProducer() override{};
59 
60  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
61 
62 private:
63  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
64 
65  std::tuple<float, float, float> calculatePtRatioRel(edm::Ptr<reco::Candidate> lep,
67  const reco::Vertex& vtx) const;
68 
69  // ----------member data ---------------------------
70 
74 };
75 
76 //
77 // constants, enums and typedefs
78 //
79 
80 //
81 // static data member definitions
82 //
83 
84 //
85 // member functions
86 //
87 
88 // ------------ method called to produce the data ------------
89 template <typename T>
92  iEvent.getByToken(srcJet_, srcJet);
94  iEvent.getByToken(srcLep_, srcLep);
96  iEvent.getByToken(srcVtx_, srcVtx);
97 
98  unsigned int nJet = srcJet->size();
99  unsigned int nLep = srcLep->size();
100 
101  std::vector<float> ptRatio(nLep, -1);
102  std::vector<float> ptRel(nLep, -1);
103  std::vector<float> jetNDauChargedMVASel(nLep, 0);
104  std::vector<reco::CandidatePtr> jetForLepJetVar(nLep, reco::CandidatePtr());
105 
106  const auto& pv = (*srcVtx)[0];
107 
108  for (unsigned int il = 0; il < nLep; il++) {
109  for (unsigned int ij = 0; ij < nJet; ij++) {
110  auto lep = srcLep->ptrAt(il);
111  auto jet = srcJet->ptrAt(ij);
112  if (matchByCommonSourceCandidatePtr(*lep, *jet)) {
113  auto res = calculatePtRatioRel(lep, jet, pv);
114  ptRatio[il] = std::get<0>(res);
115  ptRel[il] = std::get<1>(res);
116  jetNDauChargedMVASel[il] = std::get<2>(res);
117  jetForLepJetVar[il] = jet;
118  break; // take leading jet with shared source candidates
119  }
120  }
121  }
122 
123  std::unique_ptr<edm::ValueMap<float>> ptRatioV(new edm::ValueMap<float>());
124  edm::ValueMap<float>::Filler fillerRatio(*ptRatioV);
125  fillerRatio.insert(srcLep, ptRatio.begin(), ptRatio.end());
126  fillerRatio.fill();
127  iEvent.put(std::move(ptRatioV), "ptRatio");
128 
129  std::unique_ptr<edm::ValueMap<float>> ptRelV(new edm::ValueMap<float>());
130  edm::ValueMap<float>::Filler fillerRel(*ptRelV);
131  fillerRel.insert(srcLep, ptRel.begin(), ptRel.end());
132  fillerRel.fill();
133  iEvent.put(std::move(ptRelV), "ptRel");
134 
135  std::unique_ptr<edm::ValueMap<float>> jetNDauChargedMVASelV(new edm::ValueMap<float>());
136  edm::ValueMap<float>::Filler fillerNDau(*jetNDauChargedMVASelV);
137  fillerNDau.insert(srcLep, jetNDauChargedMVASel.begin(), jetNDauChargedMVASel.end());
138  fillerNDau.fill();
139  iEvent.put(std::move(jetNDauChargedMVASelV), "jetNDauChargedMVASel");
140 
141  std::unique_ptr<edm::ValueMap<reco::CandidatePtr>> jetForLepJetVarV(new edm::ValueMap<reco::CandidatePtr>());
142  edm::ValueMap<reco::CandidatePtr>::Filler fillerjetForLepJetVar(*jetForLepJetVarV);
143  fillerjetForLepJetVar.insert(srcLep, jetForLepJetVar.begin(), jetForLepJetVar.end());
144  fillerjetForLepJetVar.fill();
145  iEvent.put(std::move(jetForLepJetVarV), "jetForLepJetVar");
146 }
147 
148 template <typename T>
151  const reco::Vertex& vtx) const {
152  auto rawp4 = jet->correctedP4("Uncorrected");
153  auto lepp4 = lep->p4();
154 
155  if ((rawp4 - lepp4).R() < 1e-4)
156  return std::tuple<float, float, float>(1.0, 0.0, 0.0);
157 
158  auto l1corrFactor = jet->jecFactor("L1FastJet") / jet->jecFactor("Uncorrected");
159 
160  auto jetp4 = (rawp4 - lepp4 * (1.0 / l1corrFactor)) * (jet->pt() / rawp4.pt()) + lepp4;
161  auto ptratio = lepp4.pt() / jetp4.pt();
162  auto ptrel = lepp4.Vect().Cross((jetp4 - lepp4).Vect().Unit()).R();
163 
164  unsigned int jndau = 0;
165  for (const auto& _d : jet->daughterPtrVector()) {
166  const auto d = dynamic_cast<const pat::PackedCandidate*>(_d.get());
167  if (d->charge() == 0)
168  continue;
169  if (d->fromPV() <= 1)
170  continue;
171  if (deltaR(*d, *lep) > 0.4)
172  continue;
173  if (!(d->hasTrackDetails()))
174  continue;
175  auto tk = d->pseudoTrack();
176  if (tk.pt() > 1 && tk.hitPattern().numberOfValidHits() >= 8 && tk.hitPattern().numberOfValidPixelHits() >= 2 &&
177  tk.normalizedChi2() < 5 && fabs(tk.dxy(vtx.position())) < 0.2 && fabs(tk.dz(vtx.position())) < 17)
178  jndau++;
179  }
180 
181  return std::tuple<float, float, float>(ptratio, ptrel, float(jndau));
182 }
183 
184 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
185 template <typename T>
187  //The following says we do not know what parameters are allowed so do no validation
188  // Please change this to state exactly what you do use, even if it is no parameters
190  desc.add<edm::InputTag>("srcJet")->setComment("jet input collection");
191  desc.add<edm::InputTag>("srcLep")->setComment("lepton input collection");
192  desc.add<edm::InputTag>("srcVtx")->setComment("primary vertex input collection");
194  if (typeid(T) == typeid(pat::Muon))
195  modname += "Muon";
196  else if (typeid(T) == typeid(pat::Electron))
197  modname += "Electron";
198  modname += "JetVarProducer";
199  descriptions.add(modname, desc);
200 }
201 
204 
205 //define this as a plug-in
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::tuple< float, float, float > calculatePtRatioRel(edm::Ptr< reco::Candidate > lep, edm::Ptr< pat::Jet > jet, const reco::Vertex &vtx) const
const Point & position() const
position
Definition: Vertex.h:127
LeptonJetVarProducer< pat::Electron > ElectronJetVarProducer
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< std::vector< reco::Vertex > > srcVtx_
tuple d
Definition: ztail.py:151
int iEvent
Definition: GenABIO.cc:224
bool matchByCommonSourceCandidatePtr(const C1 &c1, const C2 &c2)
Definition: MatchingUtils.h:9
def move
Definition: eostools.py:511
ParameterDescriptionBase * add(U const &iLabel, T const &value)
LeptonJetVarProducer(const edm::ParameterSet &iConfig)
Analysis-level electron class.
Definition: Electron.h:51
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< edm::View< pat::Jet > > srcJet_
constexpr char Jet[]
Definition: modules.cc:9
long double T
Analysis-level muon class.
Definition: Muon.h:51
LeptonJetVarProducer< pat::Muon > MuonJetVarProducer
edm::EDGetTokenT< edm::View< T > > srcLep_