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.
4 
5 /* class PFRecoTauDiscriminationByFlightPathSignificance
6  * created : August 30 2010,
7  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
8  * based on H+ tau ID by Lauri Wendland
9  */
10 
18 
19 #include "TLorentzVector.h"
20 
21 using namespace reco;
22 using namespace std;
23 using namespace edm;
24 
27  public:
29  const ParameterSet& iConfig)
31  vertexAssociator_(iConfig.getParameter<ParameterSet>("qualityCuts")) {
32  flightPathSig = iConfig.getParameter<double>("flightPathSig");
33  withPVError = iConfig.getParameter<bool>("UsePVerror");
34  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
35  }
36 
38 
39  void beginEvent(const edm::Event&, const edm::EventSetup&);
40  double discriminate(const reco::PFTauRef&);
41 
42  private:
43  double threeProngFlightPathSig(const PFTauRef&);
44  double vertexSignificance(reco::Vertex&,reco::Vertex&,GlobalVector&);
45 
47 
49  double flightPathSig;
51 
53 };
54 
56  const Event& iEvent, const EventSetup& iSetup){
57 
58  vertexAssociator_.setEvent(iEvent);
59 
60  // Transient Tracks
62  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",builder);
63  transientTrackBuilder = builder.product();
64 
65 }
66 
68 
69  if(booleanOutput) return ( threeProngFlightPathSig(tau) > flightPathSig ? 1. : 0. );
70  return threeProngFlightPathSig(tau);
71 }
72 
74  const PFTauRef& tau){
75  double flightPathSignificance = 0;
76 
77  reco::VertexRef primaryVertex = vertexAssociator_.associatedVertex(*tau);
78 
79  if (primaryVertex.isNull()) {
80  edm::LogError("FlightPathSignficance") << "Could not get vertex associated"
81  << " to tau, returning -999!" << std::endl;
82  return -999;
83  }
84 
85  //Secondary vertex
86  const PFCandidateRefVector pfSignalCandidates = tau->signalPFChargedHadrCands();
87  vector<TransientTrack> transientTracks;
89  for(iTrack = pfSignalCandidates.begin(); iTrack!= pfSignalCandidates.end(); iTrack++){
90  const PFCandidate& pfCand = *(iTrack->get());
91  if(pfCand.trackRef().isNonnull()){
92  const TransientTrack transientTrack = transientTrackBuilder->build(pfCand.trackRef());
93  transientTracks.push_back(transientTrack);
94  }
95  }
96  if(transientTracks.size() > 1){
97  KalmanVertexFitter kvf(true);
98  TransientVertex tv = kvf.vertex(transientTracks);
99 
100  if(tv.isValid()){
101  GlobalVector tauDir(tau->px(),
102  tau->py(),
103  tau->pz());
104  Vertex secVer = tv;
105  // We have to un-const the PV for some reason
106  reco::Vertex primaryVertexNonConst = *primaryVertex;
107  flightPathSignificance = vertexSignificance(primaryVertexNonConst,secVer,tauDir);
108  }
109  }
110  return flightPathSignificance;
111 }
112 
114  reco::Vertex& pv, Vertex& sv,GlobalVector& direction){
115  return SecondaryVertex::computeDist3d(pv,sv,direction,withPVError).significance();
116 }
117 
119 
T getParameter(std::string const &) const
double vertexSignificance(reco::Vertex &, reco::Vertex &, GlobalVector &)
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:249
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:244
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:339
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
int iEvent
Definition: GenABIO.cc:243
bool isNull() const
Checks for null.
Definition: Ref.h:247
static Measurement1D computeDist3d(const reco::Vertex &pv, const reco::Vertex &sv, const GlobalVector &direction, bool withPVError)
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:62
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:33
bool isValid() const