Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:15

0001 /*
0002  *  See header file for a description of this class.
0003  *
0004  */
0005 
0006 #include <fmt/printf.h>
0007 
0008 #include "DQMOffline/Alignment/interface/DiMuonVertexMonitor.h"
0009 #include "DQMServices/Core/interface/DQMStore.h"
0010 #include "DataFormats/Math/interface/deltaR.h"
0011 #include "DataFormats/TrackReco/interface/Track.h"
0012 #include "DataFormats/TrackReco/interface/TrackBase.h"
0013 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0017 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
0018 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
0019 #include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
0020 #include "TrackingTools/IPTools/interface/IPTools.h"
0021 
0022 #include "TLorentzVector.h"
0023 
0024 namespace {
0025   constexpr float cmToum = 10e4;
0026   constexpr float mumass2 = 0.105658367 * 0.105658367;  //mu mass squared (GeV^2/c^4)
0027 }  // namespace
0028 
0029 DiMuonVertexMonitor::DiMuonVertexMonitor(const edm::ParameterSet& iConfig)
0030     : ttbESToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
0031       tracksToken_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("muonTracks"))),
0032       vertexToken_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))),
0033       motherName_(iConfig.getParameter<std::string>("decayMotherName")),
0034       MEFolderName_(iConfig.getParameter<std::string>("FolderName")),
0035       useClosestVertex_(iConfig.getParameter<bool>("useClosestVertex")),
0036       maxSVdist_(iConfig.getParameter<double>("maxSVdist")),
0037       CosPhi3DConfiguration_(iConfig.getParameter<edm::ParameterSet>("CosPhi3DConfig")),
0038       SVDistConfiguration_(iConfig.getParameter<edm::ParameterSet>("SVDistConfig")),
0039       SVDistSigConfiguration_(iConfig.getParameter<edm::ParameterSet>("SVDistSigConfig")),
0040       SVDist3DConfiguration_(iConfig.getParameter<edm::ParameterSet>("SVDist3DConfig")),
0041       SVDist3DSigConfiguration_(iConfig.getParameter<edm::ParameterSet>("SVDist3DSigConfig")) {
0042   if (motherName_.find('Z') != std::string::npos) {
0043     massLimits_ = std::make_pair(50., 120);
0044   } else if (motherName_.find("J/#psi") != std::string::npos) {
0045     massLimits_ = std::make_pair(2.7, 3.4);
0046   } else if (motherName_.find("#Upsilon") != std::string::npos) {
0047     massLimits_ = std::make_pair(8.9, 9.9);
0048   } else {
0049     edm::LogError("DiMuonVertexMonitor") << " unrecognized decay mother particle: " << motherName_
0050                                          << " setting the default for the Z->mm (50.,120.)" << std::endl;
0051     massLimits_ = std::make_pair(50., 120);
0052   }
0053 }
0054 
0055 void DiMuonVertexMonitor::bookHistograms(DQMStore::IBooker& iBooker, edm::Run const&, edm::EventSetup const&) {
0056   iBooker.setCurrentFolder(MEFolderName_ + "/DiMuonVertexMonitor");
0057 
0058   // clang-format off
0059   std::string ts = fmt::sprintf(";%s vertex probability;N(#mu#mu pairs)", motherName_);
0060   std::string ps = "N(#mu#mu pairs)";
0061   hSVProb_ = iBooker.book1D("VtxProb", ts, 100, 0., 1.);
0062 
0063   ts = fmt::sprintf("#chi^{2} of the %s vertex; #chi^{2} of the %s vertex; %s", motherName_, motherName_, ps);
0064   hSVChi2_ = iBooker.book1D("VtxChi2", ts, 200, 0., 200.);
0065 
0066   ts = fmt::sprintf("#chi^{2}/ndf of the %s vertex; #chi^{2}/ndf of %s vertex; %s", motherName_, motherName_, ps);
0067   hSVNormChi2_ = iBooker.book1D("VtxNormChi2", ts, 100, 0., 20.);
0068 
0069   std::string histTit = motherName_ + " #rightarrow #mu^{+}#mu^{-}";
0070   ts = fmt::sprintf("%s;PV- %sV xy distance [#mum];%s", histTit, motherName_, ps);
0071   hSVDist_ = iBooker.book1D("VtxDist", ts, 100, 0., 300.);
0072 
0073   ts = fmt::sprintf("%s;PV-%sV xy distance error [#mum];%s", histTit, motherName_, ps);
0074   hSVDistErr_ = iBooker.book1D("VtxDistErr", ts, 100, 0., 1000.);
0075 
0076   ts = fmt::sprintf("%s;PV-%sV xy distance signficance;%s", histTit, motherName_, ps);
0077   hSVDistSig_ = iBooker.book1D("VtxDistSig", ts, 100, 0., 5.);
0078 
0079   ts = fmt::sprintf("compatibility of %s vertex; compatibility of %s vertex; %s", motherName_, motherName_, ps);
0080   hSVCompatibility_ = iBooker.book1D("VtxCompatibility", ts, 100, 0., 100.);
0081 
0082   ts = fmt::sprintf("%s;PV-%sV 3D distance [#mum];%s", histTit, motherName_, ps);
0083   hSVDist3D_ = iBooker.book1D("VtxDist3D", ts, 100, 0., 300.);
0084 
0085   ts = fmt::sprintf("%s;PV-%sV 3D distance error [#mum];%s", histTit, motherName_, ps);
0086   hSVDist3DErr_ = iBooker.book1D("VtxDist3DErr", ts, 100, 0., 1000.);
0087 
0088   ts = fmt::sprintf("%s;PV-%sV 3D distance signficance;%s", histTit, motherName_, ps);
0089   hSVDist3DSig_ = iBooker.book1D("VtxDist3DSig", ts, 100, 0., 5.);
0090 
0091   ts = fmt::sprintf("3D compatibility of %s vertex;3D compatibility of %s vertex; %s", motherName_, motherName_, ps);
0092   hSVCompatibility3D_ = iBooker.book1D("VtxCompatibility3D", ts, 100, 0., 100.);
0093 
0094   hInvMass_ = iBooker.book1D("InvMass", fmt::sprintf("%s;M(#mu,#mu) [GeV];%s", histTit, ps), 70., massLimits_.first, massLimits_.second);
0095   hCosPhi_ = iBooker.book1D("CosPhi", fmt::sprintf("%s;cos(#phi_{xy});%s", histTit, ps), 50, -1., 1.);
0096   hCosPhi3D_ = iBooker.book1D("CosPhi3D", fmt::sprintf("%s;cos(#phi_{3D});%s", histTit, ps), 50, -1., 1.);
0097   hCosPhiInv_ = iBooker.book1D("CosPhiInv", fmt::sprintf("%s;inverted cos(#phi_{xy});%s", histTit, ps), 50, -1., 1.);
0098   hCosPhiInv3D_ = iBooker.book1D("CosPhiInv3D", fmt::sprintf("%s;inverted cos(#phi_{3D});%s", histTit, ps), 50, -1., 1.);
0099   hCosPhiUnbalance_ = iBooker.book1D("CosPhiUnbalance", fmt::sprintf("%s;cos(#phi_{xy}) unbalance;#Delta%s", histTit, ps), 50, -1.,1.);
0100   hCosPhi3DUnbalance_ = iBooker.book1D("CosPhi3DUnbalance", fmt::sprintf("%s;cos(#phi_{3D}) unbalance;#Delta%s", histTit, ps), 50, -1., 1.);
0101 
0102   hdxy_ = iBooker.book1D("dxy", fmt::sprintf("%s;muon track d_{xy}(PV) [#mum];muon tracks", histTit), 150, -300, 300);
0103   hdz_ = iBooker.book1D("dz", fmt::sprintf("%s;muon track d_{z}(PV) [#mum];muon tracks", histTit), 150, -300, 300);
0104   hdxyErr_ = iBooker.book1D("dxyErr", fmt::sprintf("%s;muon track err_{dxy} [#mum];muon tracks", histTit), 250, 0., 500.);
0105   hdzErr_ = iBooker.book1D("dzErr", fmt::sprintf("%s;muon track err_{dz} [#mum];muon tracks", histTit), 250, 0., 500.);
0106   hIP2d_ = iBooker.book1D("IP2d", fmt::sprintf("%s;muon track IP_{2D} [#mum];muon tracks", histTit), 150, -300, 300);
0107   hIP3d_ = iBooker.book1D("IP3d", fmt::sprintf("%s;muon track IP_{3D} [#mum];muon tracks", histTit), 150, -300, 300);
0108   hIP2dsig_ = iBooker.book1D("IP2Dsig", fmt::sprintf("%s;muon track IP_{2D} significance;muon tracks", histTit), 100, 0., 5.);
0109   hIP3dsig_ = iBooker.book1D("IP3Dsig", fmt::sprintf("%s;muon track IP_{3D} significance;muon tracks", histTit), 100, 0., 5.);
0110   // clang-format on
0111 
0112   // now book the cosphi3D plots vs kinematics
0113   iBooker.setCurrentFolder(MEFolderName_ + "/DiMuonVertexMonitor/CosPhi3DPlots");
0114   CosPhi3DPlots_.bookFromPSet(iBooker, CosPhi3DConfiguration_);
0115 
0116   // now book the PV-SV distance plots vs kinematics
0117   iBooker.setCurrentFolder(MEFolderName_ + "/DiMuonVertexMonitor/SVDistPlots");
0118   SVDistPlots_.bookFromPSet(iBooker, SVDistConfiguration_);
0119 
0120   // now book the PV-SV distance significance plots vs kinematics
0121   iBooker.setCurrentFolder(MEFolderName_ + "/DiMuonVertexMonitor/SVDistSigPlots");
0122   SVDistSigPlots_.bookFromPSet(iBooker, SVDistSigConfiguration_);
0123 
0124   // now book the PV-SV 3D distance plots vs kinematics
0125   iBooker.setCurrentFolder(MEFolderName_ + "/DiMuonVertexMonitor/SVDist3DPlots");
0126   SVDist3DPlots_.bookFromPSet(iBooker, SVDist3DConfiguration_);
0127 
0128   // now book the PV-SV 3D distance significance plots vs kinematics
0129   iBooker.setCurrentFolder(MEFolderName_ + "/DiMuonVertexMonitor/SVDist3DSigPlots");
0130   SVDist3DSigPlots_.bookFromPSet(iBooker, SVDist3DSigConfiguration_);
0131 }
0132 
0133 void DiMuonVertexMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0134   std::vector<const reco::Track*> myTracks;
0135   const auto trackHandle = iEvent.getHandle(tracksToken_);
0136   if (!trackHandle.isValid()) {
0137     edm::LogError("DiMuonVertexMonitor") << "invalid track collection encountered!";
0138     return;
0139   }
0140 
0141   for (const auto& muonTrk : *trackHandle) {
0142     myTracks.emplace_back(&muonTrk);
0143   }
0144 
0145 #ifdef EDM_ML_DEBUG
0146   for (const auto& track : myTracks) {
0147     edm::LogVerbatim("DiMuonVertexMonitor") << __PRETTY_FUNCTION__ << " pT: " << track->pt() << " GeV"
0148                                             << " , pT error: " << track->ptError() << " GeV"
0149                                             << " , eta: " << track->eta() << " , phi: " << track->phi() << std::endl;
0150   }
0151 #endif
0152 
0153   const TransientTrackBuilder* theB = &iSetup.getData(ttbESToken_);
0154   TransientVertex mumuTransientVtx;
0155   std::vector<reco::TransientTrack> tks;
0156 
0157   if (myTracks.size() != 2) {
0158     edm::LogWarning("DiMuonVertexMonitor") << "There are not enough tracks to monitor!";
0159     return;
0160   }
0161 
0162   const auto& t1 = myTracks[1]->momentum();
0163   const auto& t0 = myTracks[0]->momentum();
0164   const auto& ditrack = t1 + t0;
0165 
0166   const auto& tplus = myTracks[0]->charge() > 0 ? myTracks[0] : myTracks[1];
0167   const auto& tminus = myTracks[0]->charge() < 0 ? myTracks[0] : myTracks[1];
0168 
0169   TLorentzVector p4_tplus(tplus->px(), tplus->py(), tplus->pz(), sqrt((tplus->p() * tplus->p()) + mumass2));
0170   TLorentzVector p4_tminus(tminus->px(), tminus->py(), tminus->pz(), sqrt((tminus->p() * tminus->p()) + mumass2));
0171 
0172 #ifdef EDM_ML_DEBUG
0173   // Define a lambda function to convert TLorentzVector to a string
0174   auto tLorentzVectorToString = [](const TLorentzVector& vector) {
0175     return std::to_string(vector.Px()) + " " + std::to_string(vector.Py()) + " " + std::to_string(vector.Pz()) + " " +
0176            std::to_string(vector.E());
0177   };
0178 
0179   edm::LogVerbatim("DiMuonVertexMonitor") << "mu+" << tLorentzVectorToString(p4_tplus) << std::endl;
0180   edm::LogVerbatim("DiMuonVertexMonitor") << "mu-" << tLorentzVectorToString(p4_tminus) << std::endl;
0181 #endif
0182 
0183   const auto& Zp4 = p4_tplus + p4_tminus;
0184   float track_invMass = Zp4.M();
0185   hInvMass_->Fill(track_invMass);
0186 
0187   // creat the pair of TLorentVectors used to make the plos
0188   std::pair<TLorentzVector, TLorentzVector> tktk_p4 = std::make_pair(p4_tplus, p4_tminus);
0189 
0190   math::XYZPoint ZpT(ditrack.x(), ditrack.y(), 0);
0191   math::XYZPoint Zp(ditrack.x(), ditrack.y(), ditrack.z());
0192 
0193   for (const auto& track : myTracks) {
0194     reco::TransientTrack trajectory = theB->build(track);
0195     tks.push_back(trajectory);
0196   }
0197 
0198   KalmanVertexFitter kalman(true);
0199   mumuTransientVtx = kalman.vertex(tks);
0200 
0201   double SVProb = TMath::Prob(mumuTransientVtx.totalChiSquared(), (int)mumuTransientVtx.degreesOfFreedom());
0202   hSVProb_->Fill(SVProb);
0203   hSVChi2_->Fill(mumuTransientVtx.totalChiSquared());
0204   hSVNormChi2_->Fill(mumuTransientVtx.totalChiSquared() / (int)mumuTransientVtx.degreesOfFreedom());
0205 
0206   if (!mumuTransientVtx.isValid())
0207     return;
0208 
0209   const reco::Vertex* theClosestVertex;
0210   // get collection of reconstructed vertices from event
0211   edm::Handle<reco::VertexCollection> vertexHandle = iEvent.getHandle(vertexToken_);
0212   if (vertexHandle.isValid()) {
0213     const reco::VertexCollection* vertices = vertexHandle.product();
0214     theClosestVertex = this->findClosestVertex(mumuTransientVtx, vertices);
0215   } else {
0216     edm::LogWarning("DiMuonVertexMonitor") << "invalid vertex collection encountered Skipping event!";
0217     return;
0218   }
0219 
0220   reco::Vertex theMainVtx;
0221   if (!useClosestVertex_ || theClosestVertex == nullptr) {
0222     // if the closest vertex is not available, or explicitly not chosen
0223     theMainVtx = vertexHandle.product()->front();
0224   } else {
0225     theMainVtx = *theClosestVertex;
0226   }
0227 
0228   const math::XYZPoint theMainVtxPos(theMainVtx.position().x(), theMainVtx.position().y(), theMainVtx.position().z());
0229   const math::XYZPoint myVertex(
0230       mumuTransientVtx.position().x(), mumuTransientVtx.position().y(), mumuTransientVtx.position().z());
0231   const math::XYZPoint deltaVtx(
0232       theMainVtxPos.x() - myVertex.x(), theMainVtxPos.y() - myVertex.y(), theMainVtxPos.z() - myVertex.z());
0233 
0234 #ifdef EDM_ML_DEBUG
0235   edm::LogVerbatim("DiMuonVertexMonitor") << "mm vertex position:" << mumuTransientVtx.position().x() << ","
0236                                           << mumuTransientVtx.position().y() << "," << mumuTransientVtx.position().z();
0237 
0238   edm::LogVerbatim("DiMuonVertexMonitor") << "main vertex position:" << theMainVtx.position().x() << ","
0239                                           << theMainVtx.position().y() << "," << theMainVtx.position().z();
0240 #endif
0241 
0242   if (theMainVtx.isValid()) {
0243     // fill the impact parameter plots
0244     for (const auto& track : myTracks) {
0245       hdxy_->Fill(track->dxy(theMainVtxPos) * cmToum);
0246       hdz_->Fill(track->dz(theMainVtxPos) * cmToum);
0247       hdxyErr_->Fill(track->dxyError() * cmToum);
0248       hdzErr_->Fill(track->dzError() * cmToum);
0249 
0250       const auto& ttrk = theB->build(track);
0251       Global3DVector dir(track->px(), track->py(), track->pz());
0252       const auto& ip2d = IPTools::signedTransverseImpactParameter(ttrk, dir, theMainVtx);
0253       const auto& ip3d = IPTools::signedImpactParameter3D(ttrk, dir, theMainVtx);
0254 
0255       hIP2d_->Fill(ip2d.second.value() * cmToum);
0256       hIP3d_->Fill(ip3d.second.value() * cmToum);
0257       hIP2dsig_->Fill(ip2d.second.significance());
0258       hIP3dsig_->Fill(ip3d.second.significance());
0259     }
0260 
0261     // Z Vertex distance in the xy plane
0262     VertexDistanceXY vertTool;
0263     double distance = vertTool.distance(mumuTransientVtx, theMainVtx).value();
0264     double dist_err = vertTool.distance(mumuTransientVtx, theMainVtx).error();
0265     float compatibility = 0.;
0266 
0267     try {
0268       compatibility = vertTool.compatibility(mumuTransientVtx, theMainVtx);
0269     } catch (cms::Exception& er) {
0270       LogTrace("DiMuonVertexMonitor") << "caught std::exception " << er.what() << std::endl;
0271     }
0272 
0273     hSVCompatibility_->Fill(compatibility);
0274     hSVDist_->Fill(distance * cmToum);
0275     hSVDistErr_->Fill(dist_err * cmToum);
0276     hSVDistSig_->Fill(distance / dist_err);
0277 
0278     // Z Vertex distance in 3D
0279     VertexDistance3D vertTool3D;
0280     double distance3D = vertTool3D.distance(mumuTransientVtx, theMainVtx).value();
0281     double dist3D_err = vertTool3D.distance(mumuTransientVtx, theMainVtx).error();
0282     float compatibility3D = 0.;
0283 
0284     try {
0285       compatibility3D = vertTool3D.compatibility(mumuTransientVtx, theMainVtx);
0286     } catch (cms::Exception& er) {
0287       LogTrace("DiMuonVertexMonitor") << "caught std::exception " << er.what() << std::endl;
0288     }
0289 
0290     hSVCompatibility3D_->Fill(compatibility3D);
0291     hSVDist3D_->Fill(distance3D * cmToum);
0292     hSVDist3DErr_->Fill(dist3D_err * cmToum);
0293     hSVDist3DSig_->Fill(distance3D / dist3D_err);
0294 
0295     // creat the pair of TLorentVectors used to make the plos
0296     std::pair<TLorentzVector, TLorentzVector> tktk_p4 = std::make_pair(p4_tplus, p4_tminus);
0297 
0298     SVDistPlots_.fillPlots(distance * cmToum, tktk_p4);
0299     SVDistSigPlots_.fillPlots(distance / dist_err, tktk_p4);
0300     SVDist3DPlots_.fillPlots(distance3D * cmToum, tktk_p4);
0301     SVDist3DSigPlots_.fillPlots(distance3D / dist3D_err, tktk_p4);
0302 
0303     // cut on the PV - SV distance
0304     if (distance * cmToum < maxSVdist_) {
0305       double cosphi = (ZpT.x() * deltaVtx.x() + ZpT.y() * deltaVtx.y()) /
0306                       (sqrt(ZpT.x() * ZpT.x() + ZpT.y() * ZpT.y()) *
0307                        sqrt(deltaVtx.x() * deltaVtx.x() + deltaVtx.y() * deltaVtx.y()));
0308 
0309       double cosphi3D = (Zp.x() * deltaVtx.x() + Zp.y() * deltaVtx.y() + Zp.z() * deltaVtx.z()) /
0310                         (sqrt(Zp.x() * Zp.x() + Zp.y() * Zp.y() + Zp.z() * Zp.z()) *
0311                          sqrt(deltaVtx.x() * deltaVtx.x() + deltaVtx.y() * deltaVtx.y() + deltaVtx.z() * deltaVtx.z()));
0312 
0313       hCosPhi_->Fill(cosphi);
0314       hCosPhi3D_->Fill(cosphi3D);
0315       // inverted
0316       hCosPhiInv_->Fill(-cosphi);
0317       hCosPhiInv3D_->Fill(-cosphi3D);
0318 
0319 #ifdef EDM_ML_DEBUG
0320       edm::LogVerbatim("DiMuonVertexMonitor")
0321           << "distance " << distance * cmToum << " cosphi3D:" << cosphi3D << std::endl;
0322 #endif
0323 
0324       // unbalance
0325       hCosPhiUnbalance_->Fill(cosphi, 1.);
0326       hCosPhiUnbalance_->Fill(-cosphi, -1.);
0327       hCosPhi3DUnbalance_->Fill(cosphi3D, 1.);
0328       hCosPhi3DUnbalance_->Fill(-cosphi3D, -1.);
0329 
0330       // fill the cos(phi3D) plots
0331       CosPhi3DPlots_.fillPlots(cosphi3D, tktk_p4);
0332     }
0333   } else {
0334     edm::LogWarning("DiMuonVertexMonitor") << "hardest primary vertex in the event is not valid!";
0335   }
0336 }
0337 
0338 // compute the closest vertex to di-lepton ------------------------------------
0339 const reco::Vertex* DiMuonVertexMonitor::findClosestVertex(const TransientVertex aTransVtx,
0340                                                            const reco::VertexCollection* vertices) const {
0341   reco::Vertex* defaultVtx = nullptr;
0342 
0343   if (!aTransVtx.isValid())
0344     return defaultVtx;
0345 
0346   // find the closest vertex to the secondary vertex in 3D
0347   VertexDistance3D vertTool3D;
0348   float minD = 9999.;
0349   int closestVtxIndex = 0;
0350   int counter = 0;
0351   for (const auto& vtx : *vertices) {
0352     double dist3D = vertTool3D.distance(aTransVtx, vtx).value();
0353     if (dist3D < minD) {
0354       minD = dist3D;
0355       closestVtxIndex = counter;
0356     }
0357     counter++;
0358   }
0359 
0360   if ((*vertices).at(closestVtxIndex).isValid()) {
0361     return &(vertices->at(closestVtxIndex));
0362   } else {
0363     return defaultVtx;
0364   }
0365 }
0366 
0367 void DiMuonVertexMonitor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0368   edm::ParameterSetDescription desc;
0369   desc.add<edm::InputTag>("muonTracks", edm::InputTag("ALCARECOTkAlDiMuon"));
0370   desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices"));
0371   desc.add<std::string>("FolderName", "DiMuonVertexMonitor");
0372   desc.add<std::string>("decayMotherName", "Z");
0373   desc.add<bool>("useClosestVertex", true);
0374   desc.add<double>("maxSVdist", 50.);
0375 
0376   {
0377     edm::ParameterSetDescription psCosPhi3D;
0378     psCosPhi3D.add<std::string>("name", "CosPhi3D");
0379     psCosPhi3D.add<std::string>("title", "cos(#phi_{3D})");
0380     psCosPhi3D.add<std::string>("yUnits", "");
0381     psCosPhi3D.add<int>("NxBins", 24);
0382     psCosPhi3D.add<int>("NyBins", 50);
0383     psCosPhi3D.add<double>("ymin", -1.);
0384     psCosPhi3D.add<double>("ymax", 1.);
0385     psCosPhi3D.add<double>("maxDeltaEta", 3.7);
0386     desc.add<edm::ParameterSetDescription>("CosPhi3DConfig", psCosPhi3D);
0387   }
0388 
0389   {
0390     edm::ParameterSetDescription psSVDist;
0391     psSVDist.add<std::string>("name", "SVDist");
0392     psSVDist.add<std::string>("title", "PV-SV distance");
0393     psSVDist.add<std::string>("yUnits", "[#mum]");
0394     psSVDist.add<int>("NxBins", 24);
0395     psSVDist.add<int>("NyBins", 100);
0396     psSVDist.add<double>("ymin", 0.);
0397     psSVDist.add<double>("ymax", 300.);
0398     psSVDist.add<double>("maxDeltaEta", 3.7);
0399     desc.add<edm::ParameterSetDescription>("SVDistConfig", psSVDist);
0400   }
0401 
0402   {
0403     edm::ParameterSetDescription psSVDistSig;
0404     psSVDistSig.add<std::string>("name", "SVDistSig");
0405     psSVDistSig.add<std::string>("title", "PV-SV distance significance");
0406     psSVDistSig.add<std::string>("yUnits", "[#mum]");
0407     psSVDistSig.add<int>("NxBins", 24);
0408     psSVDistSig.add<int>("NyBins", 100);
0409     psSVDistSig.add<double>("ymin", 0.);
0410     psSVDistSig.add<double>("ymax", 5.);
0411     psSVDistSig.add<double>("maxDeltaEta", 3.7);
0412     desc.add<edm::ParameterSetDescription>("SVDistSigConfig", psSVDistSig);
0413   }
0414 
0415   {
0416     edm::ParameterSetDescription psSVDist3D;
0417     psSVDist3D.add<std::string>("name", "SVDist3D");
0418     psSVDist3D.add<std::string>("title", "PV-SV 3D distance");
0419     psSVDist3D.add<std::string>("yUnits", "[#mum]");
0420     psSVDist3D.add<int>("NxBins", 24);
0421     psSVDist3D.add<int>("NyBins", 100);
0422     psSVDist3D.add<double>("ymin", 0.);
0423     psSVDist3D.add<double>("ymax", 300.);
0424     psSVDist3D.add<double>("maxDeltaEta", 3.7);
0425     desc.add<edm::ParameterSetDescription>("SVDist3DConfig", psSVDist3D);
0426   }
0427 
0428   {
0429     edm::ParameterSetDescription psSVDist3DSig;
0430     psSVDist3DSig.add<std::string>("name", "SVDist3DSig");
0431     psSVDist3DSig.add<std::string>("title", "PV-SV 3D distance significance");
0432     psSVDist3DSig.add<std::string>("yUnits", "[#mum]");
0433     psSVDist3DSig.add<int>("NxBins", 24);
0434     psSVDist3DSig.add<int>("NyBins", 100);
0435     psSVDist3DSig.add<double>("ymin", 0.);
0436     psSVDist3DSig.add<double>("ymax", 5.);
0437     psSVDist3DSig.add<double>("maxDeltaEta", 3.7);
0438     desc.add<edm::ParameterSetDescription>("SVDist3DSigConfig", psSVDist3DSig);
0439   }
0440 
0441   descriptions.addWithDefaultLabel(desc);
0442 }
0443 
0444 DEFINE_FWK_MODULE(DiMuonVertexMonitor);