CMS 3D CMS Logo

CaloRecoTauDiscriminationByDeltaE.cc
Go to the documentation of this file.
2 
4 
7 
8 /* class CaloRecoTauDiscriminationByDeltaE
9  * created : September 23 2010,
10  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
11  * based on H+ tau ID by Lauri Wendland
12  */
13 
15 #include "TLorentzVector.h"
16 
17 using namespace reco;
18 using namespace std;
19 using namespace edm;
20 
22  public:
24  deltaEmin = iConfig.getParameter<double>("deltaEmin");
25  deltaEmax = iConfig.getParameter<double>("deltaEmax");
26  chargedPionMass = 0.139;
27  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
28  }
29 
31 
32  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
33  double discriminate(const reco::CaloTauRef&) const override;
34 
35  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
36  private:
37  double DeltaE(const CaloTauRef&) const ;
38 
40 
43 };
44 
46 }
47 
49 
50  double dE = DeltaE(tau);
51  if(booleanOutput) return ( dE > deltaEmin && dE < deltaEmax ? 1. : 0. );
52  return dE;
53 }
54 
56  double tracksE = 0;
57  reco::TrackRefVector signalTracks = tau->signalTracks();
58  for(size_t i = 0; i < signalTracks.size(); ++i){
59  TLorentzVector p4;
60  p4.SetXYZM(signalTracks[i]->px(),
61  signalTracks[i]->py(),
62  signalTracks[i]->pz(),
63  chargedPionMass);
64  tracksE += p4.E();
65  }
66  if(tau->leadTrackHCAL3x3hitsEtSum() == 0) return -1; // electron
67  return tracksE/tau->leadTrackHCAL3x3hitsEtSum() - 1.0;
68 }
69 
70 void
72  // caloRecoTauDiscriminationByDeltaE
74  desc.add<double>("deltaEmin", -0.15);
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<double>("deltaEmax", 1.0);
87  desc.add<bool>("BooleanOutput", true);
88  desc.add<edm::InputTag>("TauProducer", edm::InputTag("caloRecoTauProducer"));
89  descriptions.add("caloRecoTauDiscriminationByDeltaE", desc);
90 }
91 
93 
T getParameter(std::string const &) const
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
void beginEvent(const edm::Event &, const edm::EventSetup &) override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double discriminate(const reco::CaloTauRef &) const override
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double p4[4]
Definition: TauolaWrapper.h:92
CaloRecoTauDiscriminationByDeltaE(const ParameterSet &iConfig)
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.
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107