CMS 3D CMS Logo

CaloRecoTauDiscriminationByFlightPathSignificance.cc
Go to the documentation of this file.
3 
6 
7 /* class CaloRecoTauDiscriminationByFlightPathSignificance
8  * created : September 23 2010,
9  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
10  * based on H+ tau ID by Lauri Wendland
11  */
12 
19 
20 #include "TLorentzVector.h"
21 
22 namespace {
23 
24 using namespace reco;
25 using namespace std;
26 
27 class CaloRecoTauDiscriminationByFlightPathSignificance final : public CaloTauDiscriminationProducerBase {
28  public:
29  explicit CaloRecoTauDiscriminationByFlightPathSignificance(
30  const edm::ParameterSet& iConfig)
32  flightPathSig = iConfig.getParameter<double>("flightPathSig");
33  withPVError = iConfig.getParameter<bool>("UsePVerror");
34 
35  PVProducer = iConfig.getParameter<edm::InputTag>("PVProducer");
36 
37  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
38  }
39  ~CaloRecoTauDiscriminationByFlightPathSignificance() override{}
40  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
41  double discriminate(const reco::CaloTauRef&) const override;
42 
43  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
44  private:
45  double threeProngFlightPathSig(const CaloTauRef&) const ;
46  double vertexSignificance(reco::Vertex const&,reco::Vertex const &,GlobalVector const &) const ;
47 
48  double flightPathSig;
49  bool withPVError;
50 
52  const TransientTrackBuilder* transientTrackBuilder;
54 
55  bool booleanOutput;
56 };
57 
58 void CaloRecoTauDiscriminationByFlightPathSignificance::beginEvent(
59  const edm::Event& iEvent, const edm::EventSetup& iSetup){
60  //Primary vertex
62  iEvent.getByLabel(PVProducer, vertexHandle);
63  const edm::View<reco::Vertex>& vertexCollection(*vertexHandle);
64  primaryVertex = *(vertexCollection.begin());
65  // Transient Tracks
67  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",builder);
68  transientTrackBuilder = builder.product();
69 }
70 
71 double
72 CaloRecoTauDiscriminationByFlightPathSignificance::discriminate(
73  const CaloTauRef& tau) const {
74  if(booleanOutput)
75  return ( threeProngFlightPathSig(tau) > flightPathSig ? 1. : 0. );
76  return threeProngFlightPathSig(tau);
77 }
78 
79 double
80 CaloRecoTauDiscriminationByFlightPathSignificance::threeProngFlightPathSig(
81  const CaloTauRef& tau) const {
82  double flightPathSignificance = 0;
83  //Secondary vertex
84  reco::TrackRefVector signalTracks = tau->signalTracks();
85  vector<TransientTrack> transientTracks;
86  for(size_t i = 0; i < signalTracks.size(); ++i){
87  const TransientTrack transientTrack =
88  transientTrackBuilder->build(signalTracks[i]);
89  transientTracks.push_back(transientTrack);
90  }
91  if(transientTracks.size() > 1) {
92  KalmanVertexFitter kvf(true);
93  TransientVertex tv = kvf.vertex(transientTracks);
94  if(tv.isValid()){
95  GlobalVector tauDir(tau->px(), tau->py(), tau->pz());
96  Vertex secVer = tv;
97  flightPathSignificance = vertexSignificance(primaryVertex,secVer,tauDir);
98  }
99  }
100  return flightPathSignificance;
101 }
102 
103 double
104 CaloRecoTauDiscriminationByFlightPathSignificance::vertexSignificance(
105  reco::Vertex const & pv, Vertex const & sv,GlobalVector const & direction) const {
107  pv,sv,direction,withPVError).significance();
108 }
109 
110 }
111 
112 void
114  // caloRecoTauDiscriminationByFlightPathSignificance
116  desc.add<edm::InputTag>("CaloTauProducer", edm::InputTag("caloRecoTauProducer"));
117  desc.add<double>("flightPathSig", 1.5);
118  desc.add<edm::InputTag>("PVProducer", edm::InputTag("offlinePrimaryVertices"));
119  desc.add<bool>("BooleanOutput", true);
120 
121  edm::ParameterSetDescription pset_signalQualityCuts;
122  pset_signalQualityCuts.add<double>("maxDeltaZ", 0.4);
123  pset_signalQualityCuts.add<double>("minTrackPt", 0.5);
124  pset_signalQualityCuts.add<double>("minTrackVertexWeight", -1.0);
125  pset_signalQualityCuts.add<double>("maxTrackChi2", 100.0);
126  pset_signalQualityCuts.add<unsigned int>("minTrackPixelHits", 0);
127  pset_signalQualityCuts.add<double>("minGammaEt", 1.0);
128  pset_signalQualityCuts.add<unsigned int>("minTrackHits", 3);
129  pset_signalQualityCuts.add<double>("minNeutralHadronEt", 30.0);
130  pset_signalQualityCuts.add<double>("maxTransverseImpactParameter", 0.1);
131  pset_signalQualityCuts.addOptional<bool>("useTracksInsteadOfPFHadrons");
132 
133  edm::ParameterSetDescription pset_vxAssocQualityCuts;
134  pset_vxAssocQualityCuts.add<double>("minTrackPt", 0.5);
135  pset_vxAssocQualityCuts.add<double>("minTrackVertexWeight", -1.0);
136  pset_vxAssocQualityCuts.add<double>("maxTrackChi2", 100.0);
137  pset_vxAssocQualityCuts.add<unsigned int>("minTrackPixelHits", 0);
138  pset_vxAssocQualityCuts.add<double>("minGammaEt", 1.0);
139  pset_vxAssocQualityCuts.add<unsigned int>("minTrackHits", 3);
140  pset_vxAssocQualityCuts.add<double>("maxTransverseImpactParameter", 0.1);
141  pset_vxAssocQualityCuts.addOptional<bool>("useTracksInsteadOfPFHadrons");
142 
143  edm::ParameterSetDescription pset_isolationQualityCuts;
144  pset_isolationQualityCuts.add<double>("maxDeltaZ", 0.2);
145  pset_isolationQualityCuts.add<double>("minTrackPt", 1.0);
146  pset_isolationQualityCuts.add<double>("minTrackVertexWeight", -1.0);
147  pset_isolationQualityCuts.add<double>("maxTrackChi2", 100.0);
148  pset_isolationQualityCuts.add<unsigned int>("minTrackPixelHits", 0);
149  pset_isolationQualityCuts.add<double>("minGammaEt", 1.5);
150  pset_isolationQualityCuts.add<unsigned int>("minTrackHits", 8);
151  pset_isolationQualityCuts.add<double>("maxTransverseImpactParameter", 0.03);
152  pset_isolationQualityCuts.addOptional<bool>("useTracksInsteadOfPFHadrons");
153 
154  edm::ParameterSetDescription pset_qualityCuts;
155  pset_qualityCuts.add<edm::ParameterSetDescription>("signalQualityCuts", pset_signalQualityCuts);
156  pset_qualityCuts.add<edm::ParameterSetDescription>("vxAssocQualityCuts", pset_vxAssocQualityCuts);
157  pset_qualityCuts.add<edm::ParameterSetDescription>("isolationQualityCuts", pset_isolationQualityCuts);
158  pset_qualityCuts.add<std::string>("leadingTrkOrPFCandOption", "leadPFCand");
159  pset_qualityCuts.add<std::string>("pvFindingAlgo", "closestInDeltaZ");
160  pset_qualityCuts.add<edm::InputTag>("primaryVertexSrc", edm::InputTag("offlinePrimaryVertices"));
161  pset_qualityCuts.add<bool>("vertexTrackFiltering", false);
162  pset_qualityCuts.add<bool>("recoverLeadingTrk", false);
163 
164  desc.add<edm::ParameterSetDescription>("qualityCuts", pset_qualityCuts);
165 
166  {
168  psd0.add<std::string>("BooleanOperator", "and");
169  {
171  psd1.add<double>("cut");
172  psd1.add<edm::InputTag>("Producer");
173  psd0.addOptional<edm::ParameterSetDescription>("leadTrack", psd1);
174  }
175  desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
176  }
177 
178  desc.add<bool>("UsePVerror", true);
179  descriptions.add("caloRecoTauDiscriminationByFlightPathSignificance", desc);
180 }
181 DEFINE_FWK_MODULE(CaloRecoTauDiscriminationByFlightPathSignificance);
182 
T getParameter(std::string const &) const
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
static Measurement1D computeDist3d(const reco::Vertex &pv, const SV &sv, const GlobalVector &direction, bool withPVError)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
def pv(vc)
Definition: MetAnalyzer.py:7
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:480
double significance() const
Definition: Measurement1D.h:29
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
fixed size matrix
T get() const
Definition: EventSetup.h:71
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
T const * product() const
Definition: ESHandle.h:86
bool isValid() const