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