CMS 3D CMS Logo

PFRecoTauDiscriminationByTauPolarization.cc
Go to the documentation of this file.
4 
5 /* class PFRecoTauDiscriminationByTauPolarization
6  * created : May 26 2010,
7  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
8  */
9 
10 using namespace reco;
11 using namespace std;
12 using namespace edm;
13 
16  public:
18  const ParameterSet& iConfig)
19  :PFTauDiscriminationProducerBase(iConfig) { // retrieve quality cuts
20  rTauMin = iConfig.getParameter<double>("rtau");
21  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
22  }
23 
25 
26  void beginEvent(const Event&, const EventSetup&) override;
27  double discriminate(const PFTauRef&) const override;
28 
29  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
30  private:
32  double rTauMin;
33 };
34 
36  const Event& event, const EventSetup& eventSetup){}
37 
38 double
40 
41  double rTau = 0;
42  // rtau for PFTau has to be calculated for leading PF charged hadronic candidate
43  // calculating it from leadingTrack can (and will) give rtau > 1!
44  if(tau.isNonnull() && tau->p() > 0
45  && tau->leadChargedHadrCand().isNonnull()) {
46  rTau = tau->leadChargedHadrCand()->p()/tau->p();
47  }
48 
49  if(booleanOutput) return ( rTau > rTauMin ? 1. : 0. );
50  return rTau;
51 }
52 
53 void
55  // pfRecoTauDiscriminationByTauPolarization
57  desc.add<double>("rtau", 0.8);
58  desc.add<edm::InputTag>("PVProducer", edm::InputTag("offlinePrimaryVertices"));
59  desc.add<bool>("BooleanOutput", true);
60  desc.add<edm::InputTag>("PFTauProducer", edm::InputTag("pfRecoTauProducer"));
61 
62  {
63  edm::ParameterSetDescription pset_signalQualityCuts;
64  pset_signalQualityCuts.add<double>("maxDeltaZ", 0.4);
65  pset_signalQualityCuts.add<double>("minTrackPt", 0.5);
66  pset_signalQualityCuts.add<double>("minTrackVertexWeight", -1.0);
67  pset_signalQualityCuts.add<double>("maxTrackChi2", 100.0);
68  pset_signalQualityCuts.add<unsigned int>("minTrackPixelHits", 0);
69  pset_signalQualityCuts.add<double>("minGammaEt", 1.0);
70  pset_signalQualityCuts.add<unsigned int>("minTrackHits", 3);
71  pset_signalQualityCuts.add<double>("minNeutralHadronEt", 30.0);
72  pset_signalQualityCuts.add<double>("maxTransverseImpactParameter", 0.1);
73  pset_signalQualityCuts.addOptional<bool>("useTracksInsteadOfPFHadrons");
74 
75  edm::ParameterSetDescription pset_vxAssocQualityCuts;
76  pset_vxAssocQualityCuts.add<double>("minTrackPt", 0.5);
77  pset_vxAssocQualityCuts.add<double>("minTrackVertexWeight", -1.0);
78  pset_vxAssocQualityCuts.add<double>("maxTrackChi2", 100.0);
79  pset_vxAssocQualityCuts.add<unsigned int>("minTrackPixelHits", 0);
80  pset_vxAssocQualityCuts.add<double>("minGammaEt", 1.0);
81  pset_vxAssocQualityCuts.add<unsigned int>("minTrackHits", 3);
82  pset_vxAssocQualityCuts.add<double>("maxTransverseImpactParameter", 0.1);
83  pset_vxAssocQualityCuts.addOptional<bool>("useTracksInsteadOfPFHadrons");
84 
85  edm::ParameterSetDescription pset_isolationQualityCuts;
86  pset_isolationQualityCuts.add<double>("maxDeltaZ", 0.2);
87  pset_isolationQualityCuts.add<double>("minTrackPt", 1.0);
88  pset_isolationQualityCuts.add<double>("minTrackVertexWeight", -1.0);
89  pset_isolationQualityCuts.add<double>("maxTrackChi2", 100.0);
90  pset_isolationQualityCuts.add<unsigned int>("minTrackPixelHits", 0);
91  pset_isolationQualityCuts.add<double>("minGammaEt", 1.5);
92  pset_isolationQualityCuts.add<unsigned int>("minTrackHits", 8);
93  pset_isolationQualityCuts.add<double>("maxTransverseImpactParameter", 0.03);
94  pset_isolationQualityCuts.addOptional<bool>("useTracksInsteadOfPFHadrons");
95 
96  edm::ParameterSetDescription pset_qualityCuts;
97  pset_qualityCuts.add<edm::ParameterSetDescription>("signalQualityCuts", pset_signalQualityCuts);
98  pset_qualityCuts.add<edm::ParameterSetDescription>("vxAssocQualityCuts", pset_vxAssocQualityCuts);
99  pset_qualityCuts.add<edm::ParameterSetDescription>("isolationQualityCuts", pset_isolationQualityCuts);
100  pset_qualityCuts.add<std::string>("leadingTrkOrPFCandOption", "leadPFCand");
101  pset_qualityCuts.add<std::string>("pvFindingAlgo", "closestInDeltaZ");
102  pset_qualityCuts.add<edm::InputTag>("primaryVertexSrc", edm::InputTag("offlinePrimaryVertices"));
103  pset_qualityCuts.add<bool>("vertexTrackFiltering", false);
104  pset_qualityCuts.add<bool>("recoverLeadingTrk", false);
105 
106  desc.add<edm::ParameterSetDescription>("qualityCuts", pset_qualityCuts);
107  }
108 
109  {
111  psd0.add<std::string>("BooleanOperator", "and");
112  {
114  psd1.add<double>("cut");
115  psd1.add<edm::InputTag>("Producer");
116  psd0.addOptional<edm::ParameterSetDescription>("leadTrack", psd1);
117  }
118  desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
119  }
120 
121  descriptions.add("pfRecoTauDiscriminationByTauPolarization", desc);
122 }
123 
T getParameter(std::string const &) const
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
void beginEvent(const Event &, const EventSetup &) override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
Definition: event.py:1