CMS 3D CMS Logo

TemplatedSecondaryVertexProducer.cc
Go to the documentation of this file.
1 #include <functional>
2 #include <algorithm>
3 #include <iterator>
4 #include <cstddef>
5 #include <string>
6 #include <vector>
7 #include <map>
8 #include <set>
9 
10 #include <boost/iterator/transform_iterator.hpp>
11 
23 
28 
33 
41 
46 
52 
54 
55 #include "fastjet/JetDefinition.hh"
56 #include "fastjet/ClusterSequence.hh"
57 #include "fastjet/PseudoJet.hh"
58 
59 //
60 // constants, enums and typedefs
61 //
62 typedef boost::shared_ptr<fastjet::ClusterSequence> ClusterSequencePtr;
63 typedef boost::shared_ptr<fastjet::JetDefinition> JetDefPtr;
64 
65 
66 using namespace reco;
67 
68 namespace {
69  class VertexInfo : public fastjet::PseudoJet::UserInfoBase{
70  public:
71  VertexInfo(const int vertexIndex) :
72  m_vertexIndex(vertexIndex) { }
73 
74  inline const int vertexIndex() const { return m_vertexIndex; }
75 
76  protected:
77  int m_vertexIndex;
78  };
79 
80  template<typename T>
81  struct RefToBaseLess : public std::binary_function<edm::RefToBase<T>,
82  edm::RefToBase<T>,
83  bool> {
84  inline bool operator()(const edm::RefToBase<T> &r1,
85  const edm::RefToBase<T> &r2) const
86  {
87  return r1.id() < r2.id() ||
88  (r1.id() == r2.id() && r1.key() < r2.key());
89  }
90  };
91 }
92 
94 return GlobalVector(sv.x() - pv.x(), sv.y() - pv.y(),sv.z() - pv.z());
95 }
97 return GlobalVector(sv.vertex().x() - pv.x(), sv.vertex().y() - pv.y(),sv.vertex().z() - pv.z());
98 }
100 {return sv.position();}
102 {return sv.vertex();}
103 
104 
105 template <class IPTI,class VTX>
107  public:
108  explicit TemplatedSecondaryVertexProducer(const edm::ParameterSet &params);
110  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
111  typedef std::vector<TemplatedSecondaryVertexTagInfo<IPTI,VTX> > Product;
113  typedef typename IPTI::input_container input_container;
115  typedef typename std::vector<reco::btag::IndexedTrackData> TrackDataVector;
116  virtual void produce(edm::Event &event, const edm::EventSetup &es) override;
117 
118  private:
119  template<class CONTAINER>
120  void matchReclusteredJets(const edm::Handle<CONTAINER>& jets,
121  const std::vector<fastjet::PseudoJet>& matchedJets,
122  std::vector<int>& matchedIndices,
123  const std::string& jetType="");
124  void matchGroomedJets(const edm::Handle<edm::View<reco::Jet> >& jets,
125  const edm::Handle<edm::View<reco::Jet> >& matchedJets,
126  std::vector<int>& matchedIndices);
127  void matchSubjets(const std::vector<int>& groomedIndices,
128  const edm::Handle<edm::View<reco::Jet> >& groomedJets,
129  const edm::Handle<std::vector<IPTI> >& subjets,
130  std::vector<std::vector<int> >& matchedIndices);
131  void matchSubjets(const edm::Handle<edm::View<reco::Jet> >& fatJets,
132  const edm::Handle<std::vector<IPTI> >& subjets,
133  std::vector<std::vector<int> >& matchedIndices);
134 
135  const reco::Jet * toJet(const reco::Jet & j) { return &j; }
136  const reco::Jet * toJet(const IPTI & j) { return &(*(j.jet())); }
137 
139  CONSTRAINT_NONE = 0,
144  CONSTRAINT_PV_PRIMARIES_IN_FIT
145  };
146  static ConstraintType getConstraintType(const std::string &name);
147 
166  double rParam;
167  double jetPtMin;
174 
177 
178  void markUsedTracks(TrackDataVector & trackData, const input_container & trackRefs, const SecondaryVertex & sv,size_t idx);
179 
180  struct SVBuilder :
181  public std::unary_function<const VTX&, SecondaryVertex> {
182 
184  const GlobalVector &direction,
185  bool withPVError,
186  double minTrackWeight) :
187  pv(pv), direction(direction),
188  withPVError(withPVError),
189  minTrackWeight(minTrackWeight) {}
190  SecondaryVertex operator () (const TransientVertex &sv) const;
191 
192  SecondaryVertex operator () (const VTX &sv) const
193  { return SecondaryVertex(pv, sv, direction, withPVError); }
194 
195 
196  const Vertex &pv;
200  };
201 
202  struct SVFilter :
203  public std::unary_function<const SecondaryVertex&, bool> {
204 
206  const GlobalVector &direction) :
207  filter(filter), pv(pv), direction(direction) {}
208 
209  inline bool operator () (const SecondaryVertex &sv) const
210  { return !filter(pv, sv, direction); }
211 
213  const Vertex &pv;
215  };
216 };
217 template <class IPTI,class VTX>
220 {
221  if (name == "None")
222  return CONSTRAINT_NONE;
223  else if (name == "BeamSpot")
224  return CONSTRAINT_BEAMSPOT;
225  else if (name == "BeamSpot+PVPosition")
226  return CONSTRAINT_PV_BEAMSPOT_SIZE;
227  else if (name == "BeamSpotZ+PVErrorScaledXY")
228  return CONSTRAINT_PV_BS_Z_ERRORS_SCALED;
229  else if (name == "PVErrorScaled")
230  return CONSTRAINT_PV_ERROR_SCALED;
231  else if (name == "BeamSpot+PVTracksInFit")
232  return CONSTRAINT_PV_PRIMARIES_IN_FIT;
233  else
234  throw cms::Exception("InvalidArgument")
235  << "TemplatedSecondaryVertexProducer: ``constraint'' parameter "
236  "value \"" << name << "\" not understood."
237  << std::endl;
238 }
239 
242 {
243  if (name == "AlwaysWithGhostTrack")
245  else if (name == "SingleTracksWithGhostTrack")
247  else if (name == "RefitGhostTrackWithVertices")
249  else
250  throw cms::Exception("InvalidArgument")
251  << "TemplatedSecondaryVertexProducer: ``fitType'' "
252  "parameter value \"" << name << "\" for "
253  "GhostTrackVertexFinder settings not "
254  "understood." << std::endl;
255 }
256 
257 template <class IPTI,class VTX>
259  const edm::ParameterSet &params) :
260  sortCriterium(TrackSorting::getCriterium(params.getParameter<std::string>("trackSort"))),
261  trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")),
262  constraint(getConstraintType(params.getParameter<std::string>("constraint"))),
263  constraintScaling(1.0),
264  vtxRecoPSet(params.getParameter<edm::ParameterSet>("vertexReco")),
265  useGhostTrack(vtxRecoPSet.getParameter<std::string>("finder") == "gtvr"),
266  withPVError(params.getParameter<bool>("usePVError")),
267  minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
268  vertexFilter(params.getParameter<edm::ParameterSet>("vertexCuts")),
269  vertexSorting(params.getParameter<edm::ParameterSet>("vertexSelection"))
270 {
271  token_trackIPTagInfo = consumes<std::vector<IPTI> >(params.getParameter<edm::InputTag>("trackIPTagInfos"));
274  constraintScaling = params.getParameter<double>("pvErrorScaling");
275 
280  token_BeamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpotTag"));
281  useExternalSV = false;
282  if(params.existsAs<bool>("useExternalSV")) useExternalSV = params.getParameter<bool> ("useExternalSV");
283  if(useExternalSV) {
284  token_extSVCollection = consumes<edm::View<VTX> >(params.getParameter<edm::InputTag>("extSVCollection"));
285  extSVDeltaRToJet = params.getParameter<double>("extSVDeltaRToJet");
286  }
287  useSVClustering = ( params.existsAs<bool>("useSVClustering") ? params.getParameter<bool>("useSVClustering") : false );
288  useSVMomentum = ( params.existsAs<bool>("useSVMomentum") ? params.getParameter<bool>("useSVMomentum") : false );
289  useFatJets = ( useExternalSV && params.exists("fatJets") );
290  useGroomedFatJets = ( useExternalSV && params.exists("groomedFatJets") );
291  if( useSVClustering )
292  {
293  jetAlgorithm = params.getParameter<std::string>("jetAlgorithm");
294  rParam = params.getParameter<double>("rParam");
295  jetPtMin = 0.; // hardcoded to 0. since we simply want to recluster all input jets which already had some PtMin applied
296  ghostRescaling = ( params.existsAs<double>("ghostRescaling") ? params.getParameter<double>("ghostRescaling") : 1e-18 );
297  relPtTolerance = ( params.existsAs<double>("relPtTolerance") ? params.getParameter<double>("relPtTolerance") : 1e-03); // 0.1% relative difference in Pt should be sufficient to detect possible misconfigurations
298 
299  // set jet algorithm
300  if (jetAlgorithm=="Kt")
301  fjJetDefinition = JetDefPtr( new fastjet::JetDefinition(fastjet::kt_algorithm, rParam) );
302  else if (jetAlgorithm=="CambridgeAachen")
303  fjJetDefinition = JetDefPtr( new fastjet::JetDefinition(fastjet::cambridge_algorithm, rParam) );
304  else if (jetAlgorithm=="AntiKt")
305  fjJetDefinition = JetDefPtr( new fastjet::JetDefinition(fastjet::antikt_algorithm, rParam) );
306  else
307  throw cms::Exception("InvalidJetAlgorithm") << "Jet clustering algorithm is invalid: " << jetAlgorithm << ", use CambridgeAachen | Kt | AntiKt" << std::endl;
308  }
309  if( useFatJets )
310  token_fatJets = consumes<edm::View<reco::Jet> >(params.getParameter<edm::InputTag>("fatJets"));
311  if( useGroomedFatJets )
312  token_groomedFatJets = consumes<edm::View<reco::Jet> >(params.getParameter<edm::InputTag>("groomedFatJets"));
313  if( useFatJets && !useSVClustering )
314  rParam = params.getParameter<double>("rParam"); // will be used later as a dR cut
315 
316  produces<Product>();
317 }
318 template <class IPTI,class VTX>
320 {
321 }
322 
323 template <class IPTI,class VTX>
325  const edm::EventSetup &es)
326 {
327 // typedef std::map<TrackBaseRef, TransientTrack,
328 // RefToBaseLess<Track> > TransientTrackMap;
329  //How about good old pointers?
330  typedef std::map<const Track *, TransientTrack> TransientTrackMap;
331 
332 
334  es.get<TransientTrackRecord>().get("TransientTrackBuilder",
335  trackBuilder);
336 
338  event.getByToken(token_trackIPTagInfo, trackIPTagInfos);
339 
340  // External Sec Vertex collection (e.g. for IVF usage)
341  edm::Handle<edm::View<VTX> > extSecVertex;
342  if(useExternalSV) event.getByToken(token_extSVCollection,extSecVertex);
343 
344  edm::Handle<edm::View<reco::Jet> > fatJetsHandle;
345  edm::Handle<edm::View<reco::Jet> > groomedFatJetsHandle;
346  if( useFatJets )
347  {
348  event.getByToken(token_fatJets, fatJetsHandle);
349 
350  if( useGroomedFatJets )
351  {
352  event.getByToken(token_groomedFatJets, groomedFatJetsHandle);
353 
354  if( groomedFatJetsHandle->size() > fatJetsHandle->size() )
355  edm::LogError("TooManyGroomedJets") << "There are more groomed (" << groomedFatJetsHandle->size() << ") than original fat jets (" << fatJetsHandle->size() << "). Please check that the two jet collections belong to each other.";
356  }
357  }
358 
360  unsigned int bsCovSrc[7] = { 0, };
361  double sigmaZ = 0.0, beamWidth = 0.0;
362  switch(constraint) {
364  event.getByToken(token_BeamSpot,beamSpot);
365  bsCovSrc[3] = bsCovSrc[4] = bsCovSrc[5] = bsCovSrc[6] = 1;
366  sigmaZ = beamSpot->sigmaZ();
367  beamWidth = beamSpot->BeamWidthX();
368  break;
369 
371  event.getByToken(token_BeamSpot,beamSpot);
372  bsCovSrc[0] = bsCovSrc[1] = 2;
373  bsCovSrc[3] = bsCovSrc[4] = bsCovSrc[5] = 1;
374  sigmaZ = beamSpot->sigmaZ();
375  break;
376 
378  bsCovSrc[0] = bsCovSrc[1] = bsCovSrc[2] = 2;
379  break;
380 
381  case CONSTRAINT_BEAMSPOT:
383  event.getByToken(token_BeamSpot,beamSpot);
384  break;
385 
386  default:
387  /* nothing */;
388  }
389 
390  // ------------------------------------ SV clustering START --------------------------------------------
391  std::vector<std::vector<int> > clusteredSVs(trackIPTagInfos->size(),std::vector<int>());
392  if( useExternalSV && useSVClustering && trackIPTagInfos->size()>0 )
393  {
394  // vector of constituents for reclustering jets and "ghost" SVs
395  std::vector<fastjet::PseudoJet> fjInputs;
396  // loop over all input jets and collect all their constituents
397  if( useFatJets )
398  {
399  for(edm::View<reco::Jet>::const_iterator it = fatJetsHandle->begin(); it != fatJetsHandle->end(); ++it)
400  {
401  std::vector<edm::Ptr<reco::Candidate> > constituents = it->getJetConstituents();
402  std::vector<edm::Ptr<reco::Candidate> >::const_iterator m;
403  for( m = constituents.begin(); m != constituents.end(); ++m )
404  {
405  reco::CandidatePtr constit = *m;
406  if(constit->pt() == 0)
407  {
408  edm::LogWarning("NullTransverseMomentum") << "dropping input candidate with pt=0";
409  continue;
410  }
411  fjInputs.push_back(fastjet::PseudoJet(constit->px(),constit->py(),constit->pz(),constit->energy()));
412  }
413  }
414  }
415  else
416  {
417  for(typename std::vector<IPTI>::const_iterator it = trackIPTagInfos->begin(); it != trackIPTagInfos->end(); ++it)
418  {
419  std::vector<edm::Ptr<reco::Candidate> > constituents = it->jet()->getJetConstituents();
420  std::vector<edm::Ptr<reco::Candidate> >::const_iterator m;
421  for( m = constituents.begin(); m != constituents.end(); ++m )
422  {
423  reco::CandidatePtr constit = *m;
424  if(constit->pt() == 0)
425  {
426  edm::LogWarning("NullTransverseMomentum") << "dropping input candidate with pt=0";
427  continue;
428  }
429  fjInputs.push_back(fastjet::PseudoJet(constit->px(),constit->py(),constit->pz(),constit->energy()));
430  }
431  }
432  }
433  // insert "ghost" SVs in the vector of constituents
434  for(typename edm::View<VTX>::const_iterator it = extSecVertex->begin(); it != extSecVertex->end(); ++it)
435  {
436  const reco::Vertex &pv = *(trackIPTagInfos->front().primaryVertex());
437  GlobalVector dir = flightDirection(pv, *it);
438  dir = dir.unit();
439  fastjet::PseudoJet p(dir.x(),dir.y(),dir.z(),dir.mag()); // using SV flight direction so treating SV as massless
440  if( useSVMomentum )
441  p = fastjet::PseudoJet(it->p4().px(),it->p4().py(),it->p4().pz(),it->p4().energy());
442  p*=ghostRescaling; // rescale SV direction/momentum
443  p.set_user_info(new VertexInfo( it - extSecVertex->begin() ));
444  fjInputs.push_back(p);
445  }
446 
447  // define jet clustering sequence
448  fjClusterSeq = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs, *fjJetDefinition ) );
449  // recluster jet constituents and inserted "ghosts"
450  std::vector<fastjet::PseudoJet> inclusiveJets = fastjet::sorted_by_pt( fjClusterSeq->inclusive_jets(jetPtMin) );
451 
452  if( useFatJets )
453  {
454  if( inclusiveJets.size() < fatJetsHandle->size() )
455  edm::LogError("TooFewReclusteredJets") << "There are fewer reclustered (" << inclusiveJets.size() << ") than original fat jets (" << fatJetsHandle->size() << "). Please check that the jet algorithm and jet size match those used for the original jet collection.";
456 
457  // match reclustered and original fat jets
458  std::vector<int> reclusteredIndices;
459  matchReclusteredJets<edm::View<reco::Jet> >(fatJetsHandle,inclusiveJets,reclusteredIndices,"fat");
460 
461  // match groomed and original fat jets
462  std::vector<int> groomedIndices;
463  if( useGroomedFatJets )
464  matchGroomedJets(fatJetsHandle,groomedFatJetsHandle,groomedIndices);
465 
466  // match subjets and original fat jets
467  std::vector<std::vector<int> > subjetIndices;
468  if( useGroomedFatJets )
469  matchSubjets(groomedIndices,groomedFatJetsHandle,trackIPTagInfos,subjetIndices);
470  else
471  matchSubjets(fatJetsHandle,trackIPTagInfos,subjetIndices);
472 
473  // collect clustered SVs
474  for(size_t i=0; i<fatJetsHandle->size(); ++i)
475  {
476  if( reclusteredIndices.at(i) < 0 ) continue; // continue if matching reclustered to original jets failed
477 
478  if( fatJetsHandle->at(i).pt() == 0 ) // continue if the original jet has Pt=0
479  {
480  edm::LogWarning("NullTransverseMomentum") << "The original fat jet " << i << " has Pt=0. This is not expected so the jet will be skipped.";
481  continue;
482  }
483 
484  if( subjetIndices.at(i).size()==0 ) continue; // continue if the original jet does not have subjets assigned
485 
486  // since the "ghosts" are extremely soft, the configuration and ordering of the reclustered and original fat jets should in principle stay the same
487  if( ( std::abs( inclusiveJets.at(reclusteredIndices.at(i)).pt() - fatJetsHandle->at(i).pt() ) / fatJetsHandle->at(i).pt() ) > relPtTolerance )
488  {
489  if( fatJetsHandle->at(i).pt() < 10. ) // special handling for low-Pt jets (Pt<10 GeV)
490  edm::LogWarning("JetPtMismatchAtLowPt") << "The reclustered and original fat jet " << i << " have different Pt's (" << inclusiveJets.at(reclusteredIndices.at(i)).pt() << " vs " << fatJetsHandle->at(i).pt() << " GeV, respectively).\n"
491  << "Please check that the jet algorithm and jet size match those used for the original fat jet collection and also make sure the original fat jets are uncorrected. In addition, make sure you are not using CaloJets which are presently not supported.\n"
492  << "Since the mismatch is at low Pt, it is ignored and only a warning is issued.\n"
493  << "\nIn extremely rare instances the mismatch could be caused by a difference in the machine precision in which case make sure the original jet collection is produced and reclustering is performed in the same job.";
494  else
495  edm::LogError("JetPtMismatch") << "The reclustered and original fat jet " << i << " have different Pt's (" << inclusiveJets.at(reclusteredIndices.at(i)).pt() << " vs " << fatJetsHandle->at(i).pt() << " GeV, respectively).\n"
496  << "Please check that the jet algorithm and jet size match those used for the original fat jet collection and also make sure the original fat jets are uncorrected. In addition, make sure you are not using CaloJets which are presently not supported.\n"
497  << "\nIn extremely rare instances the mismatch could be caused by a difference in the machine precision in which case make sure the original jet collection is produced and reclustering is performed in the same job.";
498  }
499 
500  // get jet constituents
501  std::vector<fastjet::PseudoJet> constituents = inclusiveJets.at(reclusteredIndices.at(i)).constituents();
502 
503  std::vector<int> svIndices;
504  // loop over jet constituents and try to find "ghosts"
505  for(std::vector<fastjet::PseudoJet>::const_iterator it = constituents.begin(); it != constituents.end(); ++it)
506  {
507  if( !it->has_user_info() ) continue; // skip if not a "ghost"
508 
509  svIndices.push_back( it->user_info<VertexInfo>().vertexIndex() );
510  }
511 
512  // loop over clustered SVs and assign them to different subjets based on smallest dR
513  for(size_t sv=0; sv<svIndices.size(); ++sv)
514  {
515  const reco::Vertex &pv = *(trackIPTagInfos->front().primaryVertex());
516  const VTX &extSV = (*extSecVertex)[ svIndices.at(sv) ];
517  GlobalVector dir = flightDirection(pv, extSV);
518  dir = dir.unit();
519  fastjet::PseudoJet p(dir.x(),dir.y(),dir.z(),dir.mag()); // using SV flight direction so treating SV as massless
520  if( useSVMomentum )
521  p = fastjet::PseudoJet(extSV.p4().px(),extSV.p4().py(),extSV.p4().pz(),extSV.p4().energy());
522 
523  std::vector<double> dR2toSubjets;
524 
525  for(size_t sj=0; sj<subjetIndices.at(i).size(); ++sj)
526  dR2toSubjets.push_back( Geom::deltaR2( p.rapidity(), p.phi_std(), trackIPTagInfos->at(subjetIndices.at(i).at(sj)).jet()->rapidity(), trackIPTagInfos->at(subjetIndices.at(i).at(sj)).jet()->phi() ) );
527 
528  // find the closest subjet
529  int closestSubjetIdx = std::distance( dR2toSubjets.begin(), std::min_element(dR2toSubjets.begin(), dR2toSubjets.end()) );
530 
531  clusteredSVs.at(subjetIndices.at(i).at(closestSubjetIdx)).push_back( svIndices.at(sv) );
532  }
533  }
534  }
535  else
536  {
537  if( inclusiveJets.size() < trackIPTagInfos->size() )
538  edm::LogError("TooFewReclusteredJets") << "There are fewer reclustered (" << inclusiveJets.size() << ") than original jets (" << trackIPTagInfos->size() << "). Please check that the jet algorithm and jet size match those used for the original jet collection.";
539 
540  // match reclustered and original jets
541  std::vector<int> reclusteredIndices;
542  matchReclusteredJets<std::vector<IPTI> >(trackIPTagInfos,inclusiveJets,reclusteredIndices);
543 
544  // collect clustered SVs
545  for(size_t i=0; i<trackIPTagInfos->size(); ++i)
546  {
547  if( reclusteredIndices.at(i) < 0 ) continue; // continue if matching reclustered to original jets failed
548 
549  if( trackIPTagInfos->at(i).jet()->pt() == 0 ) // continue if the original jet has Pt=0
550  {
551  edm::LogWarning("NullTransverseMomentum") << "The original jet " << i << " has Pt=0. This is not expected so the jet will be skipped.";
552  continue;
553  }
554 
555  // since the "ghosts" are extremely soft, the configuration and ordering of the reclustered and original jets should in principle stay the same
556  if( ( std::abs( inclusiveJets.at(reclusteredIndices.at(i)).pt() - trackIPTagInfos->at(i).jet()->pt() ) / trackIPTagInfos->at(i).jet()->pt() ) > relPtTolerance )
557  {
558  if( trackIPTagInfos->at(i).jet()->pt() < 10. ) // special handling for low-Pt jets (Pt<10 GeV)
559  edm::LogWarning("JetPtMismatchAtLowPt") << "The reclustered and original jet " << i << " have different Pt's (" << inclusiveJets.at(reclusteredIndices.at(i)).pt() << " vs " << trackIPTagInfos->at(i).jet()->pt() << " GeV, respectively).\n"
560  << "Please check that the jet algorithm and jet size match those used for the original jet collection and also make sure the original jets are uncorrected. In addition, make sure you are not using CaloJets which are presently not supported.\n"
561  << "Since the mismatch is at low Pt, it is ignored and only a warning is issued.\n"
562  << "\nIn extremely rare instances the mismatch could be caused by a difference in the machine precision in which case make sure the original jet collection is produced and reclustering is performed in the same job.";
563  else
564  edm::LogError("JetPtMismatch") << "The reclustered and original jet " << i << " have different Pt's (" << inclusiveJets.at(reclusteredIndices.at(i)).pt() << " vs " << trackIPTagInfos->at(i).jet()->pt() << " GeV, respectively).\n"
565  << "Please check that the jet algorithm and jet size match those used for the original jet collection and also make sure the original jets are uncorrected. In addition, make sure you are not using CaloJets which are presently not supported.\n"
566  << "\nIn extremely rare instances the mismatch could be caused by a difference in the machine precision in which case make sure the original jet collection is produced and reclustering is performed in the same job.";
567  }
568 
569  // get jet constituents
570  std::vector<fastjet::PseudoJet> constituents = inclusiveJets.at(reclusteredIndices.at(i)).constituents();
571 
572  // loop over jet constituents and try to find "ghosts"
573  for(std::vector<fastjet::PseudoJet>::const_iterator it = constituents.begin(); it != constituents.end(); ++it)
574  {
575  if( !it->has_user_info() ) continue; // skip if not a "ghost"
576  // push back clustered SV indices
577  clusteredSVs.at(i).push_back( it->user_info<VertexInfo>().vertexIndex() );
578  }
579  }
580  }
581  }
582  // case where fat jets are used to associate SVs to subjets but no SV clustering is performed
583  else if( useExternalSV && !useSVClustering && trackIPTagInfos->size()>0 && useFatJets )
584  {
585  // match groomed and original fat jets
586  std::vector<int> groomedIndices;
587  if( useGroomedFatJets )
588  matchGroomedJets(fatJetsHandle,groomedFatJetsHandle,groomedIndices);
589 
590  // match subjets and original fat jets
591  std::vector<std::vector<int> > subjetIndices;
592  if( useGroomedFatJets )
593  matchSubjets(groomedIndices,groomedFatJetsHandle,trackIPTagInfos,subjetIndices);
594  else
595  matchSubjets(fatJetsHandle,trackIPTagInfos,subjetIndices);
596 
597  // loop over fat jets
598  for(size_t i=0; i<fatJetsHandle->size(); ++i)
599  {
600  if( fatJetsHandle->at(i).pt() == 0 ) // continue if the original jet has Pt=0
601  {
602  edm::LogWarning("NullTransverseMomentum") << "The original fat jet " << i << " has Pt=0. This is not expected so the jet will be skipped.";
603  continue;
604  }
605 
606  if( subjetIndices.at(i).size()==0 ) continue; // continue if the original jet does not have subjets assigned
607 
608  // loop over SVs, associate them to fat jets based on dR cone and
609  // then assign them to the closets subjet in dR
610  for(typename edm::View<VTX>::const_iterator it = extSecVertex->begin(); it != extSecVertex->end(); ++it)
611  {
612  size_t sv = ( it - extSecVertex->begin() );
613 
614  const reco::Vertex &pv = *(trackIPTagInfos->front().primaryVertex());
615  const VTX &extSV = (*extSecVertex)[sv];
616  GlobalVector dir = flightDirection(pv, extSV);
617  GlobalVector jetDir(fatJetsHandle->at(i).px(),
618  fatJetsHandle->at(i).py(),
619  fatJetsHandle->at(i).pz());
620  // skip SVs outside the dR cone
621  if( Geom::deltaR2( dir, jetDir ) > rParam*rParam ) // here using the jet clustering rParam as a dR cut
622  continue;
623 
624  dir = dir.unit();
625  fastjet::PseudoJet p(dir.x(),dir.y(),dir.z(),dir.mag()); // using SV flight direction so treating SV as massless
626  if( useSVMomentum )
627  p = fastjet::PseudoJet(extSV.p4().px(),extSV.p4().py(),extSV.p4().pz(),extSV.p4().energy());
628 
629  std::vector<double> dR2toSubjets;
630 
631  for(size_t sj=0; sj<subjetIndices.at(i).size(); ++sj)
632  dR2toSubjets.push_back( Geom::deltaR2( p.rapidity(), p.phi_std(), trackIPTagInfos->at(subjetIndices.at(i).at(sj)).jet()->rapidity(), trackIPTagInfos->at(subjetIndices.at(i).at(sj)).jet()->phi() ) );
633 
634  // find the closest subjet
635  int closestSubjetIdx = std::distance( dR2toSubjets.begin(), std::min_element(dR2toSubjets.begin(), dR2toSubjets.end()) );
636 
637  clusteredSVs.at(subjetIndices.at(i).at(closestSubjetIdx)).push_back(sv);
638  }
639  }
640  }
641  // ------------------------------------ SV clustering END ----------------------------------------------
642 
643  std::auto_ptr<ConfigurableVertexReconstructor> vertexReco;
644  std::auto_ptr<GhostTrackVertexFinder> vertexRecoGT;
645  if (useGhostTrack)
646  vertexRecoGT.reset(new GhostTrackVertexFinder(
647  vtxRecoPSet.getParameter<double>("maxFitChi2"),
648  vtxRecoPSet.getParameter<double>("mergeThreshold"),
649  vtxRecoPSet.getParameter<double>("primcut"),
650  vtxRecoPSet.getParameter<double>("seccut"),
652  else
653  vertexReco.reset(
655 
656  TransientTrackMap primariesMap;
657 
658  // result secondary vertices
659 
660  auto tagInfos = std::make_unique<Product>();
661 
662  for(typename std::vector<IPTI>::const_iterator iterJets =
663  trackIPTagInfos->begin(); iterJets != trackIPTagInfos->end();
664  ++iterJets) {
665  TrackDataVector trackData;
666 // std::cout << "Jet " << iterJets-trackIPTagInfos->begin() << std::endl;
667 
668  const Vertex &pv = *iterJets->primaryVertex();
669 
670  std::set<TransientTrack> primaries;
672  for(Vertex::trackRef_iterator iter = pv.tracks_begin();
673  iter != pv.tracks_end(); ++iter) {
674  TransientTrackMap::iterator pos =
675  primariesMap.lower_bound(iter->get());
676 
677  if (pos != primariesMap.end() &&
678  pos->first == iter->get())
679  primaries.insert(pos->second);
680  else {
682  trackBuilder->build(
683  iter->castTo<TrackRef>());
684  primariesMap.insert(pos,
685  std::make_pair(iter->get(), track));
686  primaries.insert(track);
687  }
688  }
689  }
690 
691  edm::RefToBase<Jet> jetRef = iterJets->jet();
692 
693  GlobalVector jetDir(jetRef->momentum().x(),
694  jetRef->momentum().y(),
695  jetRef->momentum().z());
696 
697  std::vector<std::size_t> indices =
698  iterJets->sortedIndexes(sortCriterium);
699 
700  input_container trackRefs = iterJets->sortedTracks(indices);
701 
702  const std::vector<reco::btag::TrackIPData> &ipData =
703  iterJets->impactParameterData();
704 
705  // build transient tracks used for vertex reconstruction
706 
707  std::vector<TransientTrack> fitTracks;
708  std::vector<GhostTrackState> gtStates;
709  std::auto_ptr<GhostTrackPrediction> gtPred;
710  if (useGhostTrack)
711  gtPred.reset(new GhostTrackPrediction(
712  *iterJets->ghostTrack()));
713 
714  for(unsigned int i = 0; i < indices.size(); i++) {
716 
717  const input_item &trackRef = trackRefs[i];
718 
719  trackData.push_back(IndexedTrackData());
720  trackData.back().first = indices[i];
721 
722  // select tracks for SV finder
723 
724  if (!trackSelector(*reco::btag::toTrack(trackRef), ipData[indices[i]], *jetRef,
726  pv.position()))) {
727  trackData.back().second.svStatus =
729  continue;
730  }
731 
732  TransientTrackMap::const_iterator pos =
733  primariesMap.find(reco::btag::toTrack((trackRef)));
734  TransientTrack fitTrack;
735  if (pos != primariesMap.end()) {
736  primaries.erase(pos->second);
737  fitTrack = pos->second;
738  } else
739  fitTrack = trackBuilder->build(trackRef);
740  fitTracks.push_back(fitTrack);
741 
742  trackData.back().second.svStatus =
744 
745  if (useGhostTrack) {
746  GhostTrackState gtState(fitTrack);
747  GlobalPoint pos =
748  ipData[indices[i]].closestToGhostTrack;
749  gtState.linearize(*gtPred, true,
750  gtPred->lambda(pos));
751  gtState.setWeight(ipData[indices[i]].ghostTrackWeight);
752  gtStates.push_back(gtState);
753  }
754  }
755 
756  std::auto_ptr<GhostTrack> ghostTrack;
757  if (useGhostTrack)
758  ghostTrack.reset(new GhostTrack(
762  GlobalVector(
763  iterJets->ghostTrack()->px(),
764  iterJets->ghostTrack()->py(),
765  iterJets->ghostTrack()->pz()),
766  0.05),
767  *gtPred, gtStates,
768  iterJets->ghostTrack()->chi2(),
769  iterJets->ghostTrack()->ndof()));
770 
771  // perform actual vertex finding
772 
773 
774  std::vector<VTX> extAssoCollection;
775  std::vector<TransientVertex> fittedSVs;
776  std::vector<SecondaryVertex> SVs;
777  if(!useExternalSV){
778  switch(constraint) {
779  case CONSTRAINT_NONE:
780  if (useGhostTrack)
781  fittedSVs = vertexRecoGT->vertices(
782  pv, *ghostTrack);
783  else
784  fittedSVs = vertexReco->vertices(fitTracks);
785  break;
786 
787  case CONSTRAINT_BEAMSPOT:
788  if (useGhostTrack)
789  fittedSVs = vertexRecoGT->vertices(
790  pv, *beamSpot, *ghostTrack);
791  else
792  fittedSVs = vertexReco->vertices(fitTracks,
793  *beamSpot);
794  break;
795 
800  for(unsigned int i = 0; i < 7; i++) {
801  unsigned int covSrc = bsCovSrc[i];
802  for(unsigned int j = 0; j < 7; j++) {
803  double v=0.0;
804  if (!covSrc || bsCovSrc[j] != covSrc)
805  v = 0.0;
806  else if (covSrc == 1)
807  v = beamSpot->covariance(i, j);
808  else if (j<3 && i<3)
809  v = pv.covariance(i, j) *
811  cov(i, j) = v;
812  }
813  }
814 
815  BeamSpot bs(pv.position(), sigmaZ,
816  beamSpot.isValid() ? beamSpot->dxdz() : 0.,
817  beamSpot.isValid() ? beamSpot->dydz() : 0.,
818  beamWidth, cov, BeamSpot::Unknown);
819 
820  if (useGhostTrack)
821  fittedSVs = vertexRecoGT->vertices(
822  pv, bs, *ghostTrack);
823  else
824  fittedSVs = vertexReco->vertices(fitTracks, bs);
825  } break;
826 
828  std::vector<TransientTrack> primaries_(
829  primaries.begin(), primaries.end());
830  if (useGhostTrack)
831  fittedSVs = vertexRecoGT->vertices(
832  pv, *beamSpot, primaries_,
833  *ghostTrack);
834  else
835  fittedSVs = vertexReco->vertices(
836  primaries_, fitTracks,
837  *beamSpot);
838  } break;
839  }
840  // build combined SV information and filter
841  SVBuilder svBuilder(pv, jetDir, withPVError, minTrackWeight);
842  std::remove_copy_if(boost::make_transform_iterator(
843  fittedSVs.begin(), svBuilder),
844  boost::make_transform_iterator(
845  fittedSVs.end(), svBuilder),
846  std::back_inserter(SVs),
847  SVFilter(vertexFilter, pv, jetDir));
848 
849  }else{
850  if( useSVClustering || useFatJets ) {
851  size_t jetIdx = ( iterJets - trackIPTagInfos->begin() );
852 
853  for(size_t iExtSv = 0; iExtSv < clusteredSVs.at(jetIdx).size(); iExtSv++){
854  const VTX & extVertex = (*extSecVertex)[ clusteredSVs.at(jetIdx).at(iExtSv) ];
855  if( extVertex.p4().M() < 0.3 )
856  continue;
857  extAssoCollection.push_back( extVertex );
858  }
859  }
860  else {
861  for(size_t iExtSv = 0; iExtSv < extSecVertex->size(); iExtSv++){
862  const VTX & extVertex = (*extSecVertex)[iExtSv];
863  if( Geom::deltaR2( ( position(extVertex) - pv.position() ), (extSVDeltaRToJet > 0) ? jetDir : -jetDir ) > extSVDeltaRToJet*extSVDeltaRToJet || extVertex.p4().M() < 0.3 )
864  continue;
865  extAssoCollection.push_back( extVertex );
866 
867  }
868  }
869  // build combined SV information and filter
870  SVBuilder svBuilder(pv, jetDir, withPVError, minTrackWeight);
871  std::remove_copy_if(boost::make_transform_iterator( extAssoCollection.begin(), svBuilder),
872  boost::make_transform_iterator(extAssoCollection.end(), svBuilder),
873  std::back_inserter(SVs),
874  SVFilter(vertexFilter, pv, jetDir));
875  }
876  // clean up now unneeded collections
877  gtPred.reset();
878  ghostTrack.reset();
879  gtStates.clear();
880  fitTracks.clear();
881  fittedSVs.clear();
882  extAssoCollection.clear();
883 
884  // sort SVs by importance
885 
886  std::vector<unsigned int> vtxIndices = vertexSorting(SVs);
887 
888  std::vector<typename TemplatedSecondaryVertexTagInfo<IPTI,VTX>::VertexData> svData;
889 
890  svData.resize(vtxIndices.size());
891  for(unsigned int idx = 0; idx < vtxIndices.size(); idx++) {
892  const SecondaryVertex &sv = SVs[vtxIndices[idx]];
893 
894  svData[idx].vertex = sv;
895  svData[idx].dist1d = sv.dist1d();
896  svData[idx].dist2d = sv.dist2d();
897  svData[idx].dist3d = sv.dist3d();
898  svData[idx].direction = flightDirection(pv,sv);
899  // mark tracks successfully used in vertex fit
900  markUsedTracks(trackData,trackRefs,sv,idx);
901  }
902 
903  // fill result into tag infos
904 
905  tagInfos->push_back(
907  trackData, svData, SVs.size(),
909  iterJets - trackIPTagInfos->begin())));
910  }
911 
912  event.put(std::move(tagInfos));
913 }
914 
915 //Need specialized template because reco::Vertex iterators are TrackBase and it is a mess to make general
916 template <>
918 {
919  for(Vertex::trackRef_iterator iter = sv.tracks_begin(); iter != sv.tracks_end(); ++iter) {
920  if (sv.trackWeight(*iter) < minTrackWeight)
921  continue;
922 
923  typename input_container::const_iterator pos =
924  std::find(trackRefs.begin(), trackRefs.end(),
925  iter->castTo<input_item>());
926 
927  if (pos == trackRefs.end() ) {
928  if(!useExternalSV)
929  throw cms::Exception("TrackNotFound")
930  << "Could not find track from secondary "
931  "vertex in original tracks."
932  << std::endl;
933  } else {
934  unsigned int index = pos - trackRefs.begin();
935  trackData[index].second.svStatus =
937  ((unsigned int)btag::TrackData::trackAssociatedToVertex + idx);
938  }
939  }
940 }
941 template <>
943 {
944  for(typename input_container::const_iterator iter = sv.daughterPtrVector().begin(); iter != sv.daughterPtrVector().end(); ++iter)
945  {
946  typename input_container::const_iterator pos =
947  std::find(trackRefs.begin(), trackRefs.end(), *iter);
948 
949  if (pos != trackRefs.end() )
950  {
951  unsigned int index = pos - trackRefs.begin();
952  trackData[index].second.svStatus =
954  ((unsigned int)btag::TrackData::trackAssociatedToVertex + idx);
955  }
956  }
957 }
958 
959 
960 template <>
963 {
964  if(sv.originalTracks().size()>0 && sv.originalTracks()[0].trackBaseRef().isNonnull())
965  return SecondaryVertex(pv, sv, direction, withPVError);
966  else
967  {
968  edm::LogError("UnexpectedInputs") << "Building from Candidates, should not happen!";
969  return SecondaryVertex(pv, sv, direction, withPVError);
970  }
971 }
972 
973 template <>
976 {
977  if(sv.originalTracks().size()>0 && sv.originalTracks()[0].trackBaseRef().isNonnull())
978  {
979  edm::LogError("UnexpectedInputs") << "Building from Tracks, should not happen!";
980  VertexCompositePtrCandidate vtxCompPtrCand;
981 
982  vtxCompPtrCand.setCovariance(sv.vertexState().error().matrix());
983  vtxCompPtrCand.setChi2AndNdof(sv.totalChiSquared(), sv.degreesOfFreedom());
984  vtxCompPtrCand.setVertex(Candidate::Point(sv.position().x(),sv.position().y(),sv.position().z()));
985 
986  return SecondaryVertex(pv, vtxCompPtrCand, direction, withPVError);
987  }
988  else
989  {
990  VertexCompositePtrCandidate vtxCompPtrCand;
991 
992  vtxCompPtrCand.setCovariance(sv.vertexState().error().matrix());
993  vtxCompPtrCand.setChi2AndNdof(sv.totalChiSquared(), sv.degreesOfFreedom());
994  vtxCompPtrCand.setVertex(Candidate::Point(sv.position().x(),sv.position().y(),sv.position().z()));
995 
997  for(std::vector<reco::TransientTrack>::const_iterator tt = sv.originalTracks().begin(); tt != sv.originalTracks().end(); ++tt)
998  {
999  if (sv.trackWeight(*tt) < minTrackWeight)
1000  continue;
1001 
1002  const CandidatePtrTransientTrack* cptt = dynamic_cast<const CandidatePtrTransientTrack*>(tt->basicTransientTrack());
1003  if ( cptt==0 )
1004  edm::LogError("DynamicCastingFailed") << "Casting of TransientTrack to CandidatePtrTransientTrack failed!";
1005  else
1006  {
1007  p4 += cptt->candidate()->p4();
1008  vtxCompPtrCand.addDaughter(cptt->candidate());
1009  }
1010  }
1011  vtxCompPtrCand.setP4(p4);
1012 
1013  return SecondaryVertex(pv, vtxCompPtrCand, direction, withPVError);
1014  }
1015 }
1016 
1017 // ------------ method that matches reclustered and original jets based on minimum dR ------------
1018 template<class IPTI,class VTX>
1019 template<class CONTAINER>
1021  const std::vector<fastjet::PseudoJet>& reclusteredJets,
1022  std::vector<int>& matchedIndices,
1023  const std::string& jetType)
1024 {
1025  std::string type = ( !jetType.empty() ? jetType + " " : jetType );
1026 
1027  std::vector<bool> matchedLocks(reclusteredJets.size(),false);
1028 
1029  for(size_t j=0; j<jets->size(); ++j)
1030  {
1031  double matchedDR2 = 1e9;
1032  int matchedIdx = -1;
1033 
1034  for(size_t rj=0; rj<reclusteredJets.size(); ++rj)
1035  {
1036  if( matchedLocks.at(rj) ) continue; // skip jets that have already been matched
1037 
1038  double tempDR2 = Geom::deltaR2( toJet(jets->at(j))->rapidity(), toJet(jets->at(j))->phi(), reclusteredJets.at(rj).rapidity(), reclusteredJets.at(rj).phi_std() );
1039  if( tempDR2 < matchedDR2 )
1040  {
1041  matchedDR2 = tempDR2;
1042  matchedIdx = rj;
1043  }
1044  }
1045 
1046  if( matchedIdx>=0 )
1047  {
1048  if ( matchedDR2 > rParam*rParam )
1049  {
1050  edm::LogError("JetMatchingFailed") << "Matched reclustered jet " << matchedIdx << " and original " << type << "jet " << j <<" are separated by dR=" << sqrt(matchedDR2) << " which is greater than the jet size R=" << rParam << ".\n"
1051  << "This is not expected so please check that the jet algorithm and jet size match those used for the original " << type << "jet collection.";
1052  }
1053  else
1054  matchedLocks.at(matchedIdx) = true;
1055  }
1056  else
1057  edm::LogError("JetMatchingFailed") << "Matching reclustered to original " << type << "jets failed. Please check that the jet algorithm and jet size match those used for the original " << type << "jet collection.";
1058 
1059  matchedIndices.push_back(matchedIdx);
1060  }
1061 }
1062 
1063 // ------------ method that matches groomed and original jets based on minimum dR ------------
1064 template<class IPTI,class VTX>
1066  const edm::Handle<edm::View<reco::Jet> >& groomedJets,
1067  std::vector<int>& matchedIndices)
1068 {
1069  std::vector<bool> jetLocks(jets->size(),false);
1070  std::vector<int> jetIndices;
1071 
1072  for(size_t gj=0; gj<groomedJets->size(); ++gj)
1073  {
1074  double matchedDR2 = 1e9;
1075  int matchedIdx = -1;
1076 
1077  if( groomedJets->at(gj).pt()>0. ) // skip pathological cases of groomed jets with Pt=0
1078  {
1079  for(size_t j=0; j<jets->size(); ++j)
1080  {
1081  if( jetLocks.at(j) ) continue; // skip jets that have already been matched
1082 
1083  double tempDR2 = Geom::deltaR2( jets->at(j).rapidity(), jets->at(j).phi(), groomedJets->at(gj).rapidity(), groomedJets->at(gj).phi() );
1084  if( tempDR2 < matchedDR2 )
1085  {
1086  matchedDR2 = tempDR2;
1087  matchedIdx = j;
1088  }
1089  }
1090  }
1091 
1092  if( matchedIdx>=0 )
1093  {
1094  if ( matchedDR2 > rParam*rParam )
1095  {
1096  edm::LogWarning("MatchedJetsFarApart") << "Matched groomed jet " << gj << " and original jet " << matchedIdx <<" are separated by dR=" << sqrt(matchedDR2) << " which is greater than the jet size R=" << rParam << ".\n"
1097  << "This is not expected so the matching of these two jets has been discarded. Please check that the two jet collections belong to each other.";
1098  matchedIdx = -1;
1099  }
1100  else
1101  jetLocks.at(matchedIdx) = true;
1102  }
1103  jetIndices.push_back(matchedIdx);
1104  }
1105 
1106  for(size_t j=0; j<jets->size(); ++j)
1107  {
1108  std::vector<int>::iterator matchedIndex = std::find( jetIndices.begin(), jetIndices.end(), j );
1109 
1110  matchedIndices.push_back( matchedIndex != jetIndices.end() ? std::distance(jetIndices.begin(),matchedIndex) : -1 );
1111  }
1112 }
1113 
1114 // ------------ method that matches subjets and original fat jets ------------
1115 template<class IPTI,class VTX>
1116 void TemplatedSecondaryVertexProducer<IPTI,VTX>::matchSubjets(const std::vector<int>& groomedIndices,
1117  const edm::Handle<edm::View<reco::Jet> >& groomedJets,
1118  const edm::Handle<std::vector<IPTI> >& subjets,
1119  std::vector<std::vector<int> >& matchedIndices)
1120 {
1121  for(size_t g=0; g<groomedIndices.size(); ++g)
1122  {
1123  std::vector<int> subjetIndices;
1124 
1125  if( groomedIndices.at(g)>=0 )
1126  {
1127  for(size_t s=0; s<groomedJets->at(groomedIndices.at(g)).numberOfDaughters(); ++s)
1128  {
1129  const edm::Ptr<reco::Candidate> & subjet = groomedJets->at(groomedIndices.at(g)).daughterPtr(s);
1130 
1131  for(size_t sj=0; sj<subjets->size(); ++sj)
1132  {
1133  const edm::RefToBase<reco::Jet> &subjetRef = subjets->at(sj).jet();
1134  if( subjet == edm::Ptr<reco::Candidate>( subjetRef.id(), subjetRef.get(), subjetRef.key() ) )
1135  {
1136  subjetIndices.push_back(sj);
1137  break;
1138  }
1139  }
1140  }
1141 
1142  if( subjetIndices.size() == 0 )
1143  edm::LogError("SubjetMatchingFailed") << "Matching subjets to original fat jets failed. Please check that the groomed fat jet and subjet collections belong to each other.";
1144 
1145  matchedIndices.push_back(subjetIndices);
1146  }
1147  else
1148  matchedIndices.push_back(subjetIndices);
1149  }
1150 }
1151 
1152 // ------------ method that matches subjets and original fat jets ------------
1153 template<class IPTI,class VTX>
1155  const edm::Handle<std::vector<IPTI> >& subjets,
1156  std::vector<std::vector<int> >& matchedIndices)
1157 {
1158  for(size_t fj=0; fj<fatJets->size(); ++fj)
1159  {
1160  std::vector<int> subjetIndices;
1161  size_t nSubjetCollections = 0;
1162  size_t nSubjets = 0;
1163 
1164  const pat::Jet * fatJet = dynamic_cast<const pat::Jet *>( fatJets->ptrAt(fj).get() );
1165 
1166  if( !fatJet )
1167  {
1168  if( fj==0 ) edm::LogError("WrongJetType") << "Wrong jet type for input fat jets. Please check that the input fat jets are of the pat::Jet type.";
1169 
1170  matchedIndices.push_back(subjetIndices);
1171  continue;
1172  }
1173  else
1174  {
1175  nSubjetCollections = fatJet->subjetCollectionNames().size();
1176 
1177  if( nSubjetCollections>0 )
1178  {
1179  for(size_t coll=0; coll<nSubjetCollections; ++coll)
1180  {
1181  const pat::JetPtrCollection & fatJetSubjets = fatJet->subjets(coll);
1182 
1183  for(size_t fjsj=0; fjsj<fatJetSubjets.size(); ++fjsj)
1184  {
1185  ++nSubjets;
1186 
1187  for(size_t sj=0; sj<subjets->size(); ++sj)
1188  {
1189  const pat::Jet * subJet = dynamic_cast<const pat::Jet *>( subjets->at(sj).jet().get() );
1190 
1191  if( !subJet )
1192  {
1193  if( fj==0 && coll==0 && fjsj==0 && sj==0 ) edm::LogError("WrongJetType") << "Wrong jet type for input subjets. Please check that the input subjets are of the pat::Jet type.";
1194 
1195  break;
1196  }
1197  else
1198  {
1199  if( subJet->originalObjectRef() == fatJetSubjets.at(fjsj)->originalObjectRef() )
1200  {
1201  subjetIndices.push_back(sj);
1202  break;
1203  }
1204  }
1205  }
1206  }
1207  }
1208 
1209  if( subjetIndices.size() == 0 && nSubjets > 0)
1210  edm::LogError("SubjetMatchingFailed") << "Matching subjets to fat jets failed. Please check that the fat jet and subjet collections belong to each other.";
1211 
1212  matchedIndices.push_back(subjetIndices);
1213  }
1214  else
1215  matchedIndices.push_back(subjetIndices);
1216  }
1217  }
1218 }
1219 
1220 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
1221 template<class IPTI,class VTX>
1223 
1225  desc.add<double>("extSVDeltaRToJet",0.3);
1226  desc.add<edm::InputTag>("beamSpotTag",edm::InputTag("offlineBeamSpot"));
1227  {
1229  vertexReco.add<double>("primcut",1.8);
1230  vertexReco.add<double>("seccut",6.0);
1231  vertexReco.add<std::string>("finder","avr");
1232  vertexReco.addOptionalNode( edm::ParameterDescription<double>("minweight",0.5, true) and
1233  edm::ParameterDescription<double>("weightthreshold",0.001, true) and
1234  edm::ParameterDescription<bool>("smoothing",false, true), true );
1235  vertexReco.addOptionalNode( edm::ParameterDescription<double>("maxFitChi2",10.0, true) and
1236  edm::ParameterDescription<double>("mergeThreshold",3.0, true) and
1237  edm::ParameterDescription<std::string>("fitType","RefitGhostTrackWithVertices", true), true );
1238  desc.add<edm::ParameterSetDescription>("vertexReco",vertexReco);
1239  }
1240  {
1242  vertexSelection.add<std::string>("sortCriterium","dist3dError");
1243  desc.add<edm::ParameterSetDescription>("vertexSelection",vertexSelection);
1244  }
1245  desc.add<std::string>("constraint","BeamSpot");
1246  desc.add<edm::InputTag>("trackIPTagInfos",edm::InputTag("impactParameterTagInfos"));
1247  {
1249  vertexCuts.add<double>("distSig3dMax",99999.9);
1250  vertexCuts.add<double>("fracPV",0.65);
1251  vertexCuts.add<double>("distVal2dMax",2.5);
1252  vertexCuts.add<bool>("useTrackWeights",true);
1253  vertexCuts.add<double>("maxDeltaRToJetAxis",0.4);
1254  {
1256  v0Filter.add<double>("k0sMassWindow",0.05);
1257  vertexCuts.add<edm::ParameterSetDescription>("v0Filter",v0Filter);
1258  }
1259  vertexCuts.add<double>("distSig2dMin",3.0);
1260  vertexCuts.add<unsigned int>("multiplicityMin",2);
1261  vertexCuts.add<double>("distVal2dMin",0.01);
1262  vertexCuts.add<double>("distSig2dMax",99999.9);
1263  vertexCuts.add<double>("distVal3dMax",99999.9);
1264  vertexCuts.add<double>("minimumTrackWeight",0.5);
1265  vertexCuts.add<double>("distVal3dMin",-99999.9);
1266  vertexCuts.add<double>("massMax",6.5);
1267  vertexCuts.add<double>("distSig3dMin",-99999.9);
1268  desc.add<edm::ParameterSetDescription>("vertexCuts",vertexCuts);
1269  }
1270  desc.add<bool>("useExternalSV",false);
1271  desc.add<double>("minimumTrackWeight",0.5);
1272  desc.add<bool>("usePVError",true);
1273  {
1275  trackSelection.add<double>("b_pT",0.3684);
1276  trackSelection.add<double>("max_pT",500);
1277  trackSelection.add<bool>("useVariableJTA",false);
1278  trackSelection.add<double>("maxDecayLen",99999.9);
1279  trackSelection.add<double>("sip3dValMin",-99999.9);
1280  trackSelection.add<double>("max_pT_dRcut",0.1);
1281  trackSelection.add<double>("a_pT",0.005263);
1282  trackSelection.add<unsigned int>("totalHitsMin",8);
1283  trackSelection.add<double>("jetDeltaRMax",0.3);
1284  trackSelection.add<double>("a_dR",-0.001053);
1285  trackSelection.add<double>("maxDistToAxis",0.2);
1286  trackSelection.add<double>("ptMin",1.0);
1287  trackSelection.add<std::string>("qualityClass","any");
1288  trackSelection.add<unsigned int>("pixelHitsMin",2);
1289  trackSelection.add<double>("sip2dValMax",99999.9);
1290  trackSelection.add<double>("max_pT_trackPTcut",3);
1291  trackSelection.add<double>("sip2dValMin",-99999.9);
1292  trackSelection.add<double>("normChi2Max",99999.9);
1293  trackSelection.add<double>("sip3dValMax",99999.9);
1294  trackSelection.add<double>("sip3dSigMin",-99999.9);
1295  trackSelection.add<double>("min_pT",120);
1296  trackSelection.add<double>("min_pT_dRcut",0.5);
1297  trackSelection.add<double>("sip2dSigMax",99999.9);
1298  trackSelection.add<double>("sip3dSigMax",99999.9);
1299  trackSelection.add<double>("sip2dSigMin",-99999.9);
1300  trackSelection.add<double>("b_dR",0.6263);
1301  desc.add<edm::ParameterSetDescription>("trackSelection",trackSelection);
1302  }
1303  desc.add<std::string>("trackSort","sip3dSig");
1304  desc.add<edm::InputTag>("extSVCollection",edm::InputTag("secondaryVertices"));
1305  desc.addOptionalNode( edm::ParameterDescription<bool>("useSVClustering",false, true) and
1306  edm::ParameterDescription<std::string>("jetAlgorithm", true) and
1307  edm::ParameterDescription<double>("rParam", true), true );
1308  desc.addOptional<bool>("useSVMomentum",false);
1309  desc.addOptional<double>("ghostRescaling",1e-18);
1310  desc.addOptional<double>("relPtTolerance",1e-03);
1311  desc.addOptional<edm::InputTag>("fatJets");
1312  desc.addOptional<edm::InputTag>("groomedFatJets");
1313  descriptions.addDefault(desc);
1314 }
1315 
1316 //define this as a plug-in
1319 
1322 
type
Definition: HCALResponse.h:21
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
reco::Vertex::Point convertPos(const GlobalPoint &p)
value_type const * get() const
Definition: RefToBase.h:234
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::View< reco::Jet > > token_fatJets
boost::shared_ptr< fastjet::JetDefinition > JetDefPtr
const math::XYZPoint & position(const reco::Vertex &sv)
reco::btag::SortCriteria getCriterium(const std::string &name)
Definition: TrackSorting.cc:11
boost::shared_ptr< fastjet::JetDefinition > JetDefPtr
virtual double pz() const =0
z coordinate of momentum vector
TemplatedSecondaryVertex< reco::Vertex > SecondaryVertex
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:81
TemplatedSecondaryVertexProducer(const edm::ParameterSet &params)
virtual void produce(edm::Event &event, const edm::EventSetup &es) override
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< reco::BeamSpot > token_BeamSpot
const AlgebraicSymMatrix33 matrix() const
double y() const
y coordinate
Definition: Vertex.h:113
Base class for all types of Jets.
Definition: Jet.h:20
float totalChiSquared() const
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
reco::TransientTrack build(const reco::Track *p) const
T y() const
Definition: PV3DBase.h:63
bool exists(std::string const &parameterName) const
checks if a parameter exists
void setWeight(double weight)
SVFilter(const VertexFilter &filter, const Vertex &pv, const GlobalVector &direction)
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:130
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
const reco::Track * toTrack(const reco::TrackBaseRef &t)
Definition: IPTagInfo.h:24
const Point & position() const
position
Definition: Vertex.h:109
ProductID id() const
Definition: RefToBase.h:242
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
TemplatedSecondaryVertexProducer< TrackIPTagInfo, reco::Vertex > SecondaryVertexProducer
void setVertex(const Point &vertex) override
set vertex
static GhostTrackVertexFinder::FitType getGhostTrackFitType(const std::string &name)
void matchReclusteredJets(const edm::Handle< CONTAINER > &jets, const std::vector< fastjet::PseudoJet > &matchedJets, std::vector< int > &matchedIndices, const std::string &jetType="")
double dydz() const
dydz slope
Definition: BeamSpot.h:84
edm::EDGetTokenT< std::vector< IPTI > > token_trackIPTagInfo
IPTI::input_container::value_type input_item
virtual double energy() const =0
energy
virtual double py() const =0
y coordinate of momentum vector
T mag() const
Definition: PV3DBase.h:67
void addDefault(ParameterSetDescription const &psetDescription)
GlobalVector flightDirection(const reco::Vertex &pv, const reco::Vertex &sv)
std::vector< reco::TransientTrack > const & originalTracks() const
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
float degreesOfFreedom() const
size_t key() const
Definition: RefToBase.h:250
T sqrt(T t)
Definition: SSEVec.h:18
double p4[4]
Definition: TauolaWrapper.h:92
GlobalPoint position() const
vector< PseudoJet > jets
std::vector< edm::Ptr< pat::Jet > > JetPtrCollection
Definition: Jet.h:77
T z() const
Definition: PV3DBase.h:64
const Point & vertex() const override
vertex position (overwritten by PF...)
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
def pv(vc)
Definition: MetAnalyzer.py:6
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double z() const
z coordinate
Definition: Vertex.h:115
ParameterDescriptionNode * addOptionalNode(ParameterDescriptionNode const &node, bool writeToCfi)
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:86
VertexSorting< SecondaryVertex > vertexSorting
std::vector< std::string > const & subjetCollectionNames() const
Subjet collection names.
Definition: Jet.h:518
static ConstraintType getConstraintType(const std::string &name)
const edm::Ptr< reco::Candidate > & originalObjectRef() const
reference to original object. Returns a null reference if not available
Definition: PATObject.h:500
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
float trackWeight(const reco::TransientTrack &track) const
CandidatePtr candidate() const override
void matchGroomedJets(const edm::Handle< edm::View< reco::Jet > > &jets, const edm::Handle< edm::View< reco::Jet > > &matchedJets, std::vector< int > &matchedIndices)
void matchSubjets(const std::vector< int > &groomedIndices, const edm::Handle< edm::View< reco::Jet > > &groomedJets, const edm::Handle< std::vector< IPTI > > &subjets, std::vector< std::vector< int > > &matchedIndices)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double dxdz() const
dxdz slope
Definition: BeamSpot.h:82
void markUsedTracks(TrackDataVector &trackData, const input_container &trackRefs, const SecondaryVertex &sv, size_t idx)
edm::EDGetTokenT< edm::View< VTX > > token_extSVCollection
Vector3DBase unit() const
Definition: Vector3DBase.h:57
SVBuilder(const reco::Vertex &pv, const GlobalVector &direction, bool withPVError, double minTrackWeight)
double x() const
x coordinate
Definition: Vertex.h:111
JetCorrectorParametersCollection coll
Definition: classes.h:10
TemplatedSecondaryVertex< VTX > SecondaryVertex
edm::EDGetTokenT< edm::View< reco::Jet > > token_groomedFatJets
virtual double pt() const =0
transverse momentum
bool linearize(const GhostTrackPrediction &pred, bool initial=false, double lambda=0.)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const T & get() const
Definition: EventSetup.h:55
double sigmaZ() const
sigma z
Definition: BeamSpot.h:80
Analysis-level calorimeter jet class.
Definition: Jet.h:80
SecondaryVertex operator()(const TransientVertex &sv) const
Error error() const
return SMatrix
Definition: Vertex.h:139
TemplatedSecondaryVertexProducer< CandIPTagInfo, reco::VertexCompositePtrCandidate > CandSecondaryVertexProducer
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
double covariance(int i, int j) const
(i,j)-th element of error matrix
Definition: BeamSpot.h:112
std::vector< reco::btag::IndexedTrackData > TrackDataVector
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
std::vector< TemplatedSecondaryVertexTagInfo< IPTI, VTX > > Product
void setCovariance(const CovarianceMatrix &m)
set covariance matrix
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
void addDaughter(const CandidatePtr &)
add a daughter via a reference
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
GlobalError error() const
Definition: VertexState.h:74
dbl *** dir
Definition: mlp_gen.cc:35
void setChi2AndNdof(double chi2, double ndof)
set chi2 and ndof
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76
virtual double px() const =0
x coordinate of momentum vector
VertexState const & vertexState() const
T x() const
Definition: PV3DBase.h:62
pat::JetPtrCollection const & subjets(unsigned int index=0) const
Access to subjet list.
void setP4(const LorentzVector &p4) final
set 4-momentum
const reco::Jet * toJet(const reco::Jet &j)
def move(src, dest)
Definition: eostools.py:510
std::pair< unsigned int, TrackData > IndexedTrackData
Definition: event.py:1
Global3DVector GlobalVector
Definition: GlobalVector.h:10