File indexing completed on 2023-10-25 10:07:37
0001
0002
0003
0004
0005
0006
0007 #include <fstream>
0008 #include <iostream>
0009 #include <memory>
0010 #include <numeric>
0011 #include <vector>
0012
0013
0014 #include "DataFormats/Common/interface/DetSetVector.h"
0015 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0016 #include "DataFormats/Common/interface/Ptr.h"
0017 #include "DataFormats/Common/interface/Ref.h"
0018 #include "DataFormats/L1TrackTrigger/interface/TTCluster.h"
0019 #include "DataFormats/L1TrackTrigger/interface/TTStub.h"
0020 #include "DataFormats/L1TrackTrigger/interface/TTTrack.h"
0021 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0022 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0024 #include "FWCore/ServiceRegistry/interface/Service.h"
0025 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0026 #include "Geometry/Records/interface/StackedTrackerGeometryRecord.h"
0027 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0028 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
0029 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0030 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
0031 #include "SimTracker/TrackTriggerAssociation/interface/TTClusterAssociationMap.h"
0032 #include "SimTracker/TrackTriggerAssociation/interface/TTStubAssociationMap.h"
0033 #include "SimTracker/TrackTriggerAssociation/interface/TTTrackAssociationMap.h"
0034
0035 #include "OuterTrackerMonitorTrackingParticles.h"
0036
0037
0038
0039
0040 OuterTrackerMonitorTrackingParticles::OuterTrackerMonitorTrackingParticles(const edm::ParameterSet &iConfig)
0041 : m_topoToken(esConsumes()), conf_(iConfig) {
0042 topFolderName_ = conf_.getParameter<std::string>("TopFolderName");
0043 trackingParticleToken_ =
0044 consumes<std::vector<TrackingParticle>>(conf_.getParameter<edm::InputTag>("trackingParticleToken"));
0045 ttStubMCTruthToken_ =
0046 consumes<TTStubAssociationMap<Ref_Phase2TrackerDigi_>>(conf_.getParameter<edm::InputTag>("MCTruthStubInputTag"));
0047 ttClusterMCTruthToken_ = consumes<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>>(
0048 conf_.getParameter<edm::InputTag>("MCTruthClusterInputTag"));
0049 ttTrackMCTruthToken_ = consumes<TTTrackAssociationMap<Ref_Phase2TrackerDigi_>>(
0050 conf_.getParameter<edm::InputTag>("MCTruthTrackInputTag"));
0051 L1Tk_minNStub = conf_.getParameter<int>("L1Tk_minNStub");
0052 L1Tk_maxChi2dof = conf_.getParameter<double>("L1Tk_maxChi2dof");
0053 TP_minNStub = conf_.getParameter<int>("TP_minNStub");
0054
0055 TP_minNLayersStub = conf_.getParameter<int>("TP_minNLayersStub");
0056 TP_minPt = conf_.getParameter<double>("TP_minPt");
0057 TP_maxEta = conf_.getParameter<double>("TP_maxEta");
0058 TP_maxVtxZ = conf_.getParameter<double>("TP_maxVtxZ");
0059 }
0060
0061 OuterTrackerMonitorTrackingParticles::~OuterTrackerMonitorTrackingParticles() = default;
0062
0063
0064
0065
0066 void OuterTrackerMonitorTrackingParticles::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0067
0068 edm::Handle<std::vector<TrackingParticle>> trackingParticleHandle;
0069 iEvent.getByToken(trackingParticleToken_, trackingParticleHandle);
0070
0071
0072 edm::Handle<TTTrackAssociationMap<Ref_Phase2TrackerDigi_>> MCTruthTTTrackHandle;
0073 iEvent.getByToken(ttTrackMCTruthToken_, MCTruthTTTrackHandle);
0074 edm::Handle<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>> MCTruthTTClusterHandle;
0075 iEvent.getByToken(ttClusterMCTruthToken_, MCTruthTTClusterHandle);
0076 edm::Handle<TTStubAssociationMap<Ref_Phase2TrackerDigi_>> MCTruthTTStubHandle;
0077 iEvent.getByToken(ttStubMCTruthToken_, MCTruthTTStubHandle);
0078
0079
0080 const TrackerTopology *const tTopo = &iSetup.getData(m_topoToken);
0081
0082
0083 int this_tp = 0;
0084 for (const auto &iterTP : *trackingParticleHandle) {
0085 edm::Ptr<TrackingParticle> tp_ptr(trackingParticleHandle, this_tp);
0086 this_tp++;
0087
0088
0089 float tmp_tp_pt = iterTP.pt();
0090 float tmp_tp_phi = iterTP.phi();
0091 float tmp_tp_eta = iterTP.eta();
0092
0093
0094 std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>>
0095 theStubRefs = MCTruthTTStubHandle->findTTStubRefs(tp_ptr);
0096
0097 int hasStubInLayer[11] = {0};
0098 for (unsigned int is = 0; is < theStubRefs.size(); is++) {
0099 DetId detid(theStubRefs.at(is)->getDetId());
0100 int layer = -1;
0101 if (detid.subdetId() == StripSubdetector::TOB)
0102 layer = static_cast<int>(tTopo->layer(detid)) - 1;
0103 else if (detid.subdetId() == StripSubdetector::TID)
0104 layer = static_cast<int>(tTopo->layer(detid)) + 5;
0105
0106
0107 if (MCTruthTTStubHandle->findTrackingParticlePtr(theStubRefs.at(is)).isNull() && hasStubInLayer[layer] < 2)
0108 hasStubInLayer[layer] = 1;
0109 else
0110 hasStubInLayer[layer] = 2;
0111 }
0112
0113 int nStubLayerTP = 0;
0114 for (int isum = 0; isum < 11; isum++) {
0115 if (hasStubInLayer[isum] >= 1)
0116 nStubLayerTP += 1;
0117 }
0118
0119 if (std::fabs(tmp_tp_eta) > TP_maxEta)
0120 continue;
0121
0122 if (tmp_tp_pt > TP_minPt && nStubLayerTP >= TP_minNLayersStub) {
0123 trackParts_Pt->Fill(tmp_tp_pt);
0124 trackParts_Eta->Fill(tmp_tp_eta);
0125 trackParts_Phi->Fill(tmp_tp_phi);
0126 }
0127
0128
0129
0130 int nStubTP = -1;
0131 if (MCTruthTTStubHandle.isValid()) {
0132 std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>>
0133 theStubRefs = MCTruthTTStubHandle->findTTStubRefs(tp_ptr);
0134 nStubTP = (int)theStubRefs.size();
0135 }
0136 if (MCTruthTTClusterHandle.isValid() && MCTruthTTClusterHandle->findTTClusterRefs(tp_ptr).empty())
0137 continue;
0138
0139 float tmp_tp_vz = iterTP.vz();
0140 float tmp_tp_vx = iterTP.vx();
0141 float tmp_tp_vy = iterTP.vy();
0142 float tmp_tp_charge = tp_ptr->charge();
0143 int tmp_tp_pdgid = iterTP.pdgId();
0144
0145
0146
0147
0148 float tmp_tp_t = tan(2.0 * atan(1.0) - 2.0 * atan(exp(-tmp_tp_eta)));
0149 float delx = -tmp_tp_vx;
0150 float dely = -tmp_tp_vy;
0151 float K = 0.01 * 0.5696 / tmp_tp_pt * tmp_tp_charge;
0152 float A = 1. / (2. * K);
0153 float tmp_tp_x0p = delx - A * sin(tmp_tp_phi);
0154 float tmp_tp_y0p = dely + A * cos(tmp_tp_phi);
0155 float tmp_tp_rp = sqrt(tmp_tp_x0p * tmp_tp_x0p + tmp_tp_y0p * tmp_tp_y0p);
0156 static double pi = 4.0 * atan(1.0);
0157 float delphi = tmp_tp_phi - atan2(-K * tmp_tp_x0p, K * tmp_tp_y0p);
0158 if (delphi < -pi)
0159 delphi += 2.0 * pi;
0160 if (delphi > pi)
0161 delphi -= 2.0 * pi;
0162
0163 float tmp_tp_VtxZ = tmp_tp_vz + tmp_tp_t * delphi / (2.0 * K);
0164 float tmp_tp_VtxR = sqrt(tmp_tp_vx * tmp_tp_vx + tmp_tp_vy * tmp_tp_vy);
0165 float tmp_tp_d0 = tmp_tp_charge * tmp_tp_rp - (1. / (2. * K));
0166
0167
0168
0169 float other_d0 = -tmp_tp_vx * sin(tmp_tp_phi) + tmp_tp_vy * cos(tmp_tp_phi);
0170 tmp_tp_d0 = tmp_tp_d0 * (-1);
0171 if (K == 0) {
0172 tmp_tp_d0 = other_d0;
0173 tmp_tp_VtxZ = tmp_tp_vz;
0174 }
0175 if (std::fabs(tmp_tp_VtxZ) > TP_maxVtxZ)
0176 continue;
0177
0178
0179 if (tmp_tp_VtxR < 1.0) {
0180 tp_pt->Fill(tmp_tp_pt);
0181 if (tmp_tp_pt <= 10)
0182 tp_pt_zoom->Fill(tmp_tp_pt);
0183 }
0184 if (tmp_tp_pt < TP_minPt)
0185 continue;
0186 tp_VtxR->Fill(tmp_tp_VtxR);
0187 if (tmp_tp_VtxR > 1.0)
0188 continue;
0189 tp_eta->Fill(tmp_tp_eta);
0190 tp_d0->Fill(tmp_tp_d0);
0191 tp_VtxZ->Fill(tmp_tp_VtxZ);
0192
0193 if (nStubTP < TP_minNStub || nStubLayerTP < TP_minNLayersStub)
0194 continue;
0195
0196
0197
0198 int tp_nMatch = 0;
0199 int i_track = -1;
0200 float i_chi2dof = 99999;
0201 if (MCTruthTTTrackHandle.isValid()) {
0202 std::vector<edm::Ptr<TTTrack<Ref_Phase2TrackerDigi_>>> matchedTracks =
0203 MCTruthTTTrackHandle->findTTTrackPtrs(tp_ptr);
0204
0205
0206
0207
0208
0209
0210 int trkCounter = 0;
0211 for (const auto &thisTrack : matchedTracks) {
0212 if (!MCTruthTTTrackHandle->isGenuine(thisTrack))
0213 continue;
0214
0215
0216
0217
0218 int tmp_trk_nstub = thisTrack->getStubRefs().size();
0219 if (tmp_trk_nstub < L1Tk_minNStub)
0220 continue;
0221 float dmatch_pt = 999;
0222 float dmatch_eta = 999;
0223 float dmatch_phi = 999;
0224 int match_id = 999;
0225
0226 edm::Ptr<TrackingParticle> my_tp = MCTruthTTTrackHandle->findTrackingParticlePtr(thisTrack);
0227 dmatch_pt = std::fabs(my_tp->p4().pt() - tmp_tp_pt);
0228 dmatch_eta = std::fabs(my_tp->p4().eta() - tmp_tp_eta);
0229 dmatch_phi = std::fabs(my_tp->p4().phi() - tmp_tp_phi);
0230 match_id = my_tp->pdgId();
0231 float tmp_trk_chi2dof = thisTrack->chi2Red();
0232
0233
0234 if (dmatch_pt < 0.1 && dmatch_eta < 0.1 && dmatch_phi < 0.1 && tmp_tp_pdgid == match_id) {
0235 tp_nMatch++;
0236 if (i_track < 0 || tmp_trk_chi2dof < i_chi2dof) {
0237 i_track = trkCounter;
0238 i_chi2dof = tmp_trk_chi2dof;
0239 }
0240 }
0241 trkCounter++;
0242 }
0243
0244 if (tp_nMatch < 1)
0245 continue;
0246
0247 float tmp_matchtrk_pt = -999;
0248 float tmp_matchtrk_eta = -999;
0249 float tmp_matchtrk_phi = -999;
0250 float tmp_matchtrk_VtxZ = -999;
0251 float tmp_matchtrk_chi2dof = -999;
0252 int tmp_matchTrk_nStub = -999;
0253 float tmp_matchtrk_d0 = -999;
0254
0255 tmp_matchtrk_pt = matchedTracks[i_track]->momentum().perp();
0256 tmp_matchtrk_eta = matchedTracks[i_track]->momentum().eta();
0257 tmp_matchtrk_phi = matchedTracks[i_track]->momentum().phi();
0258 tmp_matchtrk_VtxZ = matchedTracks[i_track]->z0();
0259 tmp_matchtrk_chi2dof = matchedTracks[i_track]->chi2Red();
0260 tmp_matchTrk_nStub = (int)matchedTracks[i_track]->getStubRefs().size();
0261
0262
0263 float tmp_matchtrk_x0 = matchedTracks[i_track]->POCA().x();
0264 float tmp_matchtrk_y0 = matchedTracks[i_track]->POCA().y();
0265 tmp_matchtrk_d0 = -tmp_matchtrk_x0 * sin(tmp_matchtrk_phi) + tmp_matchtrk_y0 * cos(tmp_matchtrk_phi);
0266
0267
0268 if (tmp_matchTrk_nStub < L1Tk_minNStub || tmp_matchtrk_chi2dof > L1Tk_maxChi2dof)
0269 continue;
0270
0271
0272 match_tp_pt->Fill(tmp_tp_pt);
0273 if (tmp_tp_pt > 0 && tmp_tp_pt <= 10)
0274 match_tp_pt_zoom->Fill(tmp_tp_pt);
0275 match_tp_eta->Fill(tmp_tp_eta);
0276 match_tp_d0->Fill(tmp_tp_d0);
0277 match_tp_VtxR->Fill(tmp_tp_VtxR);
0278 match_tp_VtxZ->Fill(tmp_tp_VtxZ);
0279
0280
0281 float pt_diff = tmp_matchtrk_pt - tmp_tp_pt;
0282 float pt_res = pt_diff / tmp_tp_pt;
0283 float eta_res = tmp_matchtrk_eta - tmp_tp_eta;
0284 float phi_res = tmp_matchtrk_phi - tmp_tp_phi;
0285 float VtxZ_res = tmp_matchtrk_VtxZ - tmp_tp_VtxZ;
0286 float d0_res = tmp_matchtrk_d0 - tmp_tp_d0;
0287
0288
0289 res_pt->Fill(pt_diff);
0290 res_ptRel->Fill(pt_res);
0291 res_eta->Fill(eta_res);
0292
0293
0294
0295 if (std::fabs(tmp_tp_eta) >= 0 && std::fabs(tmp_tp_eta) < 0.7) {
0296 reseta_eta0to0p7->Fill(eta_res);
0297 resphi_eta0to0p7->Fill(phi_res);
0298 resVtxZ_eta0to0p7->Fill(VtxZ_res);
0299 resd0_eta0to0p7->Fill(d0_res);
0300 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0301 respt_eta0to0p7_pt2to3->Fill(pt_res);
0302 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0303 respt_eta0to0p7_pt3to8->Fill(pt_res);
0304 else if (tmp_tp_pt >= 8)
0305 respt_eta0to0p7_pt8toInf->Fill(pt_res);
0306 } else if (std::fabs(tmp_tp_eta) >= 0.7 && std::fabs(tmp_tp_eta) < 1.0) {
0307 reseta_eta0p7to1->Fill(eta_res);
0308 resphi_eta0p7to1->Fill(phi_res);
0309 resVtxZ_eta0p7to1->Fill(VtxZ_res);
0310 resd0_eta0p7to1->Fill(d0_res);
0311 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0312 respt_eta0p7to1_pt2to3->Fill(pt_res);
0313 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0314 respt_eta0p7to1_pt3to8->Fill(pt_res);
0315 else if (tmp_tp_pt >= 8)
0316 respt_eta0p7to1_pt8toInf->Fill(pt_res);
0317 } else if (std::fabs(tmp_tp_eta) >= 1.0 && std::fabs(tmp_tp_eta) < 1.2) {
0318 reseta_eta1to1p2->Fill(eta_res);
0319 resphi_eta1to1p2->Fill(phi_res);
0320 resVtxZ_eta1to1p2->Fill(VtxZ_res);
0321 resd0_eta1to1p2->Fill(d0_res);
0322 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0323 respt_eta1to1p2_pt2to3->Fill(pt_res);
0324 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0325 respt_eta1to1p2_pt3to8->Fill(pt_res);
0326 else if (tmp_tp_pt >= 8)
0327 respt_eta1to1p2_pt8toInf->Fill(pt_res);
0328 } else if (std::fabs(tmp_tp_eta) >= 1.2 && std::fabs(tmp_tp_eta) < 1.6) {
0329 reseta_eta1p2to1p6->Fill(eta_res);
0330 resphi_eta1p2to1p6->Fill(phi_res);
0331 resVtxZ_eta1p2to1p6->Fill(VtxZ_res);
0332 resd0_eta1p2to1p6->Fill(d0_res);
0333 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0334 respt_eta1p2to1p6_pt2to3->Fill(pt_res);
0335 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0336 respt_eta1p2to1p6_pt3to8->Fill(pt_res);
0337 else if (tmp_tp_pt >= 8)
0338 respt_eta1p2to1p6_pt8toInf->Fill(pt_res);
0339 } else if (std::fabs(tmp_tp_eta) >= 1.6 && std::fabs(tmp_tp_eta) < 2.0) {
0340 reseta_eta1p6to2->Fill(eta_res);
0341 resphi_eta1p6to2->Fill(phi_res);
0342 resVtxZ_eta1p6to2->Fill(VtxZ_res);
0343 resd0_eta1p6to2->Fill(d0_res);
0344 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0345 respt_eta1p6to2_pt2to3->Fill(pt_res);
0346 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0347 respt_eta1p6to2_pt3to8->Fill(pt_res);
0348 else if (tmp_tp_pt >= 8)
0349 respt_eta1p6to2_pt8toInf->Fill(pt_res);
0350 } else if (std::fabs(tmp_tp_eta) >= 2.0 && std::fabs(tmp_tp_eta) <= 2.4) {
0351 reseta_eta2to2p4->Fill(eta_res);
0352 resphi_eta2to2p4->Fill(phi_res);
0353 resVtxZ_eta2to2p4->Fill(VtxZ_res);
0354 resd0_eta2to2p4->Fill(d0_res);
0355 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0356 respt_eta2to2p4_pt2to3->Fill(pt_res);
0357 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0358 respt_eta2to2p4_pt3to8->Fill(pt_res);
0359 else if (tmp_tp_pt >= 8)
0360 respt_eta2to2p4_pt8toInf->Fill(pt_res);
0361 }
0362 }
0363 }
0364 }
0365
0366
0367
0368 void OuterTrackerMonitorTrackingParticles::bookHistograms(DQMStore::IBooker &iBooker,
0369 edm::Run const &run,
0370 edm::EventSetup const &es) {
0371
0372 std::string HistoName;
0373 iBooker.setCurrentFolder(topFolderName_ + "/trackParticles");
0374
0375
0376 edm::ParameterSet psTrackParts_Pt = conf_.getParameter<edm::ParameterSet>("TH1TrackParts_Pt");
0377 HistoName = "trackParts_Pt";
0378 trackParts_Pt = iBooker.book1D(HistoName,
0379 HistoName,
0380 psTrackParts_Pt.getParameter<int32_t>("Nbinsx"),
0381 psTrackParts_Pt.getParameter<double>("xmin"),
0382 psTrackParts_Pt.getParameter<double>("xmax"));
0383 trackParts_Pt->setAxisTitle("p_{T} [GeV]", 1);
0384 trackParts_Pt->setAxisTitle("# tracking particles", 2);
0385
0386
0387 edm::ParameterSet psTrackParts_Eta = conf_.getParameter<edm::ParameterSet>("TH1TrackParts_Eta");
0388 HistoName = "trackParts_Eta";
0389 trackParts_Eta = iBooker.book1D(HistoName,
0390 HistoName,
0391 psTrackParts_Eta.getParameter<int32_t>("Nbinsx"),
0392 psTrackParts_Eta.getParameter<double>("xmin"),
0393 psTrackParts_Eta.getParameter<double>("xmax"));
0394 trackParts_Eta->setAxisTitle("#eta", 1);
0395 trackParts_Eta->setAxisTitle("# tracking particles", 2);
0396
0397
0398 edm::ParameterSet psTrackParts_Phi = conf_.getParameter<edm::ParameterSet>("TH1TrackParts_Phi");
0399 HistoName = "trackParts_Phi";
0400 trackParts_Phi = iBooker.book1D(HistoName,
0401 HistoName,
0402 psTrackParts_Phi.getParameter<int32_t>("Nbinsx"),
0403 psTrackParts_Phi.getParameter<double>("xmin"),
0404 psTrackParts_Phi.getParameter<double>("xmax"));
0405 trackParts_Phi->setAxisTitle("#phi", 1);
0406 trackParts_Phi->setAxisTitle("# tracking particles", 2);
0407
0408
0409 iBooker.setCurrentFolder(topFolderName_ + "/Tracks/Efficiency");
0410
0411 edm::ParameterSet psEffic_pt = conf_.getParameter<edm::ParameterSet>("TH1Effic_pt");
0412 HistoName = "tp_pt";
0413 tp_pt = iBooker.book1D(HistoName,
0414 HistoName,
0415 psEffic_pt.getParameter<int32_t>("Nbinsx"),
0416 psEffic_pt.getParameter<double>("xmin"),
0417 psEffic_pt.getParameter<double>("xmax"));
0418 tp_pt->setAxisTitle("p_{T} [GeV]", 1);
0419 tp_pt->setAxisTitle("# tracking particles", 2);
0420
0421
0422 HistoName = "match_tp_pt";
0423 match_tp_pt = iBooker.book1D(HistoName,
0424 HistoName,
0425 psEffic_pt.getParameter<int32_t>("Nbinsx"),
0426 psEffic_pt.getParameter<double>("xmin"),
0427 psEffic_pt.getParameter<double>("xmax"));
0428 match_tp_pt->setAxisTitle("p_{T} [GeV]", 1);
0429 match_tp_pt->setAxisTitle("# matched tracking particles", 2);
0430
0431
0432 edm::ParameterSet psEffic_pt_zoom = conf_.getParameter<edm::ParameterSet>("TH1Effic_pt_zoom");
0433 HistoName = "tp_pt_zoom";
0434 tp_pt_zoom = iBooker.book1D(HistoName,
0435 HistoName,
0436 psEffic_pt_zoom.getParameter<int32_t>("Nbinsx"),
0437 psEffic_pt_zoom.getParameter<double>("xmin"),
0438 psEffic_pt_zoom.getParameter<double>("xmax"));
0439 tp_pt_zoom->setAxisTitle("p_{T} [GeV]", 1);
0440 tp_pt_zoom->setAxisTitle("# tracking particles", 2);
0441
0442
0443 HistoName = "match_tp_pt_zoom";
0444 match_tp_pt_zoom = iBooker.book1D(HistoName,
0445 HistoName,
0446 psEffic_pt_zoom.getParameter<int32_t>("Nbinsx"),
0447 psEffic_pt_zoom.getParameter<double>("xmin"),
0448 psEffic_pt_zoom.getParameter<double>("xmax"));
0449 match_tp_pt_zoom->setAxisTitle("p_{T} [GeV]", 1);
0450 match_tp_pt_zoom->setAxisTitle("# matched tracking particles", 2);
0451
0452
0453 edm::ParameterSet psEffic_eta = conf_.getParameter<edm::ParameterSet>("TH1Effic_eta");
0454 HistoName = "tp_eta";
0455 tp_eta = iBooker.book1D(HistoName,
0456 HistoName,
0457 psEffic_eta.getParameter<int32_t>("Nbinsx"),
0458 psEffic_eta.getParameter<double>("xmin"),
0459 psEffic_eta.getParameter<double>("xmax"));
0460 tp_eta->setAxisTitle("#eta", 1);
0461 tp_eta->setAxisTitle("# tracking particles", 2);
0462
0463
0464 HistoName = "match_tp_eta";
0465 match_tp_eta = iBooker.book1D(HistoName,
0466 HistoName,
0467 psEffic_eta.getParameter<int32_t>("Nbinsx"),
0468 psEffic_eta.getParameter<double>("xmin"),
0469 psEffic_eta.getParameter<double>("xmax"));
0470 match_tp_eta->setAxisTitle("#eta", 1);
0471 match_tp_eta->setAxisTitle("# matched tracking particles", 2);
0472
0473
0474 edm::ParameterSet psEffic_d0 = conf_.getParameter<edm::ParameterSet>("TH1Effic_d0");
0475 HistoName = "tp_d0";
0476 tp_d0 = iBooker.book1D(HistoName,
0477 HistoName,
0478 psEffic_d0.getParameter<int32_t>("Nbinsx"),
0479 psEffic_d0.getParameter<double>("xmin"),
0480 psEffic_d0.getParameter<double>("xmax"));
0481 tp_d0->setAxisTitle("d_{0} [cm]", 1);
0482 tp_d0->setAxisTitle("# tracking particles", 2);
0483
0484
0485 HistoName = "match_tp_d0";
0486 match_tp_d0 = iBooker.book1D(HistoName,
0487 HistoName,
0488 psEffic_d0.getParameter<int32_t>("Nbinsx"),
0489 psEffic_d0.getParameter<double>("xmin"),
0490 psEffic_d0.getParameter<double>("xmax"));
0491 match_tp_d0->setAxisTitle("d_{0} [cm]", 1);
0492 match_tp_d0->setAxisTitle("# matched tracking particles", 2);
0493
0494
0495 edm::ParameterSet psEffic_VtxR = conf_.getParameter<edm::ParameterSet>("TH1Effic_VtxR");
0496 HistoName = "tp_VtxR";
0497 tp_VtxR = iBooker.book1D(HistoName,
0498 HistoName,
0499 psEffic_VtxR.getParameter<int32_t>("Nbinsx"),
0500 psEffic_VtxR.getParameter<double>("xmin"),
0501 psEffic_VtxR.getParameter<double>("xmax"));
0502 tp_VtxR->setAxisTitle("d_{xy} [cm]", 1);
0503 tp_VtxR->setAxisTitle("# tracking particles", 2);
0504
0505
0506 HistoName = "match_tp_VtxR";
0507 match_tp_VtxR = iBooker.book1D(HistoName,
0508 HistoName,
0509 psEffic_VtxR.getParameter<int32_t>("Nbinsx"),
0510 psEffic_VtxR.getParameter<double>("xmin"),
0511 psEffic_VtxR.getParameter<double>("xmax"));
0512 match_tp_VtxR->setAxisTitle("d_{xy} [cm]", 1);
0513 match_tp_VtxR->setAxisTitle("# matched tracking particles", 2);
0514
0515
0516 edm::ParameterSet psEffic_VtxZ = conf_.getParameter<edm::ParameterSet>("TH1Effic_VtxZ");
0517 HistoName = "tp_VtxZ";
0518 tp_VtxZ = iBooker.book1D(HistoName,
0519 HistoName,
0520 psEffic_VtxZ.getParameter<int32_t>("Nbinsx"),
0521 psEffic_VtxZ.getParameter<double>("xmin"),
0522 psEffic_VtxZ.getParameter<double>("xmax"));
0523 tp_VtxZ->setAxisTitle("z_{0} [cm]", 1);
0524 tp_VtxZ->setAxisTitle("# tracking particles", 2);
0525
0526
0527 HistoName = "match_tp_VtxZ";
0528 match_tp_VtxZ = iBooker.book1D(HistoName,
0529 HistoName,
0530 psEffic_VtxZ.getParameter<int32_t>("Nbinsx"),
0531 psEffic_VtxZ.getParameter<double>("xmin"),
0532 psEffic_VtxZ.getParameter<double>("xmax"));
0533 match_tp_VtxZ->setAxisTitle("z_{0} [cm]", 1);
0534 match_tp_VtxZ->setAxisTitle("# matched tracking particles", 2);
0535
0536
0537 iBooker.setCurrentFolder(topFolderName_ + "/Tracks/Resolution");
0538
0539 edm::ParameterSet psRes_pt = conf_.getParameter<edm::ParameterSet>("TH1Res_pt");
0540 HistoName = "res_pt";
0541 res_pt = iBooker.book1D(HistoName,
0542 HistoName,
0543 psRes_pt.getParameter<int32_t>("Nbinsx"),
0544 psRes_pt.getParameter<double>("xmin"),
0545 psRes_pt.getParameter<double>("xmax"));
0546 res_pt->setAxisTitle("p_{T} [GeV]", 1);
0547 res_pt->setAxisTitle("# tracking particles", 2);
0548
0549
0550 edm::ParameterSet psRes_eta = conf_.getParameter<edm::ParameterSet>("TH1Res_eta");
0551 HistoName = "res_eta";
0552 res_eta = iBooker.book1D(HistoName,
0553 HistoName,
0554 psRes_eta.getParameter<int32_t>("Nbinsx"),
0555 psRes_eta.getParameter<double>("xmin"),
0556 psRes_eta.getParameter<double>("xmax"));
0557 res_eta->setAxisTitle("#eta", 1);
0558 res_eta->setAxisTitle("# tracking particles", 2);
0559
0560
0561 edm::ParameterSet psRes_ptRel = conf_.getParameter<edm::ParameterSet>("TH1Res_ptRel");
0562 HistoName = "res_ptRel";
0563 res_ptRel = iBooker.book1D(HistoName,
0564 HistoName,
0565 psRes_ptRel.getParameter<int32_t>("Nbinsx"),
0566 psRes_ptRel.getParameter<double>("xmin"),
0567 psRes_ptRel.getParameter<double>("xmax"));
0568 res_ptRel->setAxisTitle("Relative p_{T} [GeV]", 1);
0569 res_ptRel->setAxisTitle("# tracking particles", 2);
0570
0571
0572
0573 HistoName = "reseta_eta0to0p7";
0574 reseta_eta0to0p7 = iBooker.book1D(HistoName,
0575 HistoName,
0576 psRes_eta.getParameter<int32_t>("Nbinsx"),
0577 psRes_eta.getParameter<double>("xmin"),
0578 psRes_eta.getParameter<double>("xmax"));
0579 reseta_eta0to0p7->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0580 reseta_eta0to0p7->setAxisTitle("# tracking particles", 2);
0581
0582
0583 HistoName = "reseta_eta0p7to1";
0584 reseta_eta0p7to1 = iBooker.book1D(HistoName,
0585 HistoName,
0586 psRes_eta.getParameter<int32_t>("Nbinsx"),
0587 psRes_eta.getParameter<double>("xmin"),
0588 psRes_eta.getParameter<double>("xmax"));
0589 reseta_eta0p7to1->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0590 reseta_eta0p7to1->setAxisTitle("# tracking particles", 2);
0591
0592
0593 HistoName = "reseta_eta1to1p2";
0594 reseta_eta1to1p2 = iBooker.book1D(HistoName,
0595 HistoName,
0596 psRes_eta.getParameter<int32_t>("Nbinsx"),
0597 psRes_eta.getParameter<double>("xmin"),
0598 psRes_eta.getParameter<double>("xmax"));
0599 reseta_eta1to1p2->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0600 reseta_eta1to1p2->setAxisTitle("# tracking particles", 2);
0601
0602
0603 HistoName = "reseta_eta1p2to1p6";
0604 reseta_eta1p2to1p6 = iBooker.book1D(HistoName,
0605 HistoName,
0606 psRes_eta.getParameter<int32_t>("Nbinsx"),
0607 psRes_eta.getParameter<double>("xmin"),
0608 psRes_eta.getParameter<double>("xmax"));
0609 reseta_eta1p2to1p6->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0610 reseta_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
0611
0612
0613 HistoName = "reseta_eta1p6to2";
0614 reseta_eta1p6to2 = iBooker.book1D(HistoName,
0615 HistoName,
0616 psRes_eta.getParameter<int32_t>("Nbinsx"),
0617 psRes_eta.getParameter<double>("xmin"),
0618 psRes_eta.getParameter<double>("xmax"));
0619 reseta_eta1p6to2->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0620 reseta_eta1p6to2->setAxisTitle("# tracking particles", 2);
0621
0622
0623 HistoName = "reseta_eta2to2p4";
0624 reseta_eta2to2p4 = iBooker.book1D(HistoName,
0625 HistoName,
0626 psRes_eta.getParameter<int32_t>("Nbinsx"),
0627 psRes_eta.getParameter<double>("xmin"),
0628 psRes_eta.getParameter<double>("xmax"));
0629 reseta_eta2to2p4->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0630 reseta_eta2to2p4->setAxisTitle("# tracking particles", 2);
0631
0632
0633
0634
0635 HistoName = "respt_eta0to0p7_pt2to3";
0636 respt_eta0to0p7_pt2to3 = iBooker.book1D(HistoName,
0637 HistoName,
0638 psRes_pt.getParameter<int32_t>("Nbinsx"),
0639 psRes_pt.getParameter<double>("xmin"),
0640 psRes_pt.getParameter<double>("xmax"));
0641 respt_eta0to0p7_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0642 respt_eta0to0p7_pt2to3->setAxisTitle("# tracking particles", 2);
0643
0644
0645 HistoName = "respt_eta0p7to1_pt2to3";
0646 respt_eta0p7to1_pt2to3 = iBooker.book1D(HistoName,
0647 HistoName,
0648 psRes_pt.getParameter<int32_t>("Nbinsx"),
0649 psRes_pt.getParameter<double>("xmin"),
0650 psRes_pt.getParameter<double>("xmax"));
0651 respt_eta0p7to1_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0652 respt_eta0p7to1_pt2to3->setAxisTitle("# tracking particles", 2);
0653
0654
0655 HistoName = "respt_eta1to1p2_pt2to3";
0656 respt_eta1to1p2_pt2to3 = iBooker.book1D(HistoName,
0657 HistoName,
0658 psRes_pt.getParameter<int32_t>("Nbinsx"),
0659 psRes_pt.getParameter<double>("xmin"),
0660 psRes_pt.getParameter<double>("xmax"));
0661 respt_eta1to1p2_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0662 respt_eta1to1p2_pt2to3->setAxisTitle("# tracking particles", 2);
0663
0664
0665 HistoName = "respt_eta1p2to1p6_pt2to3";
0666 respt_eta1p2to1p6_pt2to3 = iBooker.book1D(HistoName,
0667 HistoName,
0668 psRes_pt.getParameter<int32_t>("Nbinsx"),
0669 psRes_pt.getParameter<double>("xmin"),
0670 psRes_pt.getParameter<double>("xmax"));
0671 respt_eta1p2to1p6_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0672 respt_eta1p2to1p6_pt2to3->setAxisTitle("# tracking particles", 2);
0673
0674
0675 HistoName = "respt_eta1p6to2_pt2to3";
0676 respt_eta1p6to2_pt2to3 = iBooker.book1D(HistoName,
0677 HistoName,
0678 psRes_pt.getParameter<int32_t>("Nbinsx"),
0679 psRes_pt.getParameter<double>("xmin"),
0680 psRes_pt.getParameter<double>("xmax"));
0681 respt_eta1p6to2_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0682 respt_eta1p6to2_pt2to3->setAxisTitle("# tracking particles", 2);
0683
0684
0685 HistoName = "respt_eta2to2p4_pt2to3";
0686 respt_eta2to2p4_pt2to3 = iBooker.book1D(HistoName,
0687 HistoName,
0688 psRes_pt.getParameter<int32_t>("Nbinsx"),
0689 psRes_pt.getParameter<double>("xmin"),
0690 psRes_pt.getParameter<double>("xmax"));
0691 respt_eta2to2p4_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0692 respt_eta2to2p4_pt2to3->setAxisTitle("# tracking particles", 2);
0693
0694
0695
0696 HistoName = "respt_eta0to0p7_pt3to8";
0697 respt_eta0to0p7_pt3to8 = iBooker.book1D(HistoName,
0698 HistoName,
0699 psRes_pt.getParameter<int32_t>("Nbinsx"),
0700 psRes_pt.getParameter<double>("xmin"),
0701 psRes_pt.getParameter<double>("xmax"));
0702 respt_eta0to0p7_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0703 respt_eta0to0p7_pt3to8->setAxisTitle("# tracking particles", 2);
0704
0705
0706 HistoName = "respt_eta0p7to1_pt3to8";
0707 respt_eta0p7to1_pt3to8 = iBooker.book1D(HistoName,
0708 HistoName,
0709 psRes_pt.getParameter<int32_t>("Nbinsx"),
0710 psRes_pt.getParameter<double>("xmin"),
0711 psRes_pt.getParameter<double>("xmax"));
0712 respt_eta0p7to1_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0713 respt_eta0p7to1_pt3to8->setAxisTitle("# tracking particles", 2);
0714
0715
0716 HistoName = "respt_eta1to1p2_pt3to8";
0717 respt_eta1to1p2_pt3to8 = iBooker.book1D(HistoName,
0718 HistoName,
0719 psRes_pt.getParameter<int32_t>("Nbinsx"),
0720 psRes_pt.getParameter<double>("xmin"),
0721 psRes_pt.getParameter<double>("xmax"));
0722 respt_eta1to1p2_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0723 respt_eta1to1p2_pt3to8->setAxisTitle("# tracking particles", 2);
0724
0725
0726 HistoName = "respt_eta1p2to1p6_pt3to8";
0727 respt_eta1p2to1p6_pt3to8 = iBooker.book1D(HistoName,
0728 HistoName,
0729 psRes_pt.getParameter<int32_t>("Nbinsx"),
0730 psRes_pt.getParameter<double>("xmin"),
0731 psRes_pt.getParameter<double>("xmax"));
0732 respt_eta1p2to1p6_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0733 respt_eta1p2to1p6_pt3to8->setAxisTitle("# tracking particles", 2);
0734
0735
0736 HistoName = "respt_eta1p6to2_pt3to8";
0737 respt_eta1p6to2_pt3to8 = iBooker.book1D(HistoName,
0738 HistoName,
0739 psRes_pt.getParameter<int32_t>("Nbinsx"),
0740 psRes_pt.getParameter<double>("xmin"),
0741 psRes_pt.getParameter<double>("xmax"));
0742 respt_eta1p6to2_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0743 respt_eta1p6to2_pt3to8->setAxisTitle("# tracking particles", 2);
0744
0745
0746 HistoName = "respt_eta2to2p4_pt3to8";
0747 respt_eta2to2p4_pt3to8 = iBooker.book1D(HistoName,
0748 HistoName,
0749 psRes_pt.getParameter<int32_t>("Nbinsx"),
0750 psRes_pt.getParameter<double>("xmin"),
0751 psRes_pt.getParameter<double>("xmax"));
0752 respt_eta2to2p4_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0753 respt_eta2to2p4_pt3to8->setAxisTitle("# tracking particles", 2);
0754
0755
0756
0757 HistoName = "respt_eta0to0p7_pt8toInf";
0758 respt_eta0to0p7_pt8toInf = iBooker.book1D(HistoName,
0759 HistoName,
0760 psRes_pt.getParameter<int32_t>("Nbinsx"),
0761 psRes_pt.getParameter<double>("xmin"),
0762 psRes_pt.getParameter<double>("xmax"));
0763 respt_eta0to0p7_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0764 respt_eta0to0p7_pt8toInf->setAxisTitle("# tracking particles", 2);
0765
0766
0767 HistoName = "respt_eta0p7to1_pt8toInf";
0768 respt_eta0p7to1_pt8toInf = iBooker.book1D(HistoName,
0769 HistoName,
0770 psRes_pt.getParameter<int32_t>("Nbinsx"),
0771 psRes_pt.getParameter<double>("xmin"),
0772 psRes_pt.getParameter<double>("xmax"));
0773 respt_eta0p7to1_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0774 respt_eta0p7to1_pt8toInf->setAxisTitle("# tracking particles", 2);
0775
0776
0777 HistoName = "respt_eta1to1p2_pt8toInf";
0778 respt_eta1to1p2_pt8toInf = iBooker.book1D(HistoName,
0779 HistoName,
0780 psRes_pt.getParameter<int32_t>("Nbinsx"),
0781 psRes_pt.getParameter<double>("xmin"),
0782 psRes_pt.getParameter<double>("xmax"));
0783 respt_eta1to1p2_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0784 respt_eta1to1p2_pt8toInf->setAxisTitle("# tracking particles", 2);
0785
0786
0787 HistoName = "respt_eta1p2to1p6_pt8toInf";
0788 respt_eta1p2to1p6_pt8toInf = iBooker.book1D(HistoName,
0789 HistoName,
0790 psRes_pt.getParameter<int32_t>("Nbinsx"),
0791 psRes_pt.getParameter<double>("xmin"),
0792 psRes_pt.getParameter<double>("xmax"));
0793 respt_eta1p2to1p6_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0794 respt_eta1p2to1p6_pt8toInf->setAxisTitle("# tracking particles", 2);
0795
0796
0797 HistoName = "respt_eta1p6to2_pt8toInf";
0798 respt_eta1p6to2_pt8toInf = iBooker.book1D(HistoName,
0799 HistoName,
0800 psRes_pt.getParameter<int32_t>("Nbinsx"),
0801 psRes_pt.getParameter<double>("xmin"),
0802 psRes_pt.getParameter<double>("xmax"));
0803 respt_eta1p6to2_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0804 respt_eta1p6to2_pt8toInf->setAxisTitle("# tracking particles", 2);
0805
0806
0807 HistoName = "respt_eta2to2p4_pt8toInf";
0808 respt_eta2to2p4_pt8toInf = iBooker.book1D(HistoName,
0809 HistoName,
0810 psRes_pt.getParameter<int32_t>("Nbinsx"),
0811 psRes_pt.getParameter<double>("xmin"),
0812 psRes_pt.getParameter<double>("xmax"));
0813 respt_eta2to2p4_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0814 respt_eta2to2p4_pt8toInf->setAxisTitle("# tracking particles", 2);
0815
0816
0817
0818 edm::ParameterSet psRes_phi = conf_.getParameter<edm::ParameterSet>("TH1Res_phi");
0819 HistoName = "resphi_eta0to0p7";
0820 resphi_eta0to0p7 = iBooker.book1D(HistoName,
0821 HistoName,
0822 psRes_phi.getParameter<int32_t>("Nbinsx"),
0823 psRes_phi.getParameter<double>("xmin"),
0824 psRes_phi.getParameter<double>("xmax"));
0825 resphi_eta0to0p7->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0826 resphi_eta0to0p7->setAxisTitle("# tracking particles", 2);
0827
0828
0829 HistoName = "resphi_eta0p7to1";
0830 resphi_eta0p7to1 = iBooker.book1D(HistoName,
0831 HistoName,
0832 psRes_phi.getParameter<int32_t>("Nbinsx"),
0833 psRes_phi.getParameter<double>("xmin"),
0834 psRes_phi.getParameter<double>("xmax"));
0835 resphi_eta0p7to1->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0836 resphi_eta0p7to1->setAxisTitle("# tracking particles", 2);
0837
0838
0839 HistoName = "resphi_eta1to1p2";
0840 resphi_eta1to1p2 = iBooker.book1D(HistoName,
0841 HistoName,
0842 psRes_phi.getParameter<int32_t>("Nbinsx"),
0843 psRes_phi.getParameter<double>("xmin"),
0844 psRes_phi.getParameter<double>("xmax"));
0845 resphi_eta1to1p2->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0846 resphi_eta1to1p2->setAxisTitle("# tracking particles", 2);
0847
0848
0849 HistoName = "resphi_eta1p2to1p6";
0850 resphi_eta1p2to1p6 = iBooker.book1D(HistoName,
0851 HistoName,
0852 psRes_phi.getParameter<int32_t>("Nbinsx"),
0853 psRes_phi.getParameter<double>("xmin"),
0854 psRes_phi.getParameter<double>("xmax"));
0855 resphi_eta1p2to1p6->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0856 resphi_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
0857
0858
0859 HistoName = "resphi_eta1p6to2";
0860 resphi_eta1p6to2 = iBooker.book1D(HistoName,
0861 HistoName,
0862 psRes_phi.getParameter<int32_t>("Nbinsx"),
0863 psRes_phi.getParameter<double>("xmin"),
0864 psRes_phi.getParameter<double>("xmax"));
0865 resphi_eta1p6to2->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0866 resphi_eta1p6to2->setAxisTitle("# tracking particles", 2);
0867
0868
0869 HistoName = "resphi_eta2to2p4";
0870 resphi_eta2to2p4 = iBooker.book1D(HistoName,
0871 HistoName,
0872 psRes_phi.getParameter<int32_t>("Nbinsx"),
0873 psRes_phi.getParameter<double>("xmin"),
0874 psRes_phi.getParameter<double>("xmax"));
0875 resphi_eta2to2p4->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0876 resphi_eta2to2p4->setAxisTitle("# tracking particles", 2);
0877
0878
0879
0880 edm::ParameterSet psRes_VtxZ = conf_.getParameter<edm::ParameterSet>("TH1Res_VtxZ");
0881 HistoName = "resVtxZ_eta0to0p7";
0882 resVtxZ_eta0to0p7 = iBooker.book1D(HistoName,
0883 HistoName,
0884 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
0885 psRes_VtxZ.getParameter<double>("xmin"),
0886 psRes_VtxZ.getParameter<double>("xmax"));
0887 resVtxZ_eta0to0p7->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
0888 resVtxZ_eta0to0p7->setAxisTitle("# tracking particles", 2);
0889
0890
0891 HistoName = "resVtxZ_eta0p7to1";
0892 resVtxZ_eta0p7to1 = iBooker.book1D(HistoName,
0893 HistoName,
0894 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
0895 psRes_VtxZ.getParameter<double>("xmin"),
0896 psRes_VtxZ.getParameter<double>("xmax"));
0897 resVtxZ_eta0p7to1->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
0898 resVtxZ_eta0p7to1->setAxisTitle("# tracking particles", 2);
0899
0900
0901 HistoName = "resVtxZ_eta1to1p2";
0902 resVtxZ_eta1to1p2 = iBooker.book1D(HistoName,
0903 HistoName,
0904 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
0905 psRes_VtxZ.getParameter<double>("xmin"),
0906 psRes_VtxZ.getParameter<double>("xmax"));
0907 resVtxZ_eta1to1p2->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
0908 resVtxZ_eta1to1p2->setAxisTitle("# tracking particles", 2);
0909
0910
0911 HistoName = "resVtxZ_eta1p2to1p6";
0912 resVtxZ_eta1p2to1p6 = iBooker.book1D(HistoName,
0913 HistoName,
0914 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
0915 psRes_VtxZ.getParameter<double>("xmin"),
0916 psRes_VtxZ.getParameter<double>("xmax"));
0917 resVtxZ_eta1p2to1p6->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
0918 resVtxZ_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
0919
0920
0921 HistoName = "resVtxZ_eta1p6to2";
0922 resVtxZ_eta1p6to2 = iBooker.book1D(HistoName,
0923 HistoName,
0924 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
0925 psRes_VtxZ.getParameter<double>("xmin"),
0926 psRes_VtxZ.getParameter<double>("xmax"));
0927 resVtxZ_eta1p6to2->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
0928 resVtxZ_eta1p6to2->setAxisTitle("# tracking particles", 2);
0929
0930
0931 HistoName = "resVtxZ_eta2to2p4";
0932 resVtxZ_eta2to2p4 = iBooker.book1D(HistoName,
0933 HistoName,
0934 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
0935 psRes_VtxZ.getParameter<double>("xmin"),
0936 psRes_VtxZ.getParameter<double>("xmax"));
0937 resVtxZ_eta2to2p4->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
0938 resVtxZ_eta2to2p4->setAxisTitle("# tracking particles", 2);
0939
0940
0941
0942 edm::ParameterSet psRes_d0 = conf_.getParameter<edm::ParameterSet>("TH1Res_d0");
0943 HistoName = "resd0_eta0to0p7";
0944 resd0_eta0to0p7 = iBooker.book1D(HistoName,
0945 HistoName,
0946 psRes_d0.getParameter<int32_t>("Nbinsx"),
0947 psRes_d0.getParameter<double>("xmin"),
0948 psRes_d0.getParameter<double>("xmax"));
0949 resd0_eta0to0p7->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
0950 resd0_eta0to0p7->setAxisTitle("# tracking particles", 2);
0951
0952
0953 HistoName = "resd0_eta0p7to1";
0954 resd0_eta0p7to1 = iBooker.book1D(HistoName,
0955 HistoName,
0956 psRes_d0.getParameter<int32_t>("Nbinsx"),
0957 psRes_d0.getParameter<double>("xmin"),
0958 psRes_d0.getParameter<double>("xmax"));
0959 resd0_eta0p7to1->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
0960 resd0_eta0p7to1->setAxisTitle("# tracking particles", 2);
0961
0962
0963 HistoName = "resd0_eta1to1p2";
0964 resd0_eta1to1p2 = iBooker.book1D(HistoName,
0965 HistoName,
0966 psRes_d0.getParameter<int32_t>("Nbinsx"),
0967 psRes_d0.getParameter<double>("xmin"),
0968 psRes_d0.getParameter<double>("xmax"));
0969 resd0_eta1to1p2->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
0970 resd0_eta1to1p2->setAxisTitle("# tracking particles", 2);
0971
0972
0973 HistoName = "resd0_eta1p2to1p6";
0974 resd0_eta1p2to1p6 = iBooker.book1D(HistoName,
0975 HistoName,
0976 psRes_d0.getParameter<int32_t>("Nbinsx"),
0977 psRes_d0.getParameter<double>("xmin"),
0978 psRes_d0.getParameter<double>("xmax"));
0979 resd0_eta1p2to1p6->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
0980 resd0_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
0981
0982
0983 HistoName = "resd0_eta1p6to2";
0984 resd0_eta1p6to2 = iBooker.book1D(HistoName,
0985 HistoName,
0986 psRes_d0.getParameter<int32_t>("Nbinsx"),
0987 psRes_d0.getParameter<double>("xmin"),
0988 psRes_d0.getParameter<double>("xmax"));
0989 resd0_eta1p6to2->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
0990 resd0_eta1p6to2->setAxisTitle("# tracking particles", 2);
0991
0992
0993 HistoName = "resd0_eta2to2p4";
0994 resd0_eta2to2p4 = iBooker.book1D(HistoName,
0995 HistoName,
0996 psRes_d0.getParameter<int32_t>("Nbinsx"),
0997 psRes_d0.getParameter<double>("xmin"),
0998 psRes_d0.getParameter<double>("xmax"));
0999 resd0_eta2to2p4->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1000 resd0_eta2to2p4->setAxisTitle("# tracking particles", 2);
1001
1002 }
1003
1004 DEFINE_FWK_MODULE(OuterTrackerMonitorTrackingParticles);