CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CaloRecoTauDiscriminationByInvMass.cc
Go to the documentation of this file.
3 
4 /* class CaloRecoTauDiscriminationByInvMass
5  * created : September 23 2010,
6  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
7  * based on H+ tau ID by Lauri Wendland
8  */
9 
11 #include "TLorentzVector.h"
12 
13 namespace {
14 
15 using namespace reco;
16 using namespace std;
17 
18 class CaloRecoTauDiscriminationByInvMass final : public CaloTauDiscriminationProducerBase {
19  public:
20  explicit CaloRecoTauDiscriminationByInvMass(
21  const edm::ParameterSet& iConfig)
23  invMassMin = iConfig.getParameter<double>("invMassMin");
24  invMassMax = iConfig.getParameter<double>("invMassMax");
25  chargedPionMass = 0.139;
26  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
27  }
28 
29  ~CaloRecoTauDiscriminationByInvMass(){}
30 
31  double discriminate(const reco::CaloTauRef&) const override;
32 
33  private:
34  double threeProngInvMass(const CaloTauRef&) const ;
35  double chargedPionMass;
36  double invMassMin,invMassMax;
37  bool booleanOutput;
38 };
39 
40 double CaloRecoTauDiscriminationByInvMass::discriminate(const CaloTauRef& tau) const {
41 
42  double invMass = threeProngInvMass(tau);
43  if(booleanOutput) return (
44  invMass > invMassMin && invMass < invMassMax ? 1. : 0. );
45  return invMass;
46 }
47 
48 double CaloRecoTauDiscriminationByInvMass::threeProngInvMass(
49  const CaloTauRef& tau) const {
50  TLorentzVector sum;
51  reco::TrackRefVector signalTracks = tau->signalTracks();
52  for(size_t i = 0; i < signalTracks.size(); ++i){
53  TLorentzVector p4;
54  p4.SetXYZM(signalTracks[i]->px(),
55  signalTracks[i]->py(),
56  signalTracks[i]->pz(),
57  chargedPionMass);
58  sum += p4;
59  }
60  return sum.M();
61 }
62 
63 }
64 DEFINE_FWK_MODULE(CaloRecoTauDiscriminationByInvMass);
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
double p4[4]
Definition: TauolaWrapper.h:92
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89