CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ScoutingTestAnalyzer.cc
Go to the documentation of this file.
1 // This class is used to test the functionalities of the package
2 
4 
5 
7 
8 //------------------------------------------------------------------------------
9 // A simple constructor which takes as inoput only the name of the PF jet collection
11  :ScoutingAnalyzerBase(conf){
12  m_pfJetsCollectionTag = conf.getUntrackedParameter<edm::InputTag>("pfJetsCollectionName");
13  }
14 
15 //------------------------------------------------------------------------------
16 // Nothing to destroy: the DQM service thinks about everything
18 
19 //------------------------------------------------------------------------------
20 // Usual analyze method
22 
23 
24  edm::Handle<reco::CaloJetCollection> calojets_handle ;
25  iEvent.getByLabel(m_pfJetsCollectionTag,calojets_handle) ;
26  /* This is an example of how C++11 can simplify or lifes. The auto keyword
27  make the compiler figure out by itself which is the type of the pfjets object.
28  The qualifier const of course still apply.
29  Poor's man explaination: "compiler, make pfjets a const ref and figure out
30  for me the type"*/
31  auto const& calojets = *calojets_handle;
32 
33  // Again, C++11. A loop on a std::vector becomes as simple as this!
34  for (auto const & calojet: calojets){
35  m_jetPt->Fill(calojet.pt());
36  m_jetEtaPhi->Fill(calojet.eta(),calojet.phi());
37  }
38 
39 
40 }
41 
42 //------------------------------------------------------------------------------
43 /* Method called at the end of the Run. Ideal to finalise stuff within the
44  * DQM infrastructure, which is entirely Run based. */
46 
47  std::string collection_name = m_pfJetsCollectionTag.label();
48  /* This function is specific of this class and allows us to make a
49  * projection in one line */
50 
52  collection_name+" Jets #eta (projection)",
53  "#eta^{Jet}");
54 
56  collection_name+" Jets phi (projection)",
57  "#phi^{Jet}");
58 
59 
60 }
61 
62 //------------------------------------------------------------------------------
63 // Function to book the Monitoring Elements.
65  std::string collection_name = m_pfJetsCollectionTag.label();
66 
67  /* This method allows us to book an Histogram in one line in a completely
68  * transparent way. Take your time to put axis titles!!!!*/
69  m_jetPt = bookH1withSumw2(collection_name+"_pt",
70  collection_name+" Jet P_{T}",
71  50,0.,500.,
72  "Jet P_{T} [GeV]");
73 
74  m_jetEtaPhi = bookH2withSumw2(collection_name+"_etaphi",
75  collection_name+" #eta #phi",
76  50,-5,5,
77  50,-3.1415,+3.1415,
78  "#eta^{Jet}",
79  "#phi^{Jet}");
80 }
81 
82 //------------------------------------------------------------------------------
83 
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * bookH2withSumw2(const std::string &name, const std::string &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, const std::string &titleX="", const std::string &titleY="", Option_t *option="COLZ")
MonitorElement * bookH1withSumw2(const std::string &name, const std::string &title, int nchX, double lowX, double highX, const std::string &titleX="", const std::string &titleY="Events", Option_t *option="E1 P")
virtual void analyze(const edm::Event &, const edm::EventSetup &)
MonitorElement * m_jetEtaPhi
edm::InputTag m_pfJetsCollectionTag
void Fill(long long x)
MonitorElement * profileY(MonitorElement *me2d, const std::string &title="", const std::string &titleX="", const std::string &titleY="", Double_t minimum=-1111, Double_t maximum=-1111)
int iEvent
Definition: GenABIO.cc:243
MonitorElement * profileX(MonitorElement *me2d, const std::string &title="", const std::string &titleX="", const std::string &titleY="", Double_t minimum=-1111, Double_t maximum=-1111)
ScoutingTestAnalyzer(const edm::ParameterSet &)
MonitorElement * m_jetPt
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
tuple conf
Definition: dbtoconf.py:185
virtual void endRun(edm::Run const &, edm::EventSetup const &)
std::string const & label() const
Definition: InputTag.h:42
Definition: Run.h:36