ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/LJMet/MultivariateAnalysis/src/FirstCollisonsAnalyzer.cc
Revision: 1.1
Committed: Thu Apr 1 16:08:44 2010 UTC (15 years, 1 month ago) by msegala
Content type: text/plain
Branch: MAIN
CVS Tags: V00-03-01, ZMorph_BASE_20100408, gak040610_morphing, HEAD
Branch point for: ZMorph-V00-03-01
Log Message:
CMSSW_3_5_6

File Contents

# Content
1 // -*- C++ -*-
2
3 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
4 #include "DataFormats/PatCandidates/interface/Jet.h"
5 #include "DataFormats/PatCandidates/interface/Electron.h"
6 #include "DataFormats/PatCandidates/interface/Muon.h"
7 #include "DataFormats/PatCandidates/interface/MET.h"
8 #include "../interface/FirstCollisionsAnalyzer.h"
9 #include "LJMet/MultivariateAnalysis/interface/TopTopologicalVariables.h"
10 #include "LJMet/MultivariateAnalysis/interface/LJetsTopoVars.h"
11 #include "FWCore/Framework/interface/TriggerNames.h"
12 #include "DataFormats/Common/interface/TriggerResults.h"
13 #include "DataFormats/MuonReco/interface/MuonSelectors.h"
14 #include "DataFormats/PatCandidates/interface/TriggerPath.h"
15 #include "FWCore/ServiceRegistry/interface/Service.h"
16 #include "CommonTools/UtilAlgos/interface/TFileService.h"
17 #include "TLorentzVector.h"
18
19 using namespace std;
20 using namespace top_cafe;
21
22
23 FirstCollisonsAnalyzer::FirstCollisonsAnalyzer(const edm::ParameterSet& iConfig)
24
25 {
26
27 // init TFileService, needed for the output root file
28 edm::Service<TFileService> fs;
29 _tree = fs->make<TTree>("ttljets", "ttljets", 64000000);
30
31
32 //
33 // _____ some default values __________________________________________
34 //
35 nCaloJets_min = 1;
36 nLepton_min = 1;
37 jet_pt_min = 5.0;
38 jet_eta_max = 2.4;
39 muon_pt_min = 2.0;
40 muon_eta_max = 2.1;
41 muon_trackIso_max = 100.0;
42 muon_caloIso_max = 100.0;
43 electron_pt_min = 2.0;
44 electron_eta_max = 2.1;
45 _is_mc = true;
46
47
48
49
50 //
51 //_____ config file parameters ________________________________________
52 //
53 _jetSource = iConfig.getParameter<string>("jetSource");
54 _electronSource = iConfig.getParameter<string>("electronSource");
55 _muonSource = iConfig.getParameter<string>("muonSource");
56 _METSource = iConfig.getParameter<string>("METSource");
57 //
58 //_____ Optional cuts (better use the official PAT selection filters)__
59 //
60 if (iConfig.exists("nCaloJets_min")) nCaloJets_min = iConfig.getParameter<int>("nCaloJets_min");
61 if (iConfig.exists("nLepton_min")) nLepton_min = iConfig.getParameter<int>("nLepton_min");
62 if (iConfig.exists("jet_pt_min")) jet_pt_min = iConfig.getParameter<double>("jet_pt_min");
63 if (iConfig.exists("jet_eta_max")) jet_eta_max = iConfig.getParameter<double>("jet_eta_max");
64 if (iConfig.exists("muon_pt_min")) muon_pt_min = iConfig.getParameter<double>("muon_pt_min");
65 if (iConfig.exists("muon_eta_max")) muon_eta_max = iConfig.getParameter<double>("muon_eta_max");
66 if (iConfig.exists("muon_trackIso_max")) muon_trackIso_max = iConfig.getParameter<double>("muon_trackIso_max");
67 if (iConfig.exists("muon_caloIso_max")) muon_caloIso_max = iConfig.getParameter<double>("muon_caloIso_max");
68 if (iConfig.exists("electron_pt_min")) electron_pt_min = iConfig.getParameter<double>("electron_pt_min");
69 if (iConfig.exists("electron_eta_max")) electron_eta_max = iConfig.getParameter<double>("electron_eta_max");
70 if (iConfig.exists("electron_trackIso_max")) electron_trackIso_max = iConfig.getParameter<double>("electron_trackIso_max");
71 if (iConfig.exists("electron_caloIso_max")) electron_caloIso_max = iConfig.getParameter<double>("electron_caloIso_max");
72 if (iConfig.exists("is_mc")) _is_mc = iConfig.getParameter<bool>("is_mc");
73
74 _jet_pt = new vector<double>();
75 _jet_px = new vector<double>();
76 _jet_py = new vector<double>();
77 _jet_pz = new vector<double>();
78 _jet_eta = new vector<double>();
79 _jet_phi = new vector<double>();
80 _jet_et = new vector<double>();
81 _jet_energy = new vector<double>();
82
83
84 _muon_pt = new vector<double>();
85 _muon_pz = new vector<double>();
86 _muon_eta = new vector<double>();
87 _muon_phi = new vector<double>();
88 _muon_et = new vector<double>();
89 _muon_energy = new vector<double>();
90
91 _electron_pt = new vector<double>();
92 _electron_pz = new vector<double>();
93 _electron_eta = new vector<double>();
94 _electron_phi = new vector<double>();
95 _electron_et = new vector<double>();
96 _electron_energy = new vector<double>();
97
98
99
100
101 //
102 //_____ Create tree branches __________________________________________
103 //
104
105 //_tree -> Branch("event", &_event, "event/I" );
106 //_tree -> Branch("lumi", &ilumi, "lumi/I" );
107 //_tree -> Branch("run", &irun, "run/I" );
108
109 //
110 //___
111 //
112 b_n_jets = _tree -> Branch("n_jets", &_n_jets, "n_jets/I" );
113 b_n_met = _tree -> Branch("n_met", &_n_met, "n_met/I" );
114 b_n_muons = _tree -> Branch("n_muons", &_n_muons, "n_muons/I" );
115 b_n_electrons = _tree -> Branch("n_electrons", &_n_electrons, "n_electrons/I" );
116
117 //
118 //___MET
119 //
120 b_met_et = _tree -> Branch("met_et", &_met_et, "met_et/D" );
121 b_met_pt = _tree -> Branch("met_pt", &_met_pt, "met_pt/D" );
122 b_met_pz = _tree -> Branch("met_pz", &_met_pz, "met_pz/D" );
123 b_met_eta = _tree -> Branch("met_eta", &_met_eta, "met_eta/D" );
124 b_met_phi = _tree -> Branch("met_phi", &_met_phi, "met_phi/D" );
125 b_met_energy = _tree -> Branch("met_energy", &_met_energy, "met_energy/D" );
126
127 //
128 //___JETS
129 //
130 b_jet_pt = _tree -> Branch("jet_pt", &_jet_pt);
131 b_jet_pz = _tree -> Branch("jet_pz", &_jet_pz);
132 b_jet_eta = _tree -> Branch("jet_eta", &_jet_eta);
133 b_jet_phi = _tree -> Branch("jet_phi", &_jet_phi);
134 b_jet_et = _tree -> Branch("jet_et", &_jet_et);
135 b_jet_energy = _tree -> Branch("jet_energy", &_jet_energy);
136
137 //
138 //___MUONS
139 //
140 b_muon_pt = _tree -> Branch("v_muon_pt", &_muon_pt);
141 b_muon_pz = _tree -> Branch("v_muon_pz", &_muon_pz);
142 b_muon_eta = _tree -> Branch("v_muon_eta", &_muon_eta);
143 b_muon_phi = _tree -> Branch("v_muon_phi", &_muon_phi);
144 b_muon_et = _tree -> Branch("v_muon_et", &_muon_et);
145 b_muon_energy = _tree -> Branch("v_muon_energy", &_muon_energy);
146
147 //
148 //___ELECTRONS
149 //
150 b_electron_pt = _tree -> Branch("v_electron_pt", &_electron_pt);
151 b_electron_pz = _tree -> Branch("v_electron_pz", &_electron_pz);
152 b_electron_eta = _tree -> Branch("v_electron_eta", &_electron_eta);
153 b_electron_phi = _tree -> Branch("v_electron_phi", &_electron_phi);
154 b_electron_et = _tree -> Branch("v_electron_et", &_electron_et);
155 b_electron_energy = _tree -> Branch("v_electron_energy", &_electron_energy);
156
157
158
159 //
160 //_____ set behavior of the main event counter ________________________
161 //
162 eventCounter . setPrintCount(true);
163 eventCounter . setNewLine(true);
164 eventCounter . setMessage("Events processed: ");
165 eventCounter . setDivider(100);
166
167
168
169
170
171
172 }
173
174
175 FirstCollisonsAnalyzer::~FirstCollisonsAnalyzer()
176 {
177
178 }
179
180
181 //
182 // member functions
183 //
184
185 // ------------ method called to for each event ------------
186 void
187 FirstCollisonsAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
188 {
189
190
191 using namespace edm;
192 using namespace reco;
193 using namespace std;
194 //
195 //_____ Get objects from the event ____________________________________
196 //
197 Handle< vector< pat::Jet > > jets;
198 Handle< vector< pat::Electron > > electrons;
199 Handle< vector< pat::Muon > > muons;
200 Handle< vector< pat::MET > > METs;
201
202 iEvent . getByLabel( _jetSource, jets );
203 iEvent . getByLabel( _electronSource, electrons );
204 iEvent . getByLabel( _muonSource, muons );
205 iEvent . getByLabel( _METSource, METs );
206
207
208 _n_jets = -1;
209 _n_met = -1;
210 _n_muons = -1;
211 _n_electrons = -1;
212
213 _met_eta = -1.0;
214 _met_phi = -1.0;
215 _met_energy = -1.0;
216 _met_et = -1.0;
217 _met_pt = -1.0;
218 _met_pz = -1.0;
219
220 _jet_pt->clear();
221 _jet_pz->clear();
222 _jet_px->clear();
223 _jet_py->clear();
224 _jet_eta->clear();
225 _jet_phi->clear();
226 _jet_et->clear();
227 _jet_energy->clear();
228
229 _muon_eta->clear();
230 _muon_phi->clear();
231 _muon_energy->clear();
232 _muon_et->clear();
233 _muon_pt->clear();
234 _muon_pz->clear();
235
236 _electron_eta->clear();
237 _electron_phi->clear();
238 _electron_energy->clear();
239 _electron_et->clear();
240 _electron_pt->clear();
241 _electron_pz->clear();
242
243
244 RooGKCounter nOfGoodJets;
245 RooGKCounter nOfGoodMETs( "" );
246 RooGKCounter nOfGoodMuons;
247 RooGKCounter nOfGoodElectrons;
248
249
250
251
252
253
254 //
255 //_____ loop over jets ________________________________________________
256 //
257
258
259 int jet_ind=0;
260
261 for ( vector<pat::Jet>::const_iterator jet = jets -> begin(); jet != jets -> end(); jet++ ){
262 if (
263 ( fabs( (*jet) . eta() ) < jet_eta_max ) &&
264 ( (*jet) . pt() > jet_pt_min )
265 )
266 {
267 _jet_pt->push_back(jet->pt());
268 _jet_pz->push_back(jet->pz());
269 _jet_eta->push_back(jet->eta());
270 _jet_phi->push_back(jet->phi());
271 _jet_et->push_back(jet->et());
272 _jet_energy->push_back(jet->energy());
273 _jet_px->push_back(jet->px());
274 _jet_py->push_back(jet->py());
275
276 jet_ind++;
277 }
278 }
279
280
281
282 //
283 //_____ loop over the MET collection __________________________________
284 //
285
286 for ( vector<pat::MET>::const_iterator met = METs -> begin(); met != METs -> end(); met++){
287
288 nOfGoodMETs . count();
289 //int a = nOfGoodMETs . count();
290
291 _met_et = METs->begin()->et();
292 _met_phi = METs->begin()->phi();
293 _met_energy = METs->begin()->energy();
294 _met_eta = METs->begin()->eta();
295 _met_pt = METs->begin()->pt();
296 _met_pz = METs->begin()->pz();
297
298 }
299
300
301 //
302 //_____ loop over muons ________________________________________________
303 //
304
305
306 int muon_ind=0;
307
308 for ( vector<pat::Muon>::const_iterator muon = muons -> begin(); muon != muons -> end(); muon++ ){
309 if (
310 ( fabs( (*muon) . eta() ) < muon_eta_max ) &&
311 ( (*muon) . pt() > muon_pt_min )
312 )
313 {
314 _muon_pt->push_back(muon->pt());
315 _muon_pz->push_back(muon->pz());
316 _muon_eta->push_back(muon->eta());
317 _muon_phi->push_back(muon->phi());
318 _muon_et->push_back(muon->et());
319 _muon_energy->push_back(muon->energy());
320
321 muon_ind++;
322 }
323 }
324
325
326 //
327 //_____ loop over electrons ________________________________________________
328 //
329
330
331 int electron_ind=0;
332
333 for ( vector<pat::Electron>::const_iterator electron = electrons -> begin(); electron != electrons -> end(); electron++ ){
334 if (
335 ( fabs( (*electron) . eta() ) < electron_eta_max ) &&
336 ( (*electron) . pt() > electron_pt_min )
337 )
338 {
339 _electron_pt->push_back(electron->pt());
340 _electron_pz->push_back(electron->pz());
341 _electron_eta->push_back(electron->eta());
342 _electron_phi->push_back(electron->phi());
343 _electron_et->push_back(electron->et());
344 _electron_energy->push_back(electron->energy());
345
346 electron_ind++;
347 }
348 }
349
350
351 //
352 //_____ How many objects appear in the event ________________________________________________
353 //
354 _n_jets = nOfGoodJets.getCount();
355 _n_met = nOfGoodMETs.getCount();
356 _n_muons = nOfGoodMuons.getCount();
357 _n_electrons = nOfGoodElectrons.getCount();
358
359
360 //cout<<"jets, met, muons, electron = "<<_n_jets<<", "<<_n_met<<", "<<_n_muons<<", "<<_n_electrons<<endl;
361
362 /*
363 if (
364 (int)nOfGoodJets . getCount() >= nCaloJets_min &&
365 (int)nOfGoodMuons . getCount() >= nLepton_min &&
366 (int)nOfGoodElectrons . getCount() >= nLepton_min &&
367 nOfGoodMETs . getCount() > 0 ){
368 */
369 //
370 //____ FILL THE TREE
371 //
372
373 _tree -> Fill();
374
375
376
377 //}
378
379
380
381 }
382
383
384
385 // ------------ method called once each job just after ending the event loop ------------
386 void
387 FirstCollisonsAnalyzer::endJob() {
388
389
390 delete _jet_pt;
391 delete _jet_pz;
392 delete _jet_eta;
393 delete _jet_phi;
394 delete _jet_et;
395 delete _jet_energy;
396 //
397 delete _muon_pt;
398 delete _muon_pz;
399 delete _muon_eta;
400 delete _muon_phi;
401 delete _muon_et;
402 delete _muon_energy;
403 //
404 delete _electron_pt;
405 delete _electron_pz;
406 delete _electron_eta;
407 delete _electron_phi;
408 delete _electron_et;
409 delete _electron_energy;
410
411
412 }
413
414 //define this as a plug-in
415 //DEFINE_FWK_MODULE(FirstCollisonsAnalyzer);