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 
17 
18 #include "TLorentzVector.h"
19 
20 using namespace reco;
21 using namespace std;
22 using namespace edm;
23 
26  public:
28  qualityCuts_(iConfig.getParameter<ParameterSet>("qualityCuts")){ // retrieve quality cuts
29  flightPathSig = iConfig.getParameter<double>("flightPathSig");
30  withPVError = iConfig.getParameter<bool>("UsePVerror");
31 
32  PVProducer = iConfig.getParameter<edm::InputTag>("PVProducer");
33 
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 
48  double flightPathSig;
50 
54 
56 };
57 
59 
60  //Primary vertex
62  iEvent.getByLabel(PVProducer, vertexHandle);
63  const edm::View<reco::Vertex>& vertexCollection(*vertexHandle);
64 
65  primaryVertex = *(vertexCollection.begin());
66 
67  // Transient Tracks
69  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",builder);
70  transientTrackBuilder = builder.product();
71 
72 }
73 
75 
76  if(booleanOutput) return ( threeProngFlightPathSig(tau) > flightPathSig ? 1. : 0. );
77  return threeProngFlightPathSig(tau);
78 }
79 
81  const PFTauRef& tau){
82  double flightPathSignificance = 0;
83 
84  //Secondary vertex
85  const PFCandidateRefVector pfSignalCandidates = tau->signalPFChargedHadrCands();
86  vector<TransientTrack> transientTracks;
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  flightPathSignificance = vertexSignificance(primaryVertex,secVer,tauDir);
105  }
106  }
107  return flightPathSignificance;
108 }
109 
111  reco::Vertex& pv, Vertex& sv,GlobalVector& direction){
112  return SecondaryVertex::computeDist3d(pv,sv,direction,withPVError).significance();
113 }
114 
116 
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:243
tuple vertexCollection
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:238
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:331
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:249
int iEvent
Definition: GenABIO.cc:243
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
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
const_iterator begin() const
bool isValid() const