CMS 3D CMS Logo

CaloRecoTauDiscriminationByInvMass.cc
Go to the documentation of this file.
3 
6 
7 /* class CaloRecoTauDiscriminationByInvMass
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 
14 #include "TLorentzVector.h"
15 
16 namespace {
17 
18 using namespace reco;
19 using namespace std;
20 
21 class CaloRecoTauDiscriminationByInvMass final : public CaloTauDiscriminationProducerBase {
22  public:
23  explicit CaloRecoTauDiscriminationByInvMass(
24  const edm::ParameterSet& iConfig)
26  invMassMin = iConfig.getParameter<double>("invMassMin");
27  invMassMax = iConfig.getParameter<double>("invMassMax");
28  chargedPionMass = 0.139;
29  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
30  }
31 
32  ~CaloRecoTauDiscriminationByInvMass() override{}
33 
34  double discriminate(const reco::CaloTauRef&) const override;
35 
36  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
37  private:
38  double threeProngInvMass(const CaloTauRef&) const ;
39  double chargedPionMass;
40  double invMassMin,invMassMax;
41  bool booleanOutput;
42 };
43 
44 double CaloRecoTauDiscriminationByInvMass::discriminate(const CaloTauRef& tau) const {
45 
46  double invMass = threeProngInvMass(tau);
47  if(booleanOutput) return (
48  invMass > invMassMin && invMass < invMassMax ? 1. : 0. );
49  return invMass;
50 }
51 
52 double CaloRecoTauDiscriminationByInvMass::threeProngInvMass(
53  const CaloTauRef& tau) const {
54  TLorentzVector sum;
55  reco::TrackRefVector signalTracks = tau->signalTracks();
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  sum += p4;
63  }
64  return sum.M();
65 }
66 
67 }
68 
69 void
71  // caloRecoTauDiscriminationByInvMass
73  desc.add<double>("invMassMin", 0.0);
74  desc.add<edm::InputTag>("CaloTauProducer", edm::InputTag("caloRecoTauProducer"));
75  {
77  psd0.add<std::string>("BooleanOperator", "and");
78  {
80  psd1.add<double>("cut");
81  psd1.add<edm::InputTag>("Producer");
82  psd0.addOptional<edm::ParameterSetDescription>("leadTrack", psd1);
83  }
84  desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
85  }
86  desc.add<bool>("BooleanOutput", true);
87  desc.add<double>("invMassMax", 1.4);
88  descriptions.add("caloRecoTauDiscriminationByInvMass", desc);
89 }
90 
91 DEFINE_FWK_MODULE(CaloRecoTauDiscriminationByInvMass);
T getParameter(std::string const &) const
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double p4[4]
Definition: TauolaWrapper.h:92
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
fixed size matrix
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107