CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CandidateBoostedDoubleSecondaryVertexComputer.cc
Go to the documentation of this file.
2 
4 
16 
17 #include "fastjet/contrib/Njettiness.hh"
18 
20  beta_(parameters.getParameter<double>("beta")),
21  R0_(parameters.getParameter<double>("R0")),
22  maxSVDeltaRToJet_(parameters.getParameter<double>("maxSVDeltaRToJet")),
23  useCondDB_(parameters.getParameter<bool>("useCondDB")),
24  gbrForestLabel_(parameters.existsAs<std::string>("gbrForestLabel") ? parameters.getParameter<std::string>("gbrForestLabel") : ""),
25  weightFile_(parameters.existsAs<edm::FileInPath>("weightFile") ? parameters.getParameter<edm::FileInPath>("weightFile") : edm::FileInPath()),
26  useGBRForest_(parameters.existsAs<bool>("useGBRForest") ? parameters.getParameter<bool>("useGBRForest") : false),
27  useAdaBoost_(parameters.existsAs<bool>("useAdaBoost") ? parameters.getParameter<bool>("useAdaBoost") : false),
28  maxDistToAxis_(parameters.getParameter<edm::ParameterSet>("trackSelection").getParameter<double>("maxDistToAxis")),
29  maxDecayLen_(parameters.getParameter<edm::ParameterSet>("trackSelection").getParameter<double>("maxDecayLen")),
30  trackPairV0Filter(parameters.getParameter<edm::ParameterSet>("trackPairV0Filter")),
31  trackSelector(parameters.getParameter<edm::ParameterSet>("trackSelection"))
32 {
33  uses(0, "ipTagInfos");
34  uses(1, "svTagInfos");
35 
36  mvaID.reset(new TMVAEvaluator());
37 }
38 
40 {
41  // variable names and order need to be the same as in the training
42  std::vector<std::string> variables({"z_ratio",
43  "trackSipdSig_3","trackSipdSig_2","trackSipdSig_1","trackSipdSig_0",
44  "trackSipdSig_1_0","trackSipdSig_0_0","trackSipdSig_1_1","trackSipdSig_0_1",
45  "trackSip2dSigAboveCharm_0","trackSip2dSigAboveBottom_0","trackSip2dSigAboveBottom_1",
46  "tau0_trackEtaRel_0","tau0_trackEtaRel_1","tau0_trackEtaRel_2",
47  "tau1_trackEtaRel_0","tau1_trackEtaRel_1","tau1_trackEtaRel_2",
48  "tau_vertexMass_0","tau_vertexEnergyRatio_0","tau_vertexDeltaR_0","tau_flightDistance2dSig_0",
49  "tau_vertexMass_1","tau_vertexEnergyRatio_1","tau_flightDistance2dSig_1",
50  "jetNTracks","nSV"});
51  // book TMVA readers
52  std::vector<std::string> spectators({"massPruned", "flavour", "nbHadrons", "ptPruned", "etaPruned"});
53 
54  if (useCondDB_)
55  {
56  const GBRWrapperRcd & gbrWrapperRecord = record.getRecord<GBRWrapperRcd>();
57 
58  edm::ESHandle<GBRForest> gbrForestHandle;
59  gbrWrapperRecord.get(gbrForestLabel_.c_str(), gbrForestHandle);
60 
61  mvaID->initializeGBRForest(gbrForestHandle.product(), variables, spectators, useAdaBoost_);
62  }
63  else
64  mvaID->initialize("Color:Silent:Error", "BDT", weightFile_.fullPath(), variables, spectators, useGBRForest_, useAdaBoost_);
65 
66  // get TransientTrackBuilder
67  const TransientTrackRecord & transientTrackRcd = record.getRecord<TransientTrackRecord>();
68  transientTrackRcd.get("TransientTrackBuilder", trackBuilder);
69 }
70 
72 {
73  // get TagInfos
74  const reco::CandIPTagInfo & ipTagInfo = tagInfo.get<reco::CandIPTagInfo>(0);
76 
77  // default discriminator value
78  float value = -10.;
79 
80  // default variable values
81  float z_ratio = dummyZ_ratio;
82  float trackSip3dSig_3 = dummyTrackSip3dSig, trackSip3dSig_2 = dummyTrackSip3dSig, trackSip3dSig_1 = dummyTrackSip3dSig, trackSip3dSig_0 = dummyTrackSip3dSig;
83  float tau2_trackSip3dSig_0 = dummyTrackSip3dSig, tau1_trackSip3dSig_0 = dummyTrackSip3dSig, tau2_trackSip3dSig_1 = dummyTrackSip3dSig, tau1_trackSip3dSig_1 = dummyTrackSip3dSig;
84  float trackSip2dSigAboveCharm_0 = dummyTrackSip2dSigAbove, trackSip2dSigAboveBottom_0 = dummyTrackSip2dSigAbove, trackSip2dSigAboveBottom_1 = dummyTrackSip2dSigAbove;
85  float tau1_trackEtaRel_0 = dummyTrackEtaRel, tau1_trackEtaRel_1 = dummyTrackEtaRel, tau1_trackEtaRel_2 = dummyTrackEtaRel;
86  float tau2_trackEtaRel_0 = dummyTrackEtaRel, tau2_trackEtaRel_1 = dummyTrackEtaRel, tau2_trackEtaRel_2 = dummyTrackEtaRel;
87  float tau1_vertexMass = dummyVertexMass, tau1_vertexEnergyRatio = dummyVertexEnergyRatio, tau1_vertexDeltaR = dummyVertexDeltaR, tau1_flightDistance2dSig = dummyFlightDistance2dSig;
88  float tau2_vertexMass = dummyVertexMass, tau2_vertexEnergyRatio = dummyVertexEnergyRatio, tau2_vertexDeltaR = dummyVertexDeltaR, tau2_flightDistance2dSig = dummyFlightDistance2dSig;
89  float jetNTracks = 0, nSV = 0, tau1_nSecondaryVertices = 0, tau2_nSecondaryVertices = 0;
90 
91  // get the jet reference
92  const reco::JetBaseRef jet = svTagInfo.jet();
93 
94  std::vector<fastjet::PseudoJet> currentAxes;
95  float tau2, tau1;
96  // calculate N-subjettiness
97  calcNsubjettiness(jet, tau1, tau2, currentAxes);
98 
99  const reco::VertexRef & vertexRef = ipTagInfo.primaryVertex();
100  GlobalPoint pv(0.,0.,0.);
101  if ( ipTagInfo.primaryVertex().isNonnull() )
102  pv = GlobalPoint(vertexRef->x(),vertexRef->y(),vertexRef->z());
103 
104  const std::vector<reco::CandidatePtr> & selectedTracks = ipTagInfo.selectedTracks();
105  const std::vector<reco::btag::TrackIPData> & ipData = ipTagInfo.impactParameterData();
106  size_t trackSize = selectedTracks.size();
107 
108 
109  reco::TrackKinematics allKinematics;
110  std::vector<float> IP3Ds, IP3Ds_1, IP3Ds_2;
111  int contTrk=0;
112 
113  // loop over tracks associated to the jet
114  for (size_t itt=0; itt < trackSize; ++itt)
115  {
116  const reco::CandidatePtr trackRef = selectedTracks[itt];
117 
118  float track_PVweight = 0.;
119  setTracksPV(trackRef, vertexRef, track_PVweight);
120  if (track_PVweight>0.5) allKinematics.add(trackRef);
121 
122  const reco::btag::TrackIPData &data = ipData[itt];
123  bool isSelected = false;
124  if (trackSelector(trackRef, data, *jet, pv)) isSelected = true;
125 
126  // check if the track is from V0
127  bool isfromV0 = false, isfromV0Tight = false;
128  std::vector<reco::CandidatePtr> trackPairV0Test(2);
129 
130  trackPairV0Test[0] = trackRef;
131 
132  for (size_t jtt=0; jtt < trackSize; ++jtt)
133  {
134  if (itt == jtt) continue;
135 
136  const reco::btag::TrackIPData & pairTrackData = ipData[jtt];
137  const reco::CandidatePtr pairTrackRef = selectedTracks[jtt];
138 
139  trackPairV0Test[1] = pairTrackRef;
140 
141  if (!trackPairV0Filter(trackPairV0Test))
142  {
143  isfromV0 = true;
144 
145  if ( trackSelector(pairTrackRef, pairTrackData, *jet, pv) )
146  isfromV0Tight = true;
147  }
148 
149  if (isfromV0 && isfromV0Tight)
150  break;
151  }
152 
153  if( isSelected && !isfromV0Tight ) jetNTracks += 1.;
154 
155  reco::TransientTrack transientTrack = trackBuilder->build(trackRef);
156  GlobalVector direction(jet->px(), jet->py(), jet->pz());
157 
158  int index = 0;
159  if (currentAxes.size() > 1 && reco::deltaR2(trackRef->momentum(),currentAxes[1]) < reco::deltaR2(trackRef->momentum(),currentAxes[0]))
160  index = 1;
161  direction = GlobalVector(currentAxes[index].px(), currentAxes[index].py(), currentAxes[index].pz());
162 
163  // decay distance and track distance wrt to the closest tau axis
164  float decayLengthTau=-1;
165  float distTauAxis=-1;
166 
167  TrajectoryStateOnSurface closest = IPTools::closestApproachToJet(transientTrack.impactPointState(), *vertexRef , direction, transientTrack.field());
168  if (closest.isValid())
169  decayLengthTau = (closest.globalPosition() - RecoVertex::convertPos(vertexRef->position())).mag();
170 
171  distTauAxis = std::abs(IPTools::jetTrackDistance(transientTrack, direction, *vertexRef ).second.value());
172 
173  float IP3Dsig = ipTagInfo.impactParameterData()[itt].ip3d.significance();
174 
175  if( !isfromV0 && decayLengthTau<maxDecayLen_ && distTauAxis<maxDistToAxis_ )
176  {
177  IP3Ds.push_back( IP3Dsig<-50. ? -50. : IP3Dsig );
178  ++contTrk;
179  if (currentAxes.size() > 1)
180  {
181  if (reco::deltaR2(trackRef->momentum(),currentAxes[0]) < reco::deltaR2(trackRef->momentum(),currentAxes[1]))
182  IP3Ds_1.push_back( IP3Dsig<-50. ? -50. : IP3Dsig );
183  else
184  IP3Ds_2.push_back( IP3Dsig<-50. ? -50. : IP3Dsig );
185  }
186  else
187  IP3Ds_1.push_back( IP3Dsig<-50. ? -50. : IP3Dsig );
188  }
189  }
190 
191  std::vector<size_t> indices = ipTagInfo.sortedIndexes(reco::btag::IP2DSig);
192  bool charmThreshSet = false;
193 
195  for (size_t i=0; i<indices.size(); ++i)
196  {
197  size_t idx = indices[i];
198  const reco::btag::TrackIPData & data = ipData[idx];
199  const reco::CandidatePtr trackRef = selectedTracks[idx];
200 
201  kin.add(trackRef);
202 
203  if ( kin.vectorSum().M() > charmThreshold // charm cut
204  && !charmThreshSet )
205  {
206  trackSip2dSigAboveCharm_0 = data.ip2d.significance();
207 
208  charmThreshSet = true;
209  }
210 
211  if ( kin.vectorSum().M() > bottomThreshold ) // bottom cut
212  {
213  trackSip2dSigAboveBottom_0 = data.ip2d.significance();
214  if ( (i+1)<indices.size() ) trackSip2dSigAboveBottom_1 = (ipData[indices[i+1]]).ip2d.significance();
215 
216  break;
217  }
218  }
219 
220  float dummyTrack = -50.;
221 
222  std::sort( IP3Ds.begin(),IP3Ds.end(),std::greater<float>() );
223  std::sort( IP3Ds_1.begin(),IP3Ds_1.end(),std::greater<float>() );
224  std::sort( IP3Ds_2.begin(),IP3Ds_2.end(),std::greater<float>() );
225  int num_1 = IP3Ds_1.size();
226  int num_2 = IP3Ds_2.size();
227 
228  switch(contTrk){
229  case 0:
230 
231  trackSip3dSig_0 = dummyTrack;
232  trackSip3dSig_1 = dummyTrack;
233  trackSip3dSig_2 = dummyTrack;
234  trackSip3dSig_3 = dummyTrack;
235 
236  break;
237 
238  case 1:
239 
240  trackSip3dSig_0 = IP3Ds.at(0);
241  trackSip3dSig_1 = dummyTrack;
242  trackSip3dSig_2 = dummyTrack;
243  trackSip3dSig_3 = dummyTrack;
244 
245  break;
246 
247  case 2:
248 
249  trackSip3dSig_0 = IP3Ds.at(0);
250  trackSip3dSig_1 = IP3Ds.at(1);
251  trackSip3dSig_2 = dummyTrack;
252  trackSip3dSig_3 = dummyTrack;
253 
254  break;
255 
256  case 3:
257 
258  trackSip3dSig_0 = IP3Ds.at(0);
259  trackSip3dSig_1 = IP3Ds.at(1);
260  trackSip3dSig_2 = IP3Ds.at(2);
261  trackSip3dSig_3 = dummyTrack;
262 
263  break;
264 
265  default:
266 
267  trackSip3dSig_0 = IP3Ds.at(0);
268  trackSip3dSig_1 = IP3Ds.at(1);
269  trackSip3dSig_2 = IP3Ds.at(2);
270  trackSip3dSig_3 = IP3Ds.at(3);
271 
272  }
273 
274  switch(num_1){
275  case 0:
276 
277  tau1_trackSip3dSig_0 = dummyTrack;
278  tau1_trackSip3dSig_1 = dummyTrack;
279 
280  break;
281 
282  case 1:
283 
284  tau1_trackSip3dSig_0 = IP3Ds_1.at(0);
285  tau1_trackSip3dSig_1 = dummyTrack;
286 
287  break;
288 
289  default:
290 
291  tau1_trackSip3dSig_0 = IP3Ds_1.at(0);
292  tau1_trackSip3dSig_1 = IP3Ds_1.at(1);
293 
294  }
295 
296  switch(num_2){
297  case 0:
298 
299  tau2_trackSip3dSig_0 = dummyTrack;
300  tau2_trackSip3dSig_1 = dummyTrack;
301 
302  break;
303 
304  case 1:
305  tau2_trackSip3dSig_0 = IP3Ds_2.at(0);
306  tau2_trackSip3dSig_1 = dummyTrack;
307 
308  break;
309 
310  default:
311 
312  tau2_trackSip3dSig_0 = IP3Ds_2.at(0);
313  tau2_trackSip3dSig_1 = IP3Ds_2.at(1);
314 
315  }
316 
317  math::XYZVector jetDir = jet->momentum().Unit();
318  reco::TrackKinematics tau1Kinematics;
319  reco::TrackKinematics tau2Kinematics;
320  std::vector<float> tau1_trackEtaRels, tau2_trackEtaRels;
321 
322  std::map<double, size_t> VTXmap;
323  for (size_t vtx = 0; vtx < svTagInfo.nVertices(); ++vtx)
324  {
325  const reco::VertexCompositePtrCandidate vertex = svTagInfo.secondaryVertex(vtx);
326  // get the vertex kinematics
327  reco::TrackKinematics vertexKinematic(vertex);
328 
329  if (currentAxes.size() > 1)
330  {
331  if (reco::deltaR2(svTagInfo.flightDirection(vtx),currentAxes[1]) < reco::deltaR2(svTagInfo.flightDirection(vtx),currentAxes[0]))
332  {
333  tau2Kinematics = tau2Kinematics + vertexKinematic;
334  if( tau2_flightDistance2dSig < 0 )
335  {
336  tau2_flightDistance2dSig = svTagInfo.flightDistance(vtx,true).significance();
337  tau2_vertexDeltaR = reco::deltaR(svTagInfo.flightDirection(vtx),currentAxes[1]);
338  }
339  etaRelToTauAxis(vertex, currentAxes[1], tau2_trackEtaRels);
340  tau2_nSecondaryVertices += 1.;
341  }
342  else
343  {
344  tau1Kinematics = tau1Kinematics + vertexKinematic;
345  if( tau1_flightDistance2dSig < 0 )
346  {
347  tau1_flightDistance2dSig =svTagInfo.flightDistance(vtx,true).significance();
348  tau1_vertexDeltaR = reco::deltaR(svTagInfo.flightDirection(vtx),currentAxes[0]);
349  }
350  etaRelToTauAxis(vertex, currentAxes[0], tau1_trackEtaRels);
351  tau1_nSecondaryVertices += 1.;
352  }
353 
354  }
355  else if (currentAxes.size() > 0)
356  {
357  tau1Kinematics = tau1Kinematics + vertexKinematic;
358  if( tau1_flightDistance2dSig < 0 )
359  {
360  tau1_flightDistance2dSig =svTagInfo.flightDistance(vtx,true).significance();
361  tau1_vertexDeltaR = reco::deltaR(svTagInfo.flightDirection(vtx),currentAxes[0]);
362  }
363  etaRelToTauAxis(vertex, currentAxes[0], tau1_trackEtaRels);
364  tau1_nSecondaryVertices += 1.;
365  }
366 
367  GlobalVector flightDir = svTagInfo.flightDirection(vtx);
368  if (reco::deltaR2(flightDir, jetDir)<(maxSVDeltaRToJet_*maxSVDeltaRToJet_))
369  VTXmap[svTagInfo.flightDistance(vtx).error()]=vtx;
370  }
371  nSV = VTXmap.size();
372 
373 
374  math::XYZTLorentzVector allSum = allKinematics.weightedVectorSum() ;
375  if ( tau1_nSecondaryVertices > 0. )
376  {
377  math::XYZTLorentzVector tau1_vertexSum = tau1Kinematics.weightedVectorSum();
378  tau1_vertexEnergyRatio = tau1_vertexSum.E() / allSum.E();
379  if ( tau1_vertexEnergyRatio > 50. ) tau1_vertexEnergyRatio = 50.;
380 
381  tau1_vertexMass = tau1_vertexSum.M();
382  }
383 
384  if ( tau2_nSecondaryVertices > 0. )
385  {
386  math::XYZTLorentzVector tau2_vertexSum = tau2Kinematics.weightedVectorSum();
387  tau2_vertexEnergyRatio = tau2_vertexSum.E() / allSum.E();
388  if ( tau2_vertexEnergyRatio > 50. ) tau2_vertexEnergyRatio = 50.;
389 
390  tau2_vertexMass= tau2_vertexSum.M();
391  }
392 
393 
394  float dummyEtaRel = -1.;
395 
396  std::sort( tau1_trackEtaRels.begin(),tau1_trackEtaRels.end() );
397  std::sort( tau2_trackEtaRels.begin(),tau2_trackEtaRels.end() );
398 
399  switch(tau2_trackEtaRels.size()){
400  case 0:
401 
402  tau2_trackEtaRel_0 = dummyEtaRel;
403  tau2_trackEtaRel_1 = dummyEtaRel;
404  tau2_trackEtaRel_2 = dummyEtaRel;
405 
406  break;
407 
408  case 1:
409 
410  tau2_trackEtaRel_0 = tau2_trackEtaRels.at(0);
411  tau2_trackEtaRel_1 = dummyEtaRel;
412  tau2_trackEtaRel_2 = dummyEtaRel;
413 
414  break;
415 
416  case 2:
417 
418  tau2_trackEtaRel_0 = tau2_trackEtaRels.at(0);
419  tau2_trackEtaRel_1 = tau2_trackEtaRels.at(1);
420  tau2_trackEtaRel_2 = dummyEtaRel;
421 
422  break;
423 
424  default:
425 
426  tau2_trackEtaRel_0 = tau2_trackEtaRels.at(0);
427  tau2_trackEtaRel_1 = tau2_trackEtaRels.at(1);
428  tau2_trackEtaRel_2 = tau2_trackEtaRels.at(2);
429 
430  }
431 
432  switch(tau1_trackEtaRels.size()){
433  case 0:
434 
435  tau1_trackEtaRel_0 = dummyEtaRel;
436  tau1_trackEtaRel_1 = dummyEtaRel;
437  tau1_trackEtaRel_2 = dummyEtaRel;
438 
439  break;
440 
441  case 1:
442 
443  tau1_trackEtaRel_0 = tau1_trackEtaRels.at(0);
444  tau1_trackEtaRel_1 = dummyEtaRel;
445  tau1_trackEtaRel_2 = dummyEtaRel;
446 
447  break;
448 
449  case 2:
450 
451  tau1_trackEtaRel_0 = tau1_trackEtaRels.at(0);
452  tau1_trackEtaRel_1 = tau1_trackEtaRels.at(1);
453  tau1_trackEtaRel_2 = dummyEtaRel;
454 
455  break;
456 
457  default:
458 
459  tau1_trackEtaRel_0 = tau1_trackEtaRels.at(0);
460  tau1_trackEtaRel_1 = tau1_trackEtaRels.at(1);
461  tau1_trackEtaRel_2 = tau1_trackEtaRels.at(2);
462 
463  }
464 
465  int cont=0;
466  GlobalVector flightDir_0, flightDir_1;
467  reco::Candidate::LorentzVector SV_p4_0 , SV_p4_1;
468  double vtxMass = 0.;
469 
470  for ( std::map<double, size_t>::iterator iVtx=VTXmap.begin(); iVtx!=VTXmap.end(); ++iVtx)
471  {
472  ++cont;
473  const reco::VertexCompositePtrCandidate &vertex = svTagInfo.secondaryVertex(iVtx->second);
474  if (cont==1)
475  {
476  flightDir_0 = svTagInfo.flightDirection(iVtx->second);
477  SV_p4_0 = vertex.p4();
478  vtxMass = SV_p4_0.mass();
479 
480  if(vtxMass > 0.)
481  z_ratio = reco::deltaR(currentAxes[1],currentAxes[0])*SV_p4_0.pt()/vtxMass;
482  }
483  if (cont==2)
484  {
485  flightDir_1 = svTagInfo.flightDirection(iVtx->second);
486  SV_p4_1 = vertex.p4();
487  vtxMass = (SV_p4_1+SV_p4_0).mass();
488 
489  if(vtxMass > 0.)
490  z_ratio = reco::deltaR(flightDir_0,flightDir_1)*SV_p4_1.pt()/vtxMass;
491 
492  break;
493  }
494  }
495 
496  // when only one tau axis has SVs assigned, they are all assigned to the 1st tau axis
497  // in the special case below need to swap values
498  if( (tau1_vertexMass<0 && tau2_vertexMass>0) )
499  {
500  float temp = tau1_trackEtaRel_0;
501  tau1_trackEtaRel_0= tau2_trackEtaRel_0;
502  tau2_trackEtaRel_0= temp;
503 
504  temp = tau1_trackEtaRel_1;
505  tau1_trackEtaRel_1= tau2_trackEtaRel_1;
506  tau2_trackEtaRel_1= temp;
507 
508  temp = tau1_trackEtaRel_2;
509  tau1_trackEtaRel_2= tau2_trackEtaRel_2;
510  tau2_trackEtaRel_2= temp;
511 
512  temp = tau1_flightDistance2dSig;
513  tau1_flightDistance2dSig= tau2_flightDistance2dSig;
514  tau2_flightDistance2dSig= temp;
515 
516  tau1_vertexDeltaR= tau2_vertexDeltaR;
517 
518  temp = tau1_vertexEnergyRatio;
519  tau1_vertexEnergyRatio= tau2_vertexEnergyRatio;
520  tau2_vertexEnergyRatio= temp;
521 
522  temp = tau1_vertexMass;
523  tau1_vertexMass= tau2_vertexMass;
524  tau2_vertexMass= temp;
525  }
526 
527 
528  std::map<std::string,float> inputs;
529  inputs["z_ratio"] = z_ratio;
530  inputs["trackSipdSig_3"] = trackSip3dSig_3;
531  inputs["trackSipdSig_2"] = trackSip3dSig_2;
532  inputs["trackSipdSig_1"] = trackSip3dSig_1;
533  inputs["trackSipdSig_0"] = trackSip3dSig_0;
534  inputs["trackSipdSig_1_0"] = tau2_trackSip3dSig_0;
535  inputs["trackSipdSig_0_0"] = tau1_trackSip3dSig_0;
536  inputs["trackSipdSig_1_1"] = tau2_trackSip3dSig_1;
537  inputs["trackSipdSig_0_1"] = tau1_trackSip3dSig_1;
538  inputs["trackSip2dSigAboveCharm_0"] = trackSip2dSigAboveCharm_0;
539  inputs["trackSip2dSigAboveBottom_0"] = trackSip2dSigAboveBottom_0;
540  inputs["trackSip2dSigAboveBottom_1"] = trackSip2dSigAboveBottom_1;
541  inputs["tau1_trackEtaRel_0"] = tau2_trackEtaRel_0;
542  inputs["tau1_trackEtaRel_1"] = tau2_trackEtaRel_1;
543  inputs["tau1_trackEtaRel_2"] = tau2_trackEtaRel_2;
544  inputs["tau0_trackEtaRel_0"] = tau1_trackEtaRel_0;
545  inputs["tau0_trackEtaRel_1"] = tau1_trackEtaRel_1;
546  inputs["tau0_trackEtaRel_2"] = tau1_trackEtaRel_2;
547  inputs["tau_vertexMass_0"] = tau1_vertexMass;
548  inputs["tau_vertexEnergyRatio_0"] = tau1_vertexEnergyRatio;
549  inputs["tau_vertexDeltaR_0"] = tau1_vertexDeltaR;
550  inputs["tau_flightDistance2dSig_0"] = tau1_flightDistance2dSig;
551  inputs["tau_vertexMass_1"] = tau2_vertexMass;
552  inputs["tau_vertexEnergyRatio_1"] = tau2_vertexEnergyRatio;
553  inputs["tau_flightDistance2dSig_1"] = tau2_flightDistance2dSig;
554  inputs["jetNTracks"] = jetNTracks;
555  inputs["nSV"] = nSV;
556 
557  // evaluate the MVA
558  value = mvaID->evaluate(inputs);
559 
560  // return the final discriminator value
561  return value;
562 }
563 
564 
565 void CandidateBoostedDoubleSecondaryVertexComputer::calcNsubjettiness(const reco::JetBaseRef & jet, float & tau1, float & tau2, std::vector<fastjet::PseudoJet> & currentAxes) const
566 {
567  std::vector<fastjet::PseudoJet> fjParticles;
568 
569  // loop over jet constituents and push them in the vector of FastJet constituents
570  for(const reco::CandidatePtr & daughter : jet->daughterPtrVector())
571  {
572  if ( daughter.isNonnull() && daughter.isAvailable() )
573  {
574  const reco::Jet * subjet = dynamic_cast<const reco::Jet *>(daughter.get());
575  // if the daughter is actually a subjet
576  if( subjet && daughter->numberOfDaughters() > 1 )
577  {
578  // loop over subjet constituents and push them in the vector of FastJet constituents
579  for(size_t i=0; i<daughter->numberOfDaughters(); ++i)
580  {
581  const reco::Candidate * constit = daughter->daughter(i);
582 
583  if ( constit )
584  fjParticles.push_back( fastjet::PseudoJet( constit->px(), constit->py(), constit->pz(), constit->energy() ) );
585  else
586  edm::LogWarning("MissingJetConstituent") << "Jet constituent required for N-subjettiness computation is missing!";
587  }
588  }
589  else
590  fjParticles.push_back( fastjet::PseudoJet( daughter->px(), daughter->py(), daughter->pz(), daughter->energy() ) );
591  }
592  else
593  edm::LogWarning("MissingJetConstituent") << "Jet constituent required for N-subjettiness computation is missing!";
594  }
595 
596  // N-subjettiness calculator
597  fastjet::contrib::Njettiness njettiness(fastjet::contrib::OnePass_KT_Axes(), fastjet::contrib::NormalizedMeasure(beta_,R0_));
598 
599  // calculate N-subjettiness
600  tau1 = njettiness.getTau(1, fjParticles);
601  tau2 = njettiness.getTau(2, fjParticles);
602  currentAxes = njettiness.currentAxes();
603 }
604 
605 
606 void CandidateBoostedDoubleSecondaryVertexComputer::setTracksPVBase(const reco::TrackRef & trackRef, const reco::VertexRef & vertexRef, float & PVweight) const
607 {
608  PVweight = 0.;
609 
610  const reco::TrackBaseRef trackBaseRef( trackRef );
611 
613 
614  const reco::Vertex & vtx = *(vertexRef);
615  // loop over tracks in vertices
616  for(IT it=vtx.tracks_begin(); it!=vtx.tracks_end(); ++it)
617  {
618  const reco::TrackBaseRef & baseRef = *it;
619  // one of the tracks in the vertex is the same as the track considered in the function
620  if( baseRef == trackBaseRef )
621  {
622  PVweight = vtx.trackWeight(baseRef);
623  break;
624  }
625  }
626 }
627 
628 
629 void CandidateBoostedDoubleSecondaryVertexComputer::setTracksPV(const reco::CandidatePtr & trackRef, const reco::VertexRef & vertexRef, float & PVweight) const
630 {
631  PVweight = 0.;
632 
633  const pat::PackedCandidate * pcand = dynamic_cast<const pat::PackedCandidate *>(trackRef.get());
634 
635  if(pcand) // MiniAOD case
636  {
637  if( pcand->fromPV() == pat::PackedCandidate::PVUsedInFit )
638  {
639  PVweight = 1.;
640  }
641  }
642  else
643  {
644  const reco::PFCandidate * pfcand = dynamic_cast<const reco::PFCandidate *>(trackRef.get());
645 
646  setTracksPVBase(pfcand->trackRef(), vertexRef, PVweight);
647  }
648 }
649 
650 
652  fastjet::PseudoJet & tauAxis, std::vector<float> & tau_trackEtaRel) const
653 {
654  math::XYZVector direction(tauAxis.px(), tauAxis.py(), tauAxis.pz());
655  const std::vector<reco::CandidatePtr> & tracks = vertex.daughterPtrVector();
656 
657  for(std::vector<reco::CandidatePtr>::const_iterator track = tracks.begin(); track != tracks.end(); ++track)
658  tau_trackEtaRel.push_back(std::abs(reco::btau::etaRel(direction.Unit(), (*track)->momentum())));
659 }
reco::Vertex::Point convertPos(const GlobalPoint &p)
int i
Definition: DBlmapReader.cc:9
const math::XYZTLorentzVector & weightedVectorSum() const
virtual double energy() const =0
energy
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
const T & get(unsigned int index=0) const
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
tuple trackPairV0Filter
tuple cont
load Luminosity info ##
Definition: generateEDF.py:622
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:44
JetCorrectorParameters::Record record
Definition: classes.h:7
void add(const reco::Track &track, double weight=1.0)
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:160
Base class for all types of Jets.
Definition: Jet.h:20
double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
Measurement1D ip2d
Definition: IPTagInfo.h:31
void setTracksPV(const reco::CandidatePtr &trackRef, const reco::VertexRef &vertexRef, float &PVweight) const
virtual double pz() const =0
z coordinate of momentum vector
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
GlobalPoint globalPosition() const
const Container & selectedTracks() const
Definition: IPTagInfo.h:101
TrajectoryStateOnSurface closestApproachToJet(const TrajectoryStateOnSurface &state, const reco::Vertex &vertex, const GlobalVector &aJetDirection, const MagneticField *field)
Definition: IPTools.cc:177
const MagneticField * field() const
void etaRelToTauAxis(const reco::VertexCompositePtrCandidate &vertex, fastjet::PseudoJet &tauAxis, std::vector< float > &tau_trackEtaRel) const
const edm::Ref< VertexCollection > & primaryVertex() const
Definition: IPTagInfo.h:133
std::pair< double, Measurement1D > jetTrackDistance(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:200
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:433
U second(std::pair< T, U > const &p)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
virtual size_t numberOfDaughters() const
number of daughters
const PVAssoc fromPV(size_t ipv=0) const
void get(HolderT &iHolder) const
void uses(unsigned int id, const std::string &label)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
virtual double py() const final
y coordinate of momentum vector
float trackWeight(const TREF &r) const
returns the weight with which a Track has contributed to the vertex-fit.
Definition: Vertex.h:73
void setTracksPVBase(const reco::TrackRef &trackRef, const reco::VertexRef &vertexRef, float &PVweight) const
virtual double py() const =0
y coordinate of momentum vector
std::vector< LinkConnSpec >::const_iterator IT
std::vector< size_t > sortedIndexes(btag::SortCriteria mode=reco::btag::IP3DSig) const
Definition: IPTagInfo.h:238
virtual double pz() const final
z coordinate of momentum vector
const std::vector< btag::TrackIPData > & impactParameterData() const
Definition: IPTagInfo.h:91
virtual double px() const =0
x coordinate of momentum vector
double significance() const
Definition: Measurement1D.h:32
virtual Vector momentum() const final
spatial momentum vector
tuple tracks
Definition: testEve_cfg.py:39
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
const math::XYZTLorentzVector & vectorSum() const
void calcNsubjettiness(const reco::JetBaseRef &jet, float &tau1, float &tau2, std::vector< fastjet::PseudoJet > &currentAxes) const
const daughters & daughterPtrVector() const
references to daughtes
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
float discriminator(const TagInfoHelper &tagInfos) const override
virtual double px() const final
x coordinate of momentum vector
TrajectoryStateOnSurface impactPointState() const
volatile std::atomic< bool > shutdown_flag false
std::string fullPath() const
Definition: FileInPath.cc:184
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:39
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
Global3DVector GlobalVector
Definition: GlobalVector.h:10
tuple trackSelector
Tracks selection.
Definition: valSkim_cff.py:4