CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFRecoTauDiscriminationByFlightPathSignificance.cc
Go to the documentation of this file.
3 
4 /* class PFRecoTauDiscriminationByFlightPathSignificance
5  * created : August 30 2010,
6  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
7  * based on H+ tau ID by Lauri Wendland
8  */
9 
17 
18 #include "TLorentzVector.h"
19 
20 using namespace reco;
21 using namespace std;
22 using namespace edm;
23 
26  public:
29  flightPathSig = iConfig.getParameter<double>("flightPathSig");
30  withPVError = iConfig.getParameter<bool>("UsePVerror");
31  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
32  // edm::ConsumesCollector iC(consumesCollector());
33  vertexAssociator_ = new reco::tau::RecoTauVertexAssociator(iConfig.getParameter<ParameterSet>("qualityCuts"),consumesCollector());
34  }
35 
37 
38  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
39  double discriminate(const reco::PFTauRef&) const override;
40 
41  private:
42  double threeProngFlightPathSig(const PFTauRef&) const ;
43  double vertexSignificance(reco::Vertex const&,reco::Vertex const &,GlobalVector const&) const;
44 
46 
48  double flightPathSig;
50 
52 };
53 
55  const Event& iEvent, const EventSetup& iSetup){
56 
57  vertexAssociator_->setEvent(iEvent);
58 
59  // Transient Tracks
61  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",builder);
62  transientTrackBuilder = builder.product();
63 
64 }
65 
67 
68  if(booleanOutput) return ( threeProngFlightPathSig(tau) > flightPathSig ? 1. : 0. );
69  return threeProngFlightPathSig(tau);
70 }
71 
73  const PFTauRef& tau) const {
74  double flightPathSignificance = 0;
75 
76  reco::VertexRef primaryVertex = vertexAssociator_->associatedVertex(*tau);
77 
78  if (primaryVertex.isNull()) {
79  edm::LogError("FlightPathSignficance") << "Could not get vertex associated"
80  << " to tau, returning -999!" << std::endl;
81  return -999;
82  }
83 
84  //Secondary vertex
85  const vector<PFCandidatePtr>& pfSignalCandidates = tau->signalPFChargedHadrCands();
86  vector<TransientTrack> transientTracks;
87  vector<PFCandidatePtr>::const_iterator iTrack;
88  for(iTrack = pfSignalCandidates.begin(); iTrack!= pfSignalCandidates.end(); iTrack++){
89  const PFCandidate& pfCand = *(iTrack->get());
90  if(pfCand.trackRef().isNonnull()){
91  const TransientTrack transientTrack = transientTrackBuilder->build(pfCand.trackRef());
92  transientTracks.push_back(transientTrack);
93  }
94  }
95  if(transientTracks.size() > 1){
96  KalmanVertexFitter kvf(true);
97  TransientVertex tv = kvf.vertex(transientTracks);
98 
99  if(tv.isValid()){
100  GlobalVector tauDir(tau->px(),
101  tau->py(),
102  tau->pz());
103  Vertex secVer = tv;
104  // We have to un-const the PV for some reason
105  reco::Vertex primaryVertexNonConst = *primaryVertex;
106  flightPathSignificance = vertexSignificance(primaryVertexNonConst,secVer,tauDir);
107  }
108  }
109  return flightPathSignificance;
110 }
111 
113  reco::Vertex const & pv, Vertex const & sv,GlobalVector const & direction) const {
114  return SecondaryVertex::computeDist3d(pv,sv,direction,withPVError).significance();
115 }
116 
118 
T getParameter(std::string const &) const
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
void beginEvent(const edm::Event &, const edm::EventSetup &) override
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
static Measurement1D computeDist3d(const reco::Vertex &pv, const SV &sv, const GlobalVector &direction, bool withPVError)
double vertexSignificance(reco::Vertex const &, reco::Vertex const &, GlobalVector const &) const
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:433
int iEvent
Definition: GenABIO.cc:230
bool isNull() const
Checks for null.
Definition: Ref.h:247
DEFINE_FWK_MODULE(CosmicTrackingParticleSelector)
double significance() const
Definition: Measurement1D.h:32
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
bool isValid() const