CMS 3D CMS Logo

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

#include <GenericBenchmarkAnalyzer.h>

Inheritance diagram for GenericBenchmarkAnalyzer:
edm::EDAnalyzer GenericBenchmark edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void endJob () override
 
 GenericBenchmarkAnalyzer (const edm::ParameterSet &)
 
 ~GenericBenchmarkAnalyzer () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 
- Public Member Functions inherited from GenericBenchmark
template<typename C >
void fill (const C *RecoCollection, const C *GenCollection, bool startFromGen=false, bool PlotAgainstReco=true, bool onlyTwoJets=false, double recPt_cut=-1., double minEta_cut=-1., double maxEta_cut=-1., double deltaR_cut=-1.)
 
 GenericBenchmark ()
 
void setfile (TFile *file)
 
void setup (DQMStore *DQM=0, bool PlotAgainstReco_=true, float minDeltaEt=-100., float maxDeltaEt=50., float minDeltaPhi=-0.5, float maxDeltaPhi=0.5, bool doMetPlots=false)
 
void write (std::string Filename)
 
virtual ~GenericBenchmark ()(false)
 

Private Attributes

std::string benchmarkLabel_
 
double deltaR_cut
 
bool doMetPlots_
 
edm::InputTag inputRecoLabel_
 
edm::InputTag inputTruthLabel_
 
float maxDeltaEt_
 
float maxDeltaPhi_
 
double maxEta_cut
 
float minDeltaEt_
 
float minDeltaPhi_
 
double minEta_cut
 
edm::EDGetTokenT< edm::View< reco::Candidate > > myReco_
 
edm::EDGetTokenT< edm::View< reco::Candidate > > myTruth_
 
bool onlyTwoJets_
 
std::string outputFile_
 
bool plotAgainstRecoQuantities_
 
double recPt_cut
 
bool startFromGen_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
- Protected Attributes inherited from GenericBenchmark
PFBenchmarkAlgoalgo_
 
DQMStoredbe_
 

Detailed Description

Definition at line 17 of file GenericBenchmarkAnalyzer.h.

Constructor & Destructor Documentation

GenericBenchmarkAnalyzer::GenericBenchmarkAnalyzer ( const edm::ParameterSet iConfig)
explicit

Definition at line 37 of file GenericBenchmarkAnalyzer.cc.

References benchmarkLabel_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), and AlCaHLTBitMon_QueryRunRegistry::string.

