CMS 3D CMS Logo

PFRecoTauDiscriminationByDeltaE.cc
Go to the documentation of this file.
3 
6 
7 /* class PFRecoTauDiscriminationByDeltaE
8  * created : August 30 2010,
9  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
10  * based on H+ tau ID by Lauri Wendland
11  */
12 
13 #include "TLorentzVector.h"
14 
15 using namespace reco;
16 using namespace std;
17 using namespace edm;
18 
20  public:
22  deltaEmin = iConfig.getParameter<double>("deltaEmin");
23  deltaEmax = iConfig.getParameter<double>("deltaEmax");
24  chargedPionMass = 0.139;
25  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
26  }
27 
29 
30  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
31  double discriminate(const reco::PFTauRef&) const override;
32 
33  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
34  private:
35  double DeltaE(const PFTauRef&) const ;
36 
38 
41 };
42 
44 }
45 
47 
48  double dE = DeltaE(tau);
49  if(booleanOutput) return ( dE > deltaEmin && dE < deltaEmax ? 1. : 0. );
50  return dE;
51 }
52 
54  double tracksE = 0;
55  const std::vector<CandidatePtr>& signalTracks = tau->signalChargedHadrCands();
56  for(size_t i = 0; i < signalTracks.size(); ++i){
57  TLorentzVector p4;
58  p4.SetXYZM(signalTracks[i]->px(),
59  signalTracks[i]->py(),
60  signalTracks[i]->pz(),
61  chargedPionMass);
62  tracksE += p4.E();
63  }
64 
65  double hadrTauP = tau->momentum().r() * (1.0 - tau->emFraction());
66  if (tau->emFraction() >= 1.0) {
67  return -1.0; // electron
68  } else {
69  return tracksE / hadrTauP - 1.0;
70  }
71 }
72 
73 void
75  // pfRecoTauDiscriminationByDeltaE
77  desc.add<double>("deltaEmin", -0.15);
78  {
80  psd0.add<std::string>("BooleanOperator", "and");
81  {
83  psd1.add<double>("cut");
84  psd1.add<edm::InputTag>("Producer");
85  psd0.addOptional<edm::ParameterSetDescription>("leadTrack", psd1);
86  }
87  desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
88  }
89  desc.add<double>("deltaEmax", 1.0);
90  desc.add<bool>("BooleanOutput", true);
91  desc.add<edm::InputTag>("PFTauProducer", edm::InputTag("pfRecoTauProducer"));
92  descriptions.add("pfRecoTauDiscriminationByDeltaE", desc);
93 }
94 
T getParameter(std::string const &) const
PFRecoTauDiscriminationByDeltaE(const ParameterSet &iConfig)
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double p4[4]
Definition: TauolaWrapper.h:92
double discriminate(const reco::PFTauRef &) const override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void beginEvent(const edm::Event &, const edm::EventSetup &) override
fixed size matrix
HLT enums.