CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
AnalysisTasksAnalyzerJEC Class Reference

Example class that can be used both within FWLite and within the full framework. More...

#include "PhysicsTools/UtilAlgos/interface/AnalysisTasksAnalyzerJEC.h"

Inheritance diagram for AnalysisTasksAnalyzerJEC:
edm::BasicAnalyzer

Public Member Functions

 AnalysisTasksAnalyzerJEC (const edm::ParameterSet &cfg, TFileDirectory &fs)
 default constructor More...
 
 AnalysisTasksAnalyzerJEC (const edm::ParameterSet &cfg, TFileDirectory &fs, edm::ConsumesCollector &&iC)
 
void analyze (const edm::EventBase &event) override
 everything that needs to be done during the event loop More...
 
void beginJob () override
 everything that needs to be done before the event loop More...
 
void endJob () override
 everything that needs to be done after the event loop More...
 
 ~AnalysisTasksAnalyzerJEC () override
 default destructor More...
 
- Public Member Functions inherited from edm::BasicAnalyzer
 BasicAnalyzer (const edm::ParameterSet &cfg, TFileDirectory &fileService)
 default constructor More...
 
 BasicAnalyzer (const edm::ParameterSet &cfg, TFileDirectory &fileService, edm::ConsumesCollector &&iC)
 
virtual ~BasicAnalyzer ()
 default destructor More...
 

Private Attributes

bool help_
 
std::map< std::string, TH2 * > hists_
 histograms More...
 
std::string jecLevel_
 
unsigned int jetInEvents_
 
edm::InputTag Jets_
 input tag for mouns More...
 
edm::EDGetTokenT< std::vector< pat::Jet > > JetsToken_
 
std::string patJetCorrFactors_
 

Detailed Description

Example class that can be used both within FWLite and within the full framework.

This is an example for keeping classes that can be used both within FWLite and within the full framework. The class is derived from the BasicAnalyzer base class, which is an interface for the two wrapper classes EDAnalyzerWrapper and FWLiteAnalyzerWrapper. The latter provides basic configuration file reading and event looping equivalent to the FWLiteHistograms executable of this package. You can see the FWLiteAnalyzerWrapper class at work in the FWLiteWithBasicAnalyzer executable of this package.

Definition at line 19 of file AnalysisTasksAnalyzerJEC.h.

Constructor & Destructor Documentation

◆ AnalysisTasksAnalyzerJEC() [1/2]

AnalysisTasksAnalyzerJEC::AnalysisTasksAnalyzerJEC ( const edm::ParameterSet cfg,
TFileDirectory fs 
)

default constructor

Definition at line 8 of file AnalysisTasksAnalyzerJEC.cc.

References compareTotals::fs, hists_, and jetInEvents_.

10  Jets_(cfg.getParameter<edm::InputTag>("Jets")),
11  jecLevel_(cfg.getParameter<std::string>("jecLevel")),
12  patJetCorrFactors_(cfg.getParameter<std::string>("patJetCorrFactors")),
13  help_(cfg.getParameter<bool>("help")) {
14  hists_["Response"] = fs.make<TH2F>("Response", "response; #eta;p_{T}(reco)/p_{T}(gen)", 5, -3.3, 3.3, 100, 0.4, 1.6);
15  jetInEvents_ = 0;
16 }
edm::InputTag Jets_
input tag for mouns
BasicAnalyzer(const edm::ParameterSet &cfg, TFileDirectory &fileService)
default constructor
Definition: BasicAnalyzer.h:45
std::map< std::string, TH2 * > hists_
histograms

◆ AnalysisTasksAnalyzerJEC() [2/2]

AnalysisTasksAnalyzerJEC::AnalysisTasksAnalyzerJEC ( const edm::ParameterSet cfg,
TFileDirectory fs,
edm::ConsumesCollector &&  iC 
)

Definition at line 17 of file AnalysisTasksAnalyzerJEC.cc.

References compareTotals::fs, hists_, and jetInEvents_.

