CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
HLTExoticaSubAnalysis Class Reference

#include <HLTExoticaSubAnalysis.h>

Public Types

enum  Level { Level::GEN = 98, Level::RECO = 99 }
 

Public Member Functions

void analyze (const edm::Event &iEvent, const edm::EventSetup &iEventSetup, EVTColContainer *cols)
 
void beginJob ()
 
void beginRun (const edm::Run &iRun, const edm::EventSetup &iEventSetup)
 
 HLTExoticaSubAnalysis (const edm::ParameterSet &pset, const std::string &analysisname, edm::ConsumesCollector &&consCollector)
 Constructor. More...
 
void subAnalysisBookHistos (DQMStore::IBooker &iBooker, const edm::Run &iRun, const edm::EventSetup &iSetup)
 
 ~HLTExoticaSubAnalysis ()
 End Constructor. More...
 

Private Member Functions

void bookHist (DQMStore::IBooker &iBooker, const std::string &source, const std::string &objType, const std::string &variable)
 The internal functions to book and fill histograms. More...
 
void fillHist (const std::string &source, const std::string &objType, const std::string &variable, const float &value)
 
void getHandlesToObjects (const edm::Event &iEvent, EVTColContainer *col)
 Gets the collections themselves. More...
 
void getNamesOfObjects (const edm::ParameterSet &anpset)
 Creates the maps that map which collection should come from which label. More...
 
const std::vector< unsigned int > getObjectsType (const std::string &hltpath) const
 closes analyze method More...
 
void initSelector (const unsigned int &objtype)
 Initializes the selectors of the objects based on which object it is. More...
 
void insertCandidates (const unsigned int &objtype, const EVTColContainer *col, std::vector< reco::LeafCandidate > *matches, std::map< int, double > &theSumEt, std::map< int, std::vector< const reco::Track * > > &trkObjs)
 
void registerConsumes (edm::ConsumesCollector &consCollector)
 Registers consumption of objects. More...
 

Private Attributes

std::string _analysisname
 The name of this sub-analysis. More...
 
edm::InputTag _beamSpotLabel
 
edm::EDGetTokenT< reco::BeamSpot_bsToken
 
std::map< std::string,
MonitorElement * > 
_elements
 Structure of the MonitorElements. More...
 
std::map< unsigned int,
std::string > 
_genCut
 gen/rec objects cuts More...
 
std::map< unsigned int,
std::string > 
_genCut_leading
 gen/rec pt-leading objects cuts More...
 
StringCutObjectSelector
< reco::GenMET > * 
_genMETSelector
 
edm::InputTag _genParticleLabel
 
edm::EDGetTokenT
< reco::GenParticleCollection
_genParticleToken
 And also the tokens to get the object collections. More...
 
std::map< unsigned int,
StringCutObjectSelector
< reco::GenParticle > * > 
_genSelectorMap
 
HLTConfigProvider _hltConfig
 Interface to the HLT information. More...
 
std::set< std::string > _hltPaths
 The hlt paths found in the hltConfig. More...
 
std::vector< std::string > _hltPathsToCheck
 The hlt paths to check for. More...
 
std::string _hltProcessName
 The labels of the object collections to be used in this analysis. More...
 
StringCutObjectSelector
< l1extra::L1EtMissParticle > * 
_l1METSelector
 
unsigned int _minCandidates
 The minimum number of reco/gen candidates needed by the analysis. More...
 
std::vector< double > _parametersDxy
 
std::vector< double > _parametersEta
 Some kinematical parameters. More...
 
std::vector< double > _parametersPhi
 
std::vector< double > _parametersTurnOn
 
std::vector< double > _parametersTurnOnSumEt
 
std::vector< HLTExoticaPlotter_plotters
 The plotters: managers of each hlt path where the plots are done. More...
 
edm::ParameterSet _pset
 Internal, working copy of the PSet passed from above. More...
 
StringCutObjectSelector
< reco::CaloJet > * 
_recCaloJetSelector
 
StringCutObjectSelector
< reco::CaloMET > * 
_recCaloMETSelector
 
std::map< unsigned int,
std::string > 
_recCut
 
std::map< unsigned int,
std::string > 
_recCut_leading
 
StringCutObjectSelector
< reco::GsfElectron > * 
_recElecSelector
 
std::map< unsigned int,
edm::InputTag
_recLabels
 
StringCutObjectSelector
< reco::MET > * 
_recMETSelector
 
StringCutObjectSelector
< reco::Muon > * 
_recMuonSelector
 
StringCutObjectSelector
< reco::Track > * 
_recMuonTrkSelector
 
StringCutObjectSelector
< reco::PFJet > * 
_recPFJetSelector
 
StringCutObjectSelector
< reco::PFMET > * 
_recPFMETSelector
 
StringCutObjectSelector
< reco::PFMET > * 
_recPFMHTSelector
 
StringCutObjectSelector
< reco::PFTau > * 
_recPFTauSelector
 
StringCutObjectSelector
< reco::Photon > * 
_recPhotonSelector
 
StringCutObjectSelector
< reco::Track > * 
_recTrackSelector
 
std::map< std::string,
std::string > 
_shortpath2long
 Relation between the short and long versions of the path. More...
 
std::map< unsigned int,
edm::EDGetToken
_tokens
 
edm::InputTag _trigResultsLabel
 
edm::EDGetTokenT
< edm::TriggerResults
_trigResultsToken
 

Detailed Description

This class is the main workhorse of the package. It makes the histograms for one given analysis, taking care of all HLT paths related to that analysis.

Generate histograms for trigger efficiencies Exotica related Documentation available on the CMS TWiki: https://twiki.cern.ch/twiki/bin/view/CMS/EXOTICATriggerValidation

Author
Thiago R. Fernandez Perez Tomei Based and adapted from: J. Duarte Campderros code from HLTriggerOffline/Higgs J. Klukas, M. Vander Donckt and J. Alcaraz code from the HLTriggerOffline/Muon package.

Definition at line 66 of file HLTExoticaSubAnalysis.h.

Member Enumeration Documentation

Enumerator
GEN 
RECO 

Definition at line 68 of file HLTExoticaSubAnalysis.h.

68  {
69  GEN = 98,
70  RECO = 99
71  };

Constructor & Destructor Documentation

HLTExoticaSubAnalysis::HLTExoticaSubAnalysis ( const edm::ParameterSet pset,
const std::string &  analysisname,
edm::ConsumesCollector &&  consCollector 
)

Constructor.

Get the vector of paths to check, for this particular analysis.

Get the minimum candidates, for this particular analysis.

Definition at line 28 of file HLTExoticaSubAnalysis.cc.

References _genCut, _genCut_leading, _hltPathsToCheck, _minCandidates, _parametersDxy, _parametersEta, _parametersPhi, _parametersTurnOn, _parametersTurnOnSumEt, _pset, _recCut, _recCut_leading, _recLabels, edm::ParameterSet::exists(), edm::ParameterSet::existsAs(), getNamesOfObjects(), edm::ParameterSet::getParameter(), EVTColContainer::getTypeString(), edm::ParameterSet::getUntrackedParameter(), edm::ParameterSet::insert(), LogDebug, registerConsumes(), edm::ParameterSet::retrieve(), and AlCaHLTBitMon_QueryRunRegistry::string.

30  :
31  _pset(pset),
32  _analysisname(analysisname),
33  _minCandidates(0),
34  _hltProcessName(pset.getParameter<std::string>("hltProcessName")),
35  _genParticleLabel(pset.getParameter<std::string>("genParticleLabel")),
36  _trigResultsLabel("TriggerResults", "", _hltProcessName),
37  _beamSpotLabel(pset.getParameter<std::string>("beamSpotLabel")),
38  _parametersEta(pset.getParameter<std::vector<double> >("parametersEta")),
39  _parametersPhi(pset.getParameter<std::vector<double> >("parametersPhi")),
40  _parametersTurnOn(pset.getParameter<std::vector<double> >("parametersTurnOn")),
41  _parametersTurnOnSumEt(pset.getParameter<std::vector<double> >("parametersTurnOnSumEt")),
42  _parametersDxy(pset.getParameter<std::vector<double> >("parametersDxy")),
47  _recMETSelector(0),
50  _genMETSelector(0),
52  _l1METSelector(0),
57 {
58 
59  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::constructor()";
60 
61  // Specific parameters for this analysis
62  edm::ParameterSet anpset = pset.getParameter<edm::ParameterSet>(analysisname);
63 
64  // If this analysis has a particular set of binnings, use it.
65  // (Taken from the analysis-specific parameter set, of course)
66  // The "true" in the beginning of _pset.insert() means
67  // "overwrite the parameter if need be".
68  if (anpset.exists("parametersTurnOn")) {
69  _parametersTurnOn = anpset.getParameter<std::vector<double> >("parametersTurnOn");
70  _pset.insert(true, "parametersTurnOn", anpset.retrieve("parametersTurnOn"));
71  }
72  if (anpset.exists("parametersEta")) {
73  _parametersEta = anpset.getParameter<std::vector<double> >("parametersEta");
74  _pset.insert(true, "parametersEta", anpset.retrieve("parametersEta"));
75  }
76  if (anpset.exists("parametersPhi")) {
77  _parametersPhi = anpset.getParameter<std::vector<double> >("parametersPhi");
78  _pset.insert(true, "parametersPhi", anpset.retrieve("parametersPhi"));
79  }
80  if (anpset.exists("parametersDxy")) {
81  _parametersDxy = anpset.getParameter<std::vector<double> >("parametersDxy");
82  _pset.insert(true, "parametersDxy", anpset.retrieve("parametersDxy"));
83  }
84  if (anpset.exists("parametersTurnOnSumEt")) {
85  _parametersTurnOnSumEt = anpset.getParameter<std::vector<double> >("parametersTurnOnSumEt");
86  _pset.insert(true, "parametersTurnOnSumEt", anpset.retrieve("parametersTurnOnSumEt"));
87  }
88 
89  // Get names of objects that we may want to get from the event.
90  // Notice that genParticles are dealt with separately.
91  this->getNamesOfObjects(anpset);
92 
93  // Since now we have the names, we should register the consumption
94  // of objects.
95  this->registerConsumes(consCollector);
96 
97  // Generic objects: Initialization of basic phase space cuts.
98  for (std::map<unsigned int, edm::InputTag>::const_iterator it = _recLabels.begin();
99  it != _recLabels.end(); ++it) {
100  const std::string objStr = EVTColContainer::getTypeString(it->first);
101  _genCut[it->first] = pset.getParameter<std::string>(objStr + "_genCut");
102  _recCut[it->first] = pset.getParameter<std::string>(objStr + "_recCut");
103  auto const genCutParam = objStr + "_genCut_leading";
104  if (pset.exists(genCutParam)) {
105  _genCut_leading[it->first] = pset.getParameter<std::string>(genCutParam);
106  } else {
107  _genCut_leading[it->first] = "pt>0"; // no cut
108  }
109  auto const recCutParam = objStr + "_recCut_leading";
110  if (pset.exists(recCutParam)) {
111  _recCut_leading[it->first] = pset.getParameter<std::string>(recCutParam);
112  } else {
113  _recCut_leading[it->first] = "pt>0"; // no cut
114  }
115  }
116 
117  //--- Updating parameters if has to be modified for this particular specific analysis
118  for (std::map<unsigned int, edm::InputTag>::const_iterator it = _recLabels.begin();
119  it != _recLabels.end(); ++it) {
120  const std::string objStr = EVTColContainer::getTypeString(it->first);
121 
122  auto const genCutParam = objStr + "_genCut";
123  if(anpset.existsAs<std::string>(genCutParam,false) ) {
124  _genCut[it->first] = anpset.getUntrackedParameter<std::string>(genCutParam);
125  }
126 
127  auto const recCutParam = objStr + "_recCut";
128  if(anpset.existsAs<std::string>(recCutParam,false) ) {
129  _recCut[it->first] = anpset.getUntrackedParameter<std::string>(recCutParam);
130  }
131 
132  }
133 
135  _hltPathsToCheck = anpset.getParameter<std::vector<std::string> >("hltPathsToCheck");
137  _minCandidates = anpset.getParameter<unsigned int>("minCandidates");
138 
139 }
#define LogDebug(id)
StringCutObjectSelector< reco::Track > * _recMuonTrkSelector
T getParameter(std::string const &) const
unsigned int _minCandidates
The minimum number of reco/gen candidates needed by the analysis.
Entry const & retrieve(char const *) const
T getUntrackedParameter(std::string const &, T const &) const
void getNamesOfObjects(const edm::ParameterSet &anpset)
Creates the maps that map which collection should come from which label.
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
std::map< unsigned int, std::string > _genCut
gen/rec objects cuts
bool exists(std::string const &parameterName) const
checks if a parameter exists
StringCutObjectSelector< reco::CaloJet > * _recCaloJetSelector
StringCutObjectSelector< reco::GenMET > * _genMETSelector
std::vector< double > _parametersPhi
void insert(bool ok_to_replace, char const *, Entry const &)
StringCutObjectSelector< reco::MET > * _recMETSelector
std::map< unsigned int, std::string > _genCut_leading
gen/rec pt-leading objects cuts
StringCutObjectSelector< reco::PFJet > * _recPFJetSelector
std::map< unsigned int, edm::InputTag > _recLabels
StringCutObjectSelector< reco::GsfElectron > * _recElecSelector
std::map< unsigned int, std::string > _recCut_leading
StringCutObjectSelector< reco::Muon > * _recMuonSelector
edm::ParameterSet _pset
Internal, working copy of the PSet passed from above.
std::vector< double > _parametersEta
Some kinematical parameters.
std::map< unsigned int, std::string > _recCut
std::vector< double > _parametersDxy
StringCutObjectSelector< reco::Track > * _recTrackSelector
std::vector< std::string > _hltPathsToCheck
The hlt paths to check for.
std::string _hltProcessName
The labels of the object collections to be used in this analysis.
void registerConsumes(edm::ConsumesCollector &consCollector)
Registers consumption of objects.
StringCutObjectSelector< reco::PFMET > * _recPFMETSelector
std::string _analysisname
The name of this sub-analysis.
StringCutObjectSelector< reco::Photon > * _recPhotonSelector
StringCutObjectSelector< l1extra::L1EtMissParticle > * _l1METSelector
std::vector< double > _parametersTurnOn
StringCutObjectSelector< reco::CaloMET > * _recCaloMETSelector
std::vector< double > _parametersTurnOnSumEt
static const std::string getTypeString(const unsigned int &objtype)
Tranform types into strings.
StringCutObjectSelector< reco::PFMET > * _recPFMHTSelector
StringCutObjectSelector< reco::PFTau > * _recPFTauSelector
HLTExoticaSubAnalysis::~HLTExoticaSubAnalysis ( )

