Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408
/*
 *  See header file for a description of this class.
 *
 *  \author Suchandra Dutta , Giorgia Mila
 */

#include "DQM/TrackingMonitor/interface/TrackSplittingMonitor.h"
#include "DQMServices/Core/interface/DQMStore.h"
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h"
#include "DataFormats/TrackReco/interface/Track.h"
#include "DataFormats/TrackReco/interface/TrackFwd.h"
#include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
#include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
#include "TrackingTools/Records/interface/TransientRecHitRecord.h"
#include "TrackingTools/Records/interface/TransientTrackRecord.h"
#include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
#include "TrackingTools/TransientTrack/interface/TransientTrack.h"
#include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"

#include <string>

TrackSplittingMonitor::TrackSplittingMonitor(const edm::ParameterSet& iConfig)
    : conf_(iConfig),
      mfToken_(esConsumes()),
      tkGeomToken_(esConsumes()),
      dtGeomToken_(esConsumes()),
      cscGeomToken_(esConsumes()),
      rpcGeomToken_(esConsumes()),
      splitTracksToken_(consumes<std::vector<reco::Track> >(conf_.getParameter<edm::InputTag>("splitTrackCollection"))),
      splitMuonsToken_(mayConsume<std::vector<reco::Muon> >(conf_.getParameter<edm::InputTag>("splitMuonCollection"))),
      plotMuons_(conf_.getParameter<bool>("ifPlotMuons")),
      pixelHitsPerLeg_(conf_.getParameter<int>("pixelHitsPerLeg")),
      totalHitsPerLeg_(conf_.getParameter<int>("totalHitsPerLeg")),
      d0Cut_(conf_.getParameter<double>("d0Cut")),
      dzCut_(conf_.getParameter<double>("dzCut")),
      ptCut_(conf_.getParameter<double>("ptCut")),
      norchiCut_(conf_.getParameter<double>("norchiCut")) {}