21  Jets_(cfg.getParameter<edm::InputTag>("Jets")),
22  JetsToken_(iC.consumes<std::vector<pat::Jet> >(Jets_)),
23  jecLevel_(cfg.getParameter<std::string>("jecLevel")),
24  patJetCorrFactors_(cfg.getParameter<std::string>("patJetCorrFactors")),
25  help_(cfg.getParameter<bool>("help")) {
26  hists_["Response"] = fs.make<TH2F>("Response", "response; #eta;p_{T}(reco)/p_{T}(gen)", 5, -3.3, 3.3, 100, 0.4, 1.6);
27  jetInEvents_ = 0;
28 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::InputTag Jets_
input tag for mouns
edm::EDGetTokenT< std::vector< pat::Jet > > JetsToken_
BasicAnalyzer(const edm::ParameterSet &cfg, TFileDirectory &fileService)
default constructor
Definition: BasicAnalyzer.h:45
std::map< std::string, TH2 * > hists_
histograms

◆ ~AnalysisTasksAnalyzerJEC()

AnalysisTasksAnalyzerJEC::~AnalysisTasksAnalyzerJEC ( )
override

default destructor

deconstructor

A Response function is nicer in a TProfile:

Definition at line 30 of file AnalysisTasksAnalyzerJEC.cc.

References geometryDiff::file, hists_, and jecLevel_.

30  {
32  TFile* file = new TFile("myResponse" + TString(jecLevel_) + ".root", "RECREATE");
33  TProfile* prof = hists_["Response"]->ProfileX();
34  prof->Write();
35  file->Write();
36  file->Close();
37 }
std::map< std::string, TH2 * > hists_
histograms

Member Function Documentation

◆ analyze()

void AnalysisTasksAnalyzerJEC::analyze ( const edm::EventBase event)
overridevirtual

everything that needs to be done during the event loop

Implements edm::BasicAnalyzer.

Definition at line 39 of file AnalysisTasksAnalyzerJEC.cc.

References gather_cfg::cout, help_, hists_, jecLevel_, L1TRate_cfi::Jet, jetInEvents_, METSkim_cff::Jets, Jets_, isotrackApplyRegressor::k, and patJetCorrFactors_.

39  {
40  // define what Jet you are using; this is necessary as FWLite is not
41  // capable of reading edm::Views
42  using pat::Jet;
43 
44  // Handle to the Jet collection
46  event.getByLabel(Jets_, Jets);
47 
48  // loop Jet collection and fill histograms
49  for (std::vector<Jet>::const_iterator jet_it = Jets->begin(); jet_it != Jets->end(); ++jet_it) {
51  if (help_ == true) {
52  if (jetInEvents_ == 0) {
53  std::cout << "\n \n number of available JEC sets: " << jet_it->availableJECSets().size() << std::endl;
54  for (unsigned int k = 0; k < jet_it->availableJECSets().size(); ++k) {
55  std::cout << "\n \n available JEC Set: " << jet_it->availableJECSets()[k] << std::endl;
56  }
57  std::cout << "\n \n*** You found out which JEC Sets exist! Now correct it in your config file." << std::endl;
58  std::cout << "\n \n number of available JEC levels " << jet_it->availableJECLevels().size() << std::endl;
59  for (unsigned int k = 0; k < jet_it->availableJECLevels().size(); ++k) {
60  std::cout << "\n \n available JEC: " << jet_it->availableJECLevels(patJetCorrFactors_)[k] << std::endl;
61  }
62  std::cout << "\n \n**** You did it correctly congratulations!!!! And you found out which JEC levels are saved "
63  "within the jets. Correct this in your configuration file."
64  << std::endl;
65  }
66  }
67 
68  if (jet_it->genParticlesSize() > 0) {
69  hists_["Response"]->Fill(jet_it->correctedJet(jecLevel_).eta(),
70  jet_it->correctedJet(jecLevel_).pt() / jet_it->genParticle(0)->pt());
71  std::cout << "\n \n**** You did it correctly congratulations!!!! " << std::endl;
72  }
73  jetInEvents_ += 1;
74  }
75 }
edm::InputTag Jets_
input tag for mouns
std::map< std::string, TH2 * > hists_
histograms

◆ beginJob()

void AnalysisTasksAnalyzerJEC::beginJob ( void  )
inlineoverridevirtual

everything that needs to be done before the event loop

Implements edm::BasicAnalyzer.

Definition at line 27 of file AnalysisTasksAnalyzerJEC.h.

27 {}

◆ endJob()

void AnalysisTasksAnalyzerJEC::endJob ( void  )
inlineoverridevirtual

everything that needs to be done after the event loop

Implements edm::BasicAnalyzer.

Definition at line 29 of file AnalysisTasksAnalyzerJEC.h.

29 {}

Member Data Documentation

◆ help_

bool AnalysisTasksAnalyzerJEC::help_
private

Definition at line 39 of file AnalysisTasksAnalyzerJEC.h.

Referenced by analyze().

◆ hists_

std::map<std::string, TH2*> AnalysisTasksAnalyzerJEC::hists_
private

histograms

Definition at line 42 of file AnalysisTasksAnalyzerJEC.h.

Referenced by AnalysisTasksAnalyzerJEC(), analyze(), and ~AnalysisTasksAnalyzerJEC().

◆ jecLevel_

std::string AnalysisTasksAnalyzerJEC::jecLevel_
private

Definition at line 37 of file AnalysisTasksAnalyzerJEC.h.

Referenced by analyze(), and ~AnalysisTasksAnalyzerJEC().

◆ jetInEvents_

unsigned int AnalysisTasksAnalyzerJEC::jetInEvents_
private

Definition at line 40 of file AnalysisTasksAnalyzerJEC.h.

Referenced by AnalysisTasksAnalyzerJEC(), and analyze().

◆ Jets_

edm::InputTag AnalysisTasksAnalyzerJEC::Jets_
private

input tag for mouns

Definition at line 35 of file AnalysisTasksAnalyzerJEC.h.

Referenced by analyze().

◆ JetsToken_

edm::EDGetTokenT<std::vector<pat::Jet> > AnalysisTasksAnalyzerJEC::JetsToken_
private

Definition at line 36 of file AnalysisTasksAnalyzerJEC.h.

◆ patJetCorrFactors_

std::string AnalysisTasksAnalyzerJEC::patJetCorrFactors_
private

Definition at line 38 of file AnalysisTasksAnalyzerJEC.h.

Referenced by analyze().