End Constructor.

Definition at line 141 of file HLTExoticaSubAnalysis.cc.

References _genMETSelector, _genSelectorMap, _l1METSelector, _recCaloJetSelector, _recCaloMETSelector, _recElecSelector, _recMETSelector, _recMuonSelector, _recMuonTrkSelector, _recPFJetSelector, _recPFMETSelector, _recPFMHTSelector, _recPFTauSelector, _recPhotonSelector, _recTrackSelector, and python.multivaluedict::map().

142 {
143  for (std::map<unsigned int, StringCutObjectSelector<reco::GenParticle>* >::iterator it = _genSelectorMap.begin();
144  it != _genSelectorMap.end(); ++it) {
145  delete it->second;
146  it->second = 0;
147  }
148  delete _recMuonSelector;
149  _recMuonSelector = 0;
150  delete _recMuonTrkSelector;
152  delete _recTrackSelector;
153  _recTrackSelector = 0;
154  delete _recElecSelector;
155  _recElecSelector = 0;
156  delete _recPhotonSelector;
157  _recPhotonSelector = 0;
158  delete _recMETSelector;
159  _recMETSelector = 0;
160  delete _recPFMETSelector;
161  _recPFMETSelector = 0;
162  delete _recPFMHTSelector;
163  _recPFMHTSelector = 0;
164  delete _genMETSelector;
165  _genMETSelector = 0;
166  delete _recCaloMETSelector;
168  delete _l1METSelector;
169  _l1METSelector = 0;
170  delete _recPFTauSelector;
171  _recPFTauSelector = 0;
172  delete _recPFJetSelector;
173  _recPFJetSelector = 0;
174  delete _recCaloJetSelector;
176 }
StringCutObjectSelector< reco::Track > * _recMuonTrkSelector
StringCutObjectSelector< reco::CaloJet > * _recCaloJetSelector
StringCutObjectSelector< reco::GenMET > * _genMETSelector
StringCutObjectSelector< reco::MET > * _recMETSelector
StringCutObjectSelector< reco::PFJet > * _recPFJetSelector
StringCutObjectSelector< reco::GsfElectron > * _recElecSelector
StringCutObjectSelector< reco::Muon > * _recMuonSelector
StringCutObjectSelector< reco::Track > * _recTrackSelector
StringCutObjectSelector< reco::PFMET > * _recPFMETSelector
StringCutObjectSelector< reco::Photon > * _recPhotonSelector
std::map< unsigned int, StringCutObjectSelector< reco::GenParticle > * > _genSelectorMap
StringCutObjectSelector< l1extra::L1EtMissParticle > * _l1METSelector
StringCutObjectSelector< reco::CaloMET > * _recCaloMETSelector
StringCutObjectSelector< reco::PFMET > * _recPFMHTSelector
StringCutObjectSelector< reco::PFTau > * _recPFTauSelector

Member Function Documentation

void HLTExoticaSubAnalysis::analyze ( const edm::Event iEvent,
const edm::EventSetup iEventSetup,
EVTColContainer cols 
)

Method to fill all relevant histograms. Notice that we update the EVTColContaner to point to the collections we want.

We are going to make a fake reco::LeafCandidate, with our particleType as the pdgId. This is an alternative to the older implementation with MatchStruct.

Filling the histograms if pass the minimum amount of candidates needed by the analysis:

GEN CASE ///

Close GEN case

RECO CASE ///

Debugging.

Close RECO case

Definition at line 368 of file HLTExoticaSubAnalysis.cc.

References _analysisname, _genCut, _genCut_leading, _genSelectorMap, _minCandidates, _plotters, _recCut_leading, _recLabels, _shortpath2long, edm::HLTGlobalStatus::accept(), EVTColContainer::bs, ecal_dqm_sourceclient-live_cfg::cerr, funct::cos(), gather_cfg::cout, EVTColContainer::ELEC, eta, fillHist(), EVTColContainer::genParticles, getHandlesToObjects(), EVTColContainer::getTypeString(), i, initSelector(), insertCandidates(), j, findQualityFiles::jj, LogDebug, EVTColContainer::MUON, reco::Candidate::p4(), phi, reco::BeamSpot::position(), EnergyCorrector::pt, benchmark_cfg::select, funct::sin(), findQualityFiles::size, AlCaHLTBitMon_QueryRunRegistry::string, edm::TriggerNames::triggerIndex(), edm::Event::triggerNames(), EVTColContainer::triggerResults, trigNames, findQualityFiles::v, reco::Candidate::vertex(), reco::BeamSpot::x0(), and reco::BeamSpot::y0().

369 {
370  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::analyze()";
371 
372  if(verbose>2) std::cerr << "### Category : " << _analysisname << std::endl;
373  // Loop over _recLabels to make sure everything is alright.
374  /*
375  std::cout << "Now printing the _recLabels" << std::endl;
376  for (std::map<unsigned int, edm::InputTag>::iterator it = _recLabels.begin();
377  it != _recLabels.end(); ++it) {
378  std::cout << "Number: " << it->first << "\t" << "Label: " << it->second.label() << std::endl;
379  }
380  */
381 
382  // Initialize the collection (the ones which have not been initialiazed yet)
383  //std::cout << "Setting handles to objects..." << std::endl;
384  this->getHandlesToObjects(iEvent, cols);
385 
386  // Utility map, mapping kinds of objects (GEN, RECO) to strings ("gen","rec")
387  //std::map<Level, std::string> u2str;
388  //u2str[Level::GEN] = "gen";
389  //u2str[Level::RECO] = "rec";
390 
391  // Extract the match structure containing the gen/reco candidates (electron, muons,...). This part is common to all the SubAnalyses
392  std::vector<reco::LeafCandidate> matchesGen; matchesGen.clear();
393  std::vector<reco::LeafCandidate> matchesReco; matchesReco.clear();
394  std::map<int , double> theSumEt; // map< pdgId ; SumEt > in order to keep track of the MET type
395  std::map<int, std::vector<const reco::Track*> > trkObjs;
396 
397  // --- deal with GEN objects first.
398  // Make each good GEN object into the base cand for a MatchStruct
399  // Our definition of "good" is "passes the selector" defined in the config.py
400  // Save all the MatchStructs in the "matchesGen" vector.
401 
402  for (std::map<unsigned int, edm::InputTag>::iterator it = _recLabels.begin();
403  it != _recLabels.end(); ++it) {
404  // Here we are filling the vector of StringCutObjectSelector<reco::GenParticle>
405  // with objects constructed from the strings saved in _genCut.
406  // Initialize selectors when first event
407 
408  //std::cout << "Loop over the kinds of objects: objects of kind " << it->first << std::endl;
409 
410 
411  if (!_genSelectorMap[it->first]) {
413  }
414 
415  const std::string objTypeStr = EVTColContainer::getTypeString(it->first);
416  // genAnyMET doesn't make sense. No need their matchesGens
417  if ( TString(objTypeStr).Contains("MET") ||
418  TString(objTypeStr).Contains("MHT") ||
419  TString(objTypeStr).Contains("Jet") ) continue;
420 
421  // Now loop over the genParticles, and apply the operator() over each of them.
422  // Fancy syntax: for objects X and Y, X.operator()(Y) is the same as X(Y).
423  for (size_t i = 0; i < cols->genParticles->size(); ++i) {
424  //std::cout << "Now matchesGen.size() is " << matchesGen.size() << std::endl;
425  if (_genSelectorMap[it->first]->operator()(cols->genParticles->at(i))) {
426  const reco::Candidate* cand = &(cols->genParticles->at(i));
427  //std::cout << "Found good cand: cand->pt() = " << cand->pt() << std::endl;
428  //matchesGen.push_back(MatchStruct(cand, it->first));
431  reco::LeafCandidate v(0,cand->p4(),cand->vertex(),it->first,0,true);
432 
433  matchesGen.push_back(v);
434  }
435  }
436  }
437 
438  // Sort the matches by pT for later filling of turn-on curve
439  //std::cout << "Before sorting: matchesGen.size() = " << matchesGen.size() << std::endl;
440 
441  // GreaterByPt<reco::LeafCandidate> comparator;
442  // std::sort(matchesGen.begin(),
443  // matchesGen.end(),
444  // comparator);
445 
446  // --- same for RECO objects
447  // Extraction of the objects candidates
448  if(verbose>0) std::cout << "-- enter loop over recLabels" << std::endl;
449  for (std::map<unsigned int, edm::InputTag>::iterator it = _recLabels.begin();
450  it != _recLabels.end(); ++it) {
451  //std::cout << "Filling RECO \"matchesReco\" vector for particle kind it->first = "
452  // << it->first << ", which means " << it->second.label() << std::endl;
453  // Reco selectors (the function takes into account if it was instantiated
454  // before or not) ### Thiago ---> Then why don't we put it in the beginRun???
455  this->initSelector(it->first);
456  // -- Storing the matchesReco
457  this->insertCandidates(it->first, cols, &matchesReco, theSumEt, trkObjs);
458  if(verbose>0) std::cout << "--- " << EVTColContainer::getTypeString(it->first)
459  << " sumEt=" << theSumEt[it->first] << std::endl;
460  }
461 
462  // std::sort(matchesReco.begin(),
463  // matchesReco.end(),
464  // comparator);
465 
466  // -- Trigger Results
467  const edm::TriggerNames trigNames = iEvent.triggerNames(*(cols->triggerResults));
468 
470 
471  //for (std::map<unsigned int, std::vector<MatchStruct> >::iterator it = sourceMatchMap.begin(); it != sourceMatchMap.end(); ++it) {
472  // it->first: gen/reco it->second: HLT matches (std::vector<MatchStruct>)
473 
474  //if (it->second.size() < _minCandidates) { // FIXME: A bug is potentially here: what about the mixed channels?
475  //continue;
476  //}
477 
481  if(verbose>2) std::cerr << "### matchesGen.size() = " << matchesGen.size() << std::endl;
482  if( matchesGen.size() >= _minCandidates) { // FIXME: A bug is potentially here: what about the mixed channels?
483  // Okay, there are enough candidates. Move on!
484 
485  // Filling the gen/reco objects (eff-denominators):
486  // Just the first two different ones, if there are more
487  // The countobjects maps uints (object types, really) --> integers.
488  // Example:
489  // | uint | int |
490  // | 0 | 1 | --> 1 muon used
491  // | 1 | 2 | --> 2 electrons used
492 
493  // Initializing the count of the used objects.
494  std::map<unsigned int, int> countobjects;
495  for (std::map<unsigned int, edm::InputTag>::iterator co = _recLabels.begin();
496  co != _recLabels.end(); ++co) {
497  //countobjects->insert(std::pair<unsigned int, int>(co->first, 0));
498  countobjects.insert(std::pair<unsigned int, int>(co->first, 0));
499  }
500 
501  int counttotal = 0;
502 
503  // 3 : pt1, pt2, pt3
504  int totalobjectssize3 = 3 * countobjects.size();
505 
506 
507  bool isPassedLeadingCut = true;
508  // We will proceed only when cuts for the pt-leading are satisified.
509  for (size_t j = 0; j != matchesGen.size(); ++j) {
510  const unsigned int objType = matchesGen[j].pdgId();
511  // Cut for the pt-leading object
513  if ( !select( matchesGen[j] ) ) { // No interest case
514  isPassedLeadingCut = false; // Will skip the following matchesGen loop
515  matchesGen.clear();
516  break;
517  }
518  }
519 
520  std::vector<float> dxys; dxys.clear();
521 
522  for (size_t j = 0; ( j != matchesGen.size() ) && isPassedLeadingCut; ++j) {
523  const unsigned int objType = matchesGen[j].pdgId();
524  //std::cout << "(4) Gonna call with " << objType << std::endl;
525  const std::string objTypeStr = EVTColContainer::getTypeString(objType);
526 
527  float pt = matchesGen[j].pt();
528 
529  if (countobjects[objType] == 0) {
530  this->fillHist("gen", objTypeStr, "MaxPt1", pt);
531  ++(countobjects[objType]);
532  ++counttotal;
533  }
534  else if (countobjects[objType] == 1) {
535  this->fillHist("gen", objTypeStr, "MaxPt2", pt);
536  ++(countobjects[objType]);
537  ++counttotal;
538  }
539  else if (countobjects[objType] == 2) {
540  this->fillHist("gen", objTypeStr, "MaxPt3", pt);
541  ++(countobjects[objType]);
542  ++counttotal;
543  }
544  else {
545  // Already the minimum three objects has been filled, get out...
546  if (counttotal == totalobjectssize3) {
547  size_t max_size = matchesGen.size();
548  for ( size_t jj = j; jj < max_size; jj++ ) {
549  matchesGen.erase(matchesGen.end());
550  }
551  break;
552  }
553  }
554 
555  float eta = matchesGen[j].eta();
556  float phi = matchesGen[j].phi();
557 
558  this->fillHist("gen", objTypeStr, "Eta", eta);
559  this->fillHist("gen", objTypeStr, "Phi", phi);
560 
561  // If the target is electron or muon,
562  // we will add Dxy plots.
563  if ( objType == EVTColContainer::MUON ||
564  objType == EVTColContainer::ELEC ) {
565  const math::XYZPoint & vtx = matchesGen[j].vertex();
566  float momphi = matchesGen[j].momentum().phi();
567  float dxyGen = (-(vtx.x()-cols->bs->x0())*sin(momphi)+(vtx.y()-cols->bs->y0())*cos(momphi));
568  dxys.push_back(dxyGen);
569  this->fillHist("gen", objTypeStr, "Dxy", dxyGen);
570  }
571 
572  } // Closes loop in gen
573 
574  // Calling to the plotters analysis (where the evaluation of the different trigger paths are done)
575  //const std::string source = "gen";
576  for (std::vector<HLTExoticaPlotter>::iterator an = _plotters.begin(); an != _plotters.end(); ++an) {
577  const std::string hltPath = _shortpath2long[an->gethltpath()];
578  const bool ispassTrigger = cols->triggerResults->accept(trigNames.triggerIndex(hltPath));
579  LogDebug("ExoticaValidation") << " preparing to call the plotters analysis";
580  an->analyze(ispassTrigger, "gen", matchesGen, theSumEt, dxys);
581  LogDebug("ExoticaValidation") << " called the plotter";
582  }
583  }
584 
588  if(verbose>2) std::cerr << "### matchesReco.size() = " << matchesReco.size() << std::endl;
589 
590  {
591  if(matchesReco.size() < _minCandidates) return; // FIXME: A bug is potentially here: what about the mixed channels?
592 
593  // Okay, there are enough candidates. Move on!
594 
595  // Filling the gen/reco objects (eff-denominators):
596  // Just the first two different ones, if there are more
597  // The countobjects maps uints (object types, really) --> integers.
598  // Example:
599  // | uint | int |
600  // | 0 | 1 | --> 1 muon used
601  // | 1 | 2 | --> 2 electrons used
602  // Initializing the count of the used objects.
603  //std::map<unsigned int, int> * countobjects = new std::map<unsigned int, int>;
604  std::map<unsigned int, int> countobjects;
605  for (std::map<unsigned int, edm::InputTag>::iterator co = _recLabels.begin();
606  co != _recLabels.end(); ++co) {
607  countobjects.insert(std::pair<unsigned int, int>(co->first, 0));
608  }
609 
610  int counttotal = 0;
611 
612  // 3 : pt1, pt2, pt3
613  int totalobjectssize3 = 3 * countobjects.size();
614 
616  //std::cout << "Our RECO vector has matchesReco.size() = " << matchesReco.size() << std::endl;
617 
618  std::vector<float> dxys; dxys.clear();
619 
620  bool isPassedLeadingCut = true;
621  // We will proceed only when cuts for the pt-leading are satisified.
622  for (size_t j = 0; j != matchesReco.size(); ++j) {
623  const unsigned int objType = matchesReco[j].pdgId();
624  // Cut for the pt-leading object
626  if ( !select( matchesReco[j] ) ) { // No interest case
627  isPassedLeadingCut = false; // Will skip the following matchesReco loop
628  matchesReco.clear();
629  break;
630  }
631  }
632 
633  for (size_t j = 0; ( j != matchesReco.size() ) && isPassedLeadingCut; ++j) {
634  const unsigned int objType = matchesReco[j].pdgId();
635  //std::cout << "(4) Gonna call with " << objType << std::endl;
636  const std::string objTypeStr = EVTColContainer::getTypeString(objType);
637 
638  float pt = matchesReco[j].pt();
639 
640  if (countobjects[objType] == 0) {
641  this->fillHist("rec", objTypeStr, "MaxPt1", pt);
642  ++(countobjects[objType]);
643  ++counttotal;
644  }
645  else if (countobjects[objType] == 1) {
646  if( ! ( TString(objTypeStr).Contains("MET") || TString(objTypeStr).Contains("MHT") ) ) {
647  this->fillHist("rec", objTypeStr, "MaxPt2", pt);
648  }
649  ++(countobjects[objType]);
650  ++counttotal;
651  }
652  else if (countobjects[objType] == 2) {
653  if( ! ( TString(objTypeStr).Contains("MET") || TString(objTypeStr).Contains("MHT") ) ) {
654  this->fillHist("rec", objTypeStr, "MaxPt3", pt);
655  }
656  ++(countobjects[objType]);
657  ++counttotal;
658  }
659  else {
660  // Already the minimum three objects has been filled, get out...
661  if (counttotal == totalobjectssize3) {
662  size_t max_size = matchesReco.size();
663  for ( size_t jj = j; jj < max_size; jj++ ) {
664  matchesReco.erase(matchesReco.end());
665  }
666  break;
667  }
668  }
669 
670  float eta = matchesReco[j].eta();
671  float phi = matchesReco[j].phi();
672 
673  if ( !( TString(objTypeStr).Contains("MET") || TString(objTypeStr).Contains("MHT") ) ) {
674  this->fillHist("rec", objTypeStr, "Eta", eta);
675  this->fillHist("rec", objTypeStr, "Phi", phi);
676  }
677  else {
678  this->fillHist("rec", objTypeStr, "SumEt", theSumEt[objType]);
679  }
680 
681  if (trkObjs[objType].size()>=j+1) {
682  float dxyRec = trkObjs[objType].at(j)->dxy(cols->bs->position());
683  this->fillHist("rec", objTypeStr, "Dxy", dxyRec);
684  dxys.push_back(dxyRec);
685  }
686 
687  } // Closes loop in reco
688 
689  //LogDebug("ExoticaValidation") << " deleting countobjects";
690  //delete countobjects;
691 
692  // Calling to the plotters analysis (where the evaluation of the different trigger paths are done)
693  //const std::string source = "reco";
694  for (std::vector<HLTExoticaPlotter>::iterator an = _plotters.begin(); an != _plotters.end(); ++an) {
695  const std::string hltPath = _shortpath2long[an->gethltpath()];
696  const bool ispassTrigger = cols->triggerResults->accept(trigNames.triggerIndex(hltPath));
697  LogDebug("ExoticaValidation") << " preparing to call the plotters analysis";
698  an->analyze(ispassTrigger, "rec", matchesReco, theSumEt, dxys);
699  LogDebug("ExoticaValidation") << " called the plotter";
700  }
701  }
702 
703 }
#define LogDebug(id)
unsigned int _minCandidates
The minimum number of reco/gen candidates needed by the analysis.
int i
Definition: DBlmapReader.cc:9
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:213
std::map< unsigned int, std::string > _genCut
gen/rec objects cuts
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool accept() const
Has at least one path accepted the event?
std::vector< HLTExoticaPlotter > _plotters
The plotters: managers of each hlt path where the plots are done.
void fillHist(const std::string &source, const std::string &objType, const std::string &variable, const float &value)
const reco::GenParticleCollection * genParticles
std::map< unsigned int, std::string > _genCut_leading
gen/rec pt-leading objects cuts
std::map< unsigned int, edm::InputTag > _recLabels
std::map< unsigned int, std::string > _recCut_leading
unsigned int triggerIndex(std::string const &name) const
Definition: TriggerNames.cc:32
const reco::BeamSpot * bs
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int j
Definition: DBlmapReader.cc:9
virtual const Point & vertex() const =0
vertex position
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
void getHandlesToObjects(const edm::Event &iEvent, EVTColContainer *col)
Gets the collections themselves.
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::string _analysisname
The name of this sub-analysis.
std::map< unsigned int, StringCutObjectSelector< reco::GenParticle > * > _genSelectorMap
double y0() const
y coordinate
Definition: BeamSpot.h:66
tuple cout
Definition: gather_cfg.py:121
const Point & position() const
position
Definition: BeamSpot.h:62
void initSelector(const unsigned int &objtype)
Initializes the selectors of the objects based on which object it is.
static const std::string getTypeString(const unsigned int &objtype)
Tranform types into strings.
const edm::TriggerResults * triggerResults
tuple size
Write out results.
std::map< std::string, std::string > _shortpath2long
Relation between the short and long versions of the path.
void insertCandidates(const unsigned int &objtype, const EVTColContainer *col, std::vector< reco::LeafCandidate > *matches, std::map< int, double > &theSumEt, std::map< int, std::vector< const reco::Track * > > &trkObjs)
double x0() const
x coordinate
Definition: BeamSpot.h:64
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
void HLTExoticaSubAnalysis::beginJob ( void  )

Definition at line 179 of file HLTExoticaSubAnalysis.cc.

180 {
181 }
void HLTExoticaSubAnalysis::beginRun ( const edm::Run iRun,
const edm::EventSetup iEventSetup 
)

Construct the plotters right here. For that we need to create the _hltPaths vector.

Definition at line 263 of file HLTExoticaSubAnalysis.cc.

References _analysisname, _hltConfig, _hltPaths, _hltPathsToCheck, _hltProcessName, _plotters, _pset, _recLabels, _shortpath2long, gather_cfg::cout, newFWLiteAna::found, i, HLTConfigProvider::init(), j, LogDebug, LogTrace, cmsHarvester::path, chain::pattern, AlCaHLTBitMon_QueryRunRegistry::string, and HLTConfigProvider::triggerNames().

264 {
265 
266  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::beginRun()";
267 
270 
271  // Initialize the HLT config.
272  bool changedConfig(true);
273  if (!_hltConfig.init(iRun, iSetup, _hltProcessName, changedConfig)) {
274  edm::LogError("ExoticaValidation") << "HLTExoticaSubAnalysis::constructor(): "
275  << "Initialization of HLTConfigProvider failed!";
276  }
277 
278  // Parse the input paths to get them if they are in the table and associate
279  // them to the last filter of the path (in order to extract the objects).
280  _hltPaths.clear();
281  for (size_t i = 0; i < _hltPathsToCheck.size(); ++i) {
282  bool found = false;
283  TPRegexp pattern(_hltPathsToCheck[i]);
284 
285  // Loop over triggerNames from _hltConfig
286  for (size_t j = 0 ; j < _hltConfig.triggerNames().size(); ++j) {
287  std::string thetriggername = _hltConfig.triggerNames()[j];
288  if (TString(thetriggername).Contains(pattern)) {
289  _hltPaths.insert(thetriggername);
290  found = true;
291  }
292  if(verbose>2 && i==0)
293  std::cout << "--- TRIGGER PATH : " << thetriggername << std::endl;
294  }
295 
296  // Oh dear, the path we wanted seems to not be available
297  if (! found && verbose>2) {
298  edm::LogWarning("ExoticaValidation") << "HLTExoticaSubAnalysis::constructor(): In "
299  << _analysisname << " subfolder NOT found the path: '"
300  << _hltPathsToCheck[i] << "*'" ;
301  }
302  } // Close loop over paths to check.
303 
304  // At this point, _hltpaths contains the names of the paths to check
305  // that were found. Let's log it at trace level.
306  LogTrace("ExoticaValidation") << "SubAnalysis: " << _analysisname
307  << "\nHLT Trigger Paths found >>>";
308  for (std::set<std::string>::const_iterator iter = _hltPaths.begin();
309  iter != _hltPaths.end(); ++iter) {
310  LogTrace("ExoticaValidation") << (*iter) << "\n";
311  }
312 
313  // Initialize the plotters (analysers for each trigger path)
314  _plotters.clear();
315  for (std::set<std::string>::iterator iPath = _hltPaths.begin();
316  iPath != _hltPaths.end(); ++iPath) {
317  // Avoiding the dependence of the version number for the trigger paths
318  std::string path = * iPath;
319  std::string shortpath = path;
320  if (path.rfind("_v") < path.length()) {
321  shortpath = path.substr(0, path.rfind("_v"));
322  }
323  _shortpath2long[shortpath] = path;
324 
325  // Objects needed by the HLT path
326  // Thiago: instead of trying to decode the objects from the path,
327  // put the burden on the user to tell us which objects are needed.
328  //const std::vector<unsigned int> objsNeedHLT = this->getObjectsType(shortpath);
329  std::vector<unsigned int> objsNeedHLT;
330  for (std::map<unsigned int, edm::InputTag>::iterator it = _recLabels.begin() ;
331  it != _recLabels.end(); ++it) {
332  objsNeedHLT.push_back(it->first);
333  }
334 
335  /*std::vector<unsigned int> userInstantiate;
336  // Sanity check: the object needed by a trigger path should be
337  // introduced by the user via config python (_recLabels datamember)
338  for (std::map<unsigned int, edm::InputTag>::iterator it = _recLabels.begin() ;
339  it != _recLabels.end(); ++it) {
340  userInstantiate.push_back(it->first);
341  }
342  for (std::vector<unsigned int>::const_iterator it = objsNeedHLT.begin(); it != objsNeedHLT.end();
343  ++it) {
344  if (std::find(userInstantiate.begin(), userInstantiate.end(), *it) ==
345  userInstantiate.end()) {
346  edm::LogError("ExoticaValidation") << "In HLTExoticaSubAnalysis::beginRun, "
347  << "Incoherence found in the python configuration file!!\nThe SubAnalysis '"
348  << _analysisname << "' has been asked to evaluate the trigger path '"
349  << shortpath << "' (found it in 'hltPathsToCheck') BUT this path"
350  << " needs a '" << EVTColContainer::getTypeString(*it)
351  << "' which has not been instantiated ('recVariableLabels'"
352  << ")" ;
353  exit(-1); // This should probably throw an exception...
354  }
355  }
356  */
357  LogTrace("ExoticaValidation") << " --- " << shortpath;
358 
359  // The hlt path, the objects (electrons, muons, photons, ...)
360  // needed to evaluate the path are the argumens of the plotter
361  HLTExoticaPlotter analyzer(_pset, shortpath, objsNeedHLT);
362  _plotters.push_back(analyzer);
363  }// Okay, at this point we have prepared all the plotters.
364 
365 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
const std::vector< std::string > & triggerNames() const
names of trigger paths
std::vector< HLTExoticaPlotter > _plotters
The plotters: managers of each hlt path where the plots are done.
list pattern
Definition: chain.py:104
std::map< unsigned int, edm::InputTag > _recLabels
tuple path
else: Piece not in the list, fine.
edm::ParameterSet _pset
Internal, working copy of the PSet passed from above.
int j
Definition: DBlmapReader.cc:9
HLTConfigProvider _hltConfig
Interface to the HLT information.
std::vector< std::string > _hltPathsToCheck
The hlt paths to check for.
#define LogTrace(id)
std::string _hltProcessName
The labels of the object collections to be used in this analysis.
std::string _analysisname
The name of this sub-analysis.
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
tuple cout
Definition: gather_cfg.py:121
std::set< std::string > _hltPaths
The hlt paths found in the hltConfig.
std::map< std::string, std::string > _shortpath2long
Relation between the short and long versions of the path.
void HLTExoticaSubAnalysis::bookHist ( DQMStore::IBooker iBooker,
const std::string &  source,
const std::string &  objType,
const std::string &  variable 
)
private

The internal functions to book and fill histograms.

Definition at line 1029 of file HLTExoticaSubAnalysis.cc.

References _elements, _parametersDxy, _parametersEta, _parametersPhi, _parametersTurnOn, _parametersTurnOnSumEt, DQMStore::IBooker::book1D(), prof2calltree::edges, h, i, LogDebug, bookConverter::max, min(), mergeVDriftHistosByStation::name, source, AlCaHLTBitMon_QueryRunRegistry::string, and indexGen::title.

Referenced by subAnalysisBookHistos().

1033 {
1034  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::bookHist()";
1035  std::string sourceUpper = source;
1036  sourceUpper[0] = std::toupper(sourceUpper[0]);
1037  std::string name = source + objType + variable ;
1038  TH1F * h = 0;
1039 
1040  if (variable.find("SumEt") != std::string::npos) {
1041  std::string title = "Sum ET of " + sourceUpper + " " + objType;
1042  const size_t nBins = _parametersTurnOnSumEt.size() - 1;
1043  float * edges = new float[nBins + 1];
1044  for (size_t i = 0; i < nBins + 1; i++) {
1045  edges[i] = _parametersTurnOnSumEt[i];
1046  }
1047  h = new TH1F(name.c_str(), title.c_str(), nBins, edges);
1048  delete[] edges;
1049  }
1050 
1051  else if (variable.find("Dxy") != std::string::npos) {
1052  std::string title = "Dxy " + sourceUpper + " " + objType;
1053  int nBins = _parametersDxy[0];
1054  double min = _parametersDxy[1];
1055  double max = _parametersDxy[2];
1056  h = new TH1F(name.c_str(), title.c_str(), nBins, min, max);
1057  }
1058 
1059  else if (variable.find("MaxPt") != std::string::npos) {
1060  std::string desc = (variable == "MaxPt1") ? "Leading" : "Next-to-Leading";
1061  std::string title = "pT of " + desc + " " + sourceUpper + " " + objType;
1062  const size_t nBins = _parametersTurnOn.size() - 1;
1063  float * edges = new float[nBins + 1];
1064  for (size_t i = 0; i < nBins + 1; i++) {
1065  edges[i] = _parametersTurnOn[i];
1066  }
1067  h = new TH1F(name.c_str(), title.c_str(), nBins, edges);
1068  delete[] edges;
1069  }
1070 
1071  else {
1072  std::string symbol = (variable == "Eta") ? "#eta" : "#phi";
1073  std::string title = symbol + " of " + sourceUpper + " " + objType;
1074  std::vector<double> params = (variable == "Eta") ? _parametersEta : _parametersPhi;
1075  int nBins = (int)params[0];
1076  double min = params[1];
1077  double max = params[2];
1078  h = new TH1F(name.c_str(), title.c_str(), nBins, min, max);
1079  }
1080 
1081  h->Sumw2();
1082  // This is the trick, that takes a normal TH1F and puts it in in the DQM
1083  // machinery. Seems to be easy!
1084  // Updated to use the new iBooker machinery.
1085  _elements[name] = iBooker.book1D(name, h);
1086  delete h;
1087 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
std::vector< double > _parametersPhi
dictionary edges
std::vector< double > _parametersEta
Some kinematical parameters.
std::vector< double > _parametersDxy
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
T min(T a, T b)
Definition: MathUtil.h:58
std::vector< double > _parametersTurnOn
std::map< std::string, MonitorElement * > _elements
Structure of the MonitorElements.
std::vector< double > _parametersTurnOnSumEt
static std::string const source
Definition: EdmProvDump.cc:42
void HLTExoticaSubAnalysis::fillHist ( const std::string &  source,
const std::string &  objType,
const std::string &  variable,
const float &  value 
)
private

Definition at line 1090 of file HLTExoticaSubAnalysis.cc.

References _elements, LogDebug, mergeVDriftHistosByStation::name, source, AlCaHLTBitMon_QueryRunRegistry::string, and relativeConstraints::value.

Referenced by analyze().

1094 {
1095  std::string sourceUpper = source;
1096  sourceUpper[0] = toupper(sourceUpper[0]);
1097  std::string name = source + objType + variable ;
1098 
1099  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::fillHist() " << name << " " << value;
1100  _elements[name]->Fill(value);
1101  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::fillHist() " << name << " worked";
1102 }
#define LogDebug(id)
std::map< std::string, MonitorElement * > _elements
Structure of the MonitorElements.
static std::string const source
Definition: EdmProvDump.cc:42
void HLTExoticaSubAnalysis::getHandlesToObjects ( const edm::Event iEvent,
EVTColContainer col 
)
private

Gets the collections themselves.

Definition at line 918 of file HLTExoticaSubAnalysis.cc.

References _bsToken, _genParticleToken, _tokens, _trigResultsToken, EVTColContainer::bs, EVTColContainer::CALOJET, EVTColContainer::CALOMET, EVTColContainer::ELEC, EVTColContainer::GENMET, EVTColContainer::genParticles, edm::Event::getByToken(), EVTColContainer::isCommonInit(), edm::HandleBase::isValid(), EVTColContainer::L1MET, LogDebug, EVTColContainer::MET, EVTColContainer::MUON, EVTColContainer::MUTRK, EVTColContainer::PFJET, EVTColContainer::PFMET, EVTColContainer::PFMHT, EVTColContainer::PFTAU, EVTColContainer::PHOTON, edm::Handle< T >::product(), EVTColContainer::set(), EVTColContainer::setPFMHT(), EVTColContainer::TRACK, and EVTColContainer::triggerResults.

Referenced by analyze().

919 {
920  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::getHandlesToObjects()";
921 
922  if (! col->isCommonInit()) {
923  // Extract the trigger results (path info, pass,...)
925  iEvent.getByToken(_trigResultsToken, trigResults);
926  if (trigResults.isValid()) {
927  col->triggerResults = trigResults.product();
928  LogDebug("ExoticaValidation") << "Added handle to triggerResults";
929  }
930 
931  // Extract the genParticles
933  iEvent.getByToken(_genParticleToken, genPart);
934  if (genPart.isValid()) {
935  col->genParticles = genPart.product();
936  LogDebug("ExoticaValidation") << "Added handle to genParticles";
937  }
938 
939  // BeamSpot for dxy
941  iEvent.getByToken(_bsToken, bsHandle);
942  if (bsHandle.isValid()) {
943  col->bs = bsHandle.product();
944  }
945  }
946 
947  // Loop over the tokens and extract all other objects
948  LogDebug("ExoticaValidation") << "We have got " << _tokens.size() << "tokens";
949  for (std::map<unsigned int, edm::EDGetToken>::iterator it = _tokens.begin();
950  it != _tokens.end(); ++it) {
951  if (it->first == EVTColContainer::MUON) {
953  iEvent.getByToken(it->second, theHandle);
954  if (theHandle.isValid()) col->set(theHandle.product());
955  }
956  else if (it->first == EVTColContainer::MUTRK) {
958  iEvent.getByToken(it->second, theHandle);
959  if (theHandle.isValid()) col->set(theHandle.product());
960  }
961  else if (it->first == EVTColContainer::TRACK) {
963  iEvent.getByToken(it->second, theHandle);
964  if (theHandle.isValid()) col->set(theHandle.product());
965  }
966  else if (it->first == EVTColContainer::ELEC) {
968  iEvent.getByToken(it->second, theHandle);
969  if (theHandle.isValid()) col->set(theHandle.product());
970  }
971  else if (it->first == EVTColContainer::PHOTON) {
973  iEvent.getByToken(it->second, theHandle);
974  if (theHandle.isValid()) col->set(theHandle.product());
975  }
976  else if (it->first == EVTColContainer::MET) {
978  iEvent.getByToken(it->second, theHandle);
979  if (theHandle.isValid()) col->set(theHandle.product());
980  }
981  else if (it->first == EVTColContainer::PFMET) {
983  iEvent.getByToken(it->second, theHandle);
984  if (theHandle.isValid()) col->set(theHandle.product());
985  }
986  else if (it->first == EVTColContainer::PFMHT) {
988  iEvent.getByToken(it->second, theHandle);
989  if (theHandle.isValid()) col->setPFMHT(theHandle.product());
990  }
991  else if (it->first == EVTColContainer::GENMET) {
993  iEvent.getByToken(it->second, theHandle);
994  if (theHandle.isValid()) col->set(theHandle.product());
995  }
996  else if (it->first == EVTColContainer::CALOMET) {
998  iEvent.getByToken(it->second, theHandle);
999  if (theHandle.isValid()) col->set(theHandle.product());
1000  }
1001  else if (it->first == EVTColContainer::L1MET) {
1003  iEvent.getByToken(it->second, theHandle);
1004  if (theHandle.isValid()) col->set(theHandle.product());
1005  }
1006  else if (it->first == EVTColContainer::PFTAU) {
1008  iEvent.getByToken(it->second, theHandle);
1009  if (theHandle.isValid()) col->set(theHandle.product());
1010  }
1011  else if (it->first == EVTColContainer::PFJET) {
1013  iEvent.getByToken(it->second, theHandle);
1014  if (theHandle.isValid()) col->set(theHandle.product());
1015  }
1016  else if (it->first == EVTColContainer::CALOJET) {
1018  iEvent.getByToken(it->second, theHandle);
1019  if (theHandle.isValid()) col->set(theHandle.product());
1020  }
1021  else {
1022  edm::LogError("ExoticaValidation") << "HLTExoticaSubAnalysis::getHandlesToObjects "
1023  << " NOT IMPLEMENTED (yet) ERROR: '" << it->first << "'";
1024  }
1025  }
1026 }
#define LogDebug(id)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
const reco::GenParticleCollection * genParticles
edm::EDGetTokenT< edm::TriggerResults > _trigResultsToken
void setPFMHT(const reco::PFMETCollection *v)
void set(const reco::MuonCollection *v)
Setter: multiple overloaded function.
const reco::BeamSpot * bs
bool isValid() const
Definition: HandleBase.h:75
T const * product() const
Definition: Handle.h:81
std::map< unsigned int, edm::EDGetToken > _tokens
edm::EDGetTokenT< reco::GenParticleCollection > _genParticleToken
And also the tokens to get the object collections.
const edm::TriggerResults * triggerResults
edm::EDGetTokenT< reco::BeamSpot > _bsToken
void HLTExoticaSubAnalysis::getNamesOfObjects ( const edm::ParameterSet anpset)
private

Creates the maps that map which collection should come from which label.

Definition at line 746 of file HLTExoticaSubAnalysis.cc.

References _analysisname, _genSelectorMap, _recLabels, EVTColContainer::CALOJET, EVTColContainer::CALOMET, EVTColContainer::ELEC, edm::ParameterSet::exists(), EVTColContainer::GENMET, edm::ParameterSet::getParameter(), EVTColContainer::L1MET, LogDebug, EVTColContainer::MET, EVTColContainer::MUON, EVTColContainer::MUTRK, EVTColContainer::PFJET, EVTColContainer::PFMET, EVTColContainer::PFMHT, EVTColContainer::PFTAU, EVTColContainer::PHOTON, and EVTColContainer::TRACK.

Referenced by HLTExoticaSubAnalysis().

747 {
748  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::getNamesOfObjects()";
749 
750  if (anpset.exists("recMuonLabel")) {
751  _recLabels[EVTColContainer::MUON] = anpset.getParameter<edm::InputTag>("recMuonLabel");
753  }
754  if (anpset.exists("recMuonTrkLabel")) {
755  _recLabels[EVTColContainer::MUTRK] = anpset.getParameter<edm::InputTag>("recMuonTrkLabel");
757  }
758  if (anpset.exists("recTrackLabel")) {
759  _recLabels[EVTColContainer::TRACK] = anpset.getParameter<edm::InputTag>("recTrackLabel");
761  }
762  if (anpset.exists("recElecLabel")) {
763  _recLabels[EVTColContainer::ELEC] = anpset.getParameter<edm::InputTag>("recElecLabel");
765  }
766  if (anpset.exists("recPhotonLabel")) {
767  _recLabels[EVTColContainer::PHOTON] = anpset.getParameter<edm::InputTag>("recPhotonLabel");
769  }
770  if (anpset.exists("recMETLabel")) {
771  _recLabels[EVTColContainer::MET] = anpset.getParameter<edm::InputTag>("recMETLabel");
773  }
774  if (anpset.exists("recPFMETLabel")) {
775  _recLabels[EVTColContainer::PFMET] = anpset.getParameter<edm::InputTag>("recPFMETLabel");
777  }
778  if (anpset.exists("recPFMHTLabel")) {
779  _recLabels[EVTColContainer::PFMHT] = anpset.getParameter<edm::InputTag>("recPFMHTLabel");
781  }
782  if (anpset.exists("genMETLabel")) {
785  }
786  if (anpset.exists("recCaloMETLabel")) {
787  _recLabels[EVTColContainer::CALOMET] = anpset.getParameter<edm::InputTag>("recCaloMETLabel");
789  }
790  if (anpset.exists("hltMETLabel")) {
793  }
794  if (anpset.exists("l1METLabel")) {
797  }
798  if (anpset.exists("recPFTauLabel")) {
799  _recLabels[EVTColContainer::PFTAU] = anpset.getParameter<edm::InputTag>("recPFTauLabel");
801  }
802  if (anpset.exists("recPFJetLabel")) {
803  _recLabels[EVTColContainer::PFJET] = anpset.getParameter<edm::InputTag>("recPFJetLabel");
805  }
806  if (anpset.exists("recCaloJetLabel")) {
807  _recLabels[EVTColContainer::CALOJET] = anpset.getParameter<edm::InputTag>("recCaloJetLabel");
809  }
810 
811  if (_recLabels.size() < 1) {
812  edm::LogError("ExoticaValidation") << "HLTExoticaSubAnalysis::getNamesOfObjects, "
813  << "Not included any object (recMuonLabel, recElecLabel, ...) "
814  << "in the analysis " << _analysisname;
815  return;
816  }
817 }
#define LogDebug(id)
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::map< unsigned int, edm::InputTag > _recLabels
std::string _analysisname
The name of this sub-analysis.
std::map< unsigned int, StringCutObjectSelector< reco::GenParticle > * > _genSelectorMap
const std::vector< unsigned int > HLTExoticaSubAnalysis::getObjectsType ( const std::string &  hltpath) const
private

closes analyze method

Return the objects (muons,electrons,photons,...) needed by a HLT path. Will in general return: 0 for muon, 1 for electron, 2 for photon, 3 for PFMET, 4 for PFTau, 5 for Jet. Notice that this function is really based on a parsing of the name of the path; any incongruences there may lead to problems.

Definition at line 707 of file HLTExoticaSubAnalysis.cc.

References EVTColContainer::CALOJET, EVTColContainer::CALOMET, EVTColContainer::ELEC, EVTColContainer::GENMET, EVTColContainer::getTypeString(), i, EVTColContainer::L1MET, LogDebug, EVTColContainer::MET, EVTColContainer::MUON, EVTColContainer::MUTRK, EVTColContainer::PFJET, EVTColContainer::PFMET, EVTColContainer::PFMHT, EVTColContainer::PFTAU, EVTColContainer::PHOTON, AlCaHLTBitMon_QueryRunRegistry::string, and EVTColContainer::TRACK.

708 {
709  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::getObjectsType()";
710 
711  static const unsigned int objSize = 14;
712  static const unsigned int objtriggernames[] = {
727  };
728 
729  std::set<unsigned int> objsType;
730  // The object to deal has to be entered via the config .py
731  for (unsigned int i = 0; i < objSize; ++i) {
732  //std::cout << "(5) Gonna call with " << objtriggernames[i] << std::endl;
733  std::string objTypeStr = EVTColContainer::getTypeString(objtriggernames[i]);
734  // Check if it is needed this object for this trigger
735  if (! TString(hltPath).Contains(objTypeStr)) {
736  continue;
737  }
738 
739  objsType.insert(objtriggernames[i]);
740  }
741 
742  return std::vector<unsigned int>(objsType.begin(), objsType.end());
743 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
static const std::string getTypeString(const unsigned int &objtype)
Tranform types into strings.
void HLTExoticaSubAnalysis::initSelector ( const unsigned int &  objtype)
private

Initializes the selectors of the objects based on which object it is.

Definition at line 1105 of file HLTExoticaSubAnalysis.cc.

References _genMETSelector, _l1METSelector, _recCaloJetSelector, _recCaloMETSelector, _recCut, _recElecSelector, _recMETSelector, _recMuonSelector, _recMuonTrkSelector, _recPFJetSelector, _recPFMETSelector, _recPFMHTSelector, _recPFTauSelector, _recPhotonSelector, _recTrackSelector, EVTColContainer::CALOJET, EVTColContainer::CALOMET, EVTColContainer::ELEC, EVTColContainer::GENMET, EVTColContainer::L1MET, LogDebug, EVTColContainer::MET, EVTColContainer::MUON, EVTColContainer::MUTRK, EVTColContainer::PFJET, EVTColContainer::PFMET, EVTColContainer::PFMHT, EVTColContainer::PFTAU, EVTColContainer::PHOTON, and EVTColContainer::TRACK.

Referenced by analyze().

1106 {
1107 
1108  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::initSelector()";
1109 
1110  if (objtype == EVTColContainer::MUON && _recMuonSelector == 0) {
1112  } else if (objtype == EVTColContainer::MUTRK && _recMuonTrkSelector == 0) {
1114  } else if (objtype == EVTColContainer::TRACK && _recTrackSelector == 0) {
1116  } else if (objtype == EVTColContainer::ELEC && _recElecSelector == 0) {
1118  } else if (objtype == EVTColContainer::PHOTON && _recPhotonSelector == 0) {
1120  } else if (objtype == EVTColContainer::MET && _recMETSelector == 0) {
1122  } else if (objtype == EVTColContainer::PFMET && _recPFMETSelector == 0) {
1124  } else if (objtype == EVTColContainer::PFMHT && _recPFMHTSelector == 0) {
1126  } else if (objtype == EVTColContainer::GENMET && _genMETSelector == 0) {
1128  } else if (objtype == EVTColContainer::CALOMET && _recCaloMETSelector == 0) {
1130  } else if (objtype == EVTColContainer::L1MET && _l1METSelector == 0) {
1132  } else if (objtype == EVTColContainer::PFTAU && _recPFTauSelector == 0) {
1134  } else if (objtype == EVTColContainer::PFJET && _recPFJetSelector == 0) {
1136  } else if (objtype == EVTColContainer::CALOJET && _recCaloJetSelector == 0) {
1138  }
1139  /* else
1140  {
1141  FIXME: ERROR NOT IMPLEMENTED
1142  }*/
1143 }
#define LogDebug(id)
StringCutObjectSelector< reco::Track > * _recMuonTrkSelector
StringCutObjectSelector< reco::CaloJet > * _recCaloJetSelector
StringCutObjectSelector< reco::GenMET > * _genMETSelector
StringCutObjectSelector< reco::MET > * _recMETSelector
StringCutObjectSelector< reco::PFJet > * _recPFJetSelector
StringCutObjectSelector< reco::GsfElectron > * _recElecSelector
StringCutObjectSelector< reco::Muon > * _recMuonSelector
std::map< unsigned int, std::string > _recCut
StringCutObjectSelector< reco::Track > * _recTrackSelector
StringCutObjectSelector< reco::PFMET > * _recPFMETSelector
StringCutObjectSelector< reco::Photon > * _recPhotonSelector
StringCutObjectSelector< l1extra::L1EtMissParticle > * _l1METSelector
StringCutObjectSelector< reco::CaloMET > * _recCaloMETSelector
StringCutObjectSelector< reco::PFMET > * _recPFMHTSelector
StringCutObjectSelector< reco::PFTau > * _recPFTauSelector
void HLTExoticaSubAnalysis::insertCandidates ( const unsigned int &  objtype,
const EVTColContainer col,
std::vector< reco::LeafCandidate > *  matches,
std::map< int, double > &  theSumEt,
std::map< int, std::vector< const reco::Track * > > &  trkObjs 
)
private

This function applies the selectors initialized previously to the objects, and matches the passing objects to HLT objects.

This is a special case. Passing a PFMET* to the constructor of MatchStruct will trigger the usage of the special constructor which also sets the sumEt member.

Definition at line 1146 of file HLTExoticaSubAnalysis.cc.

References _genMETSelector, _l1METSelector, _recCaloJetSelector, _recCaloMETSelector, _recElecSelector, _recMuonSelector, _recMuonTrkSelector, _recPFJetSelector, _recPFMETSelector, _recPFMHTSelector, _recPFTauSelector, _recPhotonSelector, _recTrackSelector, EVTColContainer::CALOJET, EVTColContainer::caloJets, EVTColContainer::CALOMET, EVTColContainer::caloMETs, EVTColContainer::ELEC, EVTColContainer::electrons, EVTColContainer::GENMET, EVTColContainer::genMETs, i, EVTColContainer::L1MET, EVTColContainer::l1METs, LogDebug, visualization-live-secondInstance_cfg::m, EVTColContainer::MUON, EVTColContainer::muons, EVTColContainer::MUTRK, EVTColContainer::PFJET, EVTColContainer::pfJets, EVTColContainer::PFMET, EVTColContainer::pfMETs, EVTColContainer::PFMHT, EVTColContainer::pfMHTs, EVTColContainer::PFTAU, EVTColContainer::pfTaus, EVTColContainer::PHOTON, EVTColContainer::photons, EVTColContainer::TRACK, and EVTColContainer::tracks.

Referenced by analyze().

1148 {
1149 
1150  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::insertCandidates()";
1151 
1152  theSumEt[objType] = -1;
1153 
1154  if (objType == EVTColContainer::MUON) {
1155  for (size_t i = 0; i < cols->muons->size(); i++) {
1156  LogDebug("ExoticaValidation") << "Inserting muon " << i ;
1157  if (_recMuonSelector->operator()(cols->muons->at(i))) {
1158  reco::LeafCandidate m(0, cols->muons->at(i).p4(), cols->muons->at(i).vertex(), objType, 0, true);
1159  matches->push_back(m);
1160 
1161  // for making dxy plots
1162  trkObjs[objType].push_back(cols->muons->at(i).bestTrack());
1163  }
1164  }
1165  } else if (objType == EVTColContainer::MUTRK) {
1166  for (size_t i = 0; i < cols->tracks->size(); i++) {
1167  LogDebug("ExoticaValidation") << "Inserting muonTrack " << i ;
1168  if (_recMuonTrkSelector->operator()(cols->tracks->at(i))) {
1169  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>> mom4;
1170  ROOT::Math::XYZVector mom3 = cols->tracks->at(i).momentum();
1171  mom4.SetXYZT(mom3.x(),mom3.y(),mom3.z(),mom3.r());
1172  reco::LeafCandidate m(0, mom4, cols->tracks->at(i).vertex(), objType, 0, true);
1173  matches->push_back(m);
1174  }
1175  }
1176  } else if (objType == EVTColContainer::TRACK) {
1177  for (size_t i = 0; i < cols->tracks->size(); i++) {
1178  LogDebug("ExoticaValidation") << "Inserting Track " << i ;
1179  if (_recTrackSelector->operator()(cols->tracks->at(i))) {
1180  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>> mom4;
1181  ROOT::Math::XYZVector mom3 = cols->tracks->at(i).momentum();
1182  mom4.SetXYZT(mom3.x(),mom3.y(),mom3.z(),mom3.r());
1183  reco::LeafCandidate m(0, mom4, cols->tracks->at(i).vertex(), objType, 0, true);
1184  matches->push_back(m);
1185  }
1186  }
1187  } else if (objType == EVTColContainer::ELEC) {
1188  for (size_t i = 0; i < cols->electrons->size(); i++) {
1189  LogDebug("ExoticaValidation") << "Inserting electron " << i ;
1190  if (_recElecSelector->operator()(cols->electrons->at(i))) {
1191  reco::LeafCandidate m(0, cols->electrons->at(i).p4(), cols->electrons->at(i).vertex(), objType, 0, true);
1192  matches->push_back(m);
1193 
1194  // for making dxy plots
1195  trkObjs[objType].push_back(cols->electrons->at(i).bestTrack());
1196  }
1197  }
1198  } else if (objType == EVTColContainer::PHOTON) {
1199  for (size_t i = 0; i < cols->photons->size(); i++) {
1200  LogDebug("ExoticaValidation") << "Inserting photon " << i ;
1201  if (_recPhotonSelector->operator()(cols->photons->at(i))) {
1202  reco::LeafCandidate m(0, cols->photons->at(i).p4(), cols->photons->at(i).vertex(), objType, 0, true);
1203  matches->push_back(m);
1204  }
1205  }
1206  } else if (objType == EVTColContainer::PFMET) {
1210  for (size_t i = 0; i < cols->pfMETs->size(); i++) {
1211  LogDebug("ExoticaValidation") << "Inserting PFMET " << i ;
1212  if (_recPFMETSelector->operator()(cols->pfMETs->at(i))) {
1213  reco::LeafCandidate m(0, cols->pfMETs->at(i).p4(), cols->pfMETs->at(i).vertex(), objType, 0, true);
1214  matches->push_back(m);
1215  if(i==0) theSumEt[objType] = cols->pfMETs->at(i).sumEt();
1216  }
1217  }
1218  } else if (objType == EVTColContainer::PFMHT) {
1219  for (size_t i = 0; i < cols->pfMHTs->size(); i++) {
1220  LogDebug("ExoticaValidation") << "Inserting PFMHT " << i ;
1221  if (_recPFMHTSelector->operator()(cols->pfMHTs->at(i))) {
1222  reco::LeafCandidate m(0, cols->pfMHTs->at(i).p4(), cols->pfMHTs->at(i).vertex(), objType, 0, true);
1223  matches->push_back(m);
1224  if(i==0) theSumEt[objType] = cols->pfMHTs->at(i).sumEt();
1225  }
1226  }
1227  } else if (objType == EVTColContainer::GENMET) {
1228  for (size_t i = 0; i < cols->genMETs->size(); i++) {
1229  LogDebug("ExoticaValidation") << "Inserting GENMET " << i ;
1230  if (_genMETSelector->operator()(cols->genMETs->at(i))) {
1231  reco::LeafCandidate m(0, cols->genMETs->at(i).p4(), cols->genMETs->at(i).vertex(), objType, 0, true);
1232  matches->push_back(m);
1233  if(i==0) theSumEt[objType] = cols->genMETs->at(i).sumEt();
1234  }
1235  }
1236  } else if (objType == EVTColContainer::CALOMET) {
1237  for (size_t i = 0; i < cols->caloMETs->size(); i++) {
1238  LogDebug("ExoticaValidation") << "Inserting CALOMET " << i ;
1239  if (_recCaloMETSelector->operator()(cols->caloMETs->at(i))) {
1240  reco::LeafCandidate m(0, cols->caloMETs->at(i).p4(), cols->caloMETs->at(i).vertex(), objType, 0, true);
1241  matches->push_back(m);
1242  if(i==0) theSumEt[objType] = cols->caloMETs->at(i).sumEt();
1243  }
1244  }
1245  } else if (objType == EVTColContainer::L1MET) {
1246  for (size_t i = 0; i < cols->l1METs->size(); i++) {
1247  LogDebug("ExoticaValidation") << "Inserting L1MET " << i ;
1248  if (_l1METSelector->operator()(cols->l1METs->at(i))) {
1249  reco::LeafCandidate m(0, cols->l1METs->at(i).p4(), cols->l1METs->at(i).vertex(), objType, 0, true);
1250  matches->push_back(m);
1251  if(i==0) theSumEt[objType] = cols->l1METs->at(i).etTotal();
1252  }
1253  }
1254  } else if (objType == EVTColContainer::PFTAU) {
1255  for (size_t i = 0; i < cols->pfTaus->size(); i++) {
1256  LogDebug("ExoticaValidation") << "Inserting PFtau " << i ;
1257  if (_recPFTauSelector->operator()(cols->pfTaus->at(i))) {
1258  reco::LeafCandidate m(0, cols->pfTaus->at(i).p4(), cols->pfTaus->at(i).vertex(), objType, 0, true);
1259  matches->push_back(m);
1260  }
1261  }
1262  } else if (objType == EVTColContainer::PFJET) {
1263  for (size_t i = 0; i < cols->pfJets->size(); i++) {
1264  LogDebug("ExoticaValidation") << "Inserting jet " << i ;
1265  if (_recPFJetSelector->operator()(cols->pfJets->at(i))) {
1266  reco::LeafCandidate m(0, cols->pfJets->at(i).p4(), cols->pfJets->at(i).vertex(), objType, 0, true);
1267  matches->push_back(m);
1268  }
1269  }
1270  } else if (objType == EVTColContainer::CALOJET) {
1271  for (size_t i = 0; i < cols->caloJets->size(); i++) {
1272  LogDebug("ExoticaValidation") << "Inserting jet " << i ;
1273  if (_recCaloJetSelector->operator()(cols->caloJets->at(i))) {
1274  reco::LeafCandidate m(0, cols->caloJets->at(i).p4(), cols->caloJets->at(i).vertex(), objType, 0, true);
1275  matches->push_back(m);
1276  }
1277  }
1278  }
1279 
1280  /* else
1281  {
1282  FIXME: ERROR NOT IMPLEMENTED
1283  }*/
1284 }
#define LogDebug(id)
StringCutObjectSelector< reco::Track > * _recMuonTrkSelector
int i
Definition: DBlmapReader.cc:9
StringCutObjectSelector< reco::CaloJet > * _recCaloJetSelector
StringCutObjectSelector< reco::GenMET > * _genMETSelector
StringCutObjectSelector< reco::PFJet > * _recPFJetSelector
StringCutObjectSelector< reco::GsfElectron > * _recElecSelector
StringCutObjectSelector< reco::Muon > * _recMuonSelector
StringCutObjectSelector< reco::Track > * _recTrackSelector
Transform3DPJ::Vector XYZVector
StringCutObjectSelector< reco::PFMET > * _recPFMETSelector
StringCutObjectSelector< reco::Photon > * _recPhotonSelector
StringCutObjectSelector< l1extra::L1EtMissParticle > * _l1METSelector
StringCutObjectSelector< reco::CaloMET > * _recCaloMETSelector
StringCutObjectSelector< reco::PFMET > * _recPFMHTSelector
StringCutObjectSelector< reco::PFTau > * _recPFTauSelector
void HLTExoticaSubAnalysis::registerConsumes ( edm::ConsumesCollector consCollector)
private

Registers consumption of objects.

Definition at line 821 of file HLTExoticaSubAnalysis.cc.

References _beamSpotLabel, _bsToken, _genParticleLabel, _genParticleToken, _recLabels, _tokens, _trigResultsLabel, _trigResultsToken, EVTColContainer::CALOJET, EVTColContainer::CALOMET, edm::ConsumesCollector::consumes(), EVTColContainer::ELEC, EVTColContainer::GENMET, EVTColContainer::L1MET, LogDebug, EVTColContainer::MET, EVTColContainer::MUON, EVTColContainer::MUTRK, EVTColContainer::PFJET, EVTColContainer::PFMET, EVTColContainer::PFMHT, EVTColContainer::PFTAU, EVTColContainer::PHOTON, and EVTColContainer::TRACK.

Referenced by HLTExoticaSubAnalysis().

822 {
823  // Register that we are getting genParticles
825 
826  // Register that we are getting the trigger results
828 
829  // Register beamspot
830  _bsToken = iC.consumes<reco::BeamSpot>(_beamSpotLabel);
831 
832  // Loop over _recLabels, see what we need, and register.
833  // Then save the registered token in _tokens.
834  // Remember: _recLabels is a map<uint, edm::InputTag>
835  // Remember: _tokens is a map<uint, edm::EDGetToken>
836  LogDebug("ExoticaValidation") << "We have got " << _recLabels.size() << "recLabels";
837  for (std::map<unsigned int, edm::InputTag>::iterator it = _recLabels.begin();
838  it != _recLabels.end(); ++it) {
839  if (it->first == EVTColContainer::MUON) {
840  edm::EDGetTokenT<reco::MuonCollection> particularToken = iC.consumes<reco::MuonCollection>(it->second);
841  edm::EDGetToken token(particularToken);
842  _tokens[it->first] = token;
843  }
844  else if (it->first == EVTColContainer::MUTRK) {
845  edm::EDGetTokenT<reco::TrackCollection> particularToken = iC.consumes<reco::TrackCollection>(it->second);
846  edm::EDGetToken token(particularToken);
847  _tokens[it->first] = token;
848  }
849  else if (it->first == EVTColContainer::TRACK) {
850  edm::EDGetTokenT<reco::TrackCollection> particularToken = iC.consumes<reco::TrackCollection>(it->second);
851  edm::EDGetToken token(particularToken);
852  _tokens[it->first] = token;
853  }
854  else if (it->first == EVTColContainer::ELEC) {
855  edm::EDGetTokenT<reco::GsfElectronCollection> particularToken = iC.consumes<reco::GsfElectronCollection>(it->second);
856  edm::EDGetToken token(particularToken);
857  _tokens[it->first] = token;
858  }
859  else if (it->first == EVTColContainer::PHOTON) {
860  edm::EDGetTokenT<reco::PhotonCollection> particularToken = iC.consumes<reco::PhotonCollection>(it->second);
861  edm::EDGetToken token(particularToken);
862  _tokens[it->first] = token;
863  }
864  else if (it->first == EVTColContainer::MET) {
865  edm::EDGetTokenT<reco::METCollection> particularToken = iC.consumes<reco::METCollection>(it->second);
866  edm::EDGetToken token(particularToken);
867  _tokens[it->first] = token;
868  }
869  else if (it->first == EVTColContainer::PFMET) {
870  edm::EDGetTokenT<reco::PFMETCollection> particularToken = iC.consumes<reco::PFMETCollection>(it->second);
871  edm::EDGetToken token(particularToken);
872  _tokens[it->first] = token;
873  }
874  else if (it->first == EVTColContainer::PFMHT) {
875  edm::EDGetTokenT<reco::PFMETCollection> particularToken = iC.consumes<reco::PFMETCollection>(it->second);
876  edm::EDGetToken token(particularToken);
877  _tokens[it->first] = token;
878  }
879  else if (it->first == EVTColContainer::GENMET) {
880  edm::EDGetTokenT<reco::GenMETCollection> particularToken = iC.consumes<reco::GenMETCollection>(it->second);
881  edm::EDGetToken token(particularToken);
882  _tokens[it->first] = token;
883  }
884  else if (it->first == EVTColContainer::CALOMET) {
885  edm::EDGetTokenT<reco::CaloMETCollection> particularToken = iC.consumes<reco::CaloMETCollection>(it->second);
886  edm::EDGetToken token(particularToken);
887  _tokens[it->first] = token;
888  }
889  else if (it->first == EVTColContainer::L1MET) {
891  edm::EDGetToken token(particularToken);
892  _tokens[it->first] = token;
893  }
894  else if (it->first == EVTColContainer::PFTAU) {
895  edm::EDGetTokenT<reco::PFTauCollection> particularToken = iC.consumes<reco::PFTauCollection>(it->second);
896  edm::EDGetToken token(particularToken);
897  _tokens[it->first] = token;
898  }
899  else if (it->first == EVTColContainer::PFJET) {
900  edm::EDGetTokenT<reco::PFJetCollection> particularToken = iC.consumes<reco::PFJetCollection>(it->second);
901  edm::EDGetToken token(particularToken);
902  _tokens[it->first] = token;
903  }
904  else if (it->first == EVTColContainer::CALOJET) {
905  edm::EDGetTokenT<reco::CaloJetCollection> particularToken = iC.consumes<reco::CaloJetCollection>(it->second);
906  edm::EDGetToken token(particularToken);
907  _tokens[it->first] = token;
908  }
909  else {
910  edm::LogError("ExoticaValidation") << "HLTExoticaSubAnalysis::registerConsumes"
911  << " NOT IMPLEMENTED (yet) ERROR: '" << it->second.label() << "'";
912  }
913  }
914 
915 }
#define LogDebug(id)
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
std::vector< PFTau > PFTauCollection
collection of PFTau objects
Definition: PFTauFwd.h:9
std::vector< reco::GenMET > GenMETCollection
collection of GenMET objects
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
edm::EDGetTokenT< edm::TriggerResults > _trigResultsToken
std::vector< reco::MET > METCollection
collection of MET objects
Definition: METCollection.h:23
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
std::map< unsigned int, edm::InputTag > _recLabels
std::vector< reco::CaloMET > CaloMETCollection
collection of CaloMET objects
std::vector< Photon > PhotonCollection
collectin of Photon objects
Definition: PhotonFwd.h:9
std::map< unsigned int, edm::EDGetToken > _tokens
std::vector< PFJet > PFJetCollection
collection of PFJet objects
std::vector< reco::PFMET > PFMETCollection
collection of PFMET objects
edm::EDGetTokenT< reco::GenParticleCollection > _genParticleToken
And also the tokens to get the object collections.
std::vector< L1EtMissParticle > L1EtMissParticleCollection
edm::EDGetTokenT< reco::BeamSpot > _bsToken
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects
void HLTExoticaSubAnalysis::subAnalysisBookHistos ( DQMStore::IBooker iBooker,
const edm::Run iRun,
const edm::EventSetup iSetup 
)

Method to book all relevant histograms in the DQMStore. Uses the IBooker interface for thread safety. Intended to be called from master object.

Definition at line 190 of file HLTExoticaSubAnalysis.cc.

References _analysisname, _plotters, _recLabels, bookHist(), EVTColContainer::ELEC, EVTColContainer::getTypeString(), i, LogDebug, EVTColContainer::MUON, DQMStore::IBooker::setCurrentFolder(), source, and AlCaHLTBitMon_QueryRunRegistry::string.

193 {
194 
195  LogDebug("ExoticaValidation") << "In HLTExoticaSubAnalysis::subAnalysisBookHistos()";
196 
197  // Create the folder structure inside HLT/Exotica
198  std::string baseDir = "HLT/Exotica/" + _analysisname + "/";
199  iBooker.setCurrentFolder(baseDir);
200 
201  // Book the gen/reco analysis-dependent histograms (denominators)
202  for (std::map<unsigned int, edm::InputTag>::const_iterator it = _recLabels.begin();
203  it != _recLabels.end(); ++it) {
204  const std::string objStr = EVTColContainer::getTypeString(it->first);
205  std::vector<std::string> sources(2);
206  sources[0] = "gen";
207  sources[1] = "rec";
208 
209  for (size_t i = 0; i < sources.size(); i++) {
210  std::string source = sources[i];
211 
212  if ( source == "gen" ) {
213  if ( TString(objStr).Contains("MET") ||
214  TString(objStr).Contains("MHT") ||
215  TString(objStr).Contains("Jet") ) {
216  continue;
217  } else {
218  bookHist(iBooker, source, objStr, "MaxPt1");
219  bookHist(iBooker, source, objStr, "MaxPt2");
220  bookHist(iBooker, source, objStr, "MaxPt3");
221  bookHist(iBooker, source, objStr, "Eta");
222  bookHist(iBooker, source, objStr, "Phi");
223 
224  // If the target is electron or muon,
225  // we will add Dxy plots.
226  if ( it->first == EVTColContainer::ELEC ||
227  it->first == EVTColContainer::MUON ) {
228  bookHist(iBooker, source, objStr, "Dxy");
229  }
230  }
231  } else { // reco
232  if ( TString(objStr).Contains("MET") ||
233  TString(objStr).Contains("MHT") ) {
234  bookHist(iBooker, source, objStr, "MaxPt1");
235  bookHist(iBooker, source, objStr, "SumEt");
236  } else {
237  bookHist(iBooker, source, objStr, "MaxPt1");
238  bookHist(iBooker, source, objStr, "MaxPt2");
239  bookHist(iBooker, source, objStr, "MaxPt3");
240  bookHist(iBooker, source, objStr, "Eta");
241  bookHist(iBooker, source, objStr, "Phi");
242 
243  // If the target is electron or muon,
244  // we will add Dxy plots.
245  if ( it->first == EVTColContainer::ELEC ||
246  it->first == EVTColContainer::MUON ) {
247  bookHist(iBooker, source, objStr, "Dxy");
248  }
249  }
250  }
251 
252  }
253  } // closes loop in _recLabels
254 
255  // Call the plotterBookHistos() (which books all the path dependent histograms)
256  LogDebug("ExoticaValidation") << " number of plotters = " << _plotters.size();
257  for (std::vector<HLTExoticaPlotter>::iterator it = _plotters.begin();
258  it != _plotters.end(); ++it) {
259  it->plotterBookHistos(iBooker, iRun, iSetup);
260  }
261 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
std::vector< HLTExoticaPlotter > _plotters
The plotters: managers of each hlt path where the plots are done.
std::map< unsigned int, edm::InputTag > _recLabels
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
std::string _analysisname
The name of this sub-analysis.
static const std::string getTypeString(const unsigned int &objtype)
Tranform types into strings.
static std::string const source
Definition: EdmProvDump.cc:42
void bookHist(DQMStore::IBooker &iBooker, const std::string &source, const std::string &objType, const std::string &variable)
The internal functions to book and fill histograms.

Member Data Documentation

std::string HLTExoticaSubAnalysis::_analysisname
private

The name of this sub-analysis.

Definition at line 122 of file HLTExoticaSubAnalysis.h.

Referenced by analyze(), beginRun(), getNamesOfObjects(), and subAnalysisBookHistos().

edm::InputTag HLTExoticaSubAnalysis::_beamSpotLabel
private

Definition at line 138 of file HLTExoticaSubAnalysis.h.

Referenced by registerConsumes().

edm::EDGetTokenT<reco::BeamSpot> HLTExoticaSubAnalysis::_bsToken
private

Definition at line 143 of file HLTExoticaSubAnalysis.h.

Referenced by getHandlesToObjects(), and registerConsumes().

std::map<std::string, MonitorElement *> HLTExoticaSubAnalysis::_elements
private

Structure of the MonitorElements.

Definition at line 185 of file HLTExoticaSubAnalysis.h.

Referenced by bookHist(), and fillHist().

std::map<unsigned int, std::string> HLTExoticaSubAnalysis::_genCut
private

gen/rec objects cuts

Definition at line 154 of file HLTExoticaSubAnalysis.h.

Referenced by analyze(), and HLTExoticaSubAnalysis().

std::map<unsigned int, std::string> HLTExoticaSubAnalysis::_genCut_leading
private

gen/rec pt-leading objects cuts

Definition at line 157 of file HLTExoticaSubAnalysis.h.

Referenced by analyze(), and HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::GenMET>* HLTExoticaSubAnalysis::_genMETSelector
private

Definition at line 170 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

edm::InputTag HLTExoticaSubAnalysis::_genParticleLabel
private

Definition at line 136 of file HLTExoticaSubAnalysis.h.

Referenced by registerConsumes().

edm::EDGetTokenT<reco::GenParticleCollection> HLTExoticaSubAnalysis::_genParticleToken
private

And also the tokens to get the object collections.

Definition at line 141 of file HLTExoticaSubAnalysis.h.

Referenced by getHandlesToObjects(), and registerConsumes().

std::map<unsigned int, StringCutObjectSelector<reco::GenParticle> *> HLTExoticaSubAnalysis::_genSelectorMap
private

The concrete String selectors (use the string cuts introduced via the config python)

Definition at line 162 of file HLTExoticaSubAnalysis.h.

Referenced by analyze(), getNamesOfObjects(), and ~HLTExoticaSubAnalysis().

HLTConfigProvider HLTExoticaSubAnalysis::_hltConfig
private

Interface to the HLT information.

Definition at line 182 of file HLTExoticaSubAnalysis.h.

Referenced by beginRun().

std::set<std::string> HLTExoticaSubAnalysis::_hltPaths
private

The hlt paths found in the hltConfig.

Definition at line 130 of file HLTExoticaSubAnalysis.h.

Referenced by beginRun().

std::vector<std::string> HLTExoticaSubAnalysis::_hltPathsToCheck
private

The hlt paths to check for.

Definition at line 128 of file HLTExoticaSubAnalysis.h.

Referenced by beginRun(), and HLTExoticaSubAnalysis().

std::string HLTExoticaSubAnalysis::_hltProcessName
private

The labels of the object collections to be used in this analysis.

Definition at line 135 of file HLTExoticaSubAnalysis.h.

Referenced by beginRun().

StringCutObjectSelector<l1extra::L1EtMissParticle>* HLTExoticaSubAnalysis::_l1METSelector
private

Definition at line 172 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

unsigned int HLTExoticaSubAnalysis::_minCandidates
private

The minimum number of reco/gen candidates needed by the analysis.

Definition at line 125 of file HLTExoticaSubAnalysis.h.

Referenced by analyze(), and HLTExoticaSubAnalysis().

std::vector<double> HLTExoticaSubAnalysis::_parametersDxy
private

Definition at line 151 of file HLTExoticaSubAnalysis.h.

Referenced by bookHist(), and HLTExoticaSubAnalysis().

std::vector<double> HLTExoticaSubAnalysis::_parametersEta
private

Some kinematical parameters.

Definition at line 147 of file HLTExoticaSubAnalysis.h.

Referenced by bookHist(), and HLTExoticaSubAnalysis().

std::vector<double> HLTExoticaSubAnalysis::_parametersPhi
private

Definition at line 148 of file HLTExoticaSubAnalysis.h.

Referenced by bookHist(), and HLTExoticaSubAnalysis().

std::vector<double> HLTExoticaSubAnalysis::_parametersTurnOn
private

Definition at line 149 of file HLTExoticaSubAnalysis.h.

Referenced by bookHist(), and HLTExoticaSubAnalysis().

std::vector<double> HLTExoticaSubAnalysis::_parametersTurnOnSumEt
private

Definition at line 150 of file HLTExoticaSubAnalysis.h.

Referenced by bookHist(), and HLTExoticaSubAnalysis().

std::vector<HLTExoticaPlotter> HLTExoticaSubAnalysis::_plotters
private

The plotters: managers of each hlt path where the plots are done.

Definition at line 179 of file HLTExoticaSubAnalysis.h.

Referenced by analyze(), beginRun(), and subAnalysisBookHistos().

edm::ParameterSet HLTExoticaSubAnalysis::_pset
private

Internal, working copy of the PSet passed from above.

Definition at line 119 of file HLTExoticaSubAnalysis.h.

Referenced by beginRun(), and HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::CaloJet>* HLTExoticaSubAnalysis::_recCaloJetSelector
private

Definition at line 176 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::CaloMET>* HLTExoticaSubAnalysis::_recCaloMETSelector
private

Definition at line 171 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

std::map<unsigned int, std::string> HLTExoticaSubAnalysis::_recCut
private

Definition at line 155 of file HLTExoticaSubAnalysis.h.

Referenced by HLTExoticaSubAnalysis(), and initSelector().

std::map<unsigned int, std::string> HLTExoticaSubAnalysis::_recCut_leading
private

Definition at line 158 of file HLTExoticaSubAnalysis.h.

Referenced by analyze(), and HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::GsfElectron>* HLTExoticaSubAnalysis::_recElecSelector
private

Definition at line 166 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

std::map<unsigned int, edm::InputTag> HLTExoticaSubAnalysis::_recLabels
private
StringCutObjectSelector<reco::MET>* HLTExoticaSubAnalysis::_recMETSelector
private

Definition at line 167 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), and ~HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::Muon>* HLTExoticaSubAnalysis::_recMuonSelector
private

Definition at line 163 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::Track>* HLTExoticaSubAnalysis::_recMuonTrkSelector
private

Definition at line 164 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::PFJet>* HLTExoticaSubAnalysis::_recPFJetSelector
private

Definition at line 175 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::PFMET>* HLTExoticaSubAnalysis::_recPFMETSelector
private

Definition at line 168 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::PFMET>* HLTExoticaSubAnalysis::_recPFMHTSelector
private

Definition at line 169 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::PFTau>* HLTExoticaSubAnalysis::_recPFTauSelector
private

Definition at line 173 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::Photon>* HLTExoticaSubAnalysis::_recPhotonSelector
private

Definition at line 174 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

StringCutObjectSelector<reco::Track>* HLTExoticaSubAnalysis::_recTrackSelector
private

Definition at line 165 of file HLTExoticaSubAnalysis.h.

Referenced by initSelector(), insertCandidates(), and ~HLTExoticaSubAnalysis().

std::map<std::string, std::string> HLTExoticaSubAnalysis::_shortpath2long
private

Relation between the short and long versions of the path.

Definition at line 132 of file HLTExoticaSubAnalysis.h.

Referenced by analyze(), and beginRun().

std::map<unsigned int, edm::EDGetToken> HLTExoticaSubAnalysis::_tokens
private

Definition at line 144 of file HLTExoticaSubAnalysis.h.

Referenced by getHandlesToObjects(), and registerConsumes().

edm::InputTag HLTExoticaSubAnalysis::_trigResultsLabel
private

Definition at line 137 of file HLTExoticaSubAnalysis.h.

Referenced by registerConsumes().

edm::EDGetTokenT<edm::TriggerResults> HLTExoticaSubAnalysis::_trigResultsToken
private

Definition at line 142 of file HLTExoticaSubAnalysis.h.

Referenced by getHandlesToObjects(), and registerConsumes().