37  {
38  inputTruthLabel_ = iConfig.getParameter<edm::InputTag>("InputTruthLabel");
39  inputRecoLabel_ = iConfig.getParameter<edm::InputTag>("InputRecoLabel");
40  outputFile_ = iConfig.getUntrackedParameter<std::string>("OutputFile");
41  benchmarkLabel_ = iConfig.getParameter<std::string>("BenchmarkLabel");
42  startFromGen_ = iConfig.getParameter<bool>("StartFromGen");
43  plotAgainstRecoQuantities_ = iConfig.getParameter<bool>("PlotAgainstRecoQuantities");
44  onlyTwoJets_ = iConfig.getParameter<bool>("OnlyTwoJets");
45  recPt_cut = iConfig.getParameter<double>("recPt");
46  minEta_cut = iConfig.getParameter<double>("minEta");
47  maxEta_cut = iConfig.getParameter<double>("maxEta");
48  deltaR_cut = iConfig.getParameter<double>("deltaRMax");
49 
50  minDeltaEt_ = iConfig.getParameter<double>("minDeltaEt");
51  maxDeltaEt_ = iConfig.getParameter<double>("maxDeltaEt");
52  minDeltaPhi_ = iConfig.getParameter<double>("minDeltaPhi");
53  maxDeltaPhi_ = iConfig.getParameter<double>("maxDeltaPhi");
54  doMetPlots_ = iConfig.getParameter<bool>("doMetPlots");
55 
56  if (!outputFile_.empty())
57  edm::LogInfo("OutputInfo") << " ParticleFLow Task histograms will be saved to '" << outputFile_.c_str() << "'";
58  else
59  edm::LogInfo("OutputInfo") << " ParticleFlow Task histograms will NOT be saved";
60 
61  myTruth_ = consumes<edm::View<reco::Candidate>>(inputTruthLabel_);
62  myReco_ = consumes<edm::View<reco::Candidate>>(inputRecoLabel_);
63 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::View< reco::Candidate > > myTruth_
edm::EDGetTokenT< edm::View< reco::Candidate > > myReco_
GenericBenchmarkAnalyzer::~GenericBenchmarkAnalyzer ( )
override

Definition at line 65 of file GenericBenchmarkAnalyzer.cc.

65 {}

Member Function Documentation

void GenericBenchmarkAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 84 of file GenericBenchmarkAnalyzer.cc.

References gather_cfg::cout, lumiContext::fill, edm::Event::getByToken(), and edm::Handle< T >::product().

84  {
85  // Typedefs to use views
86  typedef edm::View<reco::Candidate> candidateCollection;
87  typedef edm::View<reco::Candidate> candidateCollection;
88 
89  const candidateCollection *truth_candidates;
90  const candidateCollection *reco_candidates;
91 
92  // ==========================================================
93  // Retrieve!
94  // ==========================================================
95 
96  {
97  // Get Truth Candidates (GenCandidates, GenJets, etc.)
99  bool isGen = iEvent.getByToken(myTruth_, truth_hnd);
100 
101  if (!isGen) {
102  std::cout << "Warning : no Gen jets in input !" << std::endl;
103  return;
104  }
105 
106  truth_candidates = truth_hnd.product();
107 
108  // Get Reco Candidates (PFlow, CaloJet, etc.)
110  bool isReco = iEvent.getByToken(myReco_, reco_hnd);
111  if (!isReco) {
112  std::cout << "Warning : no Reco jets in input !" << std::endl;
113  return;
114  }
115  reco_candidates = reco_hnd.product();
116 
117  // no longer needed with template-ized Benchmark
118  // const PFCandidateCollection *pf_candidates = reco_hnd.product();
119  // static CandidateCollection reco_storage =
120  // algo_->makeCandidateCollection(pf_candidates); reco_candidates =
121  // &reco_storage;
122  }
123  if (!truth_candidates || !reco_candidates) {
124  edm::LogInfo("OutputInfo") << " failed to retrieve data required by ParticleFlow Task";
125  edm::LogInfo("OutputInfo") << " ParticleFlow Task cannot continue...!";
126  return;
127  }
128 
129  // ==========================================================
130  // Analyze!
131  // ==========================================================
132 
133  fill(reco_candidates,
134  truth_candidates,
137  onlyTwoJets_,
138  recPt_cut,
139  minEta_cut,
140  maxEta_cut,
141  deltaR_cut);
142 }
void fill(const C *RecoCollection, const C *GenCollection, bool startFromGen=false, bool PlotAgainstReco=true, bool onlyTwoJets=false, double recPt_cut=-1., double minEta_cut=-1., double maxEta_cut=-1., double deltaR_cut=-1.)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< edm::View< reco::Candidate > > myTruth_
T const * product() const
Definition: Handle.h:74
edm::EDGetTokenT< edm::View< reco::Candidate > > myReco_
void GenericBenchmarkAnalyzer::beginJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 67 of file GenericBenchmarkAnalyzer.cc.

References benchmarkLabel_, dbe_, Utilities::operator, callgraph::path, GeneralSetup::setup(), and AlCaHLTBitMon_QueryRunRegistry::string.

67  {
68  // get ahold of back-end interface
70 
71  if (dbe_) {
72  // dbe_->setVerbose(1);
73  // string path = "PFTask/Benchmarks/" + benchmarkLabel_ + "/";
74  std::string path = "ParticleFlow/" + benchmarkLabel_ + "/";
76  path += "Reco";
77  else
78  path += "Gen";
79  dbe_->setCurrentFolder(path);
81  }
82 }
void setup(DQMStore *DQM=0, bool PlotAgainstReco_=true, float minDeltaEt=-100., float maxDeltaEt=50., float minDeltaPhi=-0.5, float maxDeltaPhi=0.5, bool doMetPlots=false)
void GenericBenchmarkAnalyzer::endJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 144 of file GenericBenchmarkAnalyzer.cc.

References dbe_.

144  {
145  // Store the DAQ Histograms
146  if (!outputFile_.empty())
147  dbe_->save(outputFile_);
148 }

Member Data Documentation

std::string GenericBenchmarkAnalyzer::benchmarkLabel_
private

Definition at line 33 of file GenericBenchmarkAnalyzer.h.

double GenericBenchmarkAnalyzer::deltaR_cut
private

Definition at line 40 of file GenericBenchmarkAnalyzer.h.

bool GenericBenchmarkAnalyzer::doMetPlots_
private

Definition at line 45 of file GenericBenchmarkAnalyzer.h.

edm::InputTag GenericBenchmarkAnalyzer::inputRecoLabel_
private

Definition at line 32 of file GenericBenchmarkAnalyzer.h.

edm::InputTag GenericBenchmarkAnalyzer::inputTruthLabel_
private

Definition at line 31 of file GenericBenchmarkAnalyzer.h.

float GenericBenchmarkAnalyzer::maxDeltaEt_
private

Definition at line 42 of file GenericBenchmarkAnalyzer.h.

float GenericBenchmarkAnalyzer::maxDeltaPhi_
private

Definition at line 44 of file GenericBenchmarkAnalyzer.h.

double GenericBenchmarkAnalyzer::maxEta_cut
private

Definition at line 39 of file GenericBenchmarkAnalyzer.h.

float GenericBenchmarkAnalyzer::minDeltaEt_
private

Definition at line 41 of file GenericBenchmarkAnalyzer.h.

float GenericBenchmarkAnalyzer::minDeltaPhi_
private

Definition at line 43 of file GenericBenchmarkAnalyzer.h.

double GenericBenchmarkAnalyzer::minEta_cut
private

Definition at line 38 of file GenericBenchmarkAnalyzer.h.

edm::EDGetTokenT<edm::View<reco::Candidate> > GenericBenchmarkAnalyzer::myReco_
private

Definition at line 30 of file GenericBenchmarkAnalyzer.h.

edm::EDGetTokenT<edm::View<reco::Candidate> > GenericBenchmarkAnalyzer::myTruth_
private

Definition at line 29 of file GenericBenchmarkAnalyzer.h.

bool GenericBenchmarkAnalyzer::onlyTwoJets_
private

Definition at line 36 of file GenericBenchmarkAnalyzer.h.

std::string GenericBenchmarkAnalyzer::outputFile_
private

Definition at line 28 of file GenericBenchmarkAnalyzer.h.

bool GenericBenchmarkAnalyzer::plotAgainstRecoQuantities_
private

Definition at line 35 of file GenericBenchmarkAnalyzer.h.

double GenericBenchmarkAnalyzer::recPt_cut
private

Definition at line 37 of file GenericBenchmarkAnalyzer.h.

bool GenericBenchmarkAnalyzer::startFromGen_
private

Definition at line 34 of file GenericBenchmarkAnalyzer.h.