File indexing completed on 2025-02-20 03:46:00
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021 #include <memory>
0022 #include <numeric>
0023 #include <vector>
0024 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0025 #include "DQMServices/Core/interface/DQMStore.h"
0026 #include "DataFormats/Common/interface/DetSetVector.h"
0027 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0028 #include "DataFormats/Common/interface/Ptr.h"
0029 #include "DataFormats/Common/interface/Ref.h"
0030 #include "DataFormats/L1TrackTrigger/interface/TTCluster.h"
0031 #include "DataFormats/L1TrackTrigger/interface/TTStub.h"
0032 #include "DataFormats/L1TrackTrigger/interface/TTTrack.h"
0033 #include "DataFormats/L1TrackTrigger/interface/TTTypes.h"
0034 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0035 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0036
0037 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0038 #include "FWCore/ServiceRegistry/interface/Service.h"
0039 #include "FWCore/Framework/interface/ESHandle.h"
0040 #include "FWCore/Framework/interface/Event.h"
0041 #include "FWCore/Framework/interface/Frameworkfwd.h"
0042 #include "FWCore/Framework/interface/MakerMacros.h"
0043 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0044 #include "FWCore/Utilities/interface/EDGetToken.h"
0045 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0046 #include "Geometry/CommonTopologies/interface/PixelGeomDetUnit.h"
0047 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0048 #include "Geometry/Records/interface/StackedTrackerGeometryRecord.h"
0049 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0050 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
0051 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0052 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
0053 #include "SimDataFormats/Associations/interface/TTClusterAssociationMap.h"
0054 #include "SimDataFormats/Associations/interface/TTStubAssociationMap.h"
0055 #include "SimDataFormats/Associations/interface/TTTrackAssociationMap.h"
0056
0057 class Phase2OTValidateReconstruction : public DQMEDAnalyzer {
0058 public:
0059 explicit Phase2OTValidateReconstruction(const edm::ParameterSet &);
0060 ~Phase2OTValidateReconstruction() override;
0061 void analyze(const edm::Event &, const edm::EventSetup &) override;
0062 void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0063 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0064
0065 MonitorElement *trackParts_Eta = nullptr;
0066 MonitorElement *trackParts_Phi = nullptr;
0067 MonitorElement *trackParts_Pt = nullptr;
0068
0069
0070 MonitorElement *tp_pt = nullptr;
0071 MonitorElement *tp_pt_zoom = nullptr;
0072 MonitorElement *tp_eta = nullptr;
0073 MonitorElement *tp_d0 = nullptr;
0074 MonitorElement *tp_VtxR = nullptr;
0075 MonitorElement *tp_VtxZ = nullptr;
0076 MonitorElement *match_tp_pt = nullptr;
0077 MonitorElement *match_tp_pt_zoom = nullptr;
0078 MonitorElement *match_tp_eta = nullptr;
0079 MonitorElement *match_tp_d0 = nullptr;
0080 MonitorElement *match_tp_VtxR = nullptr;
0081 MonitorElement *match_tp_VtxZ = nullptr;
0082
0083
0084 MonitorElement *gen_clusters_barrel = nullptr;
0085 MonitorElement *gen_clusters_zoom_barrel = nullptr;
0086 MonitorElement *gen_clusters_endcaps = nullptr;
0087 MonitorElement *gen_clusters_zoom_endcaps = nullptr;
0088 MonitorElement *gen_clusters_if_stub_barrel = nullptr;
0089 MonitorElement *gen_clusters_if_stub_zoom_barrel = nullptr;
0090 MonitorElement *gen_clusters_if_stub_endcaps = nullptr;
0091 MonitorElement *gen_clusters_if_stub_zoom_endcaps = nullptr;
0092
0093
0094 MonitorElement *res_eta = nullptr;
0095 MonitorElement *res_pt = nullptr;
0096 MonitorElement *res_ptRel = nullptr;
0097 MonitorElement *respt_eta0to0p7_pt2to3 = nullptr;
0098 MonitorElement *respt_eta0p7to1_pt2to3 = nullptr;
0099 MonitorElement *respt_eta1to1p2_pt2to3 = nullptr;
0100 MonitorElement *respt_eta1p2to1p6_pt2to3 = nullptr;
0101 MonitorElement *respt_eta1p6to2_pt2to3 = nullptr;
0102 MonitorElement *respt_eta2to2p4_pt2to3 = nullptr;
0103 MonitorElement *respt_eta0to0p7_pt3to8 = nullptr;
0104 MonitorElement *respt_eta0p7to1_pt3to8 = nullptr;
0105 MonitorElement *respt_eta1to1p2_pt3to8 = nullptr;
0106 MonitorElement *respt_eta1p2to1p6_pt3to8 = nullptr;
0107 MonitorElement *respt_eta1p6to2_pt3to8 = nullptr;
0108 MonitorElement *respt_eta2to2p4_pt3to8 = nullptr;
0109 MonitorElement *respt_eta0to0p7_pt8toInf = nullptr;
0110 MonitorElement *respt_eta0p7to1_pt8toInf = nullptr;
0111 MonitorElement *respt_eta1to1p2_pt8toInf = nullptr;
0112 MonitorElement *respt_eta1p2to1p6_pt8toInf = nullptr;
0113 MonitorElement *respt_eta1p6to2_pt8toInf = nullptr;
0114 MonitorElement *respt_eta2to2p4_pt8toInf = nullptr;
0115 MonitorElement *reseta_eta0to0p7 = nullptr;
0116 MonitorElement *reseta_eta0p7to1 = nullptr;
0117 MonitorElement *reseta_eta1to1p2 = nullptr;
0118 MonitorElement *reseta_eta1p2to1p6 = nullptr;
0119 MonitorElement *reseta_eta1p6to2 = nullptr;
0120 MonitorElement *reseta_eta2to2p4 = nullptr;
0121 MonitorElement *resphi_eta0to0p7 = nullptr;
0122 MonitorElement *resphi_eta0p7to1 = nullptr;
0123 MonitorElement *resphi_eta1to1p2 = nullptr;
0124 MonitorElement *resphi_eta1p2to1p6 = nullptr;
0125 MonitorElement *resphi_eta1p6to2 = nullptr;
0126 MonitorElement *resphi_eta2to2p4 = nullptr;
0127 MonitorElement *resVtxZ_eta0to0p7 = nullptr;
0128 MonitorElement *resVtxZ_eta0p7to1 = nullptr;
0129 MonitorElement *resVtxZ_eta1to1p2 = nullptr;
0130 MonitorElement *resVtxZ_eta1p2to1p6 = nullptr;
0131 MonitorElement *resVtxZ_eta1p6to2 = nullptr;
0132 MonitorElement *resVtxZ_eta2to2p4 = nullptr;
0133
0134
0135 MonitorElement *resd0_eta0to0p7 = nullptr;
0136 MonitorElement *resd0_eta0p7to1 = nullptr;
0137 MonitorElement *resd0_eta1to1p2 = nullptr;
0138 MonitorElement *resd0_eta1p2to1p6 = nullptr;
0139 MonitorElement *resd0_eta1p6to2 = nullptr;
0140 MonitorElement *resd0_eta2to2p4 = nullptr;
0141
0142 private:
0143 const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> m_topoToken;
0144 edm::ParameterSet conf_;
0145 edm::EDGetTokenT<std::vector<TrackingParticle>> trackingParticleToken_;
0146 edm::EDGetTokenT<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>>
0147 ttClusterMCTruthToken_;
0148 edm::EDGetTokenT<TTStubAssociationMap<Ref_Phase2TrackerDigi_>>
0149 ttStubMCTruthToken_;
0150 edm::EDGetTokenT<TTTrackAssociationMap<Ref_Phase2TrackerDigi_>>
0151 ttTrackMCTruthToken_;
0152 edm::EDGetTokenT<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>> ttStubToken_;
0153 edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> getTokenTrackerGeom_;
0154 int L1Tk_minNStub;
0155 double L1Tk_maxChi2dof;
0156 int TP_minNStub;
0157 int TP_minNLayersStub;
0158 double TP_minPt;
0159 double TP_maxEta;
0160 double TP_maxVtxZ;
0161 std::string topFolderName_;
0162 };
0163
0164
0165
0166
0167 Phase2OTValidateReconstruction::Phase2OTValidateReconstruction(const edm::ParameterSet &iConfig)
0168 : m_topoToken(esConsumes()), conf_(iConfig) {
0169 topFolderName_ = conf_.getParameter<std::string>("TopFolderName");
0170 trackingParticleToken_ =
0171 consumes<std::vector<TrackingParticle>>(conf_.getParameter<edm::InputTag>("trackingParticleToken"));
0172 ttStubMCTruthToken_ =
0173 consumes<TTStubAssociationMap<Ref_Phase2TrackerDigi_>>(conf_.getParameter<edm::InputTag>("MCTruthStubInputTag"));
0174 ttClusterMCTruthToken_ = consumes<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>>(
0175 conf_.getParameter<edm::InputTag>("MCTruthClusterInputTag"));
0176 ttTrackMCTruthToken_ = consumes<TTTrackAssociationMap<Ref_Phase2TrackerDigi_>>(
0177 conf_.getParameter<edm::InputTag>("MCTruthTrackInputTag"));
0178 ttStubToken_ = consumes<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>>(
0179 edm::InputTag("TTStubsFromPhase2TrackerDigis", "StubAccepted"));
0180 getTokenTrackerGeom_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
0181 L1Tk_minNStub = conf_.getParameter<int>("L1Tk_minNStub");
0182 L1Tk_maxChi2dof = conf_.getParameter<double>("L1Tk_maxChi2dof");
0183 TP_minNStub = conf_.getParameter<int>("TP_minNStub");
0184
0185 TP_minNLayersStub = conf_.getParameter<int>("TP_minNLayersStub");
0186 TP_minPt = conf_.getParameter<double>("TP_minPt");
0187 TP_maxEta = conf_.getParameter<double>("TP_maxEta");
0188 TP_maxVtxZ = conf_.getParameter<double>("TP_maxVtxZ");
0189 }
0190
0191 Phase2OTValidateReconstruction::~Phase2OTValidateReconstruction() = default;
0192
0193
0194
0195
0196 void Phase2OTValidateReconstruction::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0197
0198 edm::Handle<std::vector<TrackingParticle>> trackingParticleHandle;
0199 iEvent.getByToken(trackingParticleToken_, trackingParticleHandle);
0200
0201
0202 edm::Handle<TTTrackAssociationMap<Ref_Phase2TrackerDigi_>> MCTruthTTTrackHandle;
0203 iEvent.getByToken(ttTrackMCTruthToken_, MCTruthTTTrackHandle);
0204 edm::Handle<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>> MCTruthTTClusterHandle;
0205 iEvent.getByToken(ttClusterMCTruthToken_, MCTruthTTClusterHandle);
0206 edm::Handle<TTStubAssociationMap<Ref_Phase2TrackerDigi_>> MCTruthTTStubHandle;
0207 iEvent.getByToken(ttStubMCTruthToken_, MCTruthTTStubHandle);
0208 edm::Handle<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>> TTStubHandle;
0209 iEvent.getByToken(ttStubToken_, TTStubHandle);
0210
0211
0212 const TrackerTopology *const tTopo = &iSetup.getData(m_topoToken);
0213 const TrackerGeometry *theTrackerGeom = &iSetup.getData(getTokenTrackerGeom_);
0214
0215
0216 int this_tp = 0;
0217 for (const auto &iterTP : *trackingParticleHandle) {
0218 edm::Ptr<TrackingParticle> tp_ptr(trackingParticleHandle, this_tp);
0219 this_tp++;
0220
0221
0222 float tmp_tp_pt = iterTP.pt();
0223 float tmp_tp_phi = iterTP.phi();
0224 float tmp_tp_eta = iterTP.eta();
0225
0226
0227 std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>>
0228 theStubRefs = MCTruthTTStubHandle->findTTStubRefs(tp_ptr);
0229
0230 int hasStubInLayer[11] = {0};
0231 for (unsigned int is = 0; is < theStubRefs.size(); is++) {
0232 DetId detid(theStubRefs.at(is)->getDetId());
0233 int layer = -1;
0234 if (detid.subdetId() == StripSubdetector::TOB)
0235 layer = static_cast<int>(tTopo->layer(detid)) - 1;
0236 else if (detid.subdetId() == StripSubdetector::TID)
0237 layer = static_cast<int>(tTopo->layer(detid)) + 5;
0238
0239
0240 if (MCTruthTTStubHandle->findTrackingParticlePtr(theStubRefs.at(is)).isNull() && hasStubInLayer[layer] < 2)
0241 hasStubInLayer[layer] = 1;
0242 else
0243 hasStubInLayer[layer] = 2;
0244 }
0245
0246 int nStubLayerTP = 0;
0247 for (int isum = 0; isum < 11; isum++) {
0248 if (hasStubInLayer[isum] >= 1)
0249 nStubLayerTP += 1;
0250 }
0251
0252 if (std::fabs(tmp_tp_eta) > TP_maxEta)
0253 continue;
0254
0255 if (tmp_tp_pt > TP_minPt && nStubLayerTP >= TP_minNLayersStub) {
0256 trackParts_Pt->Fill(tmp_tp_pt);
0257 trackParts_Eta->Fill(tmp_tp_eta);
0258 trackParts_Phi->Fill(tmp_tp_phi);
0259 }
0260
0261
0262
0263 int nStubTP = -1;
0264 if (MCTruthTTStubHandle.isValid()) {
0265 std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>>
0266 theStubRefs = MCTruthTTStubHandle->findTTStubRefs(tp_ptr);
0267 nStubTP = (int)theStubRefs.size();
0268 }
0269 if (MCTruthTTClusterHandle.isValid() && MCTruthTTClusterHandle->findTTClusterRefs(tp_ptr).empty())
0270 continue;
0271
0272 float tmp_tp_vz = iterTP.vz();
0273 float tmp_tp_vx = iterTP.vx();
0274 float tmp_tp_vy = iterTP.vy();
0275 float tmp_tp_charge = tp_ptr->charge();
0276 int tmp_tp_pdgid = iterTP.pdgId();
0277
0278
0279
0280
0281 float tmp_tp_t = tan(2.0 * atan(1.0) - 2.0 * atan(exp(-tmp_tp_eta)));
0282 float delx = -tmp_tp_vx;
0283 float dely = -tmp_tp_vy;
0284 float K = 0.01 * 0.5696 / tmp_tp_pt * tmp_tp_charge;
0285 float A = 1. / (2. * K);
0286 float tmp_tp_x0p = delx - A * sin(tmp_tp_phi);
0287 float tmp_tp_y0p = dely + A * cos(tmp_tp_phi);
0288 float tmp_tp_rp = sqrt(tmp_tp_x0p * tmp_tp_x0p + tmp_tp_y0p * tmp_tp_y0p);
0289 static double pi = 4.0 * atan(1.0);
0290 float delphi = tmp_tp_phi - atan2(-K * tmp_tp_x0p, K * tmp_tp_y0p);
0291 if (delphi < -pi)
0292 delphi += 2.0 * pi;
0293 if (delphi > pi)
0294 delphi -= 2.0 * pi;
0295
0296 float tmp_tp_VtxZ = tmp_tp_vz + tmp_tp_t * delphi / (2.0 * K);
0297 float tmp_tp_VtxR = sqrt(tmp_tp_vx * tmp_tp_vx + tmp_tp_vy * tmp_tp_vy);
0298 float tmp_tp_d0 = tmp_tp_charge * tmp_tp_rp - (1. / (2. * K));
0299
0300
0301
0302 float other_d0 = -tmp_tp_vx * sin(tmp_tp_phi) + tmp_tp_vy * cos(tmp_tp_phi);
0303 tmp_tp_d0 = tmp_tp_d0 * (-1);
0304 if (K == 0) {
0305 tmp_tp_d0 = other_d0;
0306 tmp_tp_VtxZ = tmp_tp_vz;
0307 }
0308 if (std::fabs(tmp_tp_VtxZ) > TP_maxVtxZ)
0309 continue;
0310
0311
0312 if (tmp_tp_VtxR < 1.0) {
0313 tp_pt->Fill(tmp_tp_pt);
0314 if (tmp_tp_pt <= 10)
0315 tp_pt_zoom->Fill(tmp_tp_pt);
0316 }
0317 if (tmp_tp_pt < TP_minPt)
0318 continue;
0319 tp_VtxR->Fill(tmp_tp_VtxR);
0320 if (tmp_tp_VtxR > 1.0)
0321 continue;
0322 tp_eta->Fill(tmp_tp_eta);
0323 tp_d0->Fill(tmp_tp_d0);
0324 tp_VtxZ->Fill(tmp_tp_VtxZ);
0325
0326 if (nStubTP < TP_minNStub || nStubLayerTP < TP_minNLayersStub)
0327 continue;
0328
0329
0330 std::vector<edm::Ref<edmNew::DetSetVector<TTCluster<Ref_Phase2TrackerDigi_>>, TTCluster<Ref_Phase2TrackerDigi_>>>
0331 associatedClusters = MCTruthTTClusterHandle->findTTClusterRefs(tp_ptr);
0332
0333
0334 for (std::size_t k = 0; k < associatedClusters.size(); ++k) {
0335 auto clusA = associatedClusters[k];
0336
0337
0338 DetId clusdetid = clusA->getDetId();
0339 if (clusdetid.subdetId() != StripSubdetector::TOB && clusdetid.subdetId() != StripSubdetector::TID)
0340 continue;
0341
0342 bool isGenuine = MCTruthTTClusterHandle->isGenuine(clusA);
0343 if (!isGenuine)
0344 continue;
0345
0346 DetId detidA = tTopo->stack(clusdetid);
0347 const GeomDetUnit *detA = theTrackerGeom->idToDetUnit(clusdetid);
0348 const PixelGeomDetUnit *theGeomDetA = dynamic_cast<const PixelGeomDetUnit *>(detA);
0349 const PixelTopology *topoA = dynamic_cast<const PixelTopology *>(&(theGeomDetA->specificTopology()));
0350 GlobalPoint coordsA =
0351 theGeomDetA->surface().toGlobal(topoA->localPosition(clusA->findAverageLocalCoordinatesCentered()));
0352
0353 int isBarrel = 0;
0354 if (clusdetid.subdetId() == StripSubdetector::TOB) {
0355 isBarrel = 1;
0356 } else if (clusdetid.subdetId() == StripSubdetector::TID) {
0357 isBarrel = 0;
0358 } else {
0359 edm::LogVerbatim("Tracklet") << "WARNING -- neither TOB or TID stub, shouldn't happen...";
0360 }
0361
0362 if (isBarrel == 1) {
0363 gen_clusters_barrel->Fill(tmp_tp_pt);
0364 gen_clusters_zoom_barrel->Fill(tmp_tp_pt);
0365 } else {
0366 gen_clusters_endcaps->Fill(tmp_tp_pt);
0367 gen_clusters_zoom_endcaps->Fill(tmp_tp_pt);
0368 }
0369
0370
0371 if (TTStubHandle->find(detidA) != TTStubHandle->end()) {
0372 edmNew::DetSet<TTStub<Ref_Phase2TrackerDigi_>> stubs = (*TTStubHandle)[detidA];
0373 for (auto stubIter = stubs.begin(); stubIter != stubs.end(); ++stubIter) {
0374 auto stubRef = edmNew::makeRefTo(TTStubHandle, stubIter);
0375
0376
0377 auto clusterRefB = stubIter->clusterRef(0);
0378 auto clusterRefC = stubIter->clusterRef(1);
0379
0380
0381 DetId detIdB = stubIter->clusterRef(0)->getDetId();
0382 DetId detIdC = stubIter->clusterRef(1)->getDetId();
0383
0384 const GeomDetUnit *detB = theTrackerGeom->idToDetUnit(detIdB);
0385 const GeomDetUnit *detC = theTrackerGeom->idToDetUnit(detIdC);
0386 const PixelGeomDetUnit *theGeomDetB = dynamic_cast<const PixelGeomDetUnit *>(detB);
0387 const PixelGeomDetUnit *theGeomDetC = dynamic_cast<const PixelGeomDetUnit *>(detC);
0388 const PixelTopology *topoB = dynamic_cast<const PixelTopology *>(&(theGeomDetB->specificTopology()));
0389 const PixelTopology *topoC = dynamic_cast<const PixelTopology *>(&(theGeomDetC->specificTopology()));
0390
0391 GlobalPoint coordsB = theGeomDetB->surface().toGlobal(
0392 topoB->localPosition(stubIter->clusterRef(0)->findAverageLocalCoordinatesCentered()));
0393 GlobalPoint coordsC = theGeomDetC->surface().toGlobal(
0394 topoC->localPosition(stubIter->clusterRef(1)->findAverageLocalCoordinatesCentered()));
0395
0396 if (coordsA.x() == coordsB.x() || coordsA.x() == coordsC.x()) {
0397 edm::Ptr<TrackingParticle> stubTP =
0398 MCTruthTTStubHandle->findTrackingParticlePtr(edmNew::makeRefTo(TTStubHandle, stubIter));
0399 if (stubTP.isNull())
0400 continue;
0401 float stub_tp_pt = stubTP->pt();
0402 if (stub_tp_pt == tmp_tp_pt) {
0403 if (isBarrel == 1) {
0404 gen_clusters_if_stub_barrel->Fill(tmp_tp_pt);
0405 gen_clusters_if_stub_zoom_barrel->Fill(tmp_tp_pt);
0406 } else {
0407 gen_clusters_if_stub_endcaps->Fill(tmp_tp_pt);
0408 gen_clusters_if_stub_zoom_endcaps->Fill(tmp_tp_pt);
0409 }
0410 break;
0411 }
0412 }
0413 }
0414 }
0415 }
0416
0417
0418
0419 int tp_nMatch = 0;
0420 int i_track = -1;
0421 float i_chi2dof = 99999;
0422 if (MCTruthTTTrackHandle.isValid()) {
0423 std::vector<edm::Ptr<TTTrack<Ref_Phase2TrackerDigi_>>> matchedTracks =
0424 MCTruthTTTrackHandle->findTTTrackPtrs(tp_ptr);
0425
0426
0427
0428
0429
0430
0431 int trkCounter = 0;
0432 for (const auto &thisTrack : matchedTracks) {
0433 if (!MCTruthTTTrackHandle->isGenuine(thisTrack))
0434 continue;
0435
0436
0437
0438
0439 int tmp_trk_nstub = thisTrack->getStubRefs().size();
0440 if (tmp_trk_nstub < L1Tk_minNStub)
0441 continue;
0442 float dmatch_pt = 999;
0443 float dmatch_eta = 999;
0444 float dmatch_phi = 999;
0445 int match_id = 999;
0446
0447 edm::Ptr<TrackingParticle> my_tp = MCTruthTTTrackHandle->findTrackingParticlePtr(thisTrack);
0448 dmatch_pt = std::fabs(my_tp->p4().pt() - tmp_tp_pt);
0449 dmatch_eta = std::fabs(my_tp->p4().eta() - tmp_tp_eta);
0450 dmatch_phi = std::fabs(my_tp->p4().phi() - tmp_tp_phi);
0451 match_id = my_tp->pdgId();
0452 float tmp_trk_chi2dof = thisTrack->chi2Red();
0453
0454
0455 if (dmatch_pt < 0.1 && dmatch_eta < 0.1 && dmatch_phi < 0.1 && tmp_tp_pdgid == match_id) {
0456 tp_nMatch++;
0457 if (i_track < 0 || tmp_trk_chi2dof < i_chi2dof) {
0458 i_track = trkCounter;
0459 i_chi2dof = tmp_trk_chi2dof;
0460 }
0461 }
0462 trkCounter++;
0463 }
0464
0465 if (tp_nMatch < 1)
0466 continue;
0467
0468 float tmp_matchtrk_pt = -999;
0469 float tmp_matchtrk_eta = -999;
0470 float tmp_matchtrk_phi = -999;
0471 float tmp_matchtrk_VtxZ = -999;
0472 float tmp_matchtrk_chi2dof = -999;
0473 int tmp_matchTrk_nStub = -999;
0474 float tmp_matchtrk_d0 = -999;
0475
0476 tmp_matchtrk_pt = matchedTracks[i_track]->momentum().perp();
0477 tmp_matchtrk_eta = matchedTracks[i_track]->momentum().eta();
0478 tmp_matchtrk_phi = matchedTracks[i_track]->momentum().phi();
0479 tmp_matchtrk_VtxZ = matchedTracks[i_track]->z0();
0480 tmp_matchtrk_chi2dof = matchedTracks[i_track]->chi2Red();
0481 tmp_matchTrk_nStub = (int)matchedTracks[i_track]->getStubRefs().size();
0482
0483
0484 float tmp_matchtrk_x0 = matchedTracks[i_track]->POCA().x();
0485 float tmp_matchtrk_y0 = matchedTracks[i_track]->POCA().y();
0486 tmp_matchtrk_d0 = -tmp_matchtrk_x0 * sin(tmp_matchtrk_phi) + tmp_matchtrk_y0 * cos(tmp_matchtrk_phi);
0487
0488
0489 if (tmp_matchTrk_nStub < L1Tk_minNStub || tmp_matchtrk_chi2dof > L1Tk_maxChi2dof)
0490 continue;
0491
0492
0493 match_tp_pt->Fill(tmp_tp_pt);
0494 if (tmp_tp_pt > 0 && tmp_tp_pt <= 10)
0495 match_tp_pt_zoom->Fill(tmp_tp_pt);
0496 match_tp_eta->Fill(tmp_tp_eta);
0497 match_tp_d0->Fill(tmp_tp_d0);
0498 match_tp_VtxR->Fill(tmp_tp_VtxR);
0499 match_tp_VtxZ->Fill(tmp_tp_VtxZ);
0500
0501
0502 float pt_diff = tmp_matchtrk_pt - tmp_tp_pt;
0503 float pt_res = pt_diff / tmp_tp_pt;
0504 float eta_res = tmp_matchtrk_eta - tmp_tp_eta;
0505 float phi_res = tmp_matchtrk_phi - tmp_tp_phi;
0506 float VtxZ_res = tmp_matchtrk_VtxZ - tmp_tp_VtxZ;
0507 float d0_res = tmp_matchtrk_d0 - tmp_tp_d0;
0508
0509
0510 res_pt->Fill(pt_diff);
0511 res_ptRel->Fill(pt_res);
0512 res_eta->Fill(eta_res);
0513
0514
0515
0516 if (std::fabs(tmp_tp_eta) >= 0 && std::fabs(tmp_tp_eta) < 0.7) {
0517 reseta_eta0to0p7->Fill(eta_res);
0518 resphi_eta0to0p7->Fill(phi_res);
0519 resVtxZ_eta0to0p7->Fill(VtxZ_res);
0520 resd0_eta0to0p7->Fill(d0_res);
0521 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0522 respt_eta0to0p7_pt2to3->Fill(pt_res);
0523 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0524 respt_eta0to0p7_pt3to8->Fill(pt_res);
0525 else if (tmp_tp_pt >= 8)
0526 respt_eta0to0p7_pt8toInf->Fill(pt_res);
0527 } else if (std::fabs(tmp_tp_eta) >= 0.7 && std::fabs(tmp_tp_eta) < 1.0) {
0528 reseta_eta0p7to1->Fill(eta_res);
0529 resphi_eta0p7to1->Fill(phi_res);
0530 resVtxZ_eta0p7to1->Fill(VtxZ_res);
0531 resd0_eta0p7to1->Fill(d0_res);
0532 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0533 respt_eta0p7to1_pt2to3->Fill(pt_res);
0534 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0535 respt_eta0p7to1_pt3to8->Fill(pt_res);
0536 else if (tmp_tp_pt >= 8)
0537 respt_eta0p7to1_pt8toInf->Fill(pt_res);
0538 } else if (std::fabs(tmp_tp_eta) >= 1.0 && std::fabs(tmp_tp_eta) < 1.2) {
0539 reseta_eta1to1p2->Fill(eta_res);
0540 resphi_eta1to1p2->Fill(phi_res);
0541 resVtxZ_eta1to1p2->Fill(VtxZ_res);
0542 resd0_eta1to1p2->Fill(d0_res);
0543 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0544 respt_eta1to1p2_pt2to3->Fill(pt_res);
0545 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0546 respt_eta1to1p2_pt3to8->Fill(pt_res);
0547 else if (tmp_tp_pt >= 8)
0548 respt_eta1to1p2_pt8toInf->Fill(pt_res);
0549 } else if (std::fabs(tmp_tp_eta) >= 1.2 && std::fabs(tmp_tp_eta) < 1.6) {
0550 reseta_eta1p2to1p6->Fill(eta_res);
0551 resphi_eta1p2to1p6->Fill(phi_res);
0552 resVtxZ_eta1p2to1p6->Fill(VtxZ_res);
0553 resd0_eta1p2to1p6->Fill(d0_res);
0554 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0555 respt_eta1p2to1p6_pt2to3->Fill(pt_res);
0556 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0557 respt_eta1p2to1p6_pt3to8->Fill(pt_res);
0558 else if (tmp_tp_pt >= 8)
0559 respt_eta1p2to1p6_pt8toInf->Fill(pt_res);
0560 } else if (std::fabs(tmp_tp_eta) >= 1.6 && std::fabs(tmp_tp_eta) < 2.0) {
0561 reseta_eta1p6to2->Fill(eta_res);
0562 resphi_eta1p6to2->Fill(phi_res);
0563 resVtxZ_eta1p6to2->Fill(VtxZ_res);
0564 resd0_eta1p6to2->Fill(d0_res);
0565 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0566 respt_eta1p6to2_pt2to3->Fill(pt_res);
0567 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0568 respt_eta1p6to2_pt3to8->Fill(pt_res);
0569 else if (tmp_tp_pt >= 8)
0570 respt_eta1p6to2_pt8toInf->Fill(pt_res);
0571 } else if (std::fabs(tmp_tp_eta) >= 2.0 && std::fabs(tmp_tp_eta) <= 2.4) {
0572 reseta_eta2to2p4->Fill(eta_res);
0573 resphi_eta2to2p4->Fill(phi_res);
0574 resVtxZ_eta2to2p4->Fill(VtxZ_res);
0575 resd0_eta2to2p4->Fill(d0_res);
0576 if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0577 respt_eta2to2p4_pt2to3->Fill(pt_res);
0578 else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0579 respt_eta2to2p4_pt3to8->Fill(pt_res);
0580 else if (tmp_tp_pt >= 8)
0581 respt_eta2to2p4_pt8toInf->Fill(pt_res);
0582 }
0583 }
0584 }
0585 }
0586
0587
0588
0589 void Phase2OTValidateReconstruction::bookHistograms(DQMStore::IBooker &iBooker,
0590 edm::Run const &run,
0591 edm::EventSetup const &es) {
0592
0593 std::string HistoName;
0594 iBooker.setCurrentFolder(topFolderName_ + "/trackParticles");
0595
0596
0597 edm::ParameterSet psTrackParts_Pt = conf_.getParameter<edm::ParameterSet>("TH1TrackParts_Pt");
0598 HistoName = "trackParts_Pt";
0599 trackParts_Pt = iBooker.book1D(HistoName,
0600 HistoName,
0601 psTrackParts_Pt.getParameter<int32_t>("Nbinsx"),
0602 psTrackParts_Pt.getParameter<double>("xmin"),
0603 psTrackParts_Pt.getParameter<double>("xmax"));
0604 trackParts_Pt->setAxisTitle("p_{T} [GeV]", 1);
0605 trackParts_Pt->setAxisTitle("# tracking particles", 2);
0606
0607
0608 edm::ParameterSet psTrackParts_Eta = conf_.getParameter<edm::ParameterSet>("TH1TrackParts_Eta");
0609 HistoName = "trackParts_Eta";
0610 trackParts_Eta = iBooker.book1D(HistoName,
0611 HistoName,
0612 psTrackParts_Eta.getParameter<int32_t>("Nbinsx"),
0613 psTrackParts_Eta.getParameter<double>("xmin"),
0614 psTrackParts_Eta.getParameter<double>("xmax"));
0615 trackParts_Eta->setAxisTitle("#eta", 1);
0616 trackParts_Eta->setAxisTitle("# tracking particles", 2);
0617
0618
0619 edm::ParameterSet psTrackParts_Phi = conf_.getParameter<edm::ParameterSet>("TH1TrackParts_Phi");
0620 HistoName = "trackParts_Phi";
0621 trackParts_Phi = iBooker.book1D(HistoName,
0622 HistoName,
0623 psTrackParts_Phi.getParameter<int32_t>("Nbinsx"),
0624 psTrackParts_Phi.getParameter<double>("xmin"),
0625 psTrackParts_Phi.getParameter<double>("xmax"));
0626 trackParts_Phi->setAxisTitle("#phi", 1);
0627 trackParts_Phi->setAxisTitle("# tracking particles", 2);
0628
0629
0630 iBooker.setCurrentFolder(topFolderName_ + "/EfficiencyIngredients");
0631
0632 edm::ParameterSet psEffic_pt = conf_.getParameter<edm::ParameterSet>("TH1Effic_pt");
0633 HistoName = "tp_pt";
0634 tp_pt = iBooker.book1D(HistoName,
0635 HistoName,
0636 psEffic_pt.getParameter<int32_t>("Nbinsx"),
0637 psEffic_pt.getParameter<double>("xmin"),
0638 psEffic_pt.getParameter<double>("xmax"));
0639 tp_pt->setAxisTitle("p_{T} [GeV]", 1);
0640 tp_pt->setAxisTitle("# tracking particles", 2);
0641
0642
0643 HistoName = "match_tp_pt";
0644 match_tp_pt = iBooker.book1D(HistoName,
0645 HistoName,
0646 psEffic_pt.getParameter<int32_t>("Nbinsx"),
0647 psEffic_pt.getParameter<double>("xmin"),
0648 psEffic_pt.getParameter<double>("xmax"));
0649 match_tp_pt->setAxisTitle("p_{T} [GeV]", 1);
0650 match_tp_pt->setAxisTitle("# matched tracking particles", 2);
0651
0652
0653 HistoName = "gen_clusters_barrel";
0654 gen_clusters_barrel = iBooker.book1D(HistoName,
0655 HistoName,
0656 psEffic_pt.getParameter<int32_t>("Nbinsx"),
0657 psEffic_pt.getParameter<double>("xmin"),
0658 psEffic_pt.getParameter<double>("xmax"));
0659 gen_clusters_barrel->setAxisTitle("p_{T} [GeV]", 1);
0660 gen_clusters_barrel->setAxisTitle("# tracking particles", 2);
0661
0662
0663 HistoName = "gen_clusters_if_stub_barrel";
0664 gen_clusters_if_stub_barrel = iBooker.book1D(HistoName,
0665 HistoName,
0666 psEffic_pt.getParameter<int32_t>("Nbinsx"),
0667 psEffic_pt.getParameter<double>("xmin"),
0668 psEffic_pt.getParameter<double>("xmax"));
0669 gen_clusters_if_stub_barrel->setAxisTitle("p_{T} [GeV]", 1);
0670 gen_clusters_if_stub_barrel->setAxisTitle("# tracking particles", 2);
0671
0672
0673 HistoName = "gen_clusters_endcaps";
0674 gen_clusters_endcaps = iBooker.book1D(HistoName,
0675 HistoName,
0676 psEffic_pt.getParameter<int32_t>("Nbinsx"),
0677 psEffic_pt.getParameter<double>("xmin"),
0678 psEffic_pt.getParameter<double>("xmax"));
0679 gen_clusters_endcaps->setAxisTitle("p_{T} [GeV]", 1);
0680 gen_clusters_endcaps->setAxisTitle("# tracking particles", 2);
0681
0682
0683 HistoName = "gen_clusters_if_stub_endcaps";
0684 gen_clusters_if_stub_endcaps = iBooker.book1D(HistoName,
0685 HistoName,
0686 psEffic_pt.getParameter<int32_t>("Nbinsx"),
0687 psEffic_pt.getParameter<double>("xmin"),
0688 psEffic_pt.getParameter<double>("xmax"));
0689 gen_clusters_if_stub_endcaps->setAxisTitle("p_{T} [GeV]", 1);
0690 gen_clusters_if_stub_endcaps->setAxisTitle("# tracking particles", 2);
0691
0692
0693 edm::ParameterSet psEffic_pt_zoom = conf_.getParameter<edm::ParameterSet>("TH1Effic_pt_zoom");
0694 HistoName = "tp_pt_zoom";
0695 tp_pt_zoom = iBooker.book1D(HistoName,
0696 HistoName,
0697 psEffic_pt_zoom.getParameter<int32_t>("Nbinsx"),
0698 psEffic_pt_zoom.getParameter<double>("xmin"),
0699 psEffic_pt_zoom.getParameter<double>("xmax"));
0700 tp_pt_zoom->setAxisTitle("p_{T} [GeV]", 1);
0701 tp_pt_zoom->setAxisTitle("# tracking particles", 2);
0702
0703
0704 HistoName = "match_tp_pt_zoom";
0705 match_tp_pt_zoom = iBooker.book1D(HistoName,
0706 HistoName,
0707 psEffic_pt_zoom.getParameter<int32_t>("Nbinsx"),
0708 psEffic_pt_zoom.getParameter<double>("xmin"),
0709 psEffic_pt_zoom.getParameter<double>("xmax"));
0710 match_tp_pt_zoom->setAxisTitle("p_{T} [GeV]", 1);
0711 match_tp_pt_zoom->setAxisTitle("# matched tracking particles", 2);
0712
0713
0714 HistoName = "gen_clusters_zoom_barrel";
0715 gen_clusters_zoom_barrel = iBooker.book1D(HistoName,
0716 HistoName,
0717 psEffic_pt_zoom.getParameter<int32_t>("Nbinsx"),
0718 psEffic_pt_zoom.getParameter<double>("xmin"),
0719 psEffic_pt_zoom.getParameter<double>("xmax"));
0720 gen_clusters_zoom_barrel->setAxisTitle("p_{T} [GeV]", 1);
0721 gen_clusters_zoom_barrel->setAxisTitle("# tracking particles", 2);
0722
0723
0724 HistoName = "gen_clusters_if_stub_zoom_barrel";
0725 gen_clusters_if_stub_zoom_barrel = iBooker.book1D(HistoName,
0726 HistoName,
0727 psEffic_pt_zoom.getParameter<int32_t>("Nbinsx"),
0728 psEffic_pt_zoom.getParameter<double>("xmin"),
0729 psEffic_pt_zoom.getParameter<double>("xmax"));
0730 gen_clusters_if_stub_zoom_barrel->setAxisTitle("p_{T} [GeV]", 1);
0731 gen_clusters_if_stub_zoom_barrel->setAxisTitle("# tracking particles", 2);
0732
0733
0734 HistoName = "gen_clusters_zoom_endcaps";
0735 gen_clusters_zoom_endcaps = iBooker.book1D(HistoName,
0736 HistoName,
0737 psEffic_pt_zoom.getParameter<int32_t>("Nbinsx"),
0738 psEffic_pt_zoom.getParameter<double>("xmin"),
0739 psEffic_pt_zoom.getParameter<double>("xmax"));
0740 gen_clusters_zoom_endcaps->setAxisTitle("p_{T} [GeV]", 1);
0741 gen_clusters_zoom_endcaps->setAxisTitle("# tracking particles", 2);
0742
0743
0744 HistoName = "gen_clusters_if_stub_zoom_endcaps";
0745 gen_clusters_if_stub_zoom_endcaps = iBooker.book1D(HistoName,
0746 HistoName,
0747 psEffic_pt_zoom.getParameter<int32_t>("Nbinsx"),
0748 psEffic_pt_zoom.getParameter<double>("xmin"),
0749 psEffic_pt_zoom.getParameter<double>("xmax"));
0750 gen_clusters_if_stub_zoom_endcaps->setAxisTitle("p_{T} [GeV]", 1);
0751 gen_clusters_if_stub_zoom_endcaps->setAxisTitle("# tracking particles", 2);
0752
0753
0754 edm::ParameterSet psEffic_eta = conf_.getParameter<edm::ParameterSet>("TH1Effic_eta");
0755 HistoName = "tp_eta";
0756 tp_eta = iBooker.book1D(HistoName,
0757 HistoName,
0758 psEffic_eta.getParameter<int32_t>("Nbinsx"),
0759 psEffic_eta.getParameter<double>("xmin"),
0760 psEffic_eta.getParameter<double>("xmax"));
0761 tp_eta->setAxisTitle("#eta", 1);
0762 tp_eta->setAxisTitle("# tracking particles", 2);
0763
0764
0765 HistoName = "match_tp_eta";
0766 match_tp_eta = iBooker.book1D(HistoName,
0767 HistoName,
0768 psEffic_eta.getParameter<int32_t>("Nbinsx"),
0769 psEffic_eta.getParameter<double>("xmin"),
0770 psEffic_eta.getParameter<double>("xmax"));
0771 match_tp_eta->setAxisTitle("#eta", 1);
0772 match_tp_eta->setAxisTitle("# matched tracking particles", 2);
0773
0774
0775 edm::ParameterSet psEffic_d0 = conf_.getParameter<edm::ParameterSet>("TH1Effic_d0");
0776 HistoName = "tp_d0";
0777 tp_d0 = iBooker.book1D(HistoName,
0778 HistoName,
0779 psEffic_d0.getParameter<int32_t>("Nbinsx"),
0780 psEffic_d0.getParameter<double>("xmin"),
0781 psEffic_d0.getParameter<double>("xmax"));
0782 tp_d0->setAxisTitle("d_{0} [cm]", 1);
0783 tp_d0->setAxisTitle("# tracking particles", 2);
0784
0785
0786 HistoName = "match_tp_d0";
0787 match_tp_d0 = iBooker.book1D(HistoName,
0788 HistoName,
0789 psEffic_d0.getParameter<int32_t>("Nbinsx"),
0790 psEffic_d0.getParameter<double>("xmin"),
0791 psEffic_d0.getParameter<double>("xmax"));
0792 match_tp_d0->setAxisTitle("d_{0} [cm]", 1);
0793 match_tp_d0->setAxisTitle("# matched tracking particles", 2);
0794
0795
0796 edm::ParameterSet psEffic_VtxR = conf_.getParameter<edm::ParameterSet>("TH1Effic_VtxR");
0797 HistoName = "tp_VtxR";
0798 tp_VtxR = iBooker.book1D(HistoName,
0799 HistoName,
0800 psEffic_VtxR.getParameter<int32_t>("Nbinsx"),
0801 psEffic_VtxR.getParameter<double>("xmin"),
0802 psEffic_VtxR.getParameter<double>("xmax"));
0803 tp_VtxR->setAxisTitle("d_{xy} [cm]", 1);
0804 tp_VtxR->setAxisTitle("# tracking particles", 2);
0805
0806
0807 HistoName = "match_tp_VtxR";
0808 match_tp_VtxR = iBooker.book1D(HistoName,
0809 HistoName,
0810 psEffic_VtxR.getParameter<int32_t>("Nbinsx"),
0811 psEffic_VtxR.getParameter<double>("xmin"),
0812 psEffic_VtxR.getParameter<double>("xmax"));
0813 match_tp_VtxR->setAxisTitle("d_{xy} [cm]", 1);
0814 match_tp_VtxR->setAxisTitle("# matched tracking particles", 2);
0815
0816
0817 edm::ParameterSet psEffic_VtxZ = conf_.getParameter<edm::ParameterSet>("TH1Effic_VtxZ");
0818 HistoName = "tp_VtxZ";
0819 tp_VtxZ = iBooker.book1D(HistoName,
0820 HistoName,
0821 psEffic_VtxZ.getParameter<int32_t>("Nbinsx"),
0822 psEffic_VtxZ.getParameter<double>("xmin"),
0823 psEffic_VtxZ.getParameter<double>("xmax"));
0824 tp_VtxZ->setAxisTitle("z_{0} [cm]", 1);
0825 tp_VtxZ->setAxisTitle("# tracking particles", 2);
0826
0827
0828 HistoName = "match_tp_VtxZ";
0829 match_tp_VtxZ = iBooker.book1D(HistoName,
0830 HistoName,
0831 psEffic_VtxZ.getParameter<int32_t>("Nbinsx"),
0832 psEffic_VtxZ.getParameter<double>("xmin"),
0833 psEffic_VtxZ.getParameter<double>("xmax"));
0834 match_tp_VtxZ->setAxisTitle("z_{0} [cm]", 1);
0835 match_tp_VtxZ->setAxisTitle("# matched tracking particles", 2);
0836
0837
0838 iBooker.setCurrentFolder(topFolderName_ + "/ResolutionIngredients");
0839
0840 edm::ParameterSet psRes_pt = conf_.getParameter<edm::ParameterSet>("TH1Res_pt");
0841 HistoName = "res_pt";
0842 res_pt = iBooker.book1D(HistoName,
0843 HistoName,
0844 psRes_pt.getParameter<int32_t>("Nbinsx"),
0845 psRes_pt.getParameter<double>("xmin"),
0846 psRes_pt.getParameter<double>("xmax"));
0847 res_pt->setAxisTitle("p_{T} [GeV]", 1);
0848 res_pt->setAxisTitle("# tracking particles", 2);
0849
0850
0851 edm::ParameterSet psRes_eta = conf_.getParameter<edm::ParameterSet>("TH1Res_eta");
0852 HistoName = "res_eta";
0853 res_eta = iBooker.book1D(HistoName,
0854 HistoName,
0855 psRes_eta.getParameter<int32_t>("Nbinsx"),
0856 psRes_eta.getParameter<double>("xmin"),
0857 psRes_eta.getParameter<double>("xmax"));
0858 res_eta->setAxisTitle("#eta", 1);
0859 res_eta->setAxisTitle("# tracking particles", 2);
0860
0861
0862 edm::ParameterSet psRes_ptRel = conf_.getParameter<edm::ParameterSet>("TH1Res_ptRel");
0863 HistoName = "res_ptRel";
0864 res_ptRel = iBooker.book1D(HistoName,
0865 HistoName,
0866 psRes_ptRel.getParameter<int32_t>("Nbinsx"),
0867 psRes_ptRel.getParameter<double>("xmin"),
0868 psRes_ptRel.getParameter<double>("xmax"));
0869 res_ptRel->setAxisTitle("Relative p_{T} [GeV]", 1);
0870 res_ptRel->setAxisTitle("# tracking particles", 2);
0871
0872
0873
0874 HistoName = "reseta_eta0to0p7";
0875 reseta_eta0to0p7 = iBooker.book1D(HistoName,
0876 HistoName,
0877 psRes_eta.getParameter<int32_t>("Nbinsx"),
0878 psRes_eta.getParameter<double>("xmin"),
0879 psRes_eta.getParameter<double>("xmax"));
0880 reseta_eta0to0p7->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0881 reseta_eta0to0p7->setAxisTitle("# tracking particles", 2);
0882
0883
0884 HistoName = "reseta_eta0p7to1";
0885 reseta_eta0p7to1 = iBooker.book1D(HistoName,
0886 HistoName,
0887 psRes_eta.getParameter<int32_t>("Nbinsx"),
0888 psRes_eta.getParameter<double>("xmin"),
0889 psRes_eta.getParameter<double>("xmax"));
0890 reseta_eta0p7to1->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0891 reseta_eta0p7to1->setAxisTitle("# tracking particles", 2);
0892
0893
0894 HistoName = "reseta_eta1to1p2";
0895 reseta_eta1to1p2 = iBooker.book1D(HistoName,
0896 HistoName,
0897 psRes_eta.getParameter<int32_t>("Nbinsx"),
0898 psRes_eta.getParameter<double>("xmin"),
0899 psRes_eta.getParameter<double>("xmax"));
0900 reseta_eta1to1p2->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0901 reseta_eta1to1p2->setAxisTitle("# tracking particles", 2);
0902
0903
0904 HistoName = "reseta_eta1p2to1p6";
0905 reseta_eta1p2to1p6 = iBooker.book1D(HistoName,
0906 HistoName,
0907 psRes_eta.getParameter<int32_t>("Nbinsx"),
0908 psRes_eta.getParameter<double>("xmin"),
0909 psRes_eta.getParameter<double>("xmax"));
0910 reseta_eta1p2to1p6->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0911 reseta_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
0912
0913
0914 HistoName = "reseta_eta1p6to2";
0915 reseta_eta1p6to2 = iBooker.book1D(HistoName,
0916 HistoName,
0917 psRes_eta.getParameter<int32_t>("Nbinsx"),
0918 psRes_eta.getParameter<double>("xmin"),
0919 psRes_eta.getParameter<double>("xmax"));
0920 reseta_eta1p6to2->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0921 reseta_eta1p6to2->setAxisTitle("# tracking particles", 2);
0922
0923
0924 HistoName = "reseta_eta2to2p4";
0925 reseta_eta2to2p4 = iBooker.book1D(HistoName,
0926 HistoName,
0927 psRes_eta.getParameter<int32_t>("Nbinsx"),
0928 psRes_eta.getParameter<double>("xmin"),
0929 psRes_eta.getParameter<double>("xmax"));
0930 reseta_eta2to2p4->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0931 reseta_eta2to2p4->setAxisTitle("# tracking particles", 2);
0932
0933
0934
0935
0936 HistoName = "respt_eta0to0p7_pt2to3";
0937 respt_eta0to0p7_pt2to3 = iBooker.book1D(HistoName,
0938 HistoName,
0939 psRes_pt.getParameter<int32_t>("Nbinsx"),
0940 psRes_pt.getParameter<double>("xmin"),
0941 psRes_pt.getParameter<double>("xmax"));
0942 respt_eta0to0p7_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0943 respt_eta0to0p7_pt2to3->setAxisTitle("# tracking particles", 2);
0944
0945
0946 HistoName = "respt_eta0p7to1_pt2to3";
0947 respt_eta0p7to1_pt2to3 = iBooker.book1D(HistoName,
0948 HistoName,
0949 psRes_pt.getParameter<int32_t>("Nbinsx"),
0950 psRes_pt.getParameter<double>("xmin"),
0951 psRes_pt.getParameter<double>("xmax"));
0952 respt_eta0p7to1_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0953 respt_eta0p7to1_pt2to3->setAxisTitle("# tracking particles", 2);
0954
0955
0956 HistoName = "respt_eta1to1p2_pt2to3";
0957 respt_eta1to1p2_pt2to3 = iBooker.book1D(HistoName,
0958 HistoName,
0959 psRes_pt.getParameter<int32_t>("Nbinsx"),
0960 psRes_pt.getParameter<double>("xmin"),
0961 psRes_pt.getParameter<double>("xmax"));
0962 respt_eta1to1p2_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0963 respt_eta1to1p2_pt2to3->setAxisTitle("# tracking particles", 2);
0964
0965
0966 HistoName = "respt_eta1p2to1p6_pt2to3";
0967 respt_eta1p2to1p6_pt2to3 = iBooker.book1D(HistoName,
0968 HistoName,
0969 psRes_pt.getParameter<int32_t>("Nbinsx"),
0970 psRes_pt.getParameter<double>("xmin"),
0971 psRes_pt.getParameter<double>("xmax"));
0972 respt_eta1p2to1p6_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0973 respt_eta1p2to1p6_pt2to3->setAxisTitle("# tracking particles", 2);
0974
0975
0976 HistoName = "respt_eta1p6to2_pt2to3";
0977 respt_eta1p6to2_pt2to3 = iBooker.book1D(HistoName,
0978 HistoName,
0979 psRes_pt.getParameter<int32_t>("Nbinsx"),
0980 psRes_pt.getParameter<double>("xmin"),
0981 psRes_pt.getParameter<double>("xmax"));
0982 respt_eta1p6to2_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0983 respt_eta1p6to2_pt2to3->setAxisTitle("# tracking particles", 2);
0984
0985
0986 HistoName = "respt_eta2to2p4_pt2to3";
0987 respt_eta2to2p4_pt2to3 = iBooker.book1D(HistoName,
0988 HistoName,
0989 psRes_pt.getParameter<int32_t>("Nbinsx"),
0990 psRes_pt.getParameter<double>("xmin"),
0991 psRes_pt.getParameter<double>("xmax"));
0992 respt_eta2to2p4_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0993 respt_eta2to2p4_pt2to3->setAxisTitle("# tracking particles", 2);
0994
0995
0996
0997 HistoName = "respt_eta0to0p7_pt3to8";
0998 respt_eta0to0p7_pt3to8 = iBooker.book1D(HistoName,
0999 HistoName,
1000 psRes_pt.getParameter<int32_t>("Nbinsx"),
1001 psRes_pt.getParameter<double>("xmin"),
1002 psRes_pt.getParameter<double>("xmax"));
1003 respt_eta0to0p7_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1004 respt_eta0to0p7_pt3to8->setAxisTitle("# tracking particles", 2);
1005
1006
1007 HistoName = "respt_eta0p7to1_pt3to8";
1008 respt_eta0p7to1_pt3to8 = iBooker.book1D(HistoName,
1009 HistoName,
1010 psRes_pt.getParameter<int32_t>("Nbinsx"),
1011 psRes_pt.getParameter<double>("xmin"),
1012 psRes_pt.getParameter<double>("xmax"));
1013 respt_eta0p7to1_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1014 respt_eta0p7to1_pt3to8->setAxisTitle("# tracking particles", 2);
1015
1016
1017 HistoName = "respt_eta1to1p2_pt3to8";
1018 respt_eta1to1p2_pt3to8 = iBooker.book1D(HistoName,
1019 HistoName,
1020 psRes_pt.getParameter<int32_t>("Nbinsx"),
1021 psRes_pt.getParameter<double>("xmin"),
1022 psRes_pt.getParameter<double>("xmax"));
1023 respt_eta1to1p2_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1024 respt_eta1to1p2_pt3to8->setAxisTitle("# tracking particles", 2);
1025
1026
1027 HistoName = "respt_eta1p2to1p6_pt3to8";
1028 respt_eta1p2to1p6_pt3to8 = iBooker.book1D(HistoName,
1029 HistoName,
1030 psRes_pt.getParameter<int32_t>("Nbinsx"),
1031 psRes_pt.getParameter<double>("xmin"),
1032 psRes_pt.getParameter<double>("xmax"));
1033 respt_eta1p2to1p6_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1034 respt_eta1p2to1p6_pt3to8->setAxisTitle("# tracking particles", 2);
1035
1036
1037 HistoName = "respt_eta1p6to2_pt3to8";
1038 respt_eta1p6to2_pt3to8 = iBooker.book1D(HistoName,
1039 HistoName,
1040 psRes_pt.getParameter<int32_t>("Nbinsx"),
1041 psRes_pt.getParameter<double>("xmin"),
1042 psRes_pt.getParameter<double>("xmax"));
1043 respt_eta1p6to2_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1044 respt_eta1p6to2_pt3to8->setAxisTitle("# tracking particles", 2);
1045
1046
1047 HistoName = "respt_eta2to2p4_pt3to8";
1048 respt_eta2to2p4_pt3to8 = iBooker.book1D(HistoName,
1049 HistoName,
1050 psRes_pt.getParameter<int32_t>("Nbinsx"),
1051 psRes_pt.getParameter<double>("xmin"),
1052 psRes_pt.getParameter<double>("xmax"));
1053 respt_eta2to2p4_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1054 respt_eta2to2p4_pt3to8->setAxisTitle("# tracking particles", 2);
1055
1056
1057
1058 HistoName = "respt_eta0to0p7_pt8toInf";
1059 respt_eta0to0p7_pt8toInf = iBooker.book1D(HistoName,
1060 HistoName,
1061 psRes_pt.getParameter<int32_t>("Nbinsx"),
1062 psRes_pt.getParameter<double>("xmin"),
1063 psRes_pt.getParameter<double>("xmax"));
1064 respt_eta0to0p7_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1065 respt_eta0to0p7_pt8toInf->setAxisTitle("# tracking particles", 2);
1066
1067
1068 HistoName = "respt_eta0p7to1_pt8toInf";
1069 respt_eta0p7to1_pt8toInf = iBooker.book1D(HistoName,
1070 HistoName,
1071 psRes_pt.getParameter<int32_t>("Nbinsx"),
1072 psRes_pt.getParameter<double>("xmin"),
1073 psRes_pt.getParameter<double>("xmax"));
1074 respt_eta0p7to1_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1075 respt_eta0p7to1_pt8toInf->setAxisTitle("# tracking particles", 2);
1076
1077
1078 HistoName = "respt_eta1to1p2_pt8toInf";
1079 respt_eta1to1p2_pt8toInf = iBooker.book1D(HistoName,
1080 HistoName,
1081 psRes_pt.getParameter<int32_t>("Nbinsx"),
1082 psRes_pt.getParameter<double>("xmin"),
1083 psRes_pt.getParameter<double>("xmax"));
1084 respt_eta1to1p2_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1085 respt_eta1to1p2_pt8toInf->setAxisTitle("# tracking particles", 2);
1086
1087
1088 HistoName = "respt_eta1p2to1p6_pt8toInf";
1089 respt_eta1p2to1p6_pt8toInf = iBooker.book1D(HistoName,
1090 HistoName,
1091 psRes_pt.getParameter<int32_t>("Nbinsx"),
1092 psRes_pt.getParameter<double>("xmin"),
1093 psRes_pt.getParameter<double>("xmax"));
1094 respt_eta1p2to1p6_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1095 respt_eta1p2to1p6_pt8toInf->setAxisTitle("# tracking particles", 2);
1096
1097
1098 HistoName = "respt_eta1p6to2_pt8toInf";
1099 respt_eta1p6to2_pt8toInf = iBooker.book1D(HistoName,
1100 HistoName,
1101 psRes_pt.getParameter<int32_t>("Nbinsx"),
1102 psRes_pt.getParameter<double>("xmin"),
1103 psRes_pt.getParameter<double>("xmax"));
1104 respt_eta1p6to2_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1105 respt_eta1p6to2_pt8toInf->setAxisTitle("# tracking particles", 2);
1106
1107
1108 HistoName = "respt_eta2to2p4_pt8toInf";
1109 respt_eta2to2p4_pt8toInf = iBooker.book1D(HistoName,
1110 HistoName,
1111 psRes_pt.getParameter<int32_t>("Nbinsx"),
1112 psRes_pt.getParameter<double>("xmin"),
1113 psRes_pt.getParameter<double>("xmax"));
1114 respt_eta2to2p4_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
1115 respt_eta2to2p4_pt8toInf->setAxisTitle("# tracking particles", 2);
1116
1117
1118
1119 edm::ParameterSet psRes_phi = conf_.getParameter<edm::ParameterSet>("TH1Res_phi");
1120 HistoName = "resphi_eta0to0p7";
1121 resphi_eta0to0p7 = iBooker.book1D(HistoName,
1122 HistoName,
1123 psRes_phi.getParameter<int32_t>("Nbinsx"),
1124 psRes_phi.getParameter<double>("xmin"),
1125 psRes_phi.getParameter<double>("xmax"));
1126 resphi_eta0to0p7->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
1127 resphi_eta0to0p7->setAxisTitle("# tracking particles", 2);
1128
1129
1130 HistoName = "resphi_eta0p7to1";
1131 resphi_eta0p7to1 = iBooker.book1D(HistoName,
1132 HistoName,
1133 psRes_phi.getParameter<int32_t>("Nbinsx"),
1134 psRes_phi.getParameter<double>("xmin"),
1135 psRes_phi.getParameter<double>("xmax"));
1136 resphi_eta0p7to1->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
1137 resphi_eta0p7to1->setAxisTitle("# tracking particles", 2);
1138
1139
1140 HistoName = "resphi_eta1to1p2";
1141 resphi_eta1to1p2 = iBooker.book1D(HistoName,
1142 HistoName,
1143 psRes_phi.getParameter<int32_t>("Nbinsx"),
1144 psRes_phi.getParameter<double>("xmin"),
1145 psRes_phi.getParameter<double>("xmax"));
1146 resphi_eta1to1p2->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
1147 resphi_eta1to1p2->setAxisTitle("# tracking particles", 2);
1148
1149
1150 HistoName = "resphi_eta1p2to1p6";
1151 resphi_eta1p2to1p6 = iBooker.book1D(HistoName,
1152 HistoName,
1153 psRes_phi.getParameter<int32_t>("Nbinsx"),
1154 psRes_phi.getParameter<double>("xmin"),
1155 psRes_phi.getParameter<double>("xmax"));
1156 resphi_eta1p2to1p6->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
1157 resphi_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
1158
1159
1160 HistoName = "resphi_eta1p6to2";
1161 resphi_eta1p6to2 = iBooker.book1D(HistoName,
1162 HistoName,
1163 psRes_phi.getParameter<int32_t>("Nbinsx"),
1164 psRes_phi.getParameter<double>("xmin"),
1165 psRes_phi.getParameter<double>("xmax"));
1166 resphi_eta1p6to2->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
1167 resphi_eta1p6to2->setAxisTitle("# tracking particles", 2);
1168
1169
1170 HistoName = "resphi_eta2to2p4";
1171 resphi_eta2to2p4 = iBooker.book1D(HistoName,
1172 HistoName,
1173 psRes_phi.getParameter<int32_t>("Nbinsx"),
1174 psRes_phi.getParameter<double>("xmin"),
1175 psRes_phi.getParameter<double>("xmax"));
1176 resphi_eta2to2p4->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
1177 resphi_eta2to2p4->setAxisTitle("# tracking particles", 2);
1178
1179
1180
1181 edm::ParameterSet psRes_VtxZ = conf_.getParameter<edm::ParameterSet>("TH1Res_VtxZ");
1182 HistoName = "resVtxZ_eta0to0p7";
1183 resVtxZ_eta0to0p7 = iBooker.book1D(HistoName,
1184 HistoName,
1185 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
1186 psRes_VtxZ.getParameter<double>("xmin"),
1187 psRes_VtxZ.getParameter<double>("xmax"));
1188 resVtxZ_eta0to0p7->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
1189 resVtxZ_eta0to0p7->setAxisTitle("# tracking particles", 2);
1190
1191
1192 HistoName = "resVtxZ_eta0p7to1";
1193 resVtxZ_eta0p7to1 = iBooker.book1D(HistoName,
1194 HistoName,
1195 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
1196 psRes_VtxZ.getParameter<double>("xmin"),
1197 psRes_VtxZ.getParameter<double>("xmax"));
1198 resVtxZ_eta0p7to1->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
1199 resVtxZ_eta0p7to1->setAxisTitle("# tracking particles", 2);
1200
1201
1202 HistoName = "resVtxZ_eta1to1p2";
1203 resVtxZ_eta1to1p2 = iBooker.book1D(HistoName,
1204 HistoName,
1205 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
1206 psRes_VtxZ.getParameter<double>("xmin"),
1207 psRes_VtxZ.getParameter<double>("xmax"));
1208 resVtxZ_eta1to1p2->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
1209 resVtxZ_eta1to1p2->setAxisTitle("# tracking particles", 2);
1210
1211
1212 HistoName = "resVtxZ_eta1p2to1p6";
1213 resVtxZ_eta1p2to1p6 = iBooker.book1D(HistoName,
1214 HistoName,
1215 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
1216 psRes_VtxZ.getParameter<double>("xmin"),
1217 psRes_VtxZ.getParameter<double>("xmax"));
1218 resVtxZ_eta1p2to1p6->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
1219 resVtxZ_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
1220
1221
1222 HistoName = "resVtxZ_eta1p6to2";
1223 resVtxZ_eta1p6to2 = iBooker.book1D(HistoName,
1224 HistoName,
1225 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
1226 psRes_VtxZ.getParameter<double>("xmin"),
1227 psRes_VtxZ.getParameter<double>("xmax"));
1228 resVtxZ_eta1p6to2->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
1229 resVtxZ_eta1p6to2->setAxisTitle("# tracking particles", 2);
1230
1231
1232 HistoName = "resVtxZ_eta2to2p4";
1233 resVtxZ_eta2to2p4 = iBooker.book1D(HistoName,
1234 HistoName,
1235 psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
1236 psRes_VtxZ.getParameter<double>("xmin"),
1237 psRes_VtxZ.getParameter<double>("xmax"));
1238 resVtxZ_eta2to2p4->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
1239 resVtxZ_eta2to2p4->setAxisTitle("# tracking particles", 2);
1240
1241
1242
1243 edm::ParameterSet psRes_d0 = conf_.getParameter<edm::ParameterSet>("TH1Res_d0");
1244 HistoName = "resd0_eta0to0p7";
1245 resd0_eta0to0p7 = iBooker.book1D(HistoName,
1246 HistoName,
1247 psRes_d0.getParameter<int32_t>("Nbinsx"),
1248 psRes_d0.getParameter<double>("xmin"),
1249 psRes_d0.getParameter<double>("xmax"));
1250 resd0_eta0to0p7->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1251 resd0_eta0to0p7->setAxisTitle("# tracking particles", 2);
1252
1253
1254 HistoName = "resd0_eta0p7to1";
1255 resd0_eta0p7to1 = iBooker.book1D(HistoName,
1256 HistoName,
1257 psRes_d0.getParameter<int32_t>("Nbinsx"),
1258 psRes_d0.getParameter<double>("xmin"),
1259 psRes_d0.getParameter<double>("xmax"));
1260 resd0_eta0p7to1->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1261 resd0_eta0p7to1->setAxisTitle("# tracking particles", 2);
1262
1263
1264 HistoName = "resd0_eta1to1p2";
1265 resd0_eta1to1p2 = iBooker.book1D(HistoName,
1266 HistoName,
1267 psRes_d0.getParameter<int32_t>("Nbinsx"),
1268 psRes_d0.getParameter<double>("xmin"),
1269 psRes_d0.getParameter<double>("xmax"));
1270 resd0_eta1to1p2->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1271 resd0_eta1to1p2->setAxisTitle("# tracking particles", 2);
1272
1273
1274 HistoName = "resd0_eta1p2to1p6";
1275 resd0_eta1p2to1p6 = iBooker.book1D(HistoName,
1276 HistoName,
1277 psRes_d0.getParameter<int32_t>("Nbinsx"),
1278 psRes_d0.getParameter<double>("xmin"),
1279 psRes_d0.getParameter<double>("xmax"));
1280 resd0_eta1p2to1p6->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1281 resd0_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
1282
1283
1284 HistoName = "resd0_eta1p6to2";
1285 resd0_eta1p6to2 = iBooker.book1D(HistoName,
1286 HistoName,
1287 psRes_d0.getParameter<int32_t>("Nbinsx"),
1288 psRes_d0.getParameter<double>("xmin"),
1289 psRes_d0.getParameter<double>("xmax"));
1290 resd0_eta1p6to2->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1291 resd0_eta1p6to2->setAxisTitle("# tracking particles", 2);
1292
1293
1294 HistoName = "resd0_eta2to2p4";
1295 resd0_eta2to2p4 = iBooker.book1D(HistoName,
1296 HistoName,
1297 psRes_d0.getParameter<int32_t>("Nbinsx"),
1298 psRes_d0.getParameter<double>("xmin"),
1299 psRes_d0.getParameter<double>("xmax"));
1300 resd0_eta2to2p4->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1301 resd0_eta2to2p4->setAxisTitle("# tracking particles", 2);
1302
1303 }
1304
1305 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
1306 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
1307 void Phase2OTValidateReconstruction::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
1308
1309 edm::ParameterSetDescription desc;
1310 {
1311 edm::ParameterSetDescription psd0;
1312 psd0.add<int>("Nbinsx", 45);
1313 psd0.add<double>("xmax", 3);
1314 psd0.add<double>("xmin", -3);
1315 desc.add<edm::ParameterSetDescription>("TH1TrackParts_Eta", psd0);
1316 }
1317 {
1318 edm::ParameterSetDescription psd0;
1319 psd0.add<int>("Nbinsx", 60);
1320 psd0.add<double>("xmax", 3.141592653589793);
1321 psd0.add<double>("xmin", -3.141592653589793);
1322 desc.add<edm::ParameterSetDescription>("TH1TrackParts_Phi", psd0);
1323 }
1324 {
1325 edm::ParameterSetDescription psd0;
1326 psd0.add<int>("Nbinsx", 45);
1327 psd0.add<double>("xmax", 100);
1328 psd0.add<double>("xmin", 0);
1329 desc.add<edm::ParameterSetDescription>("TH1TrackParts_Pt", psd0);
1330 }
1331 {
1332 edm::ParameterSetDescription psd0;
1333 psd0.add<int>("Nbinsx", 200);
1334 psd0.add<double>("xmax", 0.5);
1335 psd0.add<double>("xmin", -0.5);
1336 desc.add<edm::ParameterSetDescription>("TH1Res_ptRel", psd0);
1337 }
1338 {
1339 edm::ParameterSetDescription psd0;
1340 psd0.add<int>("Nbinsx", 50);
1341 psd0.add<double>("xmax", 100);
1342 psd0.add<double>("xmin", 0);
1343 desc.add<edm::ParameterSetDescription>("TH1Effic_pt", psd0);
1344 }
1345 {
1346 edm::ParameterSetDescription psd0;
1347 psd0.add<int>("Nbinsx", 50);
1348 psd0.add<double>("xmax", 10);
1349 psd0.add<double>("xmin", 0);
1350 desc.add<edm::ParameterSetDescription>("TH1Effic_pt_zoom", psd0);
1351 }
1352 {
1353 edm::ParameterSetDescription psd0;
1354 psd0.add<int>("Nbinsx", 50);
1355 psd0.add<double>("xmax", 2.5);
1356 psd0.add<double>("xmin", -2.5);
1357 desc.add<edm::ParameterSetDescription>("TH1Effic_eta", psd0);
1358 }
1359 {
1360 edm::ParameterSetDescription psd0;
1361 psd0.add<int>("Nbinsx", 50);
1362 psd0.add<double>("xmax", 2);
1363 psd0.add<double>("xmin", -2);
1364 desc.add<edm::ParameterSetDescription>("TH1Effic_d0", psd0);
1365 }
1366 {
1367 edm::ParameterSetDescription psd0;
1368 psd0.add<int>("Nbinsx", 50);
1369 psd0.add<double>("xmax", 5);
1370 psd0.add<double>("xmin", -5);
1371 desc.add<edm::ParameterSetDescription>("TH1Effic_VtxR", psd0);
1372 }
1373 {
1374 edm::ParameterSetDescription psd0;
1375 psd0.add<int>("Nbinsx", 50);
1376 psd0.add<double>("xmax", 30);
1377 psd0.add<double>("xmin", -30);
1378 desc.add<edm::ParameterSetDescription>("TH1Effic_VtxZ", psd0);
1379 }
1380 {
1381 edm::ParameterSetDescription psd0;
1382 psd0.add<int>("Nbinsx", 100);
1383 psd0.add<double>("xmax", 0.2);
1384 psd0.add<double>("xmin", -0.2);
1385 desc.add<edm::ParameterSetDescription>("TH1Res_pt", psd0);
1386 }
1387 {
1388 edm::ParameterSetDescription psd0;
1389 psd0.add<int>("Nbinsx", 100);
1390 psd0.add<double>("xmax", 0.01);
1391 psd0.add<double>("xmin", -0.01);
1392 desc.add<edm::ParameterSetDescription>("TH1Res_eta", psd0);
1393 }
1394 {
1395 edm::ParameterSetDescription psd0;
1396 psd0.add<int>("Nbinsx", 100);
1397 psd0.add<double>("xmax", 0.01);
1398 psd0.add<double>("xmin", -0.01);
1399 desc.add<edm::ParameterSetDescription>("TH1Res_phi", psd0);
1400 }
1401 {
1402 edm::ParameterSetDescription psd0;
1403 psd0.add<int>("Nbinsx", 100);
1404 psd0.add<double>("xmax", 1.0);
1405 psd0.add<double>("xmin", -1.0);
1406 desc.add<edm::ParameterSetDescription>("TH1Res_VtxZ", psd0);
1407 }
1408 {
1409 edm::ParameterSetDescription psd0;
1410 psd0.add<int>("Nbinsx", 100);
1411 psd0.add<double>("xmax", 0.05);
1412 psd0.add<double>("xmin", -0.05);
1413 desc.add<edm::ParameterSetDescription>("TH1Res_d0", psd0);
1414 }
1415 desc.add<std::string>("TopFolderName", "TrackerPhase2OTL1TrackV");
1416 desc.add<edm::InputTag>("trackingParticleToken", edm::InputTag("mix", "MergedTrackTruth"));
1417 desc.add<edm::InputTag>("MCTruthStubInputTag", edm::InputTag("TTStubAssociatorFromPixelDigis", "StubAccepted"));
1418 desc.add<edm::InputTag>("MCTruthTrackInputTag", edm::InputTag("TTTrackAssociatorFromPixelDigis", "Level1TTTracks"));
1419 desc.add<edm::InputTag>("MCTruthClusterInputTag",
1420 edm::InputTag("TTClusterAssociatorFromPixelDigis", "ClusterInclusive"));
1421 desc.add<int>("L1Tk_minNStub", 4);
1422 desc.add<double>("L1Tk_maxChi2dof", 25.0);
1423 desc.add<int>("TP_minNStub", 4);
1424 desc.add<int>("TP_minNLayersStub", 4);
1425 desc.add<double>("TP_minPt", 1.5);
1426 desc.add<double>("TP_maxEta", 2.4);
1427 desc.add<double>("TP_maxVtxZ", 15.0);
1428 descriptions.add("Phase2OTValidateReconstruction", desc);
1429
1430
1431 }
1432
1433 DEFINE_FWK_MODULE(Phase2OTValidateReconstruction);