Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 09:43:52

0001 /*
0002  *  See header file for a description of this class.
0003  *
0004  *  \author Suchandra Dutta , Giorgia Mila
0005  */
0006 
0007 #include "DQM/TrackingMonitor/interface/TrackSplittingMonitor.h"
0008 #include "DQMServices/Core/interface/DQMStore.h"
0009 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0010 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h"
0011 #include "DataFormats/TrackReco/interface/Track.h"
0012 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0013 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
0017 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
0018 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
0019 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0020 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
0021 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0022 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0023 
0024 #include <string>
0025 
0026 TrackSplittingMonitor::TrackSplittingMonitor(const edm::ParameterSet& iConfig)
0027     : conf_(iConfig),
0028       mfToken_(esConsumes()),
0029       tkGeomToken_(esConsumes()),
0030       dtGeomToken_(esConsumes()),
0031       cscGeomToken_(esConsumes()),
0032       rpcGeomToken_(esConsumes()),
0033       splitTracksToken_(consumes<std::vector<reco::Track> >(conf_.getParameter<edm::InputTag>("splitTrackCollection"))),
0034       splitMuonsToken_(mayConsume<std::vector<reco::Muon> >(conf_.getParameter<edm::InputTag>("splitMuonCollection"))),
0035       plotMuons_(conf_.getParameter<bool>("ifPlotMuons")),
0036       pixelHitsPerLeg_(conf_.getParameter<int>("pixelHitsPerLeg")),
0037       totalHitsPerLeg_(conf_.getParameter<int>("totalHitsPerLeg")),
0038       d0Cut_(conf_.getParameter<double>("d0Cut")),
0039       dzCut_(conf_.getParameter<double>("dzCut")),
0040       ptCut_(conf_.getParameter<double>("ptCut")),
0041       norchiCut_(conf_.getParameter<double>("norchiCut")) {}
0042 
0043 void TrackSplittingMonitor::bookHistograms(DQMStore::IBooker& ibooker,
0044                                            edm::Run const& /* iRun */,
0045                                            edm::EventSetup const& /* iSetup */) {
0046   std::string MEFolderName = conf_.getParameter<std::string>("FolderName");
0047   ibooker.setCurrentFolder(MEFolderName);
0048 
0049   // bin declarations
0050   int ddxyBin = conf_.getParameter<int>("ddxyBin");
0051   double ddxyMin = conf_.getParameter<double>("ddxyMin");
0052   double ddxyMax = conf_.getParameter<double>("ddxyMax");
0053 
0054   int ddzBin = conf_.getParameter<int>("ddzBin");
0055   double ddzMin = conf_.getParameter<double>("ddzMin");
0056   double ddzMax = conf_.getParameter<double>("ddzMax");
0057 
0058   int dphiBin = conf_.getParameter<int>("dphiBin");
0059   double dphiMin = conf_.getParameter<double>("dphiMin");
0060   double dphiMax = conf_.getParameter<double>("dphiMax");
0061 
0062   int dthetaBin = conf_.getParameter<int>("dthetaBin");
0063   double dthetaMin = conf_.getParameter<double>("dthetaMin");
0064   double dthetaMax = conf_.getParameter<double>("dthetaMax");
0065 
0066   int dptBin = conf_.getParameter<int>("dptBin");
0067   double dptMin = conf_.getParameter<double>("dptMin");
0068   double dptMax = conf_.getParameter<double>("dptMax");
0069 
0070   int dcurvBin = conf_.getParameter<int>("dcurvBin");
0071   double dcurvMin = conf_.getParameter<double>("dcurvMin");
0072   double dcurvMax = conf_.getParameter<double>("dcurvMax");
0073 
0074   int normBin = conf_.getParameter<int>("normBin");
0075   double normMin = conf_.getParameter<double>("normMin");
0076   double normMax = conf_.getParameter<double>("normMax");
0077 
0078   // declare histogram
0079   ddxyAbsoluteResiduals_tracker_ =
0080       ibooker.book1D("ddxyAbsoluteResiduals_tracker", "ddxyAbsoluteResiduals_tracker", ddxyBin, ddxyMin, ddxyMax);
0081   ddzAbsoluteResiduals_tracker_ =
0082       ibooker.book1D("ddzAbsoluteResiduals_tracker", "ddzAbsoluteResiduals_tracker", ddzBin, ddzMin, ddzMax);
0083   dphiAbsoluteResiduals_tracker_ =
0084       ibooker.book1D("dphiAbsoluteResiduals_tracker", "dphiAbsoluteResiduals_tracker", dphiBin, dphiMin, dphiMax);
0085   dthetaAbsoluteResiduals_tracker_ = ibooker.book1D(
0086       "dthetaAbsoluteResiduals_tracker", "dthetaAbsoluteResiduals_tracker", dthetaBin, dthetaMin, dthetaMax);
0087   dptAbsoluteResiduals_tracker_ =
0088       ibooker.book1D("dptAbsoluteResiduals_tracker", "dptAbsoluteResiduals_tracker", dptBin, dptMin, dptMax);
0089   dcurvAbsoluteResiduals_tracker_ =
0090       ibooker.book1D("dcurvAbsoluteResiduals_tracker", "dcurvAbsoluteResiduals_tracker", dcurvBin, dcurvMin, dcurvMax);
0091 
0092   ddxyNormalizedResiduals_tracker_ =
0093       ibooker.book1D("ddxyNormalizedResiduals_tracker", "ddxyNormalizedResiduals_tracker", normBin, normMin, normMax);
0094   ddzNormalizedResiduals_tracker_ =
0095       ibooker.book1D("ddzNormalizedResiduals_tracker", "ddzNormalizedResiduals_tracker", normBin, normMin, normMax);
0096   dphiNormalizedResiduals_tracker_ =
0097       ibooker.book1D("dphiNormalizedResiduals_tracker", "dphiNormalizedResiduals_tracker", normBin, normMin, normMax);
0098   dthetaNormalizedResiduals_tracker_ = ibooker.book1D(
0099       "dthetaNormalizedResiduals_tracker", "dthetaNormalizedResiduals_tracker", normBin, normMin, normMax);
0100   dptNormalizedResiduals_tracker_ =
0101       ibooker.book1D("dptNormalizedResiduals_tracker", "dptNormalizedResiduals_tracker", normBin, normMin, normMax);
0102   dcurvNormalizedResiduals_tracker_ =
0103       ibooker.book1D("dcurvNormalizedResiduals_tracker", "dcurvNormalizedResiduals_tracker", normBin, normMin, normMax);
0104 
0105   if (plotMuons_) {
0106     ddxyAbsoluteResiduals_global_ =
0107         ibooker.book1D("ddxyAbsoluteResiduals_global", "ddxyAbsoluteResiduals_global", ddxyBin, ddxyMin, ddxyMax);
0108     ddzAbsoluteResiduals_global_ =
0109         ibooker.book1D("ddzAbsoluteResiduals_global", "ddzAbsoluteResiduals_global", ddzBin, ddzMin, ddzMax);
0110     dphiAbsoluteResiduals_global_ =
0111         ibooker.book1D("dphiAbsoluteResiduals_global", "dphiAbsoluteResiduals_global", dphiBin, dphiMin, dphiMax);
0112     dthetaAbsoluteResiduals_global_ = ibooker.book1D(
0113         "dthetaAbsoluteResiduals_global", "dthetaAbsoluteResiduals_global", dthetaBin, dthetaMin, dthetaMax);
0114     dptAbsoluteResiduals_global_ =
0115         ibooker.book1D("dptAbsoluteResiduals_global", "dptAbsoluteResiduals_global", dptBin, dptMin, dptMax);
0116     dcurvAbsoluteResiduals_global_ =
0117         ibooker.book1D("dcurvAbsoluteResiduals_global", "dcurvAbsoluteResiduals_global", dcurvBin, dcurvMin, dcurvMax);
0118 
0119     ddxyNormalizedResiduals_global_ =
0120         ibooker.book1D("ddxyNormalizedResiduals_global", "ddxyNormalizedResiduals_global", normBin, normMin, normMax);
0121     ddzNormalizedResiduals_global_ =
0122         ibooker.book1D("ddzNormalizedResiduals_global", "ddzNormalizedResiduals_global", normBin, normMin, normMax);
0123     dphiNormalizedResiduals_global_ =
0124         ibooker.book1D("dphiNormalizedResiduals_global", "dphiNormalizedResiduals_global", normBin, normMin, normMax);
0125     dthetaNormalizedResiduals_global_ = ibooker.book1D(
0126         "dthetaNormalizedResiduals_global", "dthetaNormalizedResiduals_global", normBin, normMin, normMax);
0127     dptNormalizedResiduals_global_ =
0128         ibooker.book1D("dptNormalizedResiduals_global", "dptNormalizedResiduals_global", normBin, normMin, normMax);
0129     dcurvNormalizedResiduals_global_ =
0130         ibooker.book1D("dcurvNormalizedResiduals_global", "dcurvNormalizedResiduals_global", normBin, normMin, normMax);
0131   }
0132 
0133   ddxyAbsoluteResiduals_tracker_->setAxisTitle("(#delta d_{xy})/#sqrt{2} [#mum]");
0134   ddzAbsoluteResiduals_tracker_->setAxisTitle("(#delta d_{z})/#sqrt{2} [#mum]");
0135   dphiAbsoluteResiduals_tracker_->setAxisTitle("(#delta #phi)/#sqrt{2} [mrad]");
0136   dthetaAbsoluteResiduals_tracker_->setAxisTitle("(#delta #theta)/#sqrt{2} [mrad]");
0137   dptAbsoluteResiduals_tracker_->setAxisTitle("(#delta p_{T})/#sqrt{2} [GeV]");
0138   dcurvAbsoluteResiduals_tracker_->setAxisTitle("(#delta (1/p_{T}))/#sqrt{2} [GeV^{-1}]");
0139 
0140   ddxyNormalizedResiduals_tracker_->setAxisTitle("#delta d_{xy}/#sigma(d_{xy})");
0141   ddzNormalizedResiduals_tracker_->setAxisTitle("#delta d_{z}/#sigma(d_{z})");
0142   dphiNormalizedResiduals_tracker_->setAxisTitle("#delta #phi/#sigma(d_{#phi})");
0143   dthetaNormalizedResiduals_tracker_->setAxisTitle("#delta #theta/#sigma(d_{#theta})");
0144   dptNormalizedResiduals_tracker_->setAxisTitle("#delta p_{T}/#sigma(p_{T})");
0145   dcurvNormalizedResiduals_tracker_->setAxisTitle("#delta 1/p_{T}/#sigma(1/p_{T})");
0146 
0147   if (plotMuons_) {
0148     ddxyAbsoluteResiduals_global_->setAxisTitle("(#delta d_{xy})/#sqrt{2} [#mum]");
0149     ddzAbsoluteResiduals_global_->setAxisTitle("(#delta d_{z})/#sqrt{2} [#mum]");
0150     dphiAbsoluteResiduals_global_->setAxisTitle("(#delta #phi)/#sqrt{2} [mrad]");
0151     dthetaAbsoluteResiduals_global_->setAxisTitle("(#delta #theta)/#sqrt{2} [mrad]");
0152     dptAbsoluteResiduals_global_->setAxisTitle("(#delta p_{T})/#sqrt{2} [GeV]");
0153     dcurvAbsoluteResiduals_global_->setAxisTitle("(#delta (1/p_{T}))/#sqrt{2} [GeV^{-1}]");
0154 
0155     ddxyNormalizedResiduals_global_->setAxisTitle("#delta d_{xy}/#sigma(d_{xy})");
0156     ddzNormalizedResiduals_global_->setAxisTitle("#delta d_{z}/#sigma(d_{z})");
0157     dphiNormalizedResiduals_global_->setAxisTitle("#delta #phi/#sigma(d_{#phi})");
0158     dthetaNormalizedResiduals_global_->setAxisTitle("#delta #theta/#sigma(d_{#theta})");
0159     dptNormalizedResiduals_global_->setAxisTitle("#delta p_{T}/#sigma(p_{T})");
0160     dcurvNormalizedResiduals_global_->setAxisTitle("#delta 1/p_{T}/#sigma(1/p_{T})");
0161   }
0162 }
0163 
0164 //
0165 // -- Analyse
0166 //
0167 void TrackSplittingMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0168   theMagField = &iSetup.getData(mfToken_);
0169   theGeometry = &iSetup.getData(tkGeomToken_);
0170   dtGeometry = &iSetup.getData(dtGeomToken_);
0171   cscGeometry = &iSetup.getData(cscGeomToken_);
0172   rpcGeometry = &iSetup.getData(rpcGeomToken_);
0173 
0174   edm::Handle<std::vector<reco::Track> > splitTracks = iEvent.getHandle(splitTracksToken_);
0175   if (!splitTracks.isValid())
0176     return;
0177 
0178   edm::Handle<std::vector<reco::Muon> > splitMuons;
0179   if (plotMuons_) {
0180     splitMuons = iEvent.getHandle(splitMuonsToken_);
0181   }
0182 
0183   if (splitTracks->size() == 2) {
0184     // check that there are 2 tracks in split track collection
0185     edm::LogInfo("TrackSplittingMonitor") << "Split Track size: " << splitTracks->size();
0186 
0187     // split tracks calculations
0188     reco::Track track1 = splitTracks->at(0);
0189     reco::Track track2 = splitTracks->at(1);
0190 
0191     // -------------------------- basic selection ---------------------------
0192 
0193     // hit counting
0194     // looping through the hits for track 1
0195     double nRechits1 = 0;
0196     double nRechitinBPIX1 = 0;
0197     for (auto const& iHit : track1.recHits()) {
0198       if (iHit->isValid()) {
0199         nRechits1++;
0200         int type = iHit->geographicalId().subdetId();
0201         if (type == int(PixelSubdetector::PixelBarrel)) {
0202           ++nRechitinBPIX1;
0203         }
0204       }
0205     }
0206     // looping through the hits for track 2
0207     double nRechits2 = 0;
0208     double nRechitinBPIX2 = 0;
0209     for (auto const& iHit : track2.recHits()) {
0210       if (iHit->isValid()) {
0211         nRechits2++;
0212         int type = iHit->geographicalId().subdetId();
0213         if (type == int(PixelSubdetector::PixelBarrel)) {
0214           ++nRechitinBPIX2;
0215         }
0216       }
0217     }
0218 
0219     // DCA of each track
0220     double d01 = track1.d0();
0221     double dz1 = track1.dz();
0222     double d02 = track2.d0();
0223     double dz2 = track2.dz();
0224 
0225     // pT of each track
0226     double pt1 = track1.pt();
0227     double pt2 = track2.pt();
0228 
0229     // chi2 of each track
0230     double norchi1 = track1.normalizedChi2();
0231     double norchi2 = track2.normalizedChi2();
0232 
0233     // basic selection
0234     // pixel hits and total hits
0235     if ((nRechitinBPIX1 >= pixelHitsPerLeg_) && (nRechitinBPIX2 >= pixelHitsPerLeg_) &&
0236         (nRechits1 >= totalHitsPerLeg_) && (nRechits2 >= totalHitsPerLeg_)) {
0237       // dca cut
0238       if (((std::abs(d01) < d0Cut_)) && (std::abs(d02) < d0Cut_) && (std::abs(dz1) < dzCut_) &&
0239           (std::abs(dz2) < dzCut_)) {
0240         // pt cut
0241         if ((pt1 + pt2) / 2 < ptCut_) {
0242           // chi2 cut
0243           if ((norchi1 < norchiCut_) && (norchi2 < norchiCut_)) {
0244             // passed all cuts...
0245             edm::LogInfo("TrackSplittingMonitor") << " Setected after all cuts ?";
0246 
0247             double ddxyVal = d01 - d02;
0248             double ddzVal = dz1 - dz2;
0249             double dphiVal = track1.phi() - track2.phi();
0250             double dthetaVal = track1.theta() - track2.theta();
0251             double dptVal = pt1 - pt2;
0252             double dcurvVal = (1 / pt1) - (1 / pt2);
0253 
0254             double d01ErrVal = track1.d0Error();
0255             double d02ErrVal = track2.d0Error();
0256             double dz1ErrVal = track1.dzError();
0257             double dz2ErrVal = track2.dzError();
0258             double phi1ErrVal = track1.phiError();
0259             double phi2ErrVal = track2.phiError();
0260             double theta1ErrVal = track1.thetaError();
0261             double theta2ErrVal = track2.thetaError();
0262             double pt1ErrVal = track1.ptError();
0263             double pt2ErrVal = track2.ptError();
0264 
0265             ddxyAbsoluteResiduals_tracker_->Fill(cmToUm * ddxyVal / sqrt2);
0266             ddzAbsoluteResiduals_tracker_->Fill(cmToUm * ddzVal / sqrt2);
0267             dphiAbsoluteResiduals_tracker_->Fill(radToUrad * dphiVal / sqrt2);
0268             dthetaAbsoluteResiduals_tracker_->Fill(radToUrad * dthetaVal / sqrt2);
0269             dptAbsoluteResiduals_tracker_->Fill(dptVal / sqrt2);
0270             dcurvAbsoluteResiduals_tracker_->Fill(dcurvVal / sqrt2);
0271 
0272             ddxyNormalizedResiduals_tracker_->Fill(ddxyVal / sqrt(d01ErrVal * d01ErrVal + d02ErrVal * d02ErrVal));
0273             ddzNormalizedResiduals_tracker_->Fill(ddzVal / sqrt(dz1ErrVal * dz1ErrVal + dz2ErrVal * dz2ErrVal));
0274             dphiNormalizedResiduals_tracker_->Fill(dphiVal / sqrt(phi1ErrVal * phi1ErrVal + phi2ErrVal * phi2ErrVal));
0275             dthetaNormalizedResiduals_tracker_->Fill(dthetaVal /
0276                                                      sqrt(theta1ErrVal * theta1ErrVal + theta2ErrVal * theta2ErrVal));
0277             dptNormalizedResiduals_tracker_->Fill(dptVal / sqrt(pt1ErrVal * pt1ErrVal + pt2ErrVal * pt2ErrVal));
0278             dcurvNormalizedResiduals_tracker_->Fill(
0279                 dcurvVal / sqrt(pow(pt1ErrVal, 2) / pow(pt1, 4) + pow(pt2ErrVal, 2) / pow(pt2, 4)));
0280 
0281             // if do the same for split muons
0282             if (plotMuons_ && splitMuons.isValid()) {
0283               int gmCtr = 0;
0284               bool topGlobalMuonFlag = false;
0285               bool bottomGlobalMuonFlag = false;
0286               int topGlobalMuon = -1;
0287               int bottomGlobalMuon = -1;
0288               double topGlobalMuonNorchi2 = 1e10;
0289               double bottomGlobalMuonNorchi2 = 1e10;
0290 
0291               // check if usable split global muons
0292               for (std::vector<reco::Muon>::const_iterator gmI = splitMuons->begin(); gmI != splitMuons->end(); gmI++) {
0293                 if (gmI->isTrackerMuon() && gmI->isStandAloneMuon() && gmI->isGlobalMuon()) {
0294                   reco::TrackRef trackerTrackRef1(splitTracks, 0);
0295                   reco::TrackRef trackerTrackRef2(splitTracks, 1);
0296 
0297                   if (gmI->innerTrack() == trackerTrackRef1) {
0298                     if (gmI->globalTrack()->normalizedChi2() < topGlobalMuonNorchi2) {
0299                       topGlobalMuonFlag = true;
0300                       topGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
0301                       topGlobalMuon = gmCtr;
0302                     }
0303                   }
0304                   if (gmI->innerTrack() == trackerTrackRef2) {
0305                     if (gmI->globalTrack()->normalizedChi2() < bottomGlobalMuonNorchi2) {
0306                       bottomGlobalMuonFlag = true;
0307                       bottomGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
0308                       bottomGlobalMuon = gmCtr;
0309                     }
0310                   }
0311                 }
0312                 gmCtr++;
0313               }
0314 
0315               if (bottomGlobalMuonFlag && topGlobalMuonFlag) {
0316                 reco::Muon muonTop = splitMuons->at(topGlobalMuon);
0317                 reco::Muon muonBottom = splitMuons->at(bottomGlobalMuon);
0318 
0319                 reco::TrackRef glb1 = muonTop.globalTrack();
0320                 reco::TrackRef glb2 = muonBottom.globalTrack();
0321 
0322                 double ddxyValGlb = glb1->d0() - glb2->d0();
0323                 double ddzValGlb = glb1->dz() - glb2->dz();
0324                 double dphiValGlb = glb1->phi() - glb2->phi();
0325                 double dthetaValGlb = glb1->theta() - glb2->theta();
0326                 double dptValGlb = glb1->pt() - glb2->pt();
0327                 double dcurvValGlb = (1 / glb1->pt()) - (1 / glb2->pt());
0328 
0329                 double d01ErrValGlb = glb1->d0Error();
0330                 double d02ErrValGlb = glb2->d0Error();
0331                 double dz1ErrValGlb = glb1->dzError();
0332                 double dz2ErrValGlb = glb2->dzError();
0333                 double phi1ErrValGlb = glb1->phiError();
0334                 double phi2ErrValGlb = glb2->phiError();
0335                 double theta1ErrValGlb = glb1->thetaError();
0336                 double theta2ErrValGlb = glb2->thetaError();
0337                 double pt1ErrValGlb = glb1->ptError();
0338                 double pt2ErrValGlb = glb2->ptError();
0339 
0340                 ddxyAbsoluteResiduals_global_->Fill(cmToUm * ddxyValGlb / sqrt2);
0341                 ddzAbsoluteResiduals_global_->Fill(cmToUm * ddzValGlb / sqrt2);
0342                 dphiAbsoluteResiduals_global_->Fill(radToUrad * dphiValGlb / sqrt2);
0343                 dthetaAbsoluteResiduals_global_->Fill(radToUrad * dthetaValGlb / sqrt2);
0344                 dptAbsoluteResiduals_global_->Fill(dptValGlb / sqrt2);
0345                 dcurvAbsoluteResiduals_global_->Fill(dcurvValGlb / sqrt2);
0346 
0347                 ddxyNormalizedResiduals_global_->Fill(ddxyValGlb /
0348                                                       sqrt(d01ErrValGlb * d01ErrValGlb + d02ErrValGlb * d02ErrValGlb));
0349                 ddxyNormalizedResiduals_global_->Fill(ddzValGlb /
0350                                                       sqrt(dz1ErrValGlb * dz1ErrValGlb + dz2ErrValGlb * dz2ErrValGlb));
0351                 ddxyNormalizedResiduals_global_->Fill(
0352                     dphiValGlb / sqrt(phi1ErrValGlb * phi1ErrValGlb + phi2ErrValGlb * phi2ErrValGlb));
0353                 ddxyNormalizedResiduals_global_->Fill(
0354                     dthetaValGlb / sqrt(theta1ErrValGlb * theta1ErrValGlb + theta2ErrValGlb * theta2ErrValGlb));
0355                 ddxyNormalizedResiduals_global_->Fill(dptValGlb /
0356                                                       sqrt(pt1ErrValGlb * pt1ErrValGlb + pt2ErrValGlb * pt2ErrValGlb));
0357                 ddxyNormalizedResiduals_global_->Fill(
0358                     dcurvValGlb / sqrt(pow(pt1ErrValGlb, 2) / pow(pt1, 4) + pow(pt2ErrValGlb, 2) / pow(pt2, 4)));
0359               }
0360 
0361             }  // end of split muons loop
0362           }
0363         }
0364       }
0365     }
0366   }
0367 }
0368 
0369 void TrackSplittingMonitor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0370   edm::ParameterSetDescription desc;
0371   desc.setComment(
0372       "Validates track parameters resolution by splitting cosmics tracks at the PCA and comparing the parameters of "
0373       "the two halves");
0374   desc.add<std::string>("FolderName", "TrackSplitMonitoring");
0375   desc.add<edm::InputTag>("splitTrackCollection", edm::InputTag("splittedTracksP5"));
0376   desc.add<edm::InputTag>("splitMuonCollection", edm::InputTag("splitMuons"));
0377   desc.add<bool>("ifPlotMuons", true);
0378   desc.add<int>("pixelHitsPerLeg", 1);
0379   desc.add<int>("totalHitsPerLeg", 6);
0380   desc.add<double>("d0Cut", 12.0);
0381   desc.add<double>("dzCut", 25.0);
0382   desc.add<double>("ptCut", 4.0);
0383   desc.add<double>("norchiCut", 100.0);
0384   desc.add<int>("ddxyBin", 100);
0385   desc.add<double>("ddxyMin", -200.0);
0386   desc.add<double>("ddxyMax", 200.0);
0387   desc.add<int>("ddzBin", 100);
0388   desc.add<double>("ddzMin", -400.0);
0389   desc.add<double>("ddzMax", 400.0);
0390   desc.add<int>("dphiBin", 100);
0391   desc.add<double>("dphiMin", -0.01);
0392   desc.add<double>("dphiMax", 0.01);
0393   desc.add<int>("dthetaBin", 100);
0394   desc.add<double>("dthetaMin", -0.01);
0395   desc.add<double>("dthetaMax", 0.01);
0396   desc.add<int>("dptBin", 100);
0397   desc.add<double>("dptMin", -5.0);
0398   desc.add<double>("dptMax", 5.0);
0399   desc.add<int>("dcurvBin", 100);
0400   desc.add<double>("dcurvMin", -0.005);
0401   desc.add<double>("dcurvMax", 0.005);
0402   desc.add<int>("normBin", 100);
0403   desc.add<double>("normMin", -5.);
0404   desc.add<double>("normMax", 5.);
0405   descriptions.addWithDefaultLabel(desc);
0406 }
0407 
0408 DEFINE_FWK_MODULE(TrackSplittingMonitor);