CMS 3D CMS Logo

PFTauElecRejectionBenchMarkAnalyzer.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 
4 // user include files
7 
10 
14 
16 
17 using namespace edm;
18 using namespace reco;
19 using namespace std;
20 
21 //
22 // class declaration
23 
25 public:
28 
31 
32 private:
33  void beginJob() override;
34  void analyze(const edm::Event &, const edm::EventSetup &) override;
35  void endJob() override;
36  // ----------member data ---------------------------
37 
38  string outputfile;
39  DQMStore *db;
41  double maxDeltaR;
42  double minMCPt;
43  double maxMCAbsEta;
44  double minRecoPt;
45  double maxRecoAbsEta;
48 
53 
55 };
57 
58 //
59 // constants, enums and typedefs
60 //
61 
62 //
63 // static data member definitions
64 //
65 
66 //
67 // constructors and destructor
68 //
70 
71 {
72  // now do what ever initialization is needed
73  outputfile = iConfig.getUntrackedParameter<string>("OutputFile");
74  benchmarkLabel = iConfig.getParameter<string>("BenchmarkLabel");
75  sGenParticleSource_tok_ = consumes<edm::HepMCProduct>(iConfig.getParameter<InputTag>("InputTruthLabel"));
76  maxDeltaR = iConfig.getParameter<double>("maxDeltaR");
77  minMCPt = iConfig.getParameter<double>("minMCPt");
78  maxMCAbsEta = iConfig.getParameter<double>("maxMCAbsEta");
79  minRecoPt = iConfig.getParameter<double>("minRecoPt");
80  maxRecoAbsEta = iConfig.getParameter<double>("maxRecoAbsEta");
81  pfTauProducer_tok_ = consumes<reco::PFTauCollection>(iConfig.getParameter<InputTag>("PFTauProducer"));
82  pfTauDiscriminatorByIsolationProducer_tok_ =
83  consumes<reco::PFTauDiscriminator>(iConfig.getParameter<InputTag>("PFTauDiscriminatorByIsolationProducer"));
84  pfTauDiscriminatorAgainstElectronProducer_tok_ =
85  consumes<reco::PFTauDiscriminator>(iConfig.getParameter<InputTag>("PFTauDiscriminatorAgainstElectronProducer"));
86  sGenMatchObjectLabel = iConfig.getParameter<string>("GenMatchObjectLabel");
87  applyEcalCrackCut = iConfig.getParameter<bool>("ApplyEcalCrackCut");
88 
90 
91  PFTauElecRejectionBenchmark_.setup(outputfile,
92  benchmarkLabel,
93  maxDeltaR,
94  minRecoPt,
95  maxRecoAbsEta,
96  minMCPt,
97  maxMCAbsEta,
98  sGenMatchObjectLabel,
99  applyEcalCrackCut,
100  db);
101 }
102 
104  // do anything here that needs to be done at desctruction time
105  // (e.g. close files, deallocate resources etc.)
106 }
107 
108 //
109 // member functions
110 //
111 
112 // ------------ method called to for each event ------------
114  // get gen products
115  Handle<HepMCProduct> mcevt;
116  iEvent.getByToken(sGenParticleSource_tok_, mcevt);
117 
118  // get pftau collection
119  Handle<PFTauCollection> thePFTau;
120  iEvent.getByToken(pfTauProducer_tok_, thePFTau);
121 
122  // get iso discriminator association vector
123  Handle<PFTauDiscriminator> thePFTauDiscriminatorByIsolation;
124  iEvent.getByToken(pfTauDiscriminatorByIsolationProducer_tok_, thePFTauDiscriminatorByIsolation);
125 
126  // get anti-elec discriminator association vector
127  Handle<PFTauDiscriminator> thePFTauDiscriminatorAgainstElectron;
128  iEvent.getByToken(pfTauDiscriminatorAgainstElectronProducer_tok_, thePFTauDiscriminatorAgainstElectron);
129 
130  PFTauElecRejectionBenchmark_.process(
131  mcevt, thePFTau, thePFTauDiscriminatorByIsolation, thePFTauDiscriminatorAgainstElectron);
132 }
133 
134 // ------------ method called once each job just before starting event loop
135 // ------------
137 
138 // ------------ method called once each job just after ending the event loop
139 // ------------
140 void PFTauElecRejectionBenchmarkAnalyzer::endJob() { PFTauElecRejectionBenchmark_.write(); }
141 
142 // define this as a plug-in
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< reco::PFTauCollection > pfTauProducer_tok_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
PFTauElecRejectionBenchmarkAnalyzer(const edm::ParameterSet &)
PFTauElecRejection Benchmark.
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
void beginJob()
Definition: Breakpoints.cc:14
void analyze(const edm::Event &, const edm::EventSetup &) override
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< reco::PFTauDiscriminator > pfTauDiscriminatorByIsolationProducer_tok_
edm::EDGetTokenT< reco::PFTauDiscriminator > pfTauDiscriminatorAgainstElectronProducer_tok_
fixed size matrix
HLT enums.
edm::EDGetTokenT< edm::HepMCProduct > sGenParticleSource_tok_