Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-08-09 23:47:48

0001 // Package:    SiOuterTrackerV
0002 // Class:      SiOuterTrackerV
0003 
0004 // Original Author:  Emily MacDonald
0005 
0006 // system include files
0007 #include <memory>
0008 #include <numeric>
0009 #include <vector>
0010 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0011 #include "DQMServices/Core/interface/DQMStore.h"
0012 #include "DataFormats/Common/interface/DetSetVector.h"
0013 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0014 #include "DataFormats/Common/interface/Ptr.h"
0015 #include "DataFormats/Common/interface/Ref.h"
0016 #include "DataFormats/L1TrackTrigger/interface/TTCluster.h"
0017 #include "DataFormats/L1TrackTrigger/interface/TTStub.h"
0018 #include "DataFormats/L1TrackTrigger/interface/TTTrack.h"
0019 #include "DataFormats/L1TrackTrigger/interface/TTTypes.h"
0020 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0021 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0022 
0023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0024 #include "FWCore/ServiceRegistry/interface/Service.h"
0025 #include "FWCore/Framework/interface/ESHandle.h"
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/Frameworkfwd.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/Utilities/interface/EDGetToken.h"
0031 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0032 #include "Geometry/Records/interface/StackedTrackerGeometryRecord.h"
0033 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0034 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
0035 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0036 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
0037 #include "SimDataFormats/Associations/interface/TTClusterAssociationMap.h"
0038 #include "SimDataFormats/Associations/interface/TTStubAssociationMap.h"
0039 #include "SimDataFormats/Associations/interface/TTTrackAssociationMap.h"
0040 
0041 class Phase2OTValidateTrackingParticles : public DQMEDAnalyzer {
0042 public:
0043   explicit Phase2OTValidateTrackingParticles(const edm::ParameterSet &);
0044   ~Phase2OTValidateTrackingParticles() override;
0045   void analyze(const edm::Event &, const edm::EventSetup &) override;
0046   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0047   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0048   // Tracking particle distributions
0049   MonitorElement *trackParts_Eta = nullptr;
0050   MonitorElement *trackParts_Phi = nullptr;
0051   MonitorElement *trackParts_Pt = nullptr;
0052 
0053   // pT and eta for efficiency plots
0054   MonitorElement *tp_pt = nullptr;             // denominator
0055   MonitorElement *tp_pt_zoom = nullptr;        // denominator
0056   MonitorElement *tp_eta = nullptr;            // denominator
0057   MonitorElement *tp_d0 = nullptr;             // denominator
0058   MonitorElement *tp_VtxR = nullptr;           // denominator (also known as vxy)
0059   MonitorElement *tp_VtxZ = nullptr;           // denominator
0060   MonitorElement *match_tp_pt = nullptr;       // numerator
0061   MonitorElement *match_tp_pt_zoom = nullptr;  // numerator
0062   MonitorElement *match_tp_eta = nullptr;      // numerator
0063   MonitorElement *match_tp_d0 = nullptr;       // numerator
0064   MonitorElement *match_tp_VtxR = nullptr;     // numerator (also known as vxy)
0065   MonitorElement *match_tp_VtxZ = nullptr;     // numerator
0066 
0067   // 1D intermediate resolution plots (pT and eta)
0068   MonitorElement *res_eta = nullptr;    // for all eta and pT
0069   MonitorElement *res_pt = nullptr;     // for all eta and pT
0070   MonitorElement *res_ptRel = nullptr;  // for all eta and pT (delta(pT)/pT)
0071   MonitorElement *respt_eta0to0p7_pt2to3 = nullptr;
0072   MonitorElement *respt_eta0p7to1_pt2to3 = nullptr;
0073   MonitorElement *respt_eta1to1p2_pt2to3 = nullptr;
0074   MonitorElement *respt_eta1p2to1p6_pt2to3 = nullptr;
0075   MonitorElement *respt_eta1p6to2_pt2to3 = nullptr;
0076   MonitorElement *respt_eta2to2p4_pt2to3 = nullptr;
0077   MonitorElement *respt_eta0to0p7_pt3to8 = nullptr;
0078   MonitorElement *respt_eta0p7to1_pt3to8 = nullptr;
0079   MonitorElement *respt_eta1to1p2_pt3to8 = nullptr;
0080   MonitorElement *respt_eta1p2to1p6_pt3to8 = nullptr;
0081   MonitorElement *respt_eta1p6to2_pt3to8 = nullptr;
0082   MonitorElement *respt_eta2to2p4_pt3to8 = nullptr;
0083   MonitorElement *respt_eta0to0p7_pt8toInf = nullptr;
0084   MonitorElement *respt_eta0p7to1_pt8toInf = nullptr;
0085   MonitorElement *respt_eta1to1p2_pt8toInf = nullptr;
0086   MonitorElement *respt_eta1p2to1p6_pt8toInf = nullptr;
0087   MonitorElement *respt_eta1p6to2_pt8toInf = nullptr;
0088   MonitorElement *respt_eta2to2p4_pt8toInf = nullptr;
0089   MonitorElement *reseta_eta0to0p7 = nullptr;
0090   MonitorElement *reseta_eta0p7to1 = nullptr;
0091   MonitorElement *reseta_eta1to1p2 = nullptr;
0092   MonitorElement *reseta_eta1p2to1p6 = nullptr;
0093   MonitorElement *reseta_eta1p6to2 = nullptr;
0094   MonitorElement *reseta_eta2to2p4 = nullptr;
0095   MonitorElement *resphi_eta0to0p7 = nullptr;
0096   MonitorElement *resphi_eta0p7to1 = nullptr;
0097   MonitorElement *resphi_eta1to1p2 = nullptr;
0098   MonitorElement *resphi_eta1p2to1p6 = nullptr;
0099   MonitorElement *resphi_eta1p6to2 = nullptr;
0100   MonitorElement *resphi_eta2to2p4 = nullptr;
0101   MonitorElement *resVtxZ_eta0to0p7 = nullptr;
0102   MonitorElement *resVtxZ_eta0p7to1 = nullptr;
0103   MonitorElement *resVtxZ_eta1to1p2 = nullptr;
0104   MonitorElement *resVtxZ_eta1p2to1p6 = nullptr;
0105   MonitorElement *resVtxZ_eta1p6to2 = nullptr;
0106   MonitorElement *resVtxZ_eta2to2p4 = nullptr;
0107 
0108   // For d0
0109   MonitorElement *resd0_eta0to0p7 = nullptr;
0110   MonitorElement *resd0_eta0p7to1 = nullptr;
0111   MonitorElement *resd0_eta1to1p2 = nullptr;
0112   MonitorElement *resd0_eta1p2to1p6 = nullptr;
0113   MonitorElement *resd0_eta1p6to2 = nullptr;
0114   MonitorElement *resd0_eta2to2p4 = nullptr;
0115 
0116 private:
0117   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> m_topoToken;
0118   edm::ParameterSet conf_;
0119   edm::EDGetTokenT<std::vector<TrackingParticle>> trackingParticleToken_;
0120   edm::EDGetTokenT<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>>
0121       ttClusterMCTruthToken_;  // MC truth association map for clusters
0122   edm::EDGetTokenT<TTStubAssociationMap<Ref_Phase2TrackerDigi_>>
0123       ttStubMCTruthToken_;  // MC truth association map for stubs
0124   edm::EDGetTokenT<TTTrackAssociationMap<Ref_Phase2TrackerDigi_>>
0125       ttTrackMCTruthToken_;  // MC truth association map for tracks
0126   int L1Tk_minNStub;
0127   double L1Tk_maxChi2dof;
0128   int TP_minNStub;
0129   int TP_minNLayersStub;
0130   double TP_minPt;
0131   double TP_maxEta;
0132   double TP_maxVtxZ;
0133   std::string topFolderName_;
0134 };
0135 
0136 //
0137 // constructors and destructor
0138 //
0139 Phase2OTValidateTrackingParticles::Phase2OTValidateTrackingParticles(const edm::ParameterSet &iConfig)
0140     : m_topoToken(esConsumes()), conf_(iConfig) {
0141   topFolderName_ = conf_.getParameter<std::string>("TopFolderName");
0142   trackingParticleToken_ =
0143       consumes<std::vector<TrackingParticle>>(conf_.getParameter<edm::InputTag>("trackingParticleToken"));
0144   ttStubMCTruthToken_ =
0145       consumes<TTStubAssociationMap<Ref_Phase2TrackerDigi_>>(conf_.getParameter<edm::InputTag>("MCTruthStubInputTag"));
0146   ttClusterMCTruthToken_ = consumes<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>>(
0147       conf_.getParameter<edm::InputTag>("MCTruthClusterInputTag"));
0148   ttTrackMCTruthToken_ = consumes<TTTrackAssociationMap<Ref_Phase2TrackerDigi_>>(
0149       conf_.getParameter<edm::InputTag>("MCTruthTrackInputTag"));
0150   L1Tk_minNStub = conf_.getParameter<int>("L1Tk_minNStub");         // min number of stubs in the track
0151   L1Tk_maxChi2dof = conf_.getParameter<double>("L1Tk_maxChi2dof");  // maximum chi2/dof of the track
0152   TP_minNStub = conf_.getParameter<int>("TP_minNStub");             // min number of stubs in the tracking particle to
0153   //min number of layers with stubs in the tracking particle to consider matching
0154   TP_minNLayersStub = conf_.getParameter<int>("TP_minNLayersStub");
0155   TP_minPt = conf_.getParameter<double>("TP_minPt");      // min pT to consider matching
0156   TP_maxEta = conf_.getParameter<double>("TP_maxEta");    // max eta to consider matching
0157   TP_maxVtxZ = conf_.getParameter<double>("TP_maxVtxZ");  // max vertZ (or z0) to consider matching
0158 }
0159 
0160 Phase2OTValidateTrackingParticles::~Phase2OTValidateTrackingParticles() = default;
0161 
0162 // member functions
0163 
0164 // ------------ method called for each event  ------------
0165 void Phase2OTValidateTrackingParticles::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0166   // Tracking Particles
0167   edm::Handle<std::vector<TrackingParticle>> trackingParticleHandle;
0168   iEvent.getByToken(trackingParticleToken_, trackingParticleHandle);
0169 
0170   // Truth Association Maps
0171   edm::Handle<TTTrackAssociationMap<Ref_Phase2TrackerDigi_>> MCTruthTTTrackHandle;
0172   iEvent.getByToken(ttTrackMCTruthToken_, MCTruthTTTrackHandle);
0173   edm::Handle<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>> MCTruthTTClusterHandle;
0174   iEvent.getByToken(ttClusterMCTruthToken_, MCTruthTTClusterHandle);
0175   edm::Handle<TTStubAssociationMap<Ref_Phase2TrackerDigi_>> MCTruthTTStubHandle;
0176   iEvent.getByToken(ttStubMCTruthToken_, MCTruthTTStubHandle);
0177 
0178   // Geometries
0179   const TrackerTopology *const tTopo = &iSetup.getData(m_topoToken);
0180 
0181   // Loop over tracking particles
0182   int this_tp = 0;
0183   for (const auto &iterTP : *trackingParticleHandle) {
0184     edm::Ptr<TrackingParticle> tp_ptr(trackingParticleHandle, this_tp);
0185     this_tp++;
0186 
0187     // int tmp_eventid = iterTP.eventId().event();
0188     float tmp_tp_pt = iterTP.pt();
0189     float tmp_tp_phi = iterTP.phi();
0190     float tmp_tp_eta = iterTP.eta();
0191 
0192     //Calculate nLayers variable
0193     std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>>
0194         theStubRefs = MCTruthTTStubHandle->findTTStubRefs(tp_ptr);
0195 
0196     int hasStubInLayer[11] = {0};
0197     for (unsigned int is = 0; is < theStubRefs.size(); is++) {
0198       DetId detid(theStubRefs.at(is)->getDetId());
0199       int layer = -1;
0200       if (detid.subdetId() == StripSubdetector::TOB)
0201         layer = static_cast<int>(tTopo->layer(detid)) - 1;  //fill in array as entries 0-5
0202       else if (detid.subdetId() == StripSubdetector::TID)
0203         layer = static_cast<int>(tTopo->layer(detid)) + 5;  //fill in array as entries 6-10
0204 
0205       //treat genuine stubs separately (==2 is genuine, ==1 is not)
0206       if (MCTruthTTStubHandle->findTrackingParticlePtr(theStubRefs.at(is)).isNull() && hasStubInLayer[layer] < 2)
0207         hasStubInLayer[layer] = 1;
0208       else
0209         hasStubInLayer[layer] = 2;
0210     }
0211 
0212     int nStubLayerTP = 0;
0213     for (int isum = 0; isum < 11; isum++) {
0214       if (hasStubInLayer[isum] >= 1)
0215         nStubLayerTP += 1;
0216     }
0217 
0218     if (std::fabs(tmp_tp_eta) > TP_maxEta)
0219       continue;
0220     // Fill the 1D distribution plots for tracking particles, to monitor change in stub definition
0221     if (tmp_tp_pt > TP_minPt && nStubLayerTP >= TP_minNLayersStub) {
0222       trackParts_Pt->Fill(tmp_tp_pt);
0223       trackParts_Eta->Fill(tmp_tp_eta);
0224       trackParts_Phi->Fill(tmp_tp_phi);
0225     }
0226 
0227     // if (TP_select_eventid == 0 && tmp_eventid != 0)
0228     //   continue;  //only care about tracking particles from the primary interaction for efficiency/resolution
0229     int nStubTP = -1;
0230     if (MCTruthTTStubHandle.isValid()) {
0231       std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>>
0232           theStubRefs = MCTruthTTStubHandle->findTTStubRefs(tp_ptr);
0233       nStubTP = (int)theStubRefs.size();
0234     }
0235     if (MCTruthTTClusterHandle.isValid() && MCTruthTTClusterHandle->findTTClusterRefs(tp_ptr).empty())
0236       continue;
0237 
0238     float tmp_tp_vz = iterTP.vz();
0239     float tmp_tp_vx = iterTP.vx();
0240     float tmp_tp_vy = iterTP.vy();
0241     float tmp_tp_charge = tp_ptr->charge();
0242     int tmp_tp_pdgid = iterTP.pdgId();
0243 
0244     // ----------------------------------------------------------------------------------------------
0245     // calculate d0 and VtxZ propagated back to the IP, pass if greater than max
0246     // VtxZ
0247     float tmp_tp_t = tan(2.0 * atan(1.0) - 2.0 * atan(exp(-tmp_tp_eta)));
0248     float delx = -tmp_tp_vx;
0249     float dely = -tmp_tp_vy;
0250     float K = 0.01 * 0.5696 / tmp_tp_pt * tmp_tp_charge;  // curvature correction
0251     float A = 1. / (2. * K);
0252     float tmp_tp_x0p = delx - A * sin(tmp_tp_phi);
0253     float tmp_tp_y0p = dely + A * cos(tmp_tp_phi);
0254     float tmp_tp_rp = sqrt(tmp_tp_x0p * tmp_tp_x0p + tmp_tp_y0p * tmp_tp_y0p);
0255     static double pi = 4.0 * atan(1.0);
0256     float delphi = tmp_tp_phi - atan2(-K * tmp_tp_x0p, K * tmp_tp_y0p);
0257     if (delphi < -pi)
0258       delphi += 2.0 * pi;
0259     if (delphi > pi)
0260       delphi -= 2.0 * pi;
0261 
0262     float tmp_tp_VtxZ = tmp_tp_vz + tmp_tp_t * delphi / (2.0 * K);
0263     float tmp_tp_VtxR = sqrt(tmp_tp_vx * tmp_tp_vx + tmp_tp_vy * tmp_tp_vy);
0264     float tmp_tp_d0 = tmp_tp_charge * tmp_tp_rp - (1. / (2. * K));
0265 
0266     // simpler formula for d0, in cases where the charge is zero:
0267     // https://github.com/cms-sw/cmssw/blob/master/DataFormats/TrackReco/interface/TrackBase.h
0268     float other_d0 = -tmp_tp_vx * sin(tmp_tp_phi) + tmp_tp_vy * cos(tmp_tp_phi);
0269     tmp_tp_d0 = tmp_tp_d0 * (-1);  // fix d0 sign
0270     if (K == 0) {
0271       tmp_tp_d0 = other_d0;
0272       tmp_tp_VtxZ = tmp_tp_vz;
0273     }
0274     if (std::fabs(tmp_tp_VtxZ) > TP_maxVtxZ)
0275       continue;
0276 
0277     // To make efficiency plots where the denominator has NO stub cuts
0278     if (tmp_tp_VtxR < 1.0) {
0279       tp_pt->Fill(tmp_tp_pt);  //pT effic, no cut on pT, but VtxR cut
0280       if (tmp_tp_pt <= 10)
0281         tp_pt_zoom->Fill(tmp_tp_pt);  //pT effic, no cut on pT, but VtxR cut
0282     }
0283     if (tmp_tp_pt < TP_minPt)
0284       continue;
0285     tp_VtxR->Fill(tmp_tp_VtxR);  // VtxR efficiency has no cut on VtxR
0286     if (tmp_tp_VtxR > 1.0)
0287       continue;
0288     tp_eta->Fill(tmp_tp_eta);
0289     tp_d0->Fill(tmp_tp_d0);
0290     tp_VtxZ->Fill(tmp_tp_VtxZ);
0291 
0292     if (nStubTP < TP_minNStub || nStubLayerTP < TP_minNLayersStub)
0293       continue;  //nStub cut not included in denominator of efficiency plots
0294 
0295     // ----------------------------------------------------------------------------------------------
0296     // look for L1 tracks matched to the tracking particle
0297     int tp_nMatch = 0;
0298     int i_track = -1;
0299     float i_chi2dof = 99999;
0300     if (MCTruthTTTrackHandle.isValid()) {
0301       std::vector<edm::Ptr<TTTrack<Ref_Phase2TrackerDigi_>>> matchedTracks =
0302           MCTruthTTTrackHandle->findTTTrackPtrs(tp_ptr);
0303 
0304       // ----------------------------------------------------------------------------------------------
0305       // loop over matched L1 tracks
0306       // here, "match" means tracks that can be associated to a TrackingParticle
0307       // with at least one hit of at least one of its clusters
0308       // https://twiki.cern.ch/twiki/bin/viewauth/CMS/SLHCTrackerTriggerSWTools#MC_truth_for_TTTrack
0309       int trkCounter = 0;
0310       for (const auto &thisTrack : matchedTracks) {
0311         if (!MCTruthTTTrackHandle->isGenuine(thisTrack))
0312           continue;
0313         // ----------------------------------------------------------------------------------------------
0314         // further require L1 track to be (loosely) genuine, that there is only
0315         // one TP matched to the track
0316         // + have >= L1Tk_minNStub stubs for it to be a valid match
0317         int tmp_trk_nstub = thisTrack->getStubRefs().size();
0318         if (tmp_trk_nstub < L1Tk_minNStub)
0319           continue;
0320         float dmatch_pt = 999;
0321         float dmatch_eta = 999;
0322         float dmatch_phi = 999;
0323         int match_id = 999;
0324 
0325         edm::Ptr<TrackingParticle> my_tp = MCTruthTTTrackHandle->findTrackingParticlePtr(thisTrack);
0326         dmatch_pt = std::fabs(my_tp->p4().pt() - tmp_tp_pt);
0327         dmatch_eta = std::fabs(my_tp->p4().eta() - tmp_tp_eta);
0328         dmatch_phi = std::fabs(my_tp->p4().phi() - tmp_tp_phi);
0329         match_id = my_tp->pdgId();
0330         float tmp_trk_chi2dof = thisTrack->chi2Red();
0331 
0332         // ensure that track is uniquely matched to the TP we are looking at!
0333         if (dmatch_pt < 0.1 && dmatch_eta < 0.1 && dmatch_phi < 0.1 && tmp_tp_pdgid == match_id) {
0334           tp_nMatch++;
0335           if (i_track < 0 || tmp_trk_chi2dof < i_chi2dof) {
0336             i_track = trkCounter;
0337             i_chi2dof = tmp_trk_chi2dof;
0338           }
0339         }
0340         trkCounter++;
0341       }  // end loop over matched L1 tracks
0342 
0343       if (tp_nMatch < 1)
0344         continue;
0345       // Get information on the matched tracks
0346       float tmp_matchtrk_pt = -999;
0347       float tmp_matchtrk_eta = -999;
0348       float tmp_matchtrk_phi = -999;
0349       float tmp_matchtrk_VtxZ = -999;
0350       float tmp_matchtrk_chi2dof = -999;
0351       int tmp_matchTrk_nStub = -999;
0352       float tmp_matchtrk_d0 = -999;
0353 
0354       tmp_matchtrk_pt = matchedTracks[i_track]->momentum().perp();
0355       tmp_matchtrk_eta = matchedTracks[i_track]->momentum().eta();
0356       tmp_matchtrk_phi = matchedTracks[i_track]->momentum().phi();
0357       tmp_matchtrk_VtxZ = matchedTracks[i_track]->z0();
0358       tmp_matchtrk_chi2dof = matchedTracks[i_track]->chi2Red();
0359       tmp_matchTrk_nStub = (int)matchedTracks[i_track]->getStubRefs().size();
0360 
0361       //for d0
0362       float tmp_matchtrk_x0 = matchedTracks[i_track]->POCA().x();
0363       float tmp_matchtrk_y0 = matchedTracks[i_track]->POCA().y();
0364       tmp_matchtrk_d0 = -tmp_matchtrk_x0 * sin(tmp_matchtrk_phi) + tmp_matchtrk_y0 * cos(tmp_matchtrk_phi);
0365 
0366       //Add cuts for the matched tracks, numerator
0367       if (tmp_matchTrk_nStub < L1Tk_minNStub || tmp_matchtrk_chi2dof > L1Tk_maxChi2dof)
0368         continue;
0369 
0370       // fill matched track histograms (if passes all criteria)
0371       match_tp_pt->Fill(tmp_tp_pt);
0372       if (tmp_tp_pt > 0 && tmp_tp_pt <= 10)
0373         match_tp_pt_zoom->Fill(tmp_tp_pt);
0374       match_tp_eta->Fill(tmp_tp_eta);
0375       match_tp_d0->Fill(tmp_tp_d0);
0376       match_tp_VtxR->Fill(tmp_tp_VtxR);
0377       match_tp_VtxZ->Fill(tmp_tp_VtxZ);
0378 
0379       // Eta and pT histograms for resolution
0380       float pt_diff = tmp_matchtrk_pt - tmp_tp_pt;
0381       float pt_res = pt_diff / tmp_tp_pt;
0382       float eta_res = tmp_matchtrk_eta - tmp_tp_eta;
0383       float phi_res = tmp_matchtrk_phi - tmp_tp_phi;
0384       float VtxZ_res = tmp_matchtrk_VtxZ - tmp_tp_VtxZ;
0385       float d0_res = tmp_matchtrk_d0 - tmp_tp_d0;
0386 
0387       // fill total resolution histograms
0388       res_pt->Fill(pt_diff);
0389       res_ptRel->Fill(pt_res);
0390       res_eta->Fill(eta_res);
0391 
0392       // Fill resolution plots for different abs(eta) bins:
0393       // (0, 0.7), (0.7, 1.0), (1.0, 1.2), (1.2, 1.6), (1.6, 2.0), (2.0, 2.4)
0394       if (std::fabs(tmp_tp_eta) >= 0 && std::fabs(tmp_tp_eta) < 0.7) {
0395         reseta_eta0to0p7->Fill(eta_res);
0396         resphi_eta0to0p7->Fill(phi_res);
0397         resVtxZ_eta0to0p7->Fill(VtxZ_res);
0398         resd0_eta0to0p7->Fill(d0_res);
0399         if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0400           respt_eta0to0p7_pt2to3->Fill(pt_res);
0401         else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0402           respt_eta0to0p7_pt3to8->Fill(pt_res);
0403         else if (tmp_tp_pt >= 8)
0404           respt_eta0to0p7_pt8toInf->Fill(pt_res);
0405       } else if (std::fabs(tmp_tp_eta) >= 0.7 && std::fabs(tmp_tp_eta) < 1.0) {
0406         reseta_eta0p7to1->Fill(eta_res);
0407         resphi_eta0p7to1->Fill(phi_res);
0408         resVtxZ_eta0p7to1->Fill(VtxZ_res);
0409         resd0_eta0p7to1->Fill(d0_res);
0410         if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0411           respt_eta0p7to1_pt2to3->Fill(pt_res);
0412         else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0413           respt_eta0p7to1_pt3to8->Fill(pt_res);
0414         else if (tmp_tp_pt >= 8)
0415           respt_eta0p7to1_pt8toInf->Fill(pt_res);
0416       } else if (std::fabs(tmp_tp_eta) >= 1.0 && std::fabs(tmp_tp_eta) < 1.2) {
0417         reseta_eta1to1p2->Fill(eta_res);
0418         resphi_eta1to1p2->Fill(phi_res);
0419         resVtxZ_eta1to1p2->Fill(VtxZ_res);
0420         resd0_eta1to1p2->Fill(d0_res);
0421         if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0422           respt_eta1to1p2_pt2to3->Fill(pt_res);
0423         else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0424           respt_eta1to1p2_pt3to8->Fill(pt_res);
0425         else if (tmp_tp_pt >= 8)
0426           respt_eta1to1p2_pt8toInf->Fill(pt_res);
0427       } else if (std::fabs(tmp_tp_eta) >= 1.2 && std::fabs(tmp_tp_eta) < 1.6) {
0428         reseta_eta1p2to1p6->Fill(eta_res);
0429         resphi_eta1p2to1p6->Fill(phi_res);
0430         resVtxZ_eta1p2to1p6->Fill(VtxZ_res);
0431         resd0_eta1p2to1p6->Fill(d0_res);
0432         if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0433           respt_eta1p2to1p6_pt2to3->Fill(pt_res);
0434         else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0435           respt_eta1p2to1p6_pt3to8->Fill(pt_res);
0436         else if (tmp_tp_pt >= 8)
0437           respt_eta1p2to1p6_pt8toInf->Fill(pt_res);
0438       } else if (std::fabs(tmp_tp_eta) >= 1.6 && std::fabs(tmp_tp_eta) < 2.0) {
0439         reseta_eta1p6to2->Fill(eta_res);
0440         resphi_eta1p6to2->Fill(phi_res);
0441         resVtxZ_eta1p6to2->Fill(VtxZ_res);
0442         resd0_eta1p6to2->Fill(d0_res);
0443         if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0444           respt_eta1p6to2_pt2to3->Fill(pt_res);
0445         else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0446           respt_eta1p6to2_pt3to8->Fill(pt_res);
0447         else if (tmp_tp_pt >= 8)
0448           respt_eta1p6to2_pt8toInf->Fill(pt_res);
0449       } else if (std::fabs(tmp_tp_eta) >= 2.0 && std::fabs(tmp_tp_eta) <= 2.4) {
0450         reseta_eta2to2p4->Fill(eta_res);
0451         resphi_eta2to2p4->Fill(phi_res);
0452         resVtxZ_eta2to2p4->Fill(VtxZ_res);
0453         resd0_eta2to2p4->Fill(d0_res);
0454         if (tmp_tp_pt >= 2 && tmp_tp_pt < 3)
0455           respt_eta2to2p4_pt2to3->Fill(pt_res);
0456         else if (tmp_tp_pt >= 3 && tmp_tp_pt < 8)
0457           respt_eta2to2p4_pt3to8->Fill(pt_res);
0458         else if (tmp_tp_pt >= 8)
0459           respt_eta2to2p4_pt8toInf->Fill(pt_res);
0460       }
0461     }  //if MC TTTrack handle is valid
0462   }    //end loop over tracking particles
0463 }  // end of method
0464 
0465 // ------------ method called once each job just before starting event loop
0466 // ------------
0467 void Phase2OTValidateTrackingParticles::bookHistograms(DQMStore::IBooker &iBooker,
0468                                                        edm::Run const &run,
0469                                                        edm::EventSetup const &es) {
0470   // Histogram setup and definitions
0471   std::string HistoName;
0472   iBooker.setCurrentFolder(topFolderName_ + "/trackParticles");
0473 
0474   // 1D: pT
0475   edm::ParameterSet psTrackParts_Pt = conf_.getParameter<edm::ParameterSet>("TH1TrackParts_Pt");
0476   HistoName = "trackParts_Pt";
0477   trackParts_Pt = iBooker.book1D(HistoName,
0478                                  HistoName,
0479                                  psTrackParts_Pt.getParameter<int32_t>("Nbinsx"),
0480                                  psTrackParts_Pt.getParameter<double>("xmin"),
0481                                  psTrackParts_Pt.getParameter<double>("xmax"));
0482   trackParts_Pt->setAxisTitle("p_{T} [GeV]", 1);
0483   trackParts_Pt->setAxisTitle("# tracking particles", 2);
0484 
0485   // 1D: eta
0486   edm::ParameterSet psTrackParts_Eta = conf_.getParameter<edm::ParameterSet>("TH1TrackParts_Eta");
0487   HistoName = "trackParts_Eta";
0488   trackParts_Eta = iBooker.book1D(HistoName,
0489                                   HistoName,
0490                                   psTrackParts_Eta.getParameter<int32_t>("Nbinsx"),
0491                                   psTrackParts_Eta.getParameter<double>("xmin"),
0492                                   psTrackParts_Eta.getParameter<double>("xmax"));
0493   trackParts_Eta->setAxisTitle("#eta", 1);
0494   trackParts_Eta->setAxisTitle("# tracking particles", 2);
0495 
0496   // 1D: phi
0497   edm::ParameterSet psTrackParts_Phi = conf_.getParameter<edm::ParameterSet>("TH1TrackParts_Phi");
0498   HistoName = "trackParts_Phi";
0499   trackParts_Phi = iBooker.book1D(HistoName,
0500                                   HistoName,
0501                                   psTrackParts_Phi.getParameter<int32_t>("Nbinsx"),
0502                                   psTrackParts_Phi.getParameter<double>("xmin"),
0503                                   psTrackParts_Phi.getParameter<double>("xmax"));
0504   trackParts_Phi->setAxisTitle("#phi", 1);
0505   trackParts_Phi->setAxisTitle("# tracking particles", 2);
0506 
0507   // 1D plots for efficiency
0508   iBooker.setCurrentFolder(topFolderName_ + "/Efficiency");
0509   // pT
0510   edm::ParameterSet psEffic_pt = conf_.getParameter<edm::ParameterSet>("TH1Effic_pt");
0511   HistoName = "tp_pt";
0512   tp_pt = iBooker.book1D(HistoName,
0513                          HistoName,
0514                          psEffic_pt.getParameter<int32_t>("Nbinsx"),
0515                          psEffic_pt.getParameter<double>("xmin"),
0516                          psEffic_pt.getParameter<double>("xmax"));
0517   tp_pt->setAxisTitle("p_{T} [GeV]", 1);
0518   tp_pt->setAxisTitle("# tracking particles", 2);
0519 
0520   // Matched TP's pT
0521   HistoName = "match_tp_pt";
0522   match_tp_pt = iBooker.book1D(HistoName,
0523                                HistoName,
0524                                psEffic_pt.getParameter<int32_t>("Nbinsx"),
0525                                psEffic_pt.getParameter<double>("xmin"),
0526                                psEffic_pt.getParameter<double>("xmax"));
0527   match_tp_pt->setAxisTitle("p_{T} [GeV]", 1);
0528   match_tp_pt->setAxisTitle("# matched tracking particles", 2);
0529 
0530   // pT zoom (0-10 GeV)
0531   edm::ParameterSet psEffic_pt_zoom = conf_.getParameter<edm::ParameterSet>("TH1Effic_pt_zoom");
0532   HistoName = "tp_pt_zoom";
0533   tp_pt_zoom = iBooker.book1D(HistoName,
0534                               HistoName,
0535                               psEffic_pt_zoom.getParameter<int32_t>("Nbinsx"),
0536                               psEffic_pt_zoom.getParameter<double>("xmin"),
0537                               psEffic_pt_zoom.getParameter<double>("xmax"));
0538   tp_pt_zoom->setAxisTitle("p_{T} [GeV]", 1);
0539   tp_pt_zoom->setAxisTitle("# tracking particles", 2);
0540 
0541   // Matched pT zoom (0-10 GeV)
0542   HistoName = "match_tp_pt_zoom";
0543   match_tp_pt_zoom = iBooker.book1D(HistoName,
0544                                     HistoName,
0545                                     psEffic_pt_zoom.getParameter<int32_t>("Nbinsx"),
0546                                     psEffic_pt_zoom.getParameter<double>("xmin"),
0547                                     psEffic_pt_zoom.getParameter<double>("xmax"));
0548   match_tp_pt_zoom->setAxisTitle("p_{T} [GeV]", 1);
0549   match_tp_pt_zoom->setAxisTitle("# matched tracking particles", 2);
0550 
0551   // eta
0552   edm::ParameterSet psEffic_eta = conf_.getParameter<edm::ParameterSet>("TH1Effic_eta");
0553   HistoName = "tp_eta";
0554   tp_eta = iBooker.book1D(HistoName,
0555                           HistoName,
0556                           psEffic_eta.getParameter<int32_t>("Nbinsx"),
0557                           psEffic_eta.getParameter<double>("xmin"),
0558                           psEffic_eta.getParameter<double>("xmax"));
0559   tp_eta->setAxisTitle("#eta", 1);
0560   tp_eta->setAxisTitle("# tracking particles", 2);
0561 
0562   // Matched eta
0563   HistoName = "match_tp_eta";
0564   match_tp_eta = iBooker.book1D(HistoName,
0565                                 HistoName,
0566                                 psEffic_eta.getParameter<int32_t>("Nbinsx"),
0567                                 psEffic_eta.getParameter<double>("xmin"),
0568                                 psEffic_eta.getParameter<double>("xmax"));
0569   match_tp_eta->setAxisTitle("#eta", 1);
0570   match_tp_eta->setAxisTitle("# matched tracking particles", 2);
0571 
0572   // d0
0573   edm::ParameterSet psEffic_d0 = conf_.getParameter<edm::ParameterSet>("TH1Effic_d0");
0574   HistoName = "tp_d0";
0575   tp_d0 = iBooker.book1D(HistoName,
0576                          HistoName,
0577                          psEffic_d0.getParameter<int32_t>("Nbinsx"),
0578                          psEffic_d0.getParameter<double>("xmin"),
0579                          psEffic_d0.getParameter<double>("xmax"));
0580   tp_d0->setAxisTitle("d_{0} [cm]", 1);
0581   tp_d0->setAxisTitle("# tracking particles", 2);
0582 
0583   // Matched d0
0584   HistoName = "match_tp_d0";
0585   match_tp_d0 = iBooker.book1D(HistoName,
0586                                HistoName,
0587                                psEffic_d0.getParameter<int32_t>("Nbinsx"),
0588                                psEffic_d0.getParameter<double>("xmin"),
0589                                psEffic_d0.getParameter<double>("xmax"));
0590   match_tp_d0->setAxisTitle("d_{0} [cm]", 1);
0591   match_tp_d0->setAxisTitle("# matched tracking particles", 2);
0592 
0593   // VtxR (also known as vxy)
0594   edm::ParameterSet psEffic_VtxR = conf_.getParameter<edm::ParameterSet>("TH1Effic_VtxR");
0595   HistoName = "tp_VtxR";
0596   tp_VtxR = iBooker.book1D(HistoName,
0597                            HistoName,
0598                            psEffic_VtxR.getParameter<int32_t>("Nbinsx"),
0599                            psEffic_VtxR.getParameter<double>("xmin"),
0600                            psEffic_VtxR.getParameter<double>("xmax"));
0601   tp_VtxR->setAxisTitle("d_{xy} [cm]", 1);
0602   tp_VtxR->setAxisTitle("# tracking particles", 2);
0603 
0604   // Matched VtxR
0605   HistoName = "match_tp_VtxR";
0606   match_tp_VtxR = iBooker.book1D(HistoName,
0607                                  HistoName,
0608                                  psEffic_VtxR.getParameter<int32_t>("Nbinsx"),
0609                                  psEffic_VtxR.getParameter<double>("xmin"),
0610                                  psEffic_VtxR.getParameter<double>("xmax"));
0611   match_tp_VtxR->setAxisTitle("d_{xy} [cm]", 1);
0612   match_tp_VtxR->setAxisTitle("# matched tracking particles", 2);
0613 
0614   // VtxZ
0615   edm::ParameterSet psEffic_VtxZ = conf_.getParameter<edm::ParameterSet>("TH1Effic_VtxZ");
0616   HistoName = "tp_VtxZ";
0617   tp_VtxZ = iBooker.book1D(HistoName,
0618                            HistoName,
0619                            psEffic_VtxZ.getParameter<int32_t>("Nbinsx"),
0620                            psEffic_VtxZ.getParameter<double>("xmin"),
0621                            psEffic_VtxZ.getParameter<double>("xmax"));
0622   tp_VtxZ->setAxisTitle("z_{0} [cm]", 1);
0623   tp_VtxZ->setAxisTitle("# tracking particles", 2);
0624 
0625   // Matched d0
0626   HistoName = "match_tp_VtxZ";
0627   match_tp_VtxZ = iBooker.book1D(HistoName,
0628                                  HistoName,
0629                                  psEffic_VtxZ.getParameter<int32_t>("Nbinsx"),
0630                                  psEffic_VtxZ.getParameter<double>("xmin"),
0631                                  psEffic_VtxZ.getParameter<double>("xmax"));
0632   match_tp_VtxZ->setAxisTitle("z_{0} [cm]", 1);
0633   match_tp_VtxZ->setAxisTitle("# matched tracking particles", 2);
0634 
0635   // 1D plots for resolution
0636   iBooker.setCurrentFolder(topFolderName_ + "/Resolution");
0637   // full pT
0638   edm::ParameterSet psRes_pt = conf_.getParameter<edm::ParameterSet>("TH1Res_pt");
0639   HistoName = "res_pt";
0640   res_pt = iBooker.book1D(HistoName,
0641                           HistoName,
0642                           psRes_pt.getParameter<int32_t>("Nbinsx"),
0643                           psRes_pt.getParameter<double>("xmin"),
0644                           psRes_pt.getParameter<double>("xmax"));
0645   res_pt->setAxisTitle("p_{T} [GeV]", 1);
0646   res_pt->setAxisTitle("# tracking particles", 2);
0647 
0648   // Full eta
0649   edm::ParameterSet psRes_eta = conf_.getParameter<edm::ParameterSet>("TH1Res_eta");
0650   HistoName = "res_eta";
0651   res_eta = iBooker.book1D(HistoName,
0652                            HistoName,
0653                            psRes_eta.getParameter<int32_t>("Nbinsx"),
0654                            psRes_eta.getParameter<double>("xmin"),
0655                            psRes_eta.getParameter<double>("xmax"));
0656   res_eta->setAxisTitle("#eta", 1);
0657   res_eta->setAxisTitle("# tracking particles", 2);
0658 
0659   // Relative pT
0660   edm::ParameterSet psRes_ptRel = conf_.getParameter<edm::ParameterSet>("TH1Res_ptRel");
0661   HistoName = "res_ptRel";
0662   res_ptRel = iBooker.book1D(HistoName,
0663                              HistoName,
0664                              psRes_ptRel.getParameter<int32_t>("Nbinsx"),
0665                              psRes_ptRel.getParameter<double>("xmin"),
0666                              psRes_ptRel.getParameter<double>("xmax"));
0667   res_ptRel->setAxisTitle("Relative p_{T} [GeV]", 1);
0668   res_ptRel->setAxisTitle("# tracking particles", 2);
0669 
0670   // Eta parts (for resolution)
0671   // Eta 1 (0 to 0.7)
0672   HistoName = "reseta_eta0to0p7";
0673   reseta_eta0to0p7 = iBooker.book1D(HistoName,
0674                                     HistoName,
0675                                     psRes_eta.getParameter<int32_t>("Nbinsx"),
0676                                     psRes_eta.getParameter<double>("xmin"),
0677                                     psRes_eta.getParameter<double>("xmax"));
0678   reseta_eta0to0p7->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0679   reseta_eta0to0p7->setAxisTitle("# tracking particles", 2);
0680 
0681   // Eta 2 (0.7 to 1.0)
0682   HistoName = "reseta_eta0p7to1";
0683   reseta_eta0p7to1 = iBooker.book1D(HistoName,
0684                                     HistoName,
0685                                     psRes_eta.getParameter<int32_t>("Nbinsx"),
0686                                     psRes_eta.getParameter<double>("xmin"),
0687                                     psRes_eta.getParameter<double>("xmax"));
0688   reseta_eta0p7to1->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0689   reseta_eta0p7to1->setAxisTitle("# tracking particles", 2);
0690 
0691   // Eta 3 (1.0 to 1.2)
0692   HistoName = "reseta_eta1to1p2";
0693   reseta_eta1to1p2 = iBooker.book1D(HistoName,
0694                                     HistoName,
0695                                     psRes_eta.getParameter<int32_t>("Nbinsx"),
0696                                     psRes_eta.getParameter<double>("xmin"),
0697                                     psRes_eta.getParameter<double>("xmax"));
0698   reseta_eta1to1p2->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0699   reseta_eta1to1p2->setAxisTitle("# tracking particles", 2);
0700 
0701   // Eta 4 (1.2 to 1.6)
0702   HistoName = "reseta_eta1p2to1p6";
0703   reseta_eta1p2to1p6 = iBooker.book1D(HistoName,
0704                                       HistoName,
0705                                       psRes_eta.getParameter<int32_t>("Nbinsx"),
0706                                       psRes_eta.getParameter<double>("xmin"),
0707                                       psRes_eta.getParameter<double>("xmax"));
0708   reseta_eta1p2to1p6->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0709   reseta_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
0710 
0711   // Eta 5 (1.6 to 2.0)
0712   HistoName = "reseta_eta1p6to2";
0713   reseta_eta1p6to2 = iBooker.book1D(HistoName,
0714                                     HistoName,
0715                                     psRes_eta.getParameter<int32_t>("Nbinsx"),
0716                                     psRes_eta.getParameter<double>("xmin"),
0717                                     psRes_eta.getParameter<double>("xmax"));
0718   reseta_eta1p6to2->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0719   reseta_eta1p6to2->setAxisTitle("# tracking particles", 2);
0720 
0721   // Eta 6 (2.0 to 2.4)
0722   HistoName = "reseta_eta2to2p4";
0723   reseta_eta2to2p4 = iBooker.book1D(HistoName,
0724                                     HistoName,
0725                                     psRes_eta.getParameter<int32_t>("Nbinsx"),
0726                                     psRes_eta.getParameter<double>("xmin"),
0727                                     psRes_eta.getParameter<double>("xmax"));
0728   reseta_eta2to2p4->setAxisTitle("#eta_{trk} - #eta_{tp}", 1);
0729   reseta_eta2to2p4->setAxisTitle("# tracking particles", 2);
0730 
0731   // pT parts for resolution (pT res vs eta)
0732   // pT a (2 to 3 GeV)
0733   // Eta 1 (0 to 0.7)
0734   HistoName = "respt_eta0to0p7_pt2to3";
0735   respt_eta0to0p7_pt2to3 = iBooker.book1D(HistoName,
0736                                           HistoName,
0737                                           psRes_pt.getParameter<int32_t>("Nbinsx"),
0738                                           psRes_pt.getParameter<double>("xmin"),
0739                                           psRes_pt.getParameter<double>("xmax"));
0740   respt_eta0to0p7_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0741   respt_eta0to0p7_pt2to3->setAxisTitle("# tracking particles", 2);
0742 
0743   // Eta 2 (0.7 to 1.0)
0744   HistoName = "respt_eta0p7to1_pt2to3";
0745   respt_eta0p7to1_pt2to3 = iBooker.book1D(HistoName,
0746                                           HistoName,
0747                                           psRes_pt.getParameter<int32_t>("Nbinsx"),
0748                                           psRes_pt.getParameter<double>("xmin"),
0749                                           psRes_pt.getParameter<double>("xmax"));
0750   respt_eta0p7to1_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0751   respt_eta0p7to1_pt2to3->setAxisTitle("# tracking particles", 2);
0752 
0753   // Eta 3 (1.0 to 1.2)
0754   HistoName = "respt_eta1to1p2_pt2to3";
0755   respt_eta1to1p2_pt2to3 = iBooker.book1D(HistoName,
0756                                           HistoName,
0757                                           psRes_pt.getParameter<int32_t>("Nbinsx"),
0758                                           psRes_pt.getParameter<double>("xmin"),
0759                                           psRes_pt.getParameter<double>("xmax"));
0760   respt_eta1to1p2_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0761   respt_eta1to1p2_pt2to3->setAxisTitle("# tracking particles", 2);
0762 
0763   // Eta 4 (1.2 to 1.6)
0764   HistoName = "respt_eta1p2to1p6_pt2to3";
0765   respt_eta1p2to1p6_pt2to3 = iBooker.book1D(HistoName,
0766                                             HistoName,
0767                                             psRes_pt.getParameter<int32_t>("Nbinsx"),
0768                                             psRes_pt.getParameter<double>("xmin"),
0769                                             psRes_pt.getParameter<double>("xmax"));
0770   respt_eta1p2to1p6_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0771   respt_eta1p2to1p6_pt2to3->setAxisTitle("# tracking particles", 2);
0772 
0773   // Eta 5 (1.6 to 2.0)
0774   HistoName = "respt_eta1p6to2_pt2to3";
0775   respt_eta1p6to2_pt2to3 = iBooker.book1D(HistoName,
0776                                           HistoName,
0777                                           psRes_pt.getParameter<int32_t>("Nbinsx"),
0778                                           psRes_pt.getParameter<double>("xmin"),
0779                                           psRes_pt.getParameter<double>("xmax"));
0780   respt_eta1p6to2_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0781   respt_eta1p6to2_pt2to3->setAxisTitle("# tracking particles", 2);
0782 
0783   // Eta 6 (2.0 to 2.4)
0784   HistoName = "respt_eta2to2p4_pt2to3";
0785   respt_eta2to2p4_pt2to3 = iBooker.book1D(HistoName,
0786                                           HistoName,
0787                                           psRes_pt.getParameter<int32_t>("Nbinsx"),
0788                                           psRes_pt.getParameter<double>("xmin"),
0789                                           psRes_pt.getParameter<double>("xmax"));
0790   respt_eta2to2p4_pt2to3->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0791   respt_eta2to2p4_pt2to3->setAxisTitle("# tracking particles", 2);
0792 
0793   // pT b (3 to 8 GeV)
0794   // Eta 1 (0 to 0.7)
0795   HistoName = "respt_eta0to0p7_pt3to8";
0796   respt_eta0to0p7_pt3to8 = iBooker.book1D(HistoName,
0797                                           HistoName,
0798                                           psRes_pt.getParameter<int32_t>("Nbinsx"),
0799                                           psRes_pt.getParameter<double>("xmin"),
0800                                           psRes_pt.getParameter<double>("xmax"));
0801   respt_eta0to0p7_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0802   respt_eta0to0p7_pt3to8->setAxisTitle("# tracking particles", 2);
0803 
0804   // Eta 2 (0.7 to 1.0)
0805   HistoName = "respt_eta0p7to1_pt3to8";
0806   respt_eta0p7to1_pt3to8 = iBooker.book1D(HistoName,
0807                                           HistoName,
0808                                           psRes_pt.getParameter<int32_t>("Nbinsx"),
0809                                           psRes_pt.getParameter<double>("xmin"),
0810                                           psRes_pt.getParameter<double>("xmax"));
0811   respt_eta0p7to1_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0812   respt_eta0p7to1_pt3to8->setAxisTitle("# tracking particles", 2);
0813 
0814   // Eta 3 (1.0 to 1.2)
0815   HistoName = "respt_eta1to1p2_pt3to8";
0816   respt_eta1to1p2_pt3to8 = iBooker.book1D(HistoName,
0817                                           HistoName,
0818                                           psRes_pt.getParameter<int32_t>("Nbinsx"),
0819                                           psRes_pt.getParameter<double>("xmin"),
0820                                           psRes_pt.getParameter<double>("xmax"));
0821   respt_eta1to1p2_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0822   respt_eta1to1p2_pt3to8->setAxisTitle("# tracking particles", 2);
0823 
0824   // Eta 4 (1.2 to 1.6)
0825   HistoName = "respt_eta1p2to1p6_pt3to8";
0826   respt_eta1p2to1p6_pt3to8 = iBooker.book1D(HistoName,
0827                                             HistoName,
0828                                             psRes_pt.getParameter<int32_t>("Nbinsx"),
0829                                             psRes_pt.getParameter<double>("xmin"),
0830                                             psRes_pt.getParameter<double>("xmax"));
0831   respt_eta1p2to1p6_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0832   respt_eta1p2to1p6_pt3to8->setAxisTitle("# tracking particles", 2);
0833 
0834   // Eta 5 (1.6 to 2.0)
0835   HistoName = "respt_eta1p6to2_pt3to8";
0836   respt_eta1p6to2_pt3to8 = iBooker.book1D(HistoName,
0837                                           HistoName,
0838                                           psRes_pt.getParameter<int32_t>("Nbinsx"),
0839                                           psRes_pt.getParameter<double>("xmin"),
0840                                           psRes_pt.getParameter<double>("xmax"));
0841   respt_eta1p6to2_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0842   respt_eta1p6to2_pt3to8->setAxisTitle("# tracking particles", 2);
0843 
0844   // Eta 6 (2.0 to 2.4)
0845   HistoName = "respt_eta2to2p4_pt3to8";
0846   respt_eta2to2p4_pt3to8 = iBooker.book1D(HistoName,
0847                                           HistoName,
0848                                           psRes_pt.getParameter<int32_t>("Nbinsx"),
0849                                           psRes_pt.getParameter<double>("xmin"),
0850                                           psRes_pt.getParameter<double>("xmax"));
0851   respt_eta2to2p4_pt3to8->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0852   respt_eta2to2p4_pt3to8->setAxisTitle("# tracking particles", 2);
0853 
0854   // pT c (>8 GeV)
0855   // Eta 1 (0 to 0.7)
0856   HistoName = "respt_eta0to0p7_pt8toInf";
0857   respt_eta0to0p7_pt8toInf = iBooker.book1D(HistoName,
0858                                             HistoName,
0859                                             psRes_pt.getParameter<int32_t>("Nbinsx"),
0860                                             psRes_pt.getParameter<double>("xmin"),
0861                                             psRes_pt.getParameter<double>("xmax"));
0862   respt_eta0to0p7_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0863   respt_eta0to0p7_pt8toInf->setAxisTitle("# tracking particles", 2);
0864 
0865   // Eta 2 (0.7 to 1.0)
0866   HistoName = "respt_eta0p7to1_pt8toInf";
0867   respt_eta0p7to1_pt8toInf = iBooker.book1D(HistoName,
0868                                             HistoName,
0869                                             psRes_pt.getParameter<int32_t>("Nbinsx"),
0870                                             psRes_pt.getParameter<double>("xmin"),
0871                                             psRes_pt.getParameter<double>("xmax"));
0872   respt_eta0p7to1_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0873   respt_eta0p7to1_pt8toInf->setAxisTitle("# tracking particles", 2);
0874 
0875   // Eta 3 (1.0 to 1.2)
0876   HistoName = "respt_eta1to1p2_pt8toInf";
0877   respt_eta1to1p2_pt8toInf = iBooker.book1D(HistoName,
0878                                             HistoName,
0879                                             psRes_pt.getParameter<int32_t>("Nbinsx"),
0880                                             psRes_pt.getParameter<double>("xmin"),
0881                                             psRes_pt.getParameter<double>("xmax"));
0882   respt_eta1to1p2_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0883   respt_eta1to1p2_pt8toInf->setAxisTitle("# tracking particles", 2);
0884 
0885   // Eta 4 (1.2 to 1.6)
0886   HistoName = "respt_eta1p2to1p6_pt8toInf";
0887   respt_eta1p2to1p6_pt8toInf = iBooker.book1D(HistoName,
0888                                               HistoName,
0889                                               psRes_pt.getParameter<int32_t>("Nbinsx"),
0890                                               psRes_pt.getParameter<double>("xmin"),
0891                                               psRes_pt.getParameter<double>("xmax"));
0892   respt_eta1p2to1p6_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0893   respt_eta1p2to1p6_pt8toInf->setAxisTitle("# tracking particles", 2);
0894 
0895   // Eta 5 (1.6 to 2.0)
0896   HistoName = "respt_eta1p6to2_pt8toInf";
0897   respt_eta1p6to2_pt8toInf = iBooker.book1D(HistoName,
0898                                             HistoName,
0899                                             psRes_pt.getParameter<int32_t>("Nbinsx"),
0900                                             psRes_pt.getParameter<double>("xmin"),
0901                                             psRes_pt.getParameter<double>("xmax"));
0902   respt_eta1p6to2_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0903   respt_eta1p6to2_pt8toInf->setAxisTitle("# tracking particles", 2);
0904 
0905   // Eta 6 (2.0 to 2.4)
0906   HistoName = "respt_eta2to2p4_pt8toInf";
0907   respt_eta2to2p4_pt8toInf = iBooker.book1D(HistoName,
0908                                             HistoName,
0909                                             psRes_pt.getParameter<int32_t>("Nbinsx"),
0910                                             psRes_pt.getParameter<double>("xmin"),
0911                                             psRes_pt.getParameter<double>("xmax"));
0912   respt_eta2to2p4_pt8toInf->setAxisTitle("(p_{T}(trk) - p_{T}(tp))/p_{T}(tp)", 1);
0913   respt_eta2to2p4_pt8toInf->setAxisTitle("# tracking particles", 2);
0914 
0915   // Phi parts (for resolution)
0916   // Eta 1 (0 to 0.7)
0917   edm::ParameterSet psRes_phi = conf_.getParameter<edm::ParameterSet>("TH1Res_phi");
0918   HistoName = "resphi_eta0to0p7";
0919   resphi_eta0to0p7 = iBooker.book1D(HistoName,
0920                                     HistoName,
0921                                     psRes_phi.getParameter<int32_t>("Nbinsx"),
0922                                     psRes_phi.getParameter<double>("xmin"),
0923                                     psRes_phi.getParameter<double>("xmax"));
0924   resphi_eta0to0p7->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0925   resphi_eta0to0p7->setAxisTitle("# tracking particles", 2);
0926 
0927   // Eta 2 (0.7 to 1.0)
0928   HistoName = "resphi_eta0p7to1";
0929   resphi_eta0p7to1 = iBooker.book1D(HistoName,
0930                                     HistoName,
0931                                     psRes_phi.getParameter<int32_t>("Nbinsx"),
0932                                     psRes_phi.getParameter<double>("xmin"),
0933                                     psRes_phi.getParameter<double>("xmax"));
0934   resphi_eta0p7to1->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0935   resphi_eta0p7to1->setAxisTitle("# tracking particles", 2);
0936 
0937   // Eta 3 (1.0 to 1.2)
0938   HistoName = "resphi_eta1to1p2";
0939   resphi_eta1to1p2 = iBooker.book1D(HistoName,
0940                                     HistoName,
0941                                     psRes_phi.getParameter<int32_t>("Nbinsx"),
0942                                     psRes_phi.getParameter<double>("xmin"),
0943                                     psRes_phi.getParameter<double>("xmax"));
0944   resphi_eta1to1p2->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0945   resphi_eta1to1p2->setAxisTitle("# tracking particles", 2);
0946 
0947   // Eta 4 (1.2 to 1.6)
0948   HistoName = "resphi_eta1p2to1p6";
0949   resphi_eta1p2to1p6 = iBooker.book1D(HistoName,
0950                                       HistoName,
0951                                       psRes_phi.getParameter<int32_t>("Nbinsx"),
0952                                       psRes_phi.getParameter<double>("xmin"),
0953                                       psRes_phi.getParameter<double>("xmax"));
0954   resphi_eta1p2to1p6->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0955   resphi_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
0956 
0957   // Eta 5 (1.6 to 2.0)
0958   HistoName = "resphi_eta1p6to2";
0959   resphi_eta1p6to2 = iBooker.book1D(HistoName,
0960                                     HistoName,
0961                                     psRes_phi.getParameter<int32_t>("Nbinsx"),
0962                                     psRes_phi.getParameter<double>("xmin"),
0963                                     psRes_phi.getParameter<double>("xmax"));
0964   resphi_eta1p6to2->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0965   resphi_eta1p6to2->setAxisTitle("# tracking particles", 2);
0966 
0967   // Eta 6 (2.0 to 2.4)
0968   HistoName = "resphi_eta2to2p4";
0969   resphi_eta2to2p4 = iBooker.book1D(HistoName,
0970                                     HistoName,
0971                                     psRes_phi.getParameter<int32_t>("Nbinsx"),
0972                                     psRes_phi.getParameter<double>("xmin"),
0973                                     psRes_phi.getParameter<double>("xmax"));
0974   resphi_eta2to2p4->setAxisTitle("#phi_{trk} - #phi_{tp}", 1);
0975   resphi_eta2to2p4->setAxisTitle("# tracking particles", 2);
0976 
0977   // VtxZ parts (for resolution)
0978   // Eta 1 (0 to 0.7)
0979   edm::ParameterSet psRes_VtxZ = conf_.getParameter<edm::ParameterSet>("TH1Res_VtxZ");
0980   HistoName = "resVtxZ_eta0to0p7";
0981   resVtxZ_eta0to0p7 = iBooker.book1D(HistoName,
0982                                      HistoName,
0983                                      psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
0984                                      psRes_VtxZ.getParameter<double>("xmin"),
0985                                      psRes_VtxZ.getParameter<double>("xmax"));
0986   resVtxZ_eta0to0p7->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
0987   resVtxZ_eta0to0p7->setAxisTitle("# tracking particles", 2);
0988 
0989   // Eta 2 (0.7 to 1.0)
0990   HistoName = "resVtxZ_eta0p7to1";
0991   resVtxZ_eta0p7to1 = iBooker.book1D(HistoName,
0992                                      HistoName,
0993                                      psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
0994                                      psRes_VtxZ.getParameter<double>("xmin"),
0995                                      psRes_VtxZ.getParameter<double>("xmax"));
0996   resVtxZ_eta0p7to1->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
0997   resVtxZ_eta0p7to1->setAxisTitle("# tracking particles", 2);
0998 
0999   // Eta 3 (1.0 to 1.2)
1000   HistoName = "resVtxZ_eta1to1p2";
1001   resVtxZ_eta1to1p2 = iBooker.book1D(HistoName,
1002                                      HistoName,
1003                                      psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
1004                                      psRes_VtxZ.getParameter<double>("xmin"),
1005                                      psRes_VtxZ.getParameter<double>("xmax"));
1006   resVtxZ_eta1to1p2->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
1007   resVtxZ_eta1to1p2->setAxisTitle("# tracking particles", 2);
1008 
1009   // Eta 4 (1.2 to 1.6)
1010   HistoName = "resVtxZ_eta1p2to1p6";
1011   resVtxZ_eta1p2to1p6 = iBooker.book1D(HistoName,
1012                                        HistoName,
1013                                        psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
1014                                        psRes_VtxZ.getParameter<double>("xmin"),
1015                                        psRes_VtxZ.getParameter<double>("xmax"));
1016   resVtxZ_eta1p2to1p6->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
1017   resVtxZ_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
1018 
1019   // Eta 5 (1.6 to 2.0)
1020   HistoName = "resVtxZ_eta1p6to2";
1021   resVtxZ_eta1p6to2 = iBooker.book1D(HistoName,
1022                                      HistoName,
1023                                      psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
1024                                      psRes_VtxZ.getParameter<double>("xmin"),
1025                                      psRes_VtxZ.getParameter<double>("xmax"));
1026   resVtxZ_eta1p6to2->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
1027   resVtxZ_eta1p6to2->setAxisTitle("# tracking particles", 2);
1028 
1029   // Eta 6 (2.0 to 2.4)
1030   HistoName = "resVtxZ_eta2to2p4";
1031   resVtxZ_eta2to2p4 = iBooker.book1D(HistoName,
1032                                      HistoName,
1033                                      psRes_VtxZ.getParameter<int32_t>("Nbinsx"),
1034                                      psRes_VtxZ.getParameter<double>("xmin"),
1035                                      psRes_VtxZ.getParameter<double>("xmax"));
1036   resVtxZ_eta2to2p4->setAxisTitle("VtxZ_{trk} - VtxZ_{tp} [cm]", 1);
1037   resVtxZ_eta2to2p4->setAxisTitle("# tracking particles", 2);
1038 
1039   // d0 parts (for resolution)
1040   // Eta 1 (0 to 0.7)
1041   edm::ParameterSet psRes_d0 = conf_.getParameter<edm::ParameterSet>("TH1Res_d0");
1042   HistoName = "resd0_eta0to0p7";
1043   resd0_eta0to0p7 = iBooker.book1D(HistoName,
1044                                    HistoName,
1045                                    psRes_d0.getParameter<int32_t>("Nbinsx"),
1046                                    psRes_d0.getParameter<double>("xmin"),
1047                                    psRes_d0.getParameter<double>("xmax"));
1048   resd0_eta0to0p7->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1049   resd0_eta0to0p7->setAxisTitle("# tracking particles", 2);
1050 
1051   // Eta 2 (0.7 to 1.0)
1052   HistoName = "resd0_eta0p7to1";
1053   resd0_eta0p7to1 = iBooker.book1D(HistoName,
1054                                    HistoName,
1055                                    psRes_d0.getParameter<int32_t>("Nbinsx"),
1056                                    psRes_d0.getParameter<double>("xmin"),
1057                                    psRes_d0.getParameter<double>("xmax"));
1058   resd0_eta0p7to1->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1059   resd0_eta0p7to1->setAxisTitle("# tracking particles", 2);
1060 
1061   // Eta 3 (1.0 to 1.2)
1062   HistoName = "resd0_eta1to1p2";
1063   resd0_eta1to1p2 = iBooker.book1D(HistoName,
1064                                    HistoName,
1065                                    psRes_d0.getParameter<int32_t>("Nbinsx"),
1066                                    psRes_d0.getParameter<double>("xmin"),
1067                                    psRes_d0.getParameter<double>("xmax"));
1068   resd0_eta1to1p2->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1069   resd0_eta1to1p2->setAxisTitle("# tracking particles", 2);
1070 
1071   // Eta 4 (1.2 to 1.6)
1072   HistoName = "resd0_eta1p2to1p6";
1073   resd0_eta1p2to1p6 = iBooker.book1D(HistoName,
1074                                      HistoName,
1075                                      psRes_d0.getParameter<int32_t>("Nbinsx"),
1076                                      psRes_d0.getParameter<double>("xmin"),
1077                                      psRes_d0.getParameter<double>("xmax"));
1078   resd0_eta1p2to1p6->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1079   resd0_eta1p2to1p6->setAxisTitle("# tracking particles", 2);
1080 
1081   // Eta 5 (1.6 to 2.0)
1082   HistoName = "resd0_eta1p6to2";
1083   resd0_eta1p6to2 = iBooker.book1D(HistoName,
1084                                    HistoName,
1085                                    psRes_d0.getParameter<int32_t>("Nbinsx"),
1086                                    psRes_d0.getParameter<double>("xmin"),
1087                                    psRes_d0.getParameter<double>("xmax"));
1088   resd0_eta1p6to2->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1089   resd0_eta1p6to2->setAxisTitle("# tracking particles", 2);
1090 
1091   // Eta 6 (2.0 to 2.4)
1092   HistoName = "resd0_eta2to2p4";
1093   resd0_eta2to2p4 = iBooker.book1D(HistoName,
1094                                    HistoName,
1095                                    psRes_d0.getParameter<int32_t>("Nbinsx"),
1096                                    psRes_d0.getParameter<double>("xmin"),
1097                                    psRes_d0.getParameter<double>("xmax"));
1098   resd0_eta2to2p4->setAxisTitle("d0_{trk} - d0_{tp} [cm]", 1);
1099   resd0_eta2to2p4->setAxisTitle("# tracking particles", 2);
1100 
1101 }  // end of method
1102 
1103 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
1104 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
1105 void Phase2OTValidateTrackingParticles::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
1106   // OuterTrackerMonitorTrackingParticles
1107   edm::ParameterSetDescription desc;
1108   {
1109     edm::ParameterSetDescription psd0;
1110     psd0.add<int>("Nbinsx", 45);
1111     psd0.add<double>("xmax", 3);
1112     psd0.add<double>("xmin", -3);
1113     desc.add<edm::ParameterSetDescription>("TH1TrackParts_Eta", psd0);
1114   }
1115   {
1116     edm::ParameterSetDescription psd0;
1117     psd0.add<int>("Nbinsx", 60);
1118     psd0.add<double>("xmax", 3.141592653589793);
1119     psd0.add<double>("xmin", -3.141592653589793);
1120     desc.add<edm::ParameterSetDescription>("TH1TrackParts_Phi", psd0);
1121   }
1122   {
1123     edm::ParameterSetDescription psd0;
1124     psd0.add<int>("Nbinsx", 45);
1125     psd0.add<double>("xmax", 100);
1126     psd0.add<double>("xmin", 0);
1127     desc.add<edm::ParameterSetDescription>("TH1TrackParts_Pt", psd0);
1128   }
1129   {
1130     edm::ParameterSetDescription psd0;
1131     psd0.add<int>("Nbinsx", 200);
1132     psd0.add<double>("xmax", 0.5);
1133     psd0.add<double>("xmin", -0.5);
1134     desc.add<edm::ParameterSetDescription>("TH1Res_ptRel", psd0);
1135   }
1136   {
1137     edm::ParameterSetDescription psd0;
1138     psd0.add<int>("Nbinsx", 50);
1139     psd0.add<double>("xmax", 100);
1140     psd0.add<double>("xmin", 0);
1141     desc.add<edm::ParameterSetDescription>("TH1Effic_pt", psd0);
1142   }
1143   {
1144     edm::ParameterSetDescription psd0;
1145     psd0.add<int>("Nbinsx", 50);
1146     psd0.add<double>("xmax", 10);
1147     psd0.add<double>("xmin", 0);
1148     desc.add<edm::ParameterSetDescription>("TH1Effic_pt_zoom", psd0);
1149   }
1150   {
1151     edm::ParameterSetDescription psd0;
1152     psd0.add<int>("Nbinsx", 50);
1153     psd0.add<double>("xmax", 2.5);
1154     psd0.add<double>("xmin", -2.5);
1155     desc.add<edm::ParameterSetDescription>("TH1Effic_eta", psd0);
1156   }
1157   {
1158     edm::ParameterSetDescription psd0;
1159     psd0.add<int>("Nbinsx", 50);
1160     psd0.add<double>("xmax", 2);
1161     psd0.add<double>("xmin", -2);
1162     desc.add<edm::ParameterSetDescription>("TH1Effic_d0", psd0);
1163   }
1164   {
1165     edm::ParameterSetDescription psd0;
1166     psd0.add<int>("Nbinsx", 50);
1167     psd0.add<double>("xmax", 5);
1168     psd0.add<double>("xmin", -5);
1169     desc.add<edm::ParameterSetDescription>("TH1Effic_VtxR", psd0);
1170   }
1171   {
1172     edm::ParameterSetDescription psd0;
1173     psd0.add<int>("Nbinsx", 50);
1174     psd0.add<double>("xmax", 30);
1175     psd0.add<double>("xmin", -30);
1176     desc.add<edm::ParameterSetDescription>("TH1Effic_VtxZ", psd0);
1177   }
1178   {
1179     edm::ParameterSetDescription psd0;
1180     psd0.add<int>("Nbinsx", 100);
1181     psd0.add<double>("xmax", 0.2);
1182     psd0.add<double>("xmin", -0.2);
1183     desc.add<edm::ParameterSetDescription>("TH1Res_pt", psd0);
1184   }
1185   {
1186     edm::ParameterSetDescription psd0;
1187     psd0.add<int>("Nbinsx", 100);
1188     psd0.add<double>("xmax", 0.01);
1189     psd0.add<double>("xmin", -0.01);
1190     desc.add<edm::ParameterSetDescription>("TH1Res_eta", psd0);
1191   }
1192   {
1193     edm::ParameterSetDescription psd0;
1194     psd0.add<int>("Nbinsx", 100);
1195     psd0.add<double>("xmax", 0.01);
1196     psd0.add<double>("xmin", -0.01);
1197     desc.add<edm::ParameterSetDescription>("TH1Res_phi", psd0);
1198   }
1199   {
1200     edm::ParameterSetDescription psd0;
1201     psd0.add<int>("Nbinsx", 100);
1202     psd0.add<double>("xmax", 1.0);
1203     psd0.add<double>("xmin", -1.0);
1204     desc.add<edm::ParameterSetDescription>("TH1Res_VtxZ", psd0);
1205   }
1206   {
1207     edm::ParameterSetDescription psd0;
1208     psd0.add<int>("Nbinsx", 100);
1209     psd0.add<double>("xmax", 0.05);
1210     psd0.add<double>("xmin", -0.05);
1211     desc.add<edm::ParameterSetDescription>("TH1Res_d0", psd0);
1212   }
1213   desc.add<std::string>("TopFolderName", "TrackerPhase2OTL1TrackV");
1214   desc.add<edm::InputTag>("trackingParticleToken", edm::InputTag("mix", "MergedTrackTruth"));
1215   desc.add<edm::InputTag>("MCTruthStubInputTag", edm::InputTag("TTStubAssociatorFromPixelDigis", "StubAccepted"));
1216   desc.add<edm::InputTag>("MCTruthTrackInputTag", edm::InputTag("TTTrackAssociatorFromPixelDigis", "Level1TTTracks"));
1217   desc.add<edm::InputTag>("MCTruthClusterInputTag",
1218                           edm::InputTag("TTClusterAssociatorFromPixelDigis", "ClusterAccepted"));
1219   desc.add<int>("L1Tk_minNStub", 4);
1220   desc.add<double>("L1Tk_maxChi2dof", 25.0);
1221   desc.add<int>("TP_minNStub", 4);
1222   desc.add<int>("TP_minNLayersStub", 4);
1223   desc.add<double>("TP_minPt", 2.0);
1224   desc.add<double>("TP_maxEta", 2.4);
1225   desc.add<double>("TP_maxVtxZ", 15.0);
1226   descriptions.add("Phase2OTValidateTrackingParticles", desc);
1227   // or use the following to generate the label from the module's C++ type
1228   //descriptions.addWithDefaultLabel(desc);
1229 }
1230 
1231 DEFINE_FWK_MODULE(Phase2OTValidateTrackingParticles);