CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TTbar_GenLepAnalyzer.cc
Go to the documentation of this file.
5 
7  leps_(iConfig.getParameter<edm::InputTag>("leptons"))
8 {
9  //now do what ever initialization is needed
10  dbe = 0;
12 
13  lepsToken_=consumes< edm::View<reco::Candidate> >(leps_);
14 
15 }
16 
17 
19 {
20 
21  // do anything here that needs to be done at desctruction time
22  // (e.g. close files, deallocate resources etc.)
23 
24 }
25 
26 
27 //
28 // member functions
29 //
30 
31 // ------------ method called for each event ------------
32 void
34 {
35 
36  // Handle to the Leptons collections
38  iEvent.getByToken(lepsToken_, leps);
39  if(!leps.isValid()) return;
40 
41  // loop Jet collection and fill histograms
42  int nleps = 0;
43  for(edm::View<reco::Candidate>::const_iterator lep_it=leps->begin(); lep_it!=leps->end(); ++lep_it){
44 
45  ++nleps;
46 
47  if (nleps > 0) { hists_["lepPtAll" ]->Fill( lep_it->p4().pt() );
48  hists_["lepEtaAll"]->Fill( lep_it->p4().eta() );
49  }
50  if (nleps == 1) { hists_["lepPt1" ]->Fill( lep_it->p4().pt() );
51  hists_["lepEta1"]->Fill( lep_it->p4().eta() );
52  }
53  if (nleps == 2) { hists_["lepPt2" ]->Fill( lep_it->p4().pt() );
54  hists_["lepEta2"]->Fill( lep_it->p4().eta() );
55  }
56  if (nleps == 3) { hists_["lepPt3" ]->Fill( lep_it->p4().pt() );
57  hists_["lepEta3"]->Fill( lep_it->p4().eta() );
58  }
59  if (nleps == 4) { hists_["lepPt4" ]->Fill( lep_it->p4().pt() );
60  hists_["lepEta4"]->Fill( lep_it->p4().eta() );
61  }
62  }
63 
64  hists_["lepN" ]->Fill( nleps ) ;
65 
66 
67 }
68 
69 
70 // ------------ method called once each job just before starting event loop ------------
71 void
73 {
74  if(!dbe) return;
75  dbe->setCurrentFolder("Generator/TTbar");
76  hists_["lepN" ] = dbe->book1D("TTbar_lepN" , "N" , 10, -.5, 9.5 );
77 
78  hists_["lepPtAll" ] = dbe->book1D("TTbar_lepPtAll_"+leps_.label() , "pt" , 1000, 0., 1000.);
79  hists_["lepPt1" ] = dbe->book1D("TTbar_lepPt1_"+leps_.label() , "pt" , 1000, 0., 1000.);
80  hists_["lepPt2" ] = dbe->book1D("TTbar_lepPt2_"+leps_.label() , "pt" , 1000, 0., 1000.);
81  hists_["lepPt3" ] = dbe->book1D("TTbar_lepPt3_"+leps_.label() , "pt" , 1000, 0., 1000.);
82  hists_["lepPt4" ] = dbe->book1D("TTbar_lepPt4_"+leps_.label() , "pt" , 1000, 0., 1000.);
83 
84  hists_["lepEtaAll"] = dbe->book1D("TTbar_lepEtaAll"+leps_.label(), "eta", 100, -5., 5.);
85  hists_["lepEta1" ] = dbe->book1D("TTbar_lepEta1"+leps_.label() , "eta", 100, -5., 5.);
86  hists_["lepEta2" ] = dbe->book1D("TTbar_lepEta2"+leps_.label() , "eta", 100, -5., 5.);
87  hists_["lepEta3" ] = dbe->book1D("TTbar_lepEta3"+leps_.label() , "eta", 100, -5., 5.);
88  hists_["lepEta4" ] = dbe->book1D("TTbar_lepEta4"+leps_.label() , "eta", 100, -5., 5.);
89 }
90 
91 // ------------ method called once each job just after ending the event loop ------------
92 void
94 {
95 }
96 
97 // ------------ method called when starting to processes a run ------------
98 void
100 {
101 }
102 
103 // ------------ method called when ending the processing of a run ------------
104 void
106 {
107 }
108 
109 // ------------ method called when starting to processes a luminosity block ------------
110 void
112 {
113 }
114 
115 // ------------ method called when ending the processing of a luminosity block ------------
116 void
118 {
119 }
120 
121 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
122 void
124  //The following says we do not know what parameters are allowed so do no validation
125  // Please change this to state exactly what you do use, even if it is no parameters
127  desc.setUnknown();
128  descriptions.addDefault(desc);
129 }
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
edm::EDGetTokenT< edm::View< reco::Candidate > > lepsToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
TTbar_GenLepAnalyzer(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:243
void addDefault(ParameterSetDescription const &psetDescription)
virtual void endLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &)
bool isValid() const
Definition: HandleBase.h:76
virtual void analyze(const edm::Event &, const edm::EventSetup &)
std::string const & label() const
Definition: InputTag.h:42
DQMStore * dbe
ME&#39;s &quot;container&quot;.
std::map< std::string, MonitorElement * > hists_
virtual void beginRun(edm::Run const &, edm::EventSetup const &)
virtual void endRun(edm::Run const &, edm::EventSetup const &)
virtual void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
Definition: Run.h:41