CMS 3D CMS Logo

RecoTauImpactParameterSignificancePlugin.cc
Go to the documentation of this file.
1 /*
2  * =============================================================================
3  * Filename: RecoTauImpactParameterSignificancePlugin.cc
4  *
5  * Description: Add the IP significance of the lead track w.r.t to the PV.
6  * to a PFTau.
7  * Created: 10/31/2010 13:32:14
8  *
9  * Authors: Evan K. Friis (UC Davis), evan.klose.friis@cern.ch,
10  * Simone Gennai, Ludovic Houchu
11  *
12  * =============================================================================
13  */
14 
20 
22 
26 
27 namespace reco {
28  namespace tau {
29 
31  public:
34  void operator()(PFTau& tau) const override;
35  void beginEvent() override;
36 
37  private:
41  };
42 
46  transTrackBuilderToken_(iC.esConsumes(edm::ESInputTag{"", "TransientTrackBuilder"})),
47  vertexAssociator_(pset.getParameter<edm::ParameterSet>("qualityCuts"), std::move(iC)) {}
48 
51  // Get tranisent track builder.
53  }
54 
55  namespace {
56  inline const reco::Track* getTrack(const Candidate& cand) {
57  const PFCandidate* pfCandPtr = dynamic_cast<const PFCandidate*>(&cand);
58  if (pfCandPtr) {
59  if (pfCandPtr->trackRef().isNonnull())
60  return pfCandPtr->trackRef().get();
61  else
62  return nullptr;
63  }
64 
65  const pat::PackedCandidate* packedCand = dynamic_cast<const pat::PackedCandidate*>(&cand);
66  if (packedCand && packedCand->hasTrackDetails())
67  return &packedCand->pseudoTrack();
68 
69  return nullptr;
70  }
71  } // namespace
72 
74  // Get the transient lead track
75  if (tau.leadChargedHadrCand().isNonnull()) {
76  const reco::Track* leadTrack = getTrack(*tau.leadChargedHadrCand());
77  if (leadTrack != nullptr) {
79  GlobalVector direction(tau.jetRef()->px(), tau.jetRef()->py(), tau.jetRef()->pz());
81  // Compute the significance
82  std::pair<bool, Measurement1D> ipsig = IPTools::signedImpactParameter3D(track, direction, *pv);
83  if (ipsig.first)
84  tau.setleadPFChargedHadrCandsignedSipt(ipsig.second.significance());
85  }
86  }
87  }
88 
89  } // namespace tau
90 } // namespace reco
94  "RecoTauImpactParameterSignificancePlugin");
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::EventSetup * evtSetup() const
std::pair< bool, Measurement1D > signedImpactParameter3D(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:81
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:232
void setEvent(const edm::Event &evt)
Load the vertices from the event.
reco::VertexRef associatedVertex(const Jet &jet) const
reco::TransientTrack build(const reco::Track *p) const
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transTrackBuilderToken_
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalVector
vector in glovbal coordinate system
Definition: Vector3D.h:28
static const TrackGhostTrackState * getTrack(const BasicGhostTrackState *basic)
RecoTauImpactParameterSignificancePlugin(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
fixed size matrix
HLT enums.
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:226
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:437
#define DEFINE_EDM_PLUGIN(factory, type, name)
def move(src, dest)
Definition: eostools.py:511
virtual const reco::Track & pseudoTrack() const
bool hasTrackDetails() const
Return true if a bestTrack can be extracted from this Candidate.