62 #include "TLorentzVector.h"
80 bool operator()(
const TLorentzVector&
a,
const TLorentzVector&
b)
const {
return a.Pt() >
b.Pt(); }
90 tok_pfjet_ = consumes<reco::PFJetCollection>(labelPFJet_);
91 tok_track_ = consumes<reco::TrackCollection>(labelTrack_);
92 tok_castorjet_ = consumes<reco::BasicJetCollection>(labelCastorJet_);
103 PFJetpt = bei.
book1D(
"PFJetpt",
";p_{T}(PFJet)", 100, 0.0, 100);
104 PFJeteta = bei.
book1D(
"PFJeteta",
";#eta(PFJet)", 50, -2.5, 2.5);
105 PFJetphi = bei.
book1D(
"PFJetphi",
";#phi(PFJet)", 50, -3.14, 3.14);
106 PFJetMulti = bei.
book1D(
"PFJetMulti",
";No. of PFJets", 10, -0.5, 9.5);
107 PFJetRapidity = bei.
book1D(
"PFJetRapidity",
";PFJetRapidity", 50, -6.0, 6.0);
108 CastorJetphi = bei.
book1D(
"CastorJetphi",
";#phi(CastorJet)", 50, -3.14, 3.14);
109 CastorJetMulti = bei.
book1D(
"CastorJetMulti",
";No. of CastorJets", 10, -0.5, 9.5);
110 Track_HP_Phi = bei.
book1D(
"Track_HP_Phi",
";#phi(HPtrack)", 50, -3.14, 3.14);
111 Track_HP_Eta = bei.
book1D(
"Track_HP_Eta",
";#eta(HPtrack)", 50, -2.5, 2.5);
112 Track_HP_Pt = bei.
book1D(
"Track_HP_Pt",
";p_{T}(HPtrack)", 500, 0.0, 50);
113 Track_HP_ptErr_over_pt = bei.
book1D(
"Track_HP_ptErr_over_pt",
";{p_{T}Err}/p_{T}", 100, 0, 0.1);
114 Track_HP_dzvtx_over_dzerr = bei.
book1D(
"Track_HP_dzvtx_over_dzerr",
";dZerr/dZ", 100, -10, 10);
115 Track_HP_dxyvtx_over_dxyerror = bei.
book1D(
"Track_HP_dxyvtx_over_dxyerror",
";dxyErr/dxy", 100, -10, 10);
116 NPV = bei.
book1D(
"NPV",
";NPV", 10, -0.5, 9.5);
117 PV_chi2 = bei.
book1D(
"PV_chi2",
";PV_chi2", 100, 0.0, 2.0);
118 PV_d0 = bei.
book1D(
"PV_d0",
";PV_d0", 100, -10.0, 10.0);
119 PV_numTrks = bei.
book1D(
"PV_numTrks",
";PV_numTrks", 100, -0.5, 99.5);
120 PV_sumTrks = bei.
book1D(
"PV_sumTrks",
";PV_sumTrks", 100, 0, 100);
121 h_ptsum_towards = bei.
book1D(
"h_ptsum_towards",
";h_ptsum_towards", 100, 0, 100);
122 h_ptsum_transverse = bei.
book1D(
"h_ptsum_transverse",
";h_ptsum_transverse", 100, 0, 100);
124 h_ntracks = bei.
book1D(
"h_ntracks",
";h_ntracks", 50, -0.5, 49.5);
125 h_trkptsum = bei.
book1D(
"h_trkptsum",
";h_trkptsum", 100, 0, 100);
126 h_ptsum_away = bei.
book1D(
"h_ptsum_away",
";h_ptsum_away", 100, 0, 100);
127 h_ntracks_towards = bei.
book1D(
"h_ntracks_towards",
";h_ntracks_towards", 50, -0.5, 49.5);
128 h_ntracks_transverse = bei.
book1D(
"h_ntracks_transverse",
";h_ntracks_transverse", 50, -0.5, 49.5);
129 h_ntracks_away = bei.
book1D(
"h_ntracks_away",
";h_ntracks_away", 50, -0.5, 49.5);
131 h_leadingtrkpt_ntrk_away =
132 bei.
bookProfile(
"h_leadingtrkpt_ntrk_away",
"h_leadingtrkpt_ntrk_away", 50, 0.0, 50, 0, 30,
" ");
133 h_leadingtrkpt_ntrk_towards =
134 bei.
bookProfile(
"h_leadingtrkpt_ntrk_towards",
"h_leadingtrkpt_ntrk_towards", 50, 0.0, 50, 0, 30,
" ");
135 h_leadingtrkpt_ntrk_transverse =
136 bei.
bookProfile(
"h_leadingtrkpt_ntrk_transverse",
"h_leadingtrkpt_ntrk_transverse", 50, 0.0, 50, 0, 30,
" ");
137 h_leadingtrkpt_ptsum_away =
138 bei.
bookProfile(
"h_leadingtrkpt_ptsum_away",
"h_leadingtrkpt_ptsum_away", 50, 0.0, 50, 0, 30,
" ");
139 h_leadingtrkpt_ptsum_towards =
140 bei.
bookProfile(
"h_leadingtrkpt_ptsum_towards",
"h_leadingtrkpt_ptsum_towards", 50, 0.0, 50, 0, 30,
" ");
141 h_leadingtrkpt_ptsum_transverse =
142 bei.
bookProfile(
"h_leadingtrkpt_ptsum_transverse",
"h_leadingtrkpt_ptsum_transverse", 50, 0.0, 50, 0, 30,
" ");
149 using namespace reco;
152 eventNumber_ =
iEvent.id().event();
153 lumiNumber_ =
iEvent.id().luminosityBlock();
154 bxNumber_ =
iEvent.bunchCrossing();
158 iEvent.getByToken(pvs_, privtxs);
159 double bestvz = -999.9, bestvx = -999.9, bestvy = -999.9;
160 double bestvzError = -999.9, bestvxError = -999.9, bestvyError = -999.9;
163 NPV->Fill(privtxs->size());
164 if (privtxs->begin() != privtxs->end() && !(pvtx.
isFake()) && pvtx.
position().Rho() <= 2. &&
169 bestvzError = pvtx.
zError();
170 bestvxError = pvtx.
xError();
171 bestvyError = pvtx.
yError();
173 PV_d0->Fill(
sqrt(pvtx.
x() * pvtx.
x() + pvtx.
y() * pvtx.
y()));
175 double vertex_sumTrks = 0.0;
177 vertex_sumTrks += (*iTrack)->pt();
179 PV_sumTrks->Fill(vertex_sumTrks);
188 iEvent.getByToken(tok_pfjet_, pfjetchscoll);
191 reco::PFJetCollection::const_iterator pfjetchsclus = pfchsjets->begin();
192 for (pfjetchsclus = pfchsjets->begin(); pfjetchsclus != pfchsjets->end(); ++pfjetchsclus) {
193 PFJetpt->Fill(pfjetchsclus->pt());
194 PFJeteta->Fill(pfjetchsclus->eta());
195 PFJetphi->Fill(pfjetchsclus->phi());
196 PFJetRapidity->Fill(pfjetchsclus->rapidity());
199 PFJetMulti->Fill(nPFCHSJet);
202 std::vector<Jet> recoCastorJets;
203 recoCastorJets.clear();
206 iEvent.getByToken(tok_castorjet_, castorJets);
208 for (
unsigned ijet = 0; ijet < castorJets->size(); ijet++) {
209 recoCastorJets.push_back((*castorJets)[ijet]);
211 for (
unsigned ijet = 0; ijet < recoCastorJets.size(); ijet++) {
213 CastorJetphi->Fill(recoCastorJets[ijet].phi());
215 CastorJetMulti->Fill(recoCastorJets.size());
220 iEvent.getByToken(tok_track_, itracks);
223 std::vector<TLorentzVector> T_trackRec_P4;
226 int ntracks_towards = 0;
227 int ntracks_transverse = 0;
228 int ntracks_away = 0;
231 float ptsum_towards = 0;
232 float ptsum_transverse = 0;
233 float ptsum_away = 0;
235 T_trackRec_P4.clear();
236 for (reco::TrackCollection::const_iterator iT = itracks->begin(); iT != itracks->end(); ++iT) {
237 if (iT->quality(hiPurity)) {
239 double dzvtx = iT->dz(bestvtx);
240 double dxyvtx = iT->dxy(bestvtx);
241 double dzerror =
sqrt(iT->dzError() * iT->dzError() + bestvzError * bestvzError);
242 double dxyerror =
sqrt(iT->d0Error() * iT->d0Error() + bestvxError * bestvyError);
243 if ((iT->ptError()) / iT->pt() < 0.05 && dzvtx < 3.0 && dxyvtx < 3.0) {
245 trk.SetPtEtaPhiE(iT->pt(), iT->eta(), iT->phi(), iT->p());
246 T_trackRec_P4.push_back(trk);
247 Track_HP_Eta->Fill(iT->eta());
248 Track_HP_Phi->Fill(iT->phi());
249 Track_HP_Pt->Fill(iT->pt());
250 Track_HP_ptErr_over_pt->Fill((iT->ptError()) / iT->pt());
251 Track_HP_dzvtx_over_dzerr->Fill(dzvtx / dzerror);
252 Track_HP_dxyvtx_over_dxyerror->Fill(dxyvtx / dxyerror);
257 float highest_pt_track = -999;
259 for (
unsigned int k = 0;
k < T_trackRec_P4.size();
k++) {
260 if (T_trackRec_P4.at(
k).Pt() > highest_pt_track) {
261 highest_pt_track = T_trackRec_P4.at(
k).Pt();
267 if (finalid < T_trackRec_P4.size()) {
269 for (
unsigned int itrk = 0; itrk < T_trackRec_P4.size(); itrk++) {
271 ptsum = ptsum + T_trackRec_P4.at(itrk).Pt();
272 dphi =
deltaPhi(T_trackRec_P4.at(itrk).Phi(), T_trackRec_P4.at(
index).Phi());
273 if (fabs(dphi) < 1.05) {
275 ptsum_towards = ptsum_towards + T_trackRec_P4.at(itrk).Pt();
277 if (fabs(dphi) > 1.05 && fabs(dphi) < 2.09) {
278 ++ntracks_transverse;
279 ptsum_transverse = ptsum_transverse + T_trackRec_P4.at(itrk).Pt();
281 if (fabs(dphi) > 2.09) {
283 ptsum_away = ptsum_away + T_trackRec_P4.at(itrk).Pt();
287 h_ntracks->Fill(ntracks);
288 h_trkptsum->Fill(ptsum);
289 h_ptsum_towards->Fill(ptsum_towards);
290 h_ptsum_transverse->Fill(ptsum_transverse);
291 h_ptsum_away->Fill(ptsum_away);
292 h_ntracks_towards->Fill(ntracks_towards);
293 h_ntracks_transverse->Fill(ntracks_transverse);
294 h_ntracks_away->Fill(ntracks_away);
296 if (!T_trackRec_P4.empty()) {
297 h_leadingtrkpt_ntrk_towards->Fill(T_trackRec_P4.at(
index).Pt(), ntracks_towards / 8.37);
298 h_leadingtrkpt_ntrk_transverse->Fill(T_trackRec_P4.at(
index).Pt(), ntracks_transverse / 8.37);
299 h_leadingtrkpt_ntrk_away->Fill(T_trackRec_P4.at(
index).Pt(), ntracks_away / 8.37);
300 h_leadingtrkpt_ptsum_towards->Fill(T_trackRec_P4.at(
index).Pt(), ptsum_towards / 8.37);
301 h_leadingtrkpt_ptsum_transverse->Fill(T_trackRec_P4.at(
index).Pt(), ptsum_transverse / 8.37);
302 h_leadingtrkpt_ptsum_away->Fill(T_trackRec_P4.at(
index).Pt(), ptsum_away / 8.37);