CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
B2GDQM.cc
Go to the documentation of this file.
2 
3 #include <memory>
4 
5 // DQM
8 
9 // Framework
18 
19 // Candidate handling
26 
27 // Vertex utilities
30 
31 // Other
35 
36 // Math
44 
45 // vertexing
46 
47 // Transient tracks
55 
56 //JetCorrection
58 
59 //Substructure
62 
63 // ROOT
64 #include "TLorentzVector.h"
65 
66 // STDLIB
67 #include <iostream>
68 #include <iomanip>
69 #include <stdio.h>
70 #include <string>
71 #include <sstream>
72 #include <math.h>
73 
74 using namespace edm;
75 using namespace std;
76 using namespace reco;
77 using namespace trigger;
78 
79 
80 //
81 // -- Constructor
82 //
84 
85  edm::LogInfo("B2GDQM") << " Starting B2GDQM " << "\n" ;
86 
87 
88  typedef std::vector<edm::InputTag> vtag;
89 
90  // Get parameters from configuration file
91  // Trigger
92  theTriggerResultsCollection = ps.getParameter<InputTag>("triggerResultsCollection");
93  triggerToken_ = consumes< edm::TriggerResults> ( theTriggerResultsCollection );
94 
95  // Jets
96  jetLabels_ = ps.getParameter<std::vector<edm::InputTag> >("jetLabels");
97  for ( std::vector<edm::InputTag>::const_iterator jetlabel = jetLabels_.begin(),
98  jetlabelEnd = jetLabels_.end(); jetlabel != jetlabelEnd; ++jetlabel ) {
99  jetTokens_.push_back( consumes<edm::View<reco::Jet> >( *jetlabel ) );
100  }
101  jetPtMins_ = ps.getParameter<std::vector<double> > ("jetPtMins");
102  PFJetCorService_ = ps.getParameter<std::string>("PFJetCorService");
103 
104  // MET
105  PFMETLabel_ = ps.getParameter<InputTag>("pfMETCollection");
106  PFMETToken_ = consumes<std::vector<reco::PFMET> > ( PFMETLabel_ );
107 
108 
109 }
110 
111 
112 //
113 // -- Destructor
114 //
116  edm::LogInfo("B2GDQM") << " Deleting B2GDQM " << "\n" ;
117 }
118 
119 
120 //
121 // -- Begin Job
122 //
124  nLumiSecs_ = 0;
125  nEvents_ = 0;
126 }
127 
128 
129 //
130 // -- Begin Run
131 //
132 void B2GDQM::beginRun(Run const& run, edm::EventSetup const& eSetup) {
133  edm::LogInfo ("B2GDQM") <<"[B2GDQM]: Begining of Run";
134 
135 
136  bei_ = Service<DQMStore>().operator->();
137  bei_->setCurrentFolder("Physics/B2G");
138  bookHistos(bei_);
139 
140  // passed as parameter to HLTConfigProvider::init(), not yet used
141  bool isConfigChanged = false;
142 
143  // isValidHltConfig_ used to short-circuit analyze() in case of problems
144  // const std::string hltProcessName( "HLT" );
145  const std::string hltProcessName = theTriggerResultsCollection.process();
146  isValidHltConfig_ = hltConfigProvider_.init( run, eSetup, hltProcessName, isConfigChanged );
147 
148 }
149 
150 
151 //
152 // -- Begin Luminosity Block
153 //
155  edm::EventSetup const& context) {
156  //edm::LogInfo ("B2GDQM") <<"[B2GDQM]: Begin of LS transition";
157 }
158 
159 
160 //
161 // -- Book histograms
162 //
164 
165  bei->cd();
166 
167  //--- Event Interpretation
168 
169  for ( unsigned int icoll = 0; icoll < jetLabels_.size(); ++icoll ) {
170  std::stringstream ss;
171  ss << "Physics/B2G/" << jetLabels_[icoll].label();
172  bei->setCurrentFolder(ss.str().c_str());
173  pfJet_pt .push_back( bei->book1D("pfJet_pt", "Pt of PFJet (GeV)", 50, 0.0, 1000) );
174  pfJet_y .push_back( bei->book1D("pfJet_y", "Rapidity of PFJet", 60, -6.0, 6.0) );
175  pfJet_phi .push_back( bei->book1D("pfJet_phi", "#phi of PFJet (radians)",60, -3.14159, 3.14159) );
176  pfJet_m .push_back( bei->book1D("pfJet_m", "Mass of PFJet (GeV)", 50, 0.0, 500) );
177  pfJet_chef .push_back( bei->book1D("pfJet_pfchef", "PFJetID CHEF", 50, 0.0 , 1.0));
178  pfJet_nhef .push_back( bei->book1D("pfJet_pfnhef", "PFJetID NHEF", 50, 0.0 , 1.0));
179  pfJet_cemf .push_back( bei->book1D("pfJet_pfcemf", "PFJetID CEMF", 50, 0.0 , 1.0));
180  pfJet_nemf .push_back( bei->book1D("pfJet_pfnemf", "PFJetID NEMF", 50, 0.0 , 1.0));
181 
182  boostedJet_subjetPt .push_back( bei->book1D("boostedJet_subjetPt", "Pt of subjets (GeV)", 50, 0.0 , 500));
183  boostedJet_subjetY .push_back( bei->book1D("boostedJet_subjetY", "Rapidity of subjets", 60, -6.0, 6.0));
184  boostedJet_subjetPhi .push_back( bei->book1D("boostedJet_subjetPhi", "#phi of subjets (radians)", 60, -3.14159, 3.14159));
185  boostedJet_subjetM .push_back( bei->book1D("boostedJet_subjetM", "Mass of subjets (GeV)", 50, 0.0 , 250.));
186  boostedJet_subjetN .push_back( bei->book1D("boostedJet_subjetN", "Number of subjets", 10, 0, 10));
187  boostedJet_massDrop .push_back( bei->book1D("boostedJet_massDrop", "Mass drop for W-like jets", 50, 0.0 , 1.0));
188  boostedJet_minMass .push_back( bei->book1D("boostedJet_minMass", "Minimum Mass Pairing for top-like jets", 50, 0.0 , 250.0));
189 
190  }
191 
192  bei->setCurrentFolder("Physics/B2G/MET");
193  pfMet_pt = bei->book1D("pfMet_pt", "Pf Missing p_{T}; GeV", 50, 0.0 , 500);
194  pfMet_phi = bei->book1D("pfMet_phi", "Pf Missing p_{T} #phi;#phi (radians)", 35, -3.5, 3.5 );
195 
196 
197 
198  bei->cd();
199 }
200 
201 
202 //
203 // -- Analyze
204 //
205 void B2GDQM::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup){
206 
207  analyzeEventInterpretation(iEvent, iSetup);
208 }
209 
210 
212 
213 
214  // Loop over the different types of jets,
215  // Loop over the jets in that collection,
216  // fill PF jet information as well as substructure
217  // information for boosted jets.
218  // Utilizes the CMS top-tagging algorithm and the "mass drop" W-tagger.
219  for ( unsigned int icoll = 0; icoll < jetLabels_.size(); ++icoll ) {
220 
221  edm::Handle<edm::View<reco::Jet> > pfJetCollection;
222  bool ValidPFJets = iEvent.getByToken(jetTokens_[icoll], pfJetCollection );
223  if(!ValidPFJets) continue;
224  edm::View<reco::Jet> const & pfjets = *pfJetCollection;
225 
226 
227 
228  // Jet Correction
229  int countJet = 0;
230  //const JetCorrector* pfcorrector = JetCorrector::getJetCorrector(PFJetCorService_,iSetup);
231 
233  jetEnd = pfjets.end();
234  jet != jetEnd; ++jet ) {
235  if(jet->pt()<jetPtMins_[icoll]) continue;
236  pfJet_pt [icoll]->Fill(jet->pt());
237  pfJet_y [icoll]->Fill(jet->rapidity());
238  pfJet_phi [icoll]->Fill(jet->phi());
239  pfJet_m [icoll]->Fill(jet->mass());
240 
241  // Dynamic cast the base class (reco::Jet) to the derived class (PFJet)
242  // to access the PFJet information
243  reco::PFJet const * pfjet = dynamic_cast<reco::PFJet const *>( &*jet);
244 
245  if ( pfjet != 0 ) {
246  pfJet_chef[icoll]->Fill(pfjet->chargedHadronEnergyFraction());
247  pfJet_nhef[icoll]->Fill(pfjet->neutralHadronEnergyFraction());
248  pfJet_cemf[icoll]->Fill(pfjet->chargedEmEnergyFraction());
249  pfJet_nemf[icoll]->Fill(pfjet->neutralEmEnergyFraction());
250  }
251 
252  // Dynamic cast the base class (reco::Jet) to the derived class (BasicJet)
253  // to access the substructure information
254  reco::BasicJet const * basicjet = dynamic_cast<reco::BasicJet const *>( &*jet);
255 
256  if ( basicjet != 0 ) {
257  boostedJet_subjetN[icoll]->Fill( jet->numberOfDaughters() );
258 
259  for ( unsigned int ida = 0; ida < jet->numberOfDaughters(); ++ida ) {
260  reco::Candidate const * subjet = jet->daughter(ida);
261  boostedJet_subjetPt [icoll]->Fill ( subjet->pt() );
262  boostedJet_subjetY [icoll]->Fill ( subjet->rapidity() );
263  boostedJet_subjetPhi[icoll]->Fill ( subjet->phi() );
264  boostedJet_subjetM [icoll]->Fill ( subjet->mass() );
265  }
266  // Check the various tagging algorithms
267 
268  // For top-tagging, check the minimum mass pairing
269  if ( jetLabels_[icoll].label() == "cmsTopTagPFJetsCHS") {
270  CATopJetHelper helper(171.2, 80.4);
271  reco::CATopJetProperties properties = helper(*basicjet);
272  if ( jet->numberOfDaughters() > 2 ) {
273  boostedJet_minMass[icoll]->Fill ( properties.minMass );
274  } else {
275  boostedJet_minMass[icoll]->Fill ( -1.0 );
276  }
277 
278  // For W-tagging, check the mass drop
279  } else if ( jetLabels_[icoll].label() == "ca8PFJetsCHSPruned" ) {
280  if ( jet->numberOfDaughters() > 1 ) {
281  reco::Candidate const * da0 = jet->daughter(0);
282  reco::Candidate const * da1 = jet->daughter(1);
283  if ( da0->mass() > da1->mass() ) {
284  boostedJet_massDrop[icoll]->Fill( da0->mass() / jet->mass() );
285  } else {
286  boostedJet_massDrop[icoll]->Fill( da1->mass() / jet->mass() );
287  }
288  } else {
289  boostedJet_massDrop[icoll]->Fill( -1.0 );
290  }
291 
292  } // end if collection is CA8 PFJets CHS Pruned
293 
294  } // end if basic jet != 0
295  countJet++;
296  }
297  }
298 
299 
300 
301  // PFMETs
302  edm::Handle<std::vector<reco::PFMET> > pfMETCollection;
303  bool ValidPFMET = iEvent.getByToken(PFMETToken_, pfMETCollection);
304  if(!ValidPFMET) return;
305 
306  pfMet_pt->Fill( (*pfMETCollection)[0].pt() );
307  pfMet_phi->Fill( (*pfMETCollection)[0].phi() );
308 
309 }
310 
311 // -- End Luminosity Block
312 //
314  //edm::LogInfo ("B2GDQM") <<"[B2GDQM]: End of LS transition, performing the DQM client operation";
315  nLumiSecs_++;
316  //edm::LogInfo("B2GDQM") << "============================================ "
317  //<< endl << " ===> Iteration # " << nLumiSecs_ << " " << lumiSeg.luminosityBlock()
318  //<< endl << "============================================ " << endl;
319 }
320 
321 
322 //
323 // -- End Run
324 //
325 void B2GDQM::endRun(edm::Run const& run, edm::EventSetup const& eSetup){
326 }
327 
328 
329 //
330 // -- End Job
331 //
333  //edm::LogInfo("B2GDQM") <<"[B2GDQM]: endjob called!";
334 }
335 
336 
337 
T getParameter(std::string const &) const
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
void bookHistos(DQMStore *bei)
Definition: B2GDQM.cc:163
virtual float mass() const =0
mass
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:954
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
void cd(void)
go to top directory (ie. root)
Definition: DQMStore.cc:644
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction
Definition: PFJet.h:96
virtual void beginLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context)
Definition: B2GDQM.cc:154
virtual float phi() const =0
momentum azimuthal angle
virtual double rapidity() const =0
rapidity
Jets made from CaloTowers.
Definition: BasicJet.h:20
Jets made from PFObjects.
Definition: PFJet.h:21
virtual void endRun(edm::Run const &run, edm::EventSetup const &eSetup)
Definition: B2GDQM.cc:325
void bookHistos()
Definition: Histogram.h:33
float neutralEmEnergyFraction() const
neutralEmEnergyFraction
Definition: PFJet.h:148
virtual float pt() const =0
transverse momentum
int iEvent
Definition: GenABIO.cc:230
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction
Definition: PFJet.h:100
virtual void beginJob()
Definition: B2GDQM.cc:123
float chargedEmEnergyFraction() const
chargedEmEnergyFraction
Definition: PFJet.h:140
B2GDQM(const edm::ParameterSet &ps)
Definition: B2GDQM.cc:83
virtual void beginRun(edm::Run const &run, edm::EventSetup const &eSetup)
Definition: B2GDQM.cc:132
virtual void endJob()
Definition: B2GDQM.cc:332
const_iterator begin() const
const_iterator end() const
virtual void analyzeEventInterpretation(edm::Event const &e, edm::EventSetup const &eSetup)
Definition: B2GDQM.cc:211
virtual void endLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &c)
Definition: B2GDQM.cc:313
virtual ~B2GDQM()
Definition: B2GDQM.cc:115
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:667
Definition: Run.h:41
Definition: DDAxes.h:10
virtual void analyze(edm::Event const &e, edm::EventSetup const &eSetup)
Definition: B2GDQM.cc:205
int icoll
Definition: AMPTWrapper.h:136