CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTMuonMatchAndPlot.h
Go to the documentation of this file.
1 #ifndef DQMOffline_Trigger_HLTMuonMatchAndPlot_H
2 #define DQMOffline_Trigger_HLTMuonMatchAndPlot_H
3 
20 // Base Class Headers
21 
22 //#include "FWCore/Framework/interface/Frameworkfwd.h"
27 //#include "DataFormats/Common/interface/RefToBase.h"
28 //#include "DataFormats/TrackReco/interface/Track.h"
29 //#include "DataFormats/Candidate/interface/Candidate.h"
30 
33 
43 
45 
46 #include <vector>
47 #include "TFile.h"
48 #include "TNtuple.h"
49 #include "TString.h"
50 #include "TPRegexp.h"
51 
52 
53 
56 
58 
59 const double NOMATCH = 999.;
60 const std::string EFFICIENCY_SUFFIXES[2] = {"denom", "numer"};
61 
62 
65 
67 
68  public:
69 
72  std::vector<std::string>);
73 
74  // Analyzer Methods
75  void beginRun(const edm::Run &, const edm::EventSetup &);
76  void analyze(const edm::Event &, const edm::EventSetup &);
77  void endRun(const edm::Run &, const edm::EventSetup &);
78 
79  // Helper Methods
80  void fillEdges(size_t & nBins, float * & edges, std::vector<double> binning);
81  template <class T> void
82  fillMapFromPSet(std::map<std::string, T> &, edm::ParameterSet, std::string);
83  template <class T1, class T2> std::vector<size_t>
84  matchByDeltaR(const std::vector<T1> &, const std::vector<T2> &,
85  const double maxDeltaR = NOMATCH);
86 
87  private:
88 
89  // Internal Methods
93  const reco::MuonCollection &,
94  const reco::BeamSpot &,
95  bool,
97  double, double);
100  const trigger::TriggerEvent &,
101  const edm::ParameterSet &);
102 
103  // Input from Configuration File
106  std::vector<std::string> requiredTriggers_;
107  std::map<std::string, edm::InputTag> inputTags_;
108  std::map<std::string, std::vector<double> > binParams_;
109  std::map<std::string, double> plotCuts_;
112 
113  // Member Variables
115  unsigned int cutMinPt_;
117  std::vector<std::string> moduleLabels_;
119  std::map<std::string, MonitorElement *> hists_;
120 
121  // Selectors
124 
126  double targetZ0Cut_;
127  double targetD0Cut_;
129  double probeZ0Cut_;
130  double probeD0Cut_;
131 };
132 
133 #endif
reco::MuonCollection selectedMuons(const reco::MuonCollection &, const reco::BeamSpot &, bool, const StringCutObjectSelector< reco::Muon > &, double, double)
StringCutObjectSelector< reco::Muon > probeMuonSelector_
The single EDProduct to be saved for each event (AOD case)
Definition: TriggerEvent.h:27
const std::string EFFICIENCY_SUFFIXES[2]
edm::ParameterSet targetParams_
std::vector< std::string > moduleLabels_
void fillEdges(size_t &nBins, float *&edges, std::vector< double > binning)
void book1D(std::string, std::string, std::string)
std::map< std::string, std::vector< double > > binParams_
std::map< std::string, MonitorElement * > hists_
dictionary edges
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
std::vector< std::string > requiredTriggers_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
const double NOMATCH
void analyze(const edm::Event &, const edm::EventSetup &)
void beginRun(const edm::Run &, const edm::EventSetup &)
std::vector< size_t > matchByDeltaR(const std::vector< T1 > &, const std::vector< T2 > &, const double maxDeltaR=NOMATCH)
HLTMuonMatchAndPlot(const edm::ParameterSet &, std::string, std::vector< std::string >)
Constructor.
void book2D(std::string, std::string, std::string, std::string)
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:83
StringCutObjectSelector< reco::Muon > targetMuonSelector_
std::map< std::string, edm::InputTag > inputTags_
trigger::TriggerObjectCollection selectedTriggerObjects(const trigger::TriggerObjectCollection &, const trigger::TriggerEvent &, const edm::ParameterSet &)
edm::ParameterSet probeParams_
void fillMapFromPSet(std::map< std::string, T > &, edm::ParameterSet, std::string)
void endRun(const edm::Run &, const edm::EventSetup &)
std::map< std::string, double > plotCuts_
Definition: Run.h:36
math::PtEtaPhiELorentzVectorF LorentzVector