void TrackSplittingMonitor::bookHistograms(DQMStore::IBooker& ibooker,
                                           edm::Run const& /* iRun */,
                                           edm::EventSetup const& /* iSetup */) {
  std::string MEFolderName = conf_.getParameter<std::string>("FolderName");
  ibooker.setCurrentFolder(MEFolderName);

  // bin declarations
  int ddxyBin = conf_.getParameter<int>("ddxyBin");
  double ddxyMin = conf_.getParameter<double>("ddxyMin");
  double ddxyMax = conf_.getParameter<double>("ddxyMax");

  int ddzBin = conf_.getParameter<int>("ddzBin");
  double ddzMin = conf_.getParameter<double>("ddzMin");
  double ddzMax = conf_.getParameter<double>("ddzMax");

  int dphiBin = conf_.getParameter<int>("dphiBin");
  double dphiMin = conf_.getParameter<double>("dphiMin");
  double dphiMax = conf_.getParameter<double>("dphiMax");

  int dthetaBin = conf_.getParameter<int>("dthetaBin");
  double dthetaMin = conf_.getParameter<double>("dthetaMin");
  double dthetaMax = conf_.getParameter<double>("dthetaMax");

  int dptBin = conf_.getParameter<int>("dptBin");
  double dptMin = conf_.getParameter<double>("dptMin");
  double dptMax = conf_.getParameter<double>("dptMax");

  int dcurvBin = conf_.getParameter<int>("dcurvBin");
  double dcurvMin = conf_.getParameter<double>("dcurvMin");
  double dcurvMax = conf_.getParameter<double>("dcurvMax");

  int normBin = conf_.getParameter<int>("normBin");
  double normMin = conf_.getParameter<double>("normMin");
  double normMax = conf_.getParameter<double>("normMax");

  // declare histogram
  ddxyAbsoluteResiduals_tracker_ =
      ibooker.book1D("ddxyAbsoluteResiduals_tracker", "ddxyAbsoluteResiduals_tracker", ddxyBin, ddxyMin, ddxyMax);
  ddzAbsoluteResiduals_tracker_ =
      ibooker.book1D("ddzAbsoluteResiduals_tracker", "ddzAbsoluteResiduals_tracker", ddzBin, ddzMin, ddzMax);
  dphiAbsoluteResiduals_tracker_ =
      ibooker.book1D("dphiAbsoluteResiduals_tracker", "dphiAbsoluteResiduals_tracker", dphiBin, dphiMin, dphiMax);
  dthetaAbsoluteResiduals_tracker_ = ibooker.book1D(
      "dthetaAbsoluteResiduals_tracker", "dthetaAbsoluteResiduals_tracker", dthetaBin, dthetaMin, dthetaMax);
  dptAbsoluteResiduals_tracker_ =
      ibooker.book1D("dptAbsoluteResiduals_tracker", "dptAbsoluteResiduals_tracker", dptBin, dptMin, dptMax);
  dcurvAbsoluteResiduals_tracker_ =
      ibooker.book1D("dcurvAbsoluteResiduals_tracker", "dcurvAbsoluteResiduals_tracker", dcurvBin, dcurvMin, dcurvMax);

  ddxyNormalizedResiduals_tracker_ =
      ibooker.book1D("ddxyNormalizedResiduals_tracker", "ddxyNormalizedResiduals_tracker", normBin, normMin, normMax);
  ddzNormalizedResiduals_tracker_ =
      ibooker.book1D("ddzNormalizedResiduals_tracker", "ddzNormalizedResiduals_tracker", normBin, normMin, normMax);
  dphiNormalizedResiduals_tracker_ =
      ibooker.book1D("dphiNormalizedResiduals_tracker", "dphiNormalizedResiduals_tracker", normBin, normMin, normMax);
  dthetaNormalizedResiduals_tracker_ = ibooker.book1D(
      "dthetaNormalizedResiduals_tracker", "dthetaNormalizedResiduals_tracker", normBin, normMin, normMax);
  dptNormalizedResiduals_tracker_ =
      ibooker.book1D("dptNormalizedResiduals_tracker", "dptNormalizedResiduals_tracker", normBin, normMin, normMax);
  dcurvNormalizedResiduals_tracker_ =
      ibooker.book1D("dcurvNormalizedResiduals_tracker", "dcurvNormalizedResiduals_tracker", normBin, normMin, normMax);

  if (plotMuons_) {
    ddxyAbsoluteResiduals_global_ =
        ibooker.book1D("ddxyAbsoluteResiduals_global", "ddxyAbsoluteResiduals_global", ddxyBin, ddxyMin, ddxyMax);
    ddzAbsoluteResiduals_global_ =
        ibooker.book1D("ddzAbsoluteResiduals_global", "ddzAbsoluteResiduals_global", ddzBin, ddzMin, ddzMax);
    dphiAbsoluteResiduals_global_ =
        ibooker.book1D("dphiAbsoluteResiduals_global", "dphiAbsoluteResiduals_global", dphiBin, dphiMin, dphiMax);
    dthetaAbsoluteResiduals_global_ = ibooker.book1D(
        "dthetaAbsoluteResiduals_global", "dthetaAbsoluteResiduals_global", dthetaBin, dthetaMin, dthetaMax);
    dptAbsoluteResiduals_global_ =
        ibooker.book1D("dptAbsoluteResiduals_global", "dptAbsoluteResiduals_global", dptBin, dptMin, dptMax);
    dcurvAbsoluteResiduals_global_ =
        ibooker.book1D("dcurvAbsoluteResiduals_global", "dcurvAbsoluteResiduals_global", dcurvBin, dcurvMin, dcurvMax);

    ddxyNormalizedResiduals_global_ =
        ibooker.book1D("ddxyNormalizedResiduals_global", "ddxyNormalizedResiduals_global", normBin, normMin, normMax);
    ddzNormalizedResiduals_global_ =
        ibooker.book1D("ddzNormalizedResiduals_global", "ddzNormalizedResiduals_global", normBin, normMin, normMax);
    dphiNormalizedResiduals_global_ =
        ibooker.book1D("dphiNormalizedResiduals_global", "dphiNormalizedResiduals_global", normBin, normMin, normMax);
    dthetaNormalizedResiduals_global_ = ibooker.book1D(
        "dthetaNormalizedResiduals_global", "dthetaNormalizedResiduals_global", normBin, normMin, normMax);
    dptNormalizedResiduals_global_ =
        ibooker.book1D("dptNormalizedResiduals_global", "dptNormalizedResiduals_global", normBin, normMin, normMax);
    dcurvNormalizedResiduals_global_ =
        ibooker.book1D("dcurvNormalizedResiduals_global", "dcurvNormalizedResiduals_global", normBin, normMin, normMax);
  }

  ddxyAbsoluteResiduals_tracker_->setAxisTitle("(#delta d_{xy})/#sqrt{2} [#mum]");
  ddzAbsoluteResiduals_tracker_->setAxisTitle("(#delta d_{z})/#sqrt{2} [#mum]");
  dphiAbsoluteResiduals_tracker_->setAxisTitle("(#delta #phi)/#sqrt{2} [mrad]");
  dthetaAbsoluteResiduals_tracker_->setAxisTitle("(#delta #theta)/#sqrt{2} [mrad]");
  dptAbsoluteResiduals_tracker_->setAxisTitle("(#delta p_{T})/#sqrt{2} [GeV]");
  dcurvAbsoluteResiduals_tracker_->setAxisTitle("(#delta (1/p_{T}))/#sqrt{2} [GeV^{-1}]");

  ddxyNormalizedResiduals_tracker_->setAxisTitle("#delta d_{xy}/#sigma(d_{xy})");
  ddzNormalizedResiduals_tracker_->setAxisTitle("#delta d_{z}/#sigma(d_{z})");
  dphiNormalizedResiduals_tracker_->setAxisTitle("#delta #phi/#sigma(d_{#phi})");
  dthetaNormalizedResiduals_tracker_->setAxisTitle("#delta #theta/#sigma(d_{#theta})");
  dptNormalizedResiduals_tracker_->setAxisTitle("#delta p_{T}/#sigma(p_{T})");
  dcurvNormalizedResiduals_tracker_->setAxisTitle("#delta 1/p_{T}/#sigma(1/p_{T})");

  if (plotMuons_) {
    ddxyAbsoluteResiduals_global_->setAxisTitle("(#delta d_{xy})/#sqrt{2} [#mum]");
    ddzAbsoluteResiduals_global_->setAxisTitle("(#delta d_{z})/#sqrt{2} [#mum]");
    dphiAbsoluteResiduals_global_->setAxisTitle("(#delta #phi)/#sqrt{2} [mrad]");
    dthetaAbsoluteResiduals_global_->setAxisTitle("(#delta #theta)/#sqrt{2} [mrad]");
    dptAbsoluteResiduals_global_->setAxisTitle("(#delta p_{T})/#sqrt{2} [GeV]");
    dcurvAbsoluteResiduals_global_->setAxisTitle("(#delta (1/p_{T}))/#sqrt{2} [GeV^{-1}]");

    ddxyNormalizedResiduals_global_->setAxisTitle("#delta d_{xy}/#sigma(d_{xy})");
    ddzNormalizedResiduals_global_->setAxisTitle("#delta d_{z}/#sigma(d_{z})");
    dphiNormalizedResiduals_global_->setAxisTitle("#delta #phi/#sigma(d_{#phi})");
    dthetaNormalizedResiduals_global_->setAxisTitle("#delta #theta/#sigma(d_{#theta})");
    dptNormalizedResiduals_global_->setAxisTitle("#delta p_{T}/#sigma(p_{T})");
    dcurvNormalizedResiduals_global_->setAxisTitle("#delta 1/p_{T}/#sigma(1/p_{T})");
  }
}

//
// -- Analyse
//
void TrackSplittingMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
  theMagField = &iSetup.getData(mfToken_);
  theGeometry = &iSetup.getData(tkGeomToken_);
  dtGeometry = &iSetup.getData(dtGeomToken_);
  cscGeometry = &iSetup.getData(cscGeomToken_);
  rpcGeometry = &iSetup.getData(rpcGeomToken_);

  edm::Handle<std::vector<reco::Track> > splitTracks = iEvent.getHandle(splitTracksToken_);
  if (!splitTracks.isValid())
    return;

  edm::Handle<std::vector<reco::Muon> > splitMuons;
  if (plotMuons_) {
    splitMuons = iEvent.getHandle(splitMuonsToken_);
  }

  if (splitTracks->size() == 2) {
    // check that there are 2 tracks in split track collection
    edm::LogInfo("TrackSplittingMonitor") << "Split Track size: " << splitTracks->size();

    // split tracks calculations
    reco::Track track1 = splitTracks->at(0);
    reco::Track track2 = splitTracks->at(1);

    // -------------------------- basic selection ---------------------------

    // hit counting
    // looping through the hits for track 1
    double nRechits1 = 0;
    double nRechitinBPIX1 = 0;
    for (auto const& iHit : track1.recHits()) {
      if (iHit->isValid()) {
        nRechits1++;
        int type = iHit->geographicalId().subdetId();
        if (type == int(PixelSubdetector::PixelBarrel)) {
          ++nRechitinBPIX1;
        }
      }
    }
    // looping through the hits for track 2
    double nRechits2 = 0;
    double nRechitinBPIX2 = 0;
    for (auto const& iHit : track2.recHits()) {
      if (iHit->isValid()) {
        nRechits2++;
        int type = iHit->geographicalId().subdetId();
        if (type == int(PixelSubdetector::PixelBarrel)) {
          ++nRechitinBPIX2;
        }
      }
    }

    // DCA of each track
    double d01 = track1.d0();
    double dz1 = track1.dz();
    double d02 = track2.d0();
    double dz2 = track2.dz();

    // pT of each track
    double pt1 = track1.pt();
    double pt2 = track2.pt();

    // chi2 of each track
    double norchi1 = track1.normalizedChi2();
    double norchi2 = track2.normalizedChi2();

    // basic selection
    // pixel hits and total hits
    if ((nRechitinBPIX1 >= pixelHitsPerLeg_) && (nRechitinBPIX2 >= pixelHitsPerLeg_) &&
        (nRechits1 >= totalHitsPerLeg_) && (nRechits2 >= totalHitsPerLeg_)) {
      // dca cut
      if (((std::abs(d01) < d0Cut_)) && (std::abs(d02) < d0Cut_) && (std::abs(dz1) < dzCut_) &&
          (std::abs(dz2) < dzCut_)) {
        // pt cut
        if ((pt1 + pt2) / 2 < ptCut_) {
          // chi2 cut
          if ((norchi1 < norchiCut_) && (norchi2 < norchiCut_)) {
            // passed all cuts...
            edm::LogInfo("TrackSplittingMonitor") << " Setected after all cuts ?";

            double ddxyVal = d01 - d02;
            double ddzVal = dz1 - dz2;
            double dphiVal = track1.phi() - track2.phi();
            double dthetaVal = track1.theta() - track2.theta();
            double dptVal = pt1 - pt2;
            double dcurvVal = (1 / pt1) - (1 / pt2);

            double d01ErrVal = track1.d0Error();
            double d02ErrVal = track2.d0Error();
            double dz1ErrVal = track1.dzError();
            double dz2ErrVal = track2.dzError();
            double phi1ErrVal = track1.phiError();
            double phi2ErrVal = track2.phiError();
            double theta1ErrVal = track1.thetaError();
            double theta2ErrVal = track2.thetaError();
            double pt1ErrVal = track1.ptError();
            double pt2ErrVal = track2.ptError();

            ddxyAbsoluteResiduals_tracker_->Fill(cmToUm * ddxyVal / sqrt2);
            ddzAbsoluteResiduals_tracker_->Fill(cmToUm * ddzVal / sqrt2);
            dphiAbsoluteResiduals_tracker_->Fill(radToUrad * dphiVal / sqrt2);
            dthetaAbsoluteResiduals_tracker_->Fill(radToUrad * dthetaVal / sqrt2);
            dptAbsoluteResiduals_tracker_->Fill(dptVal / sqrt2);
            dcurvAbsoluteResiduals_tracker_->Fill(dcurvVal / sqrt2);

            ddxyNormalizedResiduals_tracker_->Fill(ddxyVal / sqrt(d01ErrVal * d01ErrVal + d02ErrVal * d02ErrVal));
            ddzNormalizedResiduals_tracker_->Fill(ddzVal / sqrt(dz1ErrVal * dz1ErrVal + dz2ErrVal * dz2ErrVal));
            dphiNormalizedResiduals_tracker_->Fill(dphiVal / sqrt(phi1ErrVal * phi1ErrVal + phi2ErrVal * phi2ErrVal));
            dthetaNormalizedResiduals_tracker_->Fill(dthetaVal /
                                                     sqrt(theta1ErrVal * theta1ErrVal + theta2ErrVal * theta2ErrVal));
            dptNormalizedResiduals_tracker_->Fill(dptVal / sqrt(pt1ErrVal * pt1ErrVal + pt2ErrVal * pt2ErrVal));
            dcurvNormalizedResiduals_tracker_->Fill(
                dcurvVal / sqrt(pow(pt1ErrVal, 2) / pow(pt1, 4) + pow(pt2ErrVal, 2) / pow(pt2, 4)));

            // if do the same for split muons
            if (plotMuons_ && splitMuons.isValid()) {
              int gmCtr = 0;
              bool topGlobalMuonFlag = false;
              bool bottomGlobalMuonFlag = false;
              int topGlobalMuon = -1;
              int bottomGlobalMuon = -1;
              double topGlobalMuonNorchi2 = 1e10;
              double bottomGlobalMuonNorchi2 = 1e10;

              // check if usable split global muons
              for (std::vector<reco::Muon>::const_iterator gmI = splitMuons->begin(); gmI != splitMuons->end(); gmI++) {
                if (gmI->isTrackerMuon() && gmI->isStandAloneMuon() && gmI->isGlobalMuon()) {
                  reco::TrackRef trackerTrackRef1(splitTracks, 0);
                  reco::TrackRef trackerTrackRef2(splitTracks, 1);

                  if (gmI->innerTrack() == trackerTrackRef1) {
                    if (gmI->globalTrack()->normalizedChi2() < topGlobalMuonNorchi2) {
                      topGlobalMuonFlag = true;
                      topGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
                      topGlobalMuon = gmCtr;
                    }
                  }
                  if (gmI->innerTrack() == trackerTrackRef2) {
                    if (gmI->globalTrack()->normalizedChi2() < bottomGlobalMuonNorchi2) {
                      bottomGlobalMuonFlag = true;
                      bottomGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
                      bottomGlobalMuon = gmCtr;
                    }
                  }
                }
                gmCtr++;
              }

              if (bottomGlobalMuonFlag && topGlobalMuonFlag) {
                reco::Muon muonTop = splitMuons->at(topGlobalMuon);
                reco::Muon muonBottom = splitMuons->at(bottomGlobalMuon);

                reco::TrackRef glb1 = muonTop.globalTrack();
                reco::TrackRef glb2 = muonBottom.globalTrack();

                double ddxyValGlb = glb1->d0() - glb2->d0();
                double ddzValGlb = glb1->dz() - glb2->dz();
                double dphiValGlb = glb1->phi() - glb2->phi();
                double dthetaValGlb = glb1->theta() - glb2->theta();
                double dptValGlb = glb1->pt() - glb2->pt();
                double dcurvValGlb = (1 / glb1->pt()) - (1 / glb2->pt());

                double d01ErrValGlb = glb1->d0Error();
                double d02ErrValGlb = glb2->d0Error();
                double dz1ErrValGlb = glb1->dzError();
                double dz2ErrValGlb = glb2->dzError();
                double phi1ErrValGlb = glb1->phiError();
                double phi2ErrValGlb = glb2->phiError();
                double theta1ErrValGlb = glb1->thetaError();
                double theta2ErrValGlb = glb2->thetaError();
                double pt1ErrValGlb = glb1->ptError();
                double pt2ErrValGlb = glb2->ptError();

                ddxyAbsoluteResiduals_global_->Fill(cmToUm * ddxyValGlb / sqrt2);
                ddzAbsoluteResiduals_global_->Fill(cmToUm * ddzValGlb / sqrt2);
                dphiAbsoluteResiduals_global_->Fill(radToUrad * dphiValGlb / sqrt2);
                dthetaAbsoluteResiduals_global_->Fill(radToUrad * dthetaValGlb / sqrt2);
                dptAbsoluteResiduals_global_->Fill(dptValGlb / sqrt2);
                dcurvAbsoluteResiduals_global_->Fill(dcurvValGlb / sqrt2);

                ddxyNormalizedResiduals_global_->Fill(ddxyValGlb /
                                                      sqrt(d01ErrValGlb * d01ErrValGlb + d02ErrValGlb * d02ErrValGlb));
                ddxyNormalizedResiduals_global_->Fill(ddzValGlb /
                                                      sqrt(dz1ErrValGlb * dz1ErrValGlb + dz2ErrValGlb * dz2ErrValGlb));
                ddxyNormalizedResiduals_global_->Fill(
                    dphiValGlb / sqrt(phi1ErrValGlb * phi1ErrValGlb + phi2ErrValGlb * phi2ErrValGlb));
                ddxyNormalizedResiduals_global_->Fill(
                    dthetaValGlb / sqrt(theta1ErrValGlb * theta1ErrValGlb + theta2ErrValGlb * theta2ErrValGlb));
                ddxyNormalizedResiduals_global_->Fill(dptValGlb /
                                                      sqrt(pt1ErrValGlb * pt1ErrValGlb + pt2ErrValGlb * pt2ErrValGlb));
                ddxyNormalizedResiduals_global_->Fill(
                    dcurvValGlb / sqrt(pow(pt1ErrValGlb, 2) / pow(pt1, 4) + pow(pt2ErrValGlb, 2) / pow(pt2, 4)));
              }

            }  // end of split muons loop
          }
        }
      }
    }
  }
}

void TrackSplittingMonitor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
  edm::ParameterSetDescription desc;
  desc.setComment(
      "Validates track parameters resolution by splitting cosmics tracks at the PCA and comparing the parameters of "
      "the two halves");
  desc.add<std::string>("FolderName", "TrackSplitMonitoring");
  desc.add<edm::InputTag>("splitTrackCollection", edm::InputTag("splittedTracksP5"));
  desc.add<edm::InputTag>("splitMuonCollection", edm::InputTag("splitMuons"));
  desc.add<bool>("ifPlotMuons", true);
  desc.add<int>("pixelHitsPerLeg", 1);
  desc.add<int>("totalHitsPerLeg", 6);
  desc.add<double>("d0Cut", 12.0);
  desc.add<double>("dzCut", 25.0);
  desc.add<double>("ptCut", 4.0);
  desc.add<double>("norchiCut", 100.0);
  desc.add<int>("ddxyBin", 100);
  desc.add<double>("ddxyMin", -200.0);
  desc.add<double>("ddxyMax", 200.0);
  desc.add<int>("ddzBin", 100);
  desc.add<double>("ddzMin", -400.0);
  desc.add<double>("ddzMax", 400.0);
  desc.add<int>("dphiBin", 100);
  desc.add<double>("dphiMin", -0.01);
  desc.add<double>("dphiMax", 0.01);
  desc.add<int>("dthetaBin", 100);
  desc.add<double>("dthetaMin", -0.01);
  desc.add<double>("dthetaMax", 0.01);
  desc.add<int>("dptBin", 100);
  desc.add<double>("dptMin", -5.0);
  desc.add<double>("dptMax", 5.0);
  desc.add<int>("dcurvBin", 100);
  desc.add<double>("dcurvMin", -0.005);
  desc.add<double>("dcurvMax", 0.005);
  desc.add<int>("normBin", 100);
  desc.add<double>("normMin", -5.);
  desc.add<double>("normMax", 5.);
  descriptions.addWithDefaultLabel(desc);
}

DEFINE_FWK_MODULE(TrackSplittingMonitor);