CMS 3D CMS Logo

DPFIsolation.cc
Go to the documentation of this file.
1 /*
2  * \class DPFIsolation
3  *
4  * Deep ParticleFlow tau isolation using Deep NN.
5  *
6  * \author Owen Colegrove, UCSB
7  */
8 
10 
11 namespace {
12 inline int getPFCandidateIndex(const edm::Handle<pat::PackedCandidateCollection>& pfcands,
13  const reco::CandidatePtr& cptr)
14 {
15  for(unsigned int i = 0; i < pfcands->size(); ++i) {
16  if(reco::CandidatePtr(pfcands,i) == cptr)
17  return i;
18  }
19  return -1;
20 }
21 } // anonymous namespace
22 
24 public:
25  static const OutputCollection& GetOutputs()
26  {
27  const size_t tau_index = 0;
28  static const OutputCollection outputs_ = { { "VSall", Output({tau_index}, {}) } };
29  return outputs_;
30  };
31 
32  static unsigned getNumberOfParticles(unsigned graphVersion)
33  {
34  static const std::map<unsigned, unsigned> nparticles { { 0, 60 }, { 1, 36 } };
35  return nparticles.at(graphVersion);
36  }
37 
38  static unsigned GetNumberOfFeatures(unsigned graphVersion)
39  {
40  static const std::map<unsigned, unsigned> nfeatures { { 0, 47 }, { 1, 51 } };
41  return nfeatures.at(graphVersion);
42  }
43 
45  {
47  desc.add<edm::InputTag>("pfcands", edm::InputTag("packedPFCandidates"));
48  desc.add<edm::InputTag>("taus", edm::InputTag("slimmedTaus"));
49  desc.add<edm::InputTag>("vertices", edm::InputTag("offlineSlimmedPrimaryVertices"));
50  desc.add<std::string>("graph_file", "RecoTauTag/TrainingFiles/data/DPFTauId/DPFIsolation_2017v0_quantized.pb");
51  desc.add<unsigned>("version", 0);
52  desc.add<bool>("mem_mapped", false);
53 
55  descWP.add<std::string>("VVVLoose", "0");
56  descWP.add<std::string>("VVLoose", "0");
57  descWP.add<std::string>("VLoose", "0");
58  descWP.add<std::string>("Loose", "0");
59  descWP.add<std::string>("Medium", "0");
60  descWP.add<std::string>("Tight", "0");
61  descWP.add<std::string>("VTight", "0");
62  descWP.add<std::string>("VVTight", "0");
63  descWP.add<std::string>("VVVTight", "0");
64  desc.add<edm::ParameterSetDescription>("VSallWP", descWP);
65  descriptions.add("DPFTau2016v0", desc);
66  }
67 
69  DeepTauBase(cfg, GetOutputs(), cache),
70  pfcand_token(consumes<pat::PackedCandidateCollection>(cfg.getParameter<edm::InputTag>("pfcands"))),
71  vtx_token(consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("vertices"))),
72  graphVersion(cfg.getParameter<unsigned>("version"))
73  {
74  const auto& shape = cache_->getGraph().node(0).attr().at("shape").shape();
75 
76  if(!(graphVersion == 1 || graphVersion == 0 ))
77  throw cms::Exception("DPFIsolation") << "unknown version of the graph_ file.";
78 
79  if(!(shape.dim(1).size() == getNumberOfParticles(graphVersion) && shape.dim(2).size() == GetNumberOfFeatures(graphVersion)))
80  throw cms::Exception("DPFIsolation") << "number of inputs does not match the expected inputs for the given version";
81 
82  }
83 
84 private:
85  tensorflow::Tensor getPredictions(edm::Event& event, const edm::EventSetup& es,
87  {
89  event.getByToken(pfcand_token, pfcands);
90 
92  event.getByToken(vtx_token, vertices);
93 
94  tensorflow::Tensor tensor(tensorflow::DT_FLOAT, {1,
95  static_cast<int>(getNumberOfParticles(graphVersion)), static_cast<int>(GetNumberOfFeatures(graphVersion))});
96 
97  tensorflow::Tensor predictions(tensorflow::DT_FLOAT, { static_cast<int>(taus->size()), 1});
98 
99  std::vector<tensorflow::Tensor> outputs_;
100 
101  float pfCandPt, pfCandPz, pfCandPtRel, pfCandPzRel, pfCandDr, pfCandDEta, pfCandDPhi, pfCandEta, pfCandDz,
102  pfCandDzErr, pfCandD0, pfCandD0D0, pfCandD0Dz, pfCandD0Dphi, pfCandPuppiWeight,
103  pfCandPixHits, pfCandHits, pfCandLostInnerHits, pfCandPdgID, pfCandCharge, pfCandFromPV,
104  pfCandVtxQuality, pfCandHighPurityTrk, pfCandTauIndMatch, pfCandDzSig, pfCandD0Sig, pfCandD0Err,
105  pfCandPtRelPtRel, pfCandDzDz, pfCandDVx_1, pfCandDVy_1, pfCandDVz_1, pfCandD_1;
106  float pvx = !vertices->empty() ? (*vertices)[0].x() : -1;
107  float pvy = !vertices->empty() ? (*vertices)[0].y() : -1;
108  float pvz = !vertices->empty() ? (*vertices)[0].z() : -1;
109 
110  bool pfCandIsBarrel;
111 
112  // These variables define ranges further used for standardization
113  static constexpr float pfCandPt_max = 500.f;
114  static constexpr float pfCandPz_max = 1000.f;
115  static constexpr float pfCandPtRel_max = 1.f;
116  static constexpr float pfCandPzRel_max = 100.f;
117  static constexpr float pfCandPtRelPtRel_max = 1.f;
118  static constexpr float pfCandD0_max = 5.f;
119  static constexpr float pfCandDz_max = 5.f;
120  static constexpr float pfCandDVx_y_z_1_max = 0.05f;
121  static constexpr float pfCandD_1_max = 0.1f;
122  static constexpr float pfCandD0_z_Err_max = 1.f;
123  static constexpr float pfCandDzSig_max = 3.f;
124  static constexpr float pfCandD0Sig_max = 1.f;
125  static constexpr float pfCandDr_max = 0.5f;
126  static constexpr float pfCandEta_max = 2.75f;
127  static constexpr float pfCandDEta_max = 0.5f;
128  static constexpr float pfCandDPhi_max = 0.5f;
129  static constexpr float pfCandPixHits_max = 7.f;
130  static constexpr float pfCandHits_max = 30.f;
131 
132  for(size_t tau_index = 0; tau_index < taus->size(); tau_index++) {
133  pat::Tau tau = taus->at(tau_index);
134  bool isGoodTau = false;
135  const float lepRecoPt = tau.pt();
136  const float lepRecoPz = std::abs(tau.pz());
137  const float lepRecoEta = tau.eta();
138  const float lepRecoPhi = tau.phi();
139 
140  if(lepRecoPt >= 30 && std::abs(lepRecoEta) < 2.3 && tau.isTauIDAvailable("againstMuonLoose3") &&
141  tau.isTauIDAvailable("againstElectronVLooseMVA6")) {
142  isGoodTau = (tau.tauID("againstElectronVLooseMVA6") && tau.tauID("againstMuonLoose3") );
143  }
144 
145  if (!isGoodTau) {
146  predictions.matrix<float>()(tau_index, 0) = -1;
147  continue;
148  }
149 
150  std::vector<unsigned int> signalCandidateInds;
151 
152  for(const auto c : tau.signalCands())
153  signalCandidateInds.push_back(getPFCandidateIndex(pfcands,c));
154 
155  // Use of setZero results in warnings in eigen library during compilation.
156  //tensor.flat<float>().setZero();
158  for(unsigned input_idx = 0; input_idx < n_inputs; ++input_idx)
159  tensor.flat<float>()(input_idx) = 0;
160 
161  unsigned int iPF = 0;
162  const unsigned max_iPF = getNumberOfParticles(graphVersion);
163 
164  std::vector<unsigned int> sorted_inds(pfcands->size());
165  std::size_t n = 0;
166  std::generate(std::begin(sorted_inds), std::end(sorted_inds), [&]{ return n++; });
167 
168  std::sort(std::begin(sorted_inds), std::end(sorted_inds),
169  [&](int i1, int i2) { return pfcands->at(i1).pt() > pfcands->at(i2).pt(); } );
170 
171  for(size_t pf_index = 0; pf_index < pfcands->size() && iPF < max_iPF; pf_index++) {
172  pat::PackedCandidate p = pfcands->at(sorted_inds.at(pf_index));
173  float deltaR_tau_p = deltaR(p.p4(),tau.p4());
174 
175  if (p.pt() < 0.5) continue;
176  if (p.fromPV() < 0) continue;
177  if (deltaR_tau_p > 0.5) continue;
178  if (p.fromPV() < 1 && p.charge() != 0) continue;
179  pfCandPt = p.pt();
180  pfCandPtRel = p.pt()/lepRecoPt;
181 
182  pfCandDr = deltaR_tau_p;
183  pfCandDEta = std::abs(lepRecoEta - p.eta());
184  pfCandDPhi = std::abs(deltaPhi(lepRecoPhi, p.phi()));
185  pfCandEta = p.eta();
186  pfCandIsBarrel = (std::abs(pfCandEta) < 1.4);
187  pfCandPz = std::abs(std::sinh(pfCandEta)*pfCandPt);
188  pfCandPzRel = pfCandPz/lepRecoPz;
189  pfCandPdgID = std::abs(p.pdgId());
190  pfCandCharge = p.charge();
191  pfCandDVx_1 = p.vx() - pvx;
192  pfCandDVy_1 = p.vy() - pvy;
193  pfCandDVz_1 = p.vz() - pvz;
194 
195  pfCandD_1 = std::sqrt(pfCandDVx_1*pfCandDVx_1 + pfCandDVy_1*pfCandDVy_1 + pfCandDVz_1*pfCandDVz_1);
196 
197  if (pfCandCharge != 0 and p.hasTrackDetails()){
198  pfCandDz = p.dz();
199  pfCandDzErr = p.dzError();
200  pfCandDzSig = (std::abs(p.dz()) + 0.000001)/(p.dzError() + 0.00001);
201  pfCandD0 = p.dxy();
202  pfCandD0Err = p.dxyError();
203  pfCandD0Sig = (std::abs(p.dxy()) + 0.000001)/ (p.dxyError() + 0.00001);
204  pfCandPixHits = p.numberOfPixelHits();
205  pfCandHits = p.numberOfHits();
206  pfCandLostInnerHits = p.lostInnerHits();
207  } else {
208  float disp = 1;
209  int psudorand = p.pt()*1000000;
210  if (psudorand%2 == 0) disp = -1;
211  pfCandDz = 5*disp;
212  pfCandDzErr = 0;
213  pfCandDzSig = 0;
214  pfCandD0 = 5*disp;
215  pfCandD0Err = 0;
216  pfCandD0Sig = 0;
217  pfCandPixHits = 0;
218  pfCandHits = 0;
219  pfCandLostInnerHits = 2.;
220  pfCandDVx_1 = 1;
221  pfCandDVy_1 = 1;
222  pfCandDVz_1 = 1;
223  pfCandD_1 = 1;
224  }
225 
226  pfCandPuppiWeight = p.puppiWeight();
227  pfCandFromPV = p.fromPV();
228  pfCandVtxQuality = p.pvAssociationQuality();
229  pfCandHighPurityTrk = p.trackHighPurity();
230  float pfCandTauIndMatch_temp = 0;
231 
232  for (auto i : signalCandidateInds) {
233  if (i == sorted_inds.at(pf_index)) pfCandTauIndMatch_temp = 1;
234  }
235 
236  pfCandTauIndMatch = pfCandTauIndMatch_temp;
237  pfCandPtRelPtRel = pfCandPtRel*pfCandPtRel;
238  pfCandPt = std::min(pfCandPt, pfCandPt_max);
239  pfCandPt = pfCandPt/pfCandPt_max;
240 
241  pfCandPz = std::min(pfCandPz, pfCandPz_max);
242  pfCandPz = pfCandPz/pfCandPz_max;
243 
244  pfCandPtRel = std::min(pfCandPtRel, pfCandPtRel_max);
245  pfCandPzRel = std::min(pfCandPzRel, pfCandPzRel_max);
246  pfCandPzRel = pfCandPzRel/pfCandPzRel_max;
247  pfCandDr = pfCandDr/pfCandDr_max;
248  pfCandEta = pfCandEta/pfCandEta_max;
249  pfCandDEta = pfCandDEta/pfCandDEta_max;
250  pfCandDPhi = pfCandDPhi/pfCandDPhi_max;
251  pfCandPixHits = pfCandPixHits/pfCandPixHits_max;
252  pfCandHits = pfCandHits/pfCandHits_max;
253 
254  pfCandPtRelPtRel = std::min(pfCandPtRelPtRel, pfCandPtRelPtRel_max);
255 
256  pfCandD0 = std::clamp(pfCandD0, -pfCandD0_max, pfCandD0_max);
257  pfCandD0 = pfCandD0/pfCandD0_max;
258 
259  pfCandDz = std::clamp(pfCandDz, -pfCandDz_max, pfCandDz_max);
260  pfCandDz = pfCandDz/pfCandDz_max;
261 
262  pfCandD0Err = std::min(pfCandD0Err, pfCandD0_z_Err_max);
263  pfCandDzErr = std::min(pfCandDzErr, pfCandD0_z_Err_max);
264  pfCandDzSig = std::min(pfCandDzSig, pfCandDzSig_max);
265  pfCandDzSig = pfCandDzSig/pfCandDzSig_max;
266 
267  pfCandD0Sig = std::min(pfCandD0Sig, pfCandD0Sig_max);
268  pfCandD0D0 = pfCandD0*pfCandD0;
269  pfCandDzDz = pfCandDz*pfCandDz;
270  pfCandD0Dz = pfCandD0*pfCandDz;
271  pfCandD0Dphi = pfCandD0*pfCandDPhi;
272 
273  pfCandDVx_1 = std::clamp(pfCandDVx_1, -pfCandDVx_y_z_1_max, pfCandDVx_y_z_1_max);
274  pfCandDVx_1 = pfCandDVx_1/pfCandDVx_y_z_1_max;
275 
276  pfCandDVy_1 = std::clamp(pfCandDVy_1, -pfCandDVx_y_z_1_max, pfCandDVx_y_z_1_max);
277  pfCandDVy_1 = pfCandDVy_1/pfCandDVx_y_z_1_max;
278 
279  pfCandDVz_1 = std::clamp(pfCandDVz_1, -pfCandDVx_y_z_1_max, pfCandDVx_y_z_1_max);
280  pfCandDVz_1 = pfCandDVz_1/pfCandDVx_y_z_1_max;
281 
282  pfCandD_1 = std::clamp(pfCandD_1, -pfCandD_1_max, pfCandD_1_max);
283  pfCandD_1 = pfCandD_1/ pfCandD_1_max;
284 
285  if (graphVersion == 0) {
286  tensor.tensor<float,3>()( 0, 60-1-iPF, 0) = pfCandPt;
287  tensor.tensor<float,3>()( 0, 60-1-iPF, 1) = pfCandPz;
288  tensor.tensor<float,3>()( 0, 60-1-iPF, 2) = pfCandPtRel;
289  tensor.tensor<float,3>()( 0, 60-1-iPF, 3) = pfCandPzRel;
290  tensor.tensor<float,3>()( 0, 60-1-iPF, 4) = pfCandDr;
291  tensor.tensor<float,3>()( 0, 60-1-iPF, 5) = pfCandDEta;
292  tensor.tensor<float,3>()( 0, 60-1-iPF, 6) = pfCandDPhi;
293  tensor.tensor<float,3>()( 0, 60-1-iPF, 7) = pfCandEta;
294  tensor.tensor<float,3>()( 0, 60-1-iPF, 8) = pfCandDz;
295  tensor.tensor<float,3>()( 0, 60-1-iPF, 9) = pfCandDzSig;
296  tensor.tensor<float,3>()( 0, 60-1-iPF, 10) = pfCandD0;
297  tensor.tensor<float,3>()( 0, 60-1-iPF, 11) = pfCandD0Sig;
298  tensor.tensor<float,3>()( 0, 60-1-iPF, 12) = pfCandDzErr;
299  tensor.tensor<float,3>()( 0, 60-1-iPF, 13) = pfCandD0Err;
300  tensor.tensor<float,3>()( 0, 60-1-iPF, 14) = pfCandD0D0;
301  tensor.tensor<float,3>()( 0, 60-1-iPF, 15) = pfCandCharge==0;
302  tensor.tensor<float,3>()( 0, 60-1-iPF, 16) = pfCandCharge==1;
303  tensor.tensor<float,3>()( 0, 60-1-iPF, 17) = pfCandCharge==-1;
304  tensor.tensor<float,3>()( 0, 60-1-iPF, 18) = pfCandPdgID>22;
305  tensor.tensor<float,3>()( 0, 60-1-iPF, 19) = pfCandPdgID==22;
306  tensor.tensor<float,3>()( 0, 60-1-iPF, 20) = pfCandDzDz;
307  tensor.tensor<float,3>()( 0, 60-1-iPF, 21) = pfCandD0Dz;
308  tensor.tensor<float,3>()( 0, 60-1-iPF, 22) = pfCandD0Dphi;
309  tensor.tensor<float,3>()( 0, 60-1-iPF, 23) = pfCandPtRelPtRel;
310  tensor.tensor<float,3>()( 0, 60-1-iPF, 24) = pfCandPixHits;
311  tensor.tensor<float,3>()( 0, 60-1-iPF, 25) = pfCandHits;
312  tensor.tensor<float,3>()( 0, 60-1-iPF, 26) = pfCandLostInnerHits==-1;
313  tensor.tensor<float,3>()( 0, 60-1-iPF, 27) = pfCandLostInnerHits==0;
314  tensor.tensor<float,3>()( 0, 60-1-iPF, 28) = pfCandLostInnerHits==1;
315  tensor.tensor<float,3>()( 0, 60-1-iPF, 29) = pfCandLostInnerHits==2;
316  tensor.tensor<float,3>()( 0, 60-1-iPF, 30) = pfCandPuppiWeight;
317  tensor.tensor<float,3>()( 0, 60-1-iPF, 31) = (pfCandVtxQuality == 1);
318  tensor.tensor<float,3>()( 0, 60-1-iPF, 32) = (pfCandVtxQuality == 5);
319  tensor.tensor<float,3>()( 0, 60-1-iPF, 33) = (pfCandVtxQuality == 6);
320  tensor.tensor<float,3>()( 0, 60-1-iPF, 34) = (pfCandVtxQuality == 7);
321  tensor.tensor<float,3>()( 0, 60-1-iPF, 35) = (pfCandFromPV == 1);
322  tensor.tensor<float,3>()( 0, 60-1-iPF, 36) = (pfCandFromPV == 2);
323  tensor.tensor<float,3>()( 0, 60-1-iPF, 37) = (pfCandFromPV == 3);
324  tensor.tensor<float,3>()( 0, 60-1-iPF, 38) = pfCandIsBarrel;
325  tensor.tensor<float,3>()( 0, 60-1-iPF, 39) = pfCandHighPurityTrk;
326  tensor.tensor<float,3>()( 0, 60-1-iPF, 40) = pfCandPdgID==1;
327  tensor.tensor<float,3>()( 0, 60-1-iPF, 41) = pfCandPdgID==2;
328  tensor.tensor<float,3>()( 0, 60-1-iPF, 42) = pfCandPdgID==11;
329  tensor.tensor<float,3>()( 0, 60-1-iPF, 43) = pfCandPdgID==13;
330  tensor.tensor<float,3>()( 0, 60-1-iPF, 44) = pfCandPdgID==130;
331  tensor.tensor<float,3>()( 0, 60-1-iPF, 45) = pfCandPdgID==211;
332  tensor.tensor<float,3>()( 0, 60-1-iPF, 46) = pfCandTauIndMatch;
333  }
334 
335  if (graphVersion == 1) {
336  tensor.tensor<float,3>()( 0, 36-1-iPF, 0) = pfCandPt;
337  tensor.tensor<float,3>()( 0, 36-1-iPF, 1) = pfCandPz;
338  tensor.tensor<float,3>()( 0, 36-1-iPF, 2) = pfCandPtRel;
339  tensor.tensor<float,3>()( 0, 36-1-iPF, 3) = pfCandPzRel;
340  tensor.tensor<float,3>()( 0, 36-1-iPF, 4) = pfCandDr;
341  tensor.tensor<float,3>()( 0, 36-1-iPF, 5) = pfCandDEta;
342  tensor.tensor<float,3>()( 0, 36-1-iPF, 6) = pfCandDPhi;
343  tensor.tensor<float,3>()( 0, 36-1-iPF, 7) = pfCandEta;
344  tensor.tensor<float,3>()( 0, 36-1-iPF, 8) = pfCandDz;
345  tensor.tensor<float,3>()( 0, 36-1-iPF, 9) = pfCandDzSig;
346  tensor.tensor<float,3>()( 0, 36-1-iPF, 10) = pfCandD0;
347  tensor.tensor<float,3>()( 0, 36-1-iPF, 11) = pfCandD0Sig;
348  tensor.tensor<float,3>()( 0, 36-1-iPF, 12) = pfCandDzErr;
349  tensor.tensor<float,3>()( 0, 36-1-iPF, 13) = pfCandD0Err;
350  tensor.tensor<float,3>()( 0, 36-1-iPF, 14) = pfCandD0D0;
351  tensor.tensor<float,3>()( 0, 36-1-iPF, 15) = pfCandCharge==0;
352  tensor.tensor<float,3>()( 0, 36-1-iPF, 16) = pfCandCharge==1;
353  tensor.tensor<float,3>()( 0, 36-1-iPF, 17) = pfCandCharge==-1;
354  tensor.tensor<float,3>()( 0, 36-1-iPF, 18) = pfCandPdgID>22;
355  tensor.tensor<float,3>()( 0, 36-1-iPF, 19) = pfCandPdgID==22;
356  tensor.tensor<float,3>()( 0, 36-1-iPF, 20) = pfCandDVx_1;
357  tensor.tensor<float,3>()( 0, 36-1-iPF, 21) = pfCandDVy_1;
358  tensor.tensor<float,3>()( 0, 36-1-iPF, 22) = pfCandDVz_1;
359  tensor.tensor<float,3>()( 0, 36-1-iPF, 23) = pfCandD_1;
360  tensor.tensor<float,3>()( 0, 36-1-iPF, 24) = pfCandDzDz;
361  tensor.tensor<float,3>()( 0, 36-1-iPF, 25) = pfCandD0Dz;
362  tensor.tensor<float,3>()( 0, 36-1-iPF, 26) = pfCandD0Dphi;
363  tensor.tensor<float,3>()( 0, 36-1-iPF, 27) = pfCandPtRelPtRel;
364  tensor.tensor<float,3>()( 0, 36-1-iPF, 28) = pfCandPixHits;
365  tensor.tensor<float,3>()( 0, 36-1-iPF, 29) = pfCandHits;
366  tensor.tensor<float,3>()( 0, 36-1-iPF, 30) = pfCandLostInnerHits==-1;
367  tensor.tensor<float,3>()( 0, 36-1-iPF, 31) = pfCandLostInnerHits==0;
368  tensor.tensor<float,3>()( 0, 36-1-iPF, 32) = pfCandLostInnerHits==1;
369  tensor.tensor<float,3>()( 0, 36-1-iPF, 33) = pfCandLostInnerHits==2;
370  tensor.tensor<float,3>()( 0, 36-1-iPF, 34) = pfCandPuppiWeight;
371  tensor.tensor<float,3>()( 0, 36-1-iPF, 35) = (pfCandVtxQuality == 1);
372  tensor.tensor<float,3>()( 0, 36-1-iPF, 36) = (pfCandVtxQuality == 5);
373  tensor.tensor<float,3>()( 0, 36-1-iPF, 37) = (pfCandVtxQuality == 6);
374  tensor.tensor<float,3>()( 0, 36-1-iPF, 38) = (pfCandVtxQuality == 7);
375  tensor.tensor<float,3>()( 0, 36-1-iPF, 39) = (pfCandFromPV == 1);
376  tensor.tensor<float,3>()( 0, 36-1-iPF, 40) = (pfCandFromPV == 2);
377  tensor.tensor<float,3>()( 0, 36-1-iPF, 41) = (pfCandFromPV == 3);
378  tensor.tensor<float,3>()( 0, 36-1-iPF, 42) = pfCandIsBarrel;
379  tensor.tensor<float,3>()( 0, 36-1-iPF, 43) = pfCandHighPurityTrk;
380  tensor.tensor<float,3>()( 0, 36-1-iPF, 44) = pfCandPdgID==1;
381  tensor.tensor<float,3>()( 0, 36-1-iPF, 45) = pfCandPdgID==2;
382  tensor.tensor<float,3>()( 0, 36-1-iPF, 46) = pfCandPdgID==11;
383  tensor.tensor<float,3>()( 0, 36-1-iPF, 47) = pfCandPdgID==13;
384  tensor.tensor<float,3>()( 0, 36-1-iPF, 48) = pfCandPdgID==130;
385  tensor.tensor<float,3>()( 0, 36-1-iPF, 49) = pfCandPdgID==211;
386  tensor.tensor<float,3>()( 0, 36-1-iPF, 50) = pfCandTauIndMatch;
387  }
388  iPF++;
389  }
390  tensorflow::run(&(cache_->getSession()), { { "input_1", tensor } }, { "output_node0" }, {}, &outputs_);
391  predictions.matrix<float>()(tau_index, 0) = outputs_[0].flat<float>()(0);
392  }
393  return predictions;
394  }
395 
396 private:
399  unsigned graphVersion;
400 };
401 
float puppiWeight() const
Set both weights at once (with option for only full PUPPI)
virtual float dz(size_t ipv=0) const
dz with respect to the PV[ipv]
edm::EDGetTokenT< pat::PackedCandidateCollection > pfcand_token
float dxyError() const override
uncertainty on dxy
double eta() const final
momentum pseudorapidity
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: DPFIsolation.cc:44
int numberOfHits() const
int pdgId() const override
PDG identifier.
static unsigned getNumberOfParticles(unsigned graphVersion)
Definition: DPFIsolation.cc:32
const tensorflow::GraphDef & getGraph() const
Definition: DeepTauBase.h:50
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
tensorflow::Session & getSession() const
Definition: DeepTauBase.h:49
const DeepTauCache * cache_
Definition: DeepTauBase.h:101
std::vector< pat::PackedCandidate > PackedCandidateCollection
const LorentzVector & p4() const override
four-momentum Lorentz vecto r
float dzError() const override
uncertainty on dz
double pt() const final
transverse momentum
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
float tauID(const std::string &name) const
reco::CandidatePtrVector signalCands() const
double vz() const override
z coordinate of vertex position
std::map< std::string, Output > OutputCollection
Definition: DeepTauBase.h:82
#define constexpr
Definition: HeavyIon.h:7
def generate(map_blobs=False, class_name=None)
Definition: models.py:187
OutputCollection outputs_
Definition: DeepTauBase.h:100
int charge() const override
electric charge
double vy() const override
y coordinate of vertex position
static unsigned GetNumberOfFeatures(unsigned graphVersion)
Definition: DPFIsolation.cc:38
const PVAssociationQuality pvAssociationQuality() const
double pz() const final
z coordinate of momentum vector
const PVAssoc fromPV(size_t ipv=0) const
T sqrt(T t)
Definition: SSEVec.h:18
bool trackHighPurity() const
true if the track had the highPurity quality bit
LostInnerHits lostInnerHits() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
double pt() const override
transverse momentum
DeepTauBase(const edm::ParameterSet &cfg, const OutputCollection &outputs, const DeepTauCache *cache)
Definition: DeepTauBase.cc:79
#define end
Definition: vmac.h:39
T min(T a, T b)
Definition: MathUtil.h:58
ParameterDescriptionBase * add(U const &iLabel, T const &value)
int numberOfPixelHits() const
Analysis-level tau class.
Definition: Tau.h:55
tensorflow::Tensor getPredictions(edm::Event &event, const edm::EventSetup &es, edm::Handle< TauCollection > taus) override
Definition: DPFIsolation.cc:85
def cache(function)
bool hasTrackDetails() const
Return true if a bestTrack can be extracted from this Candidate.
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
bool isTauIDAvailable(const std::string &name) const
Returns true if a specific ID is available in this pat::Tau.
double eta() const override
momentum pseudorapidity
double phi() const override
momentum azimuthal angle
edm::EDGetTokenT< reco::VertexCollection > vtx_token
DPFIsolation(const edm::ParameterSet &cfg, const deep_tau::DeepTauCache *cache)
Definition: DPFIsolation.cc:68
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
#define begin
Definition: vmac.h:32
HLT enums.
#define Output(cl)
Definition: vmac.h:194
unsigned graphVersion
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, const std::vector< std::string > &targetNodes, std::vector< Tensor > *outputs)
Definition: TensorFlow.cc:210
static const OutputCollection & GetOutputs()
Definition: DPFIsolation.cc:25
virtual float dxy() const
dxy with respect to the PV ref
double phi() const final
momentum azimuthal angle
double vx() const override
x coordinate of vertex position
Definition: event.py:1