Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:57:02

0001 // -*- C++ -*-
0002 //
0003 // Package:    Alignment/OfflineValidation
0004 // Class:      DiElectronVertexValidation
0005 //
0006 /**\class DiElectronVertexValidation DiElectronVertexValidation.cc Alignment/OfflineValidation/plugins/DiElectronVertexValidation.cc
0007 
0008  Description: [one line class summary]
0009 
0010  Implementation:
0011      [Notes on implementation]
0012 */
0013 //
0014 // Original Author:  Marco Musich
0015 //         Created:  Thu, 13 May 2021 10:24:07 GMT
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0025 
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028 
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/Utilities/interface/InputTag.h"
0031 
0032 // electrons
0033 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
0034 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0035 #include "DataFormats/EgammaCandidates/interface/GsfElectronCore.h"
0036 #include "DataFormats/ParticleFlowReco/interface/GsfPFRecTrack.h"
0037 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0038 
0039 // tracks
0040 #include "DataFormats/TrackReco/interface/Track.h"
0041 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0042 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
0043 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
0044 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0045 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0046 
0047 // vertices
0048 #include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
0049 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
0050 
0051 // TFileService
0052 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0053 #include "FWCore/ServiceRegistry/interface/Service.h"
0054 
0055 // utilities
0056 #include "DataFormats/Math/interface/deltaR.h"
0057 #include "Alignment/OfflineValidation/interface/DiLeptonVertexHelpers.h"
0058 
0059 // ROOT
0060 #include "TLorentzVector.h"
0061 #include "TH1F.h"
0062 #include "TH2F.h"
0063 
0064 //
0065 // class declaration
0066 //
0067 class DiElectronVertexValidation : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0068 public:
0069   explicit DiElectronVertexValidation(const edm::ParameterSet&);
0070   ~DiElectronVertexValidation() override;
0071 
0072   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0073 
0074 private:
0075   void beginJob() override;
0076   void analyze(const edm::Event&, const edm::EventSetup&) override;
0077   void endJob() override;
0078   bool passLooseSelection(const reco::GsfElectron& electron);
0079 
0080   // ----------member data ---------------------------
0081   DiLeptonHelp::Counts myCounts;
0082   std::vector<double> pTthresholds_;
0083   float maxSVdist_;
0084 
0085   // plot configurations
0086 
0087   edm::ParameterSet CosPhiConfiguration_;
0088   edm::ParameterSet CosPhi3DConfiguration_;
0089   edm::ParameterSet VtxProbConfiguration_;
0090   edm::ParameterSet VtxDistConfiguration_;
0091   edm::ParameterSet VtxDist3DConfiguration_;
0092   edm::ParameterSet VtxDistSigConfiguration_;
0093   edm::ParameterSet VtxDist3DSigConfiguration_;
0094   edm::ParameterSet DiMuMassConfiguration_;
0095 
0096   const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> ttbESToken_;
0097   edm::EDGetTokenT<reco::GsfElectronCollection>
0098       electronsToken_;  //used to select what electrons to read from configuration
0099   edm::EDGetTokenT<reco::GsfTrackCollection>
0100       gsfTracksToken_;                                    //used to select what tracks to read from configuration file
0101   edm::EDGetTokenT<reco::VertexCollection> vertexToken_;  //used to select what vertices to read from configuration file
0102 
0103   TH1F* hSVProb_;
0104   TH1F* hSVDist_;
0105   TH1F* hGSFMult_;
0106   TH1F* hGSFMultAftPt_;
0107   TH1F* hGSF0Pt_;
0108   TH1F* hGSF0Eta_;
0109   TH1F* hGSF1Pt_;
0110   TH1F* hGSF1Eta_;
0111   TH1F* hSVDistSig_;
0112   TH1F* hSVDist3D_;
0113   TH1F* hSVDist3DSig_;
0114   TH1F* hCosPhi_;
0115   TH1F* hCosPhi3D_;
0116   TH1F* hTrackInvMass_;
0117   TH1F* hInvMass_;
0118   TH1I* hClosestVtxIndex_;
0119   TH1F* hCutFlow_;
0120 
0121   // 2D maps
0122 
0123   DiLeptonHelp::PlotsVsKinematics CosPhiPlots = DiLeptonHelp::PlotsVsKinematics(DiLeptonHelp::EE);
0124   DiLeptonHelp::PlotsVsKinematics CosPhi3DPlots = DiLeptonHelp::PlotsVsKinematics(DiLeptonHelp::EE);
0125   DiLeptonHelp::PlotsVsKinematics VtxProbPlots = DiLeptonHelp::PlotsVsKinematics(DiLeptonHelp::EE);
0126   DiLeptonHelp::PlotsVsKinematics VtxDistPlots = DiLeptonHelp::PlotsVsKinematics(DiLeptonHelp::EE);
0127   DiLeptonHelp::PlotsVsKinematics VtxDist3DPlots = DiLeptonHelp::PlotsVsKinematics(DiLeptonHelp::EE);
0128   DiLeptonHelp::PlotsVsKinematics VtxDistSigPlots = DiLeptonHelp::PlotsVsKinematics(DiLeptonHelp::EE);
0129   DiLeptonHelp::PlotsVsKinematics VtxDist3DSigPlots = DiLeptonHelp::PlotsVsKinematics(DiLeptonHelp::EE);
0130   DiLeptonHelp::PlotsVsKinematics ZMassPlots = DiLeptonHelp::PlotsVsKinematics(DiLeptonHelp::EE);
0131 };
0132 
0133 //
0134 // constants, enums and typedefs
0135 //
0136 
0137 static constexpr float cmToum = 10e4;
0138 static constexpr float emass2 = 0.0005109990615 * 0.0005109990615;  //electron mass squared  (GeV^2/c^4)
0139 
0140 //
0141 // constructors and destructor
0142 //
0143 DiElectronVertexValidation::DiElectronVertexValidation(const edm::ParameterSet& iConfig)
0144     : pTthresholds_(iConfig.getParameter<std::vector<double>>("pTThresholds")),
0145       maxSVdist_(iConfig.getParameter<double>("maxSVdist")),
0146       CosPhiConfiguration_(iConfig.getParameter<edm::ParameterSet>("CosPhiConfig")),
0147       CosPhi3DConfiguration_(iConfig.getParameter<edm::ParameterSet>("CosPhi3DConfig")),
0148       VtxProbConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxProbConfig")),
0149       VtxDistConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDistConfig")),
0150       VtxDist3DConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDist3DConfig")),
0151       VtxDistSigConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDistSigConfig")),
0152       VtxDist3DSigConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDist3DSigConfig")),
0153       DiMuMassConfiguration_(iConfig.getParameter<edm::ParameterSet>("DiMuMassConfig")),
0154       ttbESToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
0155       electronsToken_(consumes<reco::GsfElectronCollection>(iConfig.getParameter<edm::InputTag>("electrons"))),
0156       gsfTracksToken_(consumes<reco::GsfTrackCollection>(iConfig.getParameter<edm::InputTag>("gsfTracks"))),
0157       vertexToken_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))) {
0158   usesResource(TFileService::kSharedResource);
0159 
0160   // sort the vector of thresholds
0161   std::sort(pTthresholds_.begin(), pTthresholds_.end(), [](const double& lhs, const double& rhs) { return lhs > rhs; });
0162 
0163   edm::LogInfo("DiElectronVertexValidation") << __FUNCTION__;
0164   for (const auto& thr : pTthresholds_) {
0165     edm::LogInfo("DiElectronVertexValidation") << " Threshold: " << thr << " ";
0166   }
0167   edm::LogInfo("DiElectronVertexValidation") << "Max SV distance: " << maxSVdist_ << " ";
0168 }
0169 
0170 DiElectronVertexValidation::~DiElectronVertexValidation() = default;
0171 
0172 //
0173 // member functions
0174 //
0175 
0176 // ------------ method called for each event  ------------
0177 void DiElectronVertexValidation::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0178   using namespace edm;
0179 
0180   myCounts.eventsTotal++;
0181 
0182   std::vector<const reco::GsfElectron*> myGoodGsfElectrons;
0183 
0184   int totGsfCounter = 0;
0185   for (const auto& gsfEle : iEvent.get(electronsToken_)) {
0186     totGsfCounter++;
0187     if (gsfEle.pt() > pTthresholds_[1] && passLooseSelection(gsfEle)) {
0188       myGoodGsfElectrons.emplace_back(&gsfEle);
0189     }
0190   }
0191 
0192   hGSFMult_->Fill(totGsfCounter);
0193 
0194   std::sort(myGoodGsfElectrons.begin(),
0195             myGoodGsfElectrons.end(),
0196             [](const reco::GsfElectron*& lhs, const reco::GsfElectron*& rhs) { return lhs->pt() > rhs->pt(); });
0197 
0198   hGSFMultAftPt_->Fill(myGoodGsfElectrons.size());
0199 
0200   // reject if there's no Z
0201   if (myGoodGsfElectrons.size() < 2)
0202     return;
0203 
0204   myCounts.eventsAfterMult++;
0205 
0206   if ((myGoodGsfElectrons[0]->pt()) < pTthresholds_[0] || (myGoodGsfElectrons[1]->pt() < pTthresholds_[1]))
0207     return;
0208 
0209   myCounts.eventsAfterPt++;
0210 
0211   if (myGoodGsfElectrons[0]->charge() * myGoodGsfElectrons[1]->charge() > 0)
0212     return;
0213 
0214   if (std::max(std::abs(myGoodGsfElectrons[0]->eta()), std::abs(myGoodGsfElectrons[1]->eta())) > 2.4)
0215     return;
0216 
0217   myCounts.eventsAfterEta++;
0218 
0219   const auto& ele1 = myGoodGsfElectrons[1]->p4();
0220   const auto& ele0 = myGoodGsfElectrons[0]->p4();
0221   const auto& mother = ele1 + ele0;
0222 
0223   float invMass = mother.M();
0224   hInvMass_->Fill(invMass);
0225 
0226   // just copy the top two muons
0227   std::vector<const reco::GsfElectron*> theZElectronVector;
0228   theZElectronVector.reserve(2);
0229   theZElectronVector.emplace_back(myGoodGsfElectrons[1]);
0230   theZElectronVector.emplace_back(myGoodGsfElectrons[0]);
0231 
0232   // do the matching of Z muons with inner tracks
0233 
0234   std::vector<const reco::GsfTrack*> myGoodGsfTracks;
0235 
0236   for (const auto& electron : theZElectronVector) {
0237     float minD = 1000.;
0238     const reco::GsfTrack* theMatch = nullptr;
0239     for (const auto& track : iEvent.get(gsfTracksToken_)) {
0240       float D = ::deltaR(electron->gsfTrack()->eta(), electron->gsfTrack()->phi(), track.eta(), track.phi());
0241       if (D < minD) {
0242         minD = D;
0243         theMatch = &track;
0244       }
0245     }
0246     myGoodGsfTracks.emplace_back(theMatch);
0247   }
0248 
0249   hGSF0Pt_->Fill(myGoodGsfTracks[0]->pt());
0250   hGSF0Eta_->Fill(myGoodGsfTracks[0]->eta());
0251   hGSF1Pt_->Fill(myGoodGsfTracks[1]->pt());
0252   hGSF1Eta_->Fill(myGoodGsfTracks[1]->eta());
0253 
0254   const TransientTrackBuilder* theB = &iSetup.getData(ttbESToken_);
0255   TransientVertex aTransVtx;
0256   std::vector<reco::TransientTrack> tks;
0257 
0258   std::vector<const reco::GsfTrack*> myTracks;
0259   myTracks.emplace_back(myGoodGsfTracks[0]);
0260   myTracks.emplace_back(myGoodGsfTracks[1]);
0261 
0262   if (myTracks.size() != 2)
0263     return;
0264 
0265   const auto& e1 = myTracks[1]->momentum();
0266   const auto& e0 = myTracks[0]->momentum();
0267   const auto& ditrack = e1 + e0;
0268 
0269   const auto& tplus = myTracks[0]->charge() > 0 ? myTracks[0] : myTracks[1];
0270   const auto& tminus = myTracks[0]->charge() < 0 ? myTracks[0] : myTracks[1];
0271 
0272   TLorentzVector p4_tplus(tplus->px(), tplus->py(), tplus->pz(), sqrt((tplus->p() * tplus->p()) + emass2));
0273   TLorentzVector p4_tminus(tminus->px(), tminus->py(), tminus->pz(), sqrt((tminus->p() * tminus->p()) + emass2));
0274 
0275   // creat the pair of TLorentVectors used to make the plos
0276   std::pair<TLorentzVector, TLorentzVector> tktk_p4 = std::make_pair(p4_tplus, p4_tminus);
0277 
0278   const auto& Zp4 = p4_tplus + p4_tminus;
0279   float track_invMass = Zp4.M();
0280   hTrackInvMass_->Fill(track_invMass);
0281 
0282   // fill the z->mm mass plots
0283   ZMassPlots.fillPlots(track_invMass, tktk_p4);
0284 
0285   math::XYZPoint ZpT(ditrack.x(), ditrack.y(), 0);
0286   math::XYZPoint Zp(ditrack.x(), ditrack.y(), ditrack.z());
0287 
0288   for (const auto& track : myTracks) {
0289     reco::TransientTrack trajectory = theB->build(track);
0290     tks.push_back(trajectory);
0291   }
0292 
0293   KalmanVertexFitter kalman(true);
0294   aTransVtx = kalman.vertex(tks);
0295 
0296   double SVProb = TMath::Prob(aTransVtx.totalChiSquared(), (int)aTransVtx.degreesOfFreedom());
0297   hSVProb_->Fill(SVProb);
0298 
0299   if (!aTransVtx.isValid())
0300     return;
0301 
0302   myCounts.eventsAfterVtx++;
0303 
0304   // fill the VtxProb plots
0305   VtxProbPlots.fillPlots(SVProb, tktk_p4);
0306 
0307   // get collection of reconstructed vertices from event
0308   edm::Handle<reco::VertexCollection> vertexHandle = iEvent.getHandle(vertexToken_);
0309 
0310   math::XYZPoint mainVtx(0, 0, 0);
0311   reco::Vertex TheMainVtx;  // = vertexHandle.product()->front();
0312 
0313   VertexDistanceXY vertTool;
0314   VertexDistance3D vertTool3D;
0315 
0316   if (vertexHandle.isValid()) {
0317     const reco::VertexCollection* vertices = vertexHandle.product();
0318     float minD = 9999.;
0319     int closestVtxIndex = 0;
0320     int counter = 0;
0321     for (const auto& vtx : *vertices) {
0322       double dist3D = vertTool3D.distance(aTransVtx, vtx).value();
0323       if (dist3D < minD) {
0324         minD = dist3D;
0325         closestVtxIndex = counter;
0326       }
0327       counter++;
0328     }
0329     if ((*vertices).at(closestVtxIndex).isValid()) {
0330       hClosestVtxIndex_->Fill(closestVtxIndex);
0331       TheMainVtx = (*vertices).at(closestVtxIndex);
0332       mainVtx.SetXYZ(TheMainVtx.position().x(), TheMainVtx.position().y(), TheMainVtx.position().z());
0333     }
0334   }
0335 
0336   const math::XYZPoint myVertex(aTransVtx.position().x(), aTransVtx.position().y(), aTransVtx.position().z());
0337   const math::XYZPoint deltaVtx(mainVtx.x() - myVertex.x(), mainVtx.y() - myVertex.y(), mainVtx.z() - myVertex.z());
0338 
0339   if (TheMainVtx.isValid()) {
0340     // Z Vertex distance in the xy plane
0341     double distance = vertTool.distance(aTransVtx, TheMainVtx).value();
0342     double dist_err = vertTool.distance(aTransVtx, TheMainVtx).error();
0343 
0344     hSVDist_->Fill(distance * cmToum);
0345     hSVDistSig_->Fill(distance / dist_err);
0346 
0347     // fill the VtxDist plots
0348     VtxDistPlots.fillPlots(distance * cmToum, tktk_p4);
0349 
0350     // fill the VtxDisSig plots
0351     VtxDistSigPlots.fillPlots(distance / dist_err, tktk_p4);
0352 
0353     // Z Vertex distance in 3D
0354     double distance3D = vertTool3D.distance(aTransVtx, TheMainVtx).value();
0355     double dist3D_err = vertTool3D.distance(aTransVtx, TheMainVtx).error();
0356 
0357     hSVDist3D_->Fill(distance3D * cmToum);
0358     hSVDist3DSig_->Fill(distance3D / dist3D_err);
0359 
0360     // fill the VtxDist3D plots
0361     VtxDist3DPlots.fillPlots(distance3D * cmToum, tktk_p4);
0362 
0363     // fill the VtxDisSig plots
0364     VtxDist3DSigPlots.fillPlots(distance3D / dist3D_err, tktk_p4);
0365 
0366     // cut on the PV - SV distance
0367     if (distance * cmToum < maxSVdist_) {
0368       myCounts.eventsAfterDist++;
0369 
0370       double cosphi = (ZpT.x() * deltaVtx.x() + ZpT.y() * deltaVtx.y()) /
0371                       (sqrt(ZpT.x() * ZpT.x() + ZpT.y() * ZpT.y()) *
0372                        sqrt(deltaVtx.x() * deltaVtx.x() + deltaVtx.y() * deltaVtx.y()));
0373 
0374       double cosphi3D = (Zp.x() * deltaVtx.x() + Zp.y() * deltaVtx.y() + Zp.z() * deltaVtx.z()) /
0375                         (sqrt(Zp.x() * Zp.x() + Zp.y() * Zp.y() + Zp.z() * Zp.z()) *
0376                          sqrt(deltaVtx.x() * deltaVtx.x() + deltaVtx.y() * deltaVtx.y() + deltaVtx.z() * deltaVtx.z()));
0377 
0378       hCosPhi_->Fill(cosphi);
0379       hCosPhi3D_->Fill(cosphi3D);
0380 
0381       // fill the cosphi plots
0382       CosPhiPlots.fillPlots(cosphi, tktk_p4);
0383 
0384       // fill the VtxDisSig plots
0385       CosPhi3DPlots.fillPlots(cosphi3D, tktk_p4);
0386     }
0387   }
0388 }
0389 
0390 bool DiElectronVertexValidation::passLooseSelection(const reco::GsfElectron& el) {
0391   float dEtaln = fabs(el.deltaEtaSuperClusterTrackAtVtx());
0392   float dPhiln = fabs(el.deltaPhiSuperClusterTrackAtVtx());
0393   float sigmaletaleta = el.full5x5_sigmaIetaIeta();
0394   float hem = el.hadronicOverEm();
0395   double resol = fabs((1 / el.ecalEnergy()) - (el.eSuperClusterOverP() / el.ecalEnergy()));
0396   double mHits = el.gsfTrack()->hitPattern().numberOfAllHits(reco::HitPattern::MISSING_INNER_HITS);
0397   bool barrel = (fabs(el.superCluster()->eta()) <= 1.479);
0398   bool endcap = (!barrel && fabs(el.superCluster()->eta()) < 2.5);
0399 
0400   // loose electron ID
0401 
0402   if (barrel && dEtaln < 0.00477 && dPhiln < 0.222 && sigmaletaleta < 0.011 && hem < 0.298 && resol < 0.241 &&
0403       mHits <= 1)
0404     return true;
0405   if (endcap && dEtaln < 0.00868 && dPhiln < 0.213 && sigmaletaleta < 0.0314 && hem < 0.101 && resol < 0.14 &&
0406       mHits <= 1)
0407     return true;
0408 
0409   return false;
0410 }
0411 
0412 // ------------ method called once each job just before starting event loop  ------------
0413 void DiElectronVertexValidation::beginJob() {
0414   // please remove this method if not needed
0415   edm::Service<TFileService> fs;
0416 
0417   // clang-format off
0418   TH1F::SetDefaultSumw2(kTRUE);
0419   
0420   hGSFMult_= fs->make<TH1F>("GSFMult", ";# gsf tracks;N. events", 20, 0., 20.);
0421   hGSFMultAftPt_= fs->make<TH1F>("GSFMultAftPt", ";# gsf tracks;N. events", 20, 0., 20.);
0422   hGSF0Pt_=  fs->make<TH1F>("GSF0Pt", ";leading GSF track p_{T};N. GSF tracks", 100, 0., 100.);
0423   hGSF0Eta_= fs->make<TH1F>("GSF0Eta", ";leading GSF track #eta;N. GSF tracks", 50, -2.5, 2.5);
0424   hGSF1Pt_=  fs->make<TH1F>("GSF1Pt", ";sub-leading GSF track p_{T};N. GSF tracks", 100, 0., 100.);
0425   hGSF1Eta_= fs->make<TH1F>("GSF1Eta", ";sub-leading GSF track #eta;N. GSF tracks", 50, -2.5, 2.5);
0426 
0427   hSVProb_ = fs->make<TH1F>("VtxProb", ";ZV vertex probability;N(e^{+}e^{-} pairs)", 100, 0., 1.);
0428 
0429   hSVDist_ = fs->make<TH1F>("VtxDist", ";PV-ZV xy distance [#mum];N(e^{+}e^{-} pairs)", 100, 0., 1000.);
0430   hSVDistSig_ = fs->make<TH1F>("VtxDistSig", ";PV-ZV xy distance signficance;N(e^{+}e^{-} pairs)", 100, 0., 5.);
0431 
0432   hSVDist3D_ = fs->make<TH1F>("VtxDist3D", ";PV-ZV 3D distance [#mum];N(e^{+}e^{-} pairs)", 100, 0., 1000.);
0433   hSVDist3DSig_ = fs->make<TH1F>("VtxDist3DSig", ";PV-ZV 3D distance signficance;N(e^{+}e^{-} pairs)", 100, 0., 5.);
0434 
0435   hCosPhi_ = fs->make<TH1F>("CosPhi", ";cos(#phi_{xy});N(ee pairs)", 50, -1., 1.);
0436   hCosPhi3D_ = fs->make<TH1F>("CosPhi3D", ";cos(#phi_{3D});N(ee pairs)", 50, -1., 1.);
0437   hTrackInvMass_ = fs->make<TH1F>("TkTkInvMass", ";M(tk,tk) [GeV];N(tk tk pairs)", 70., 50., 120.);
0438   hInvMass_ = fs->make<TH1F>("InvMass", ";M(#mu#mu) [GeV];N(#mu#mu pairs)", 70., 50., 120.);
0439 
0440   hClosestVtxIndex_ = fs->make<TH1I>("ClosestVtxIndex", ";closest vertex index;N(tk tk pairs)", 20, -0.5, 19.5);
0441 
0442   // 2D Maps
0443 
0444   TFileDirectory dirCosPhi = fs->mkdir("CosPhiPlots");
0445   CosPhiPlots.bookFromPSet(dirCosPhi, CosPhiConfiguration_);
0446 
0447   TFileDirectory dirCosPhi3D = fs->mkdir("CosPhi3DPlots");
0448   CosPhi3DPlots.bookFromPSet(dirCosPhi3D, CosPhi3DConfiguration_);
0449 
0450   TFileDirectory dirVtxProb = fs->mkdir("VtxProbPlots");
0451   VtxProbPlots.bookFromPSet(dirVtxProb, VtxProbConfiguration_);
0452 
0453   TFileDirectory dirVtxDist = fs->mkdir("VtxDistPlots");
0454   VtxDistPlots.bookFromPSet(dirVtxDist, VtxDistConfiguration_);
0455 
0456   TFileDirectory dirVtxDist3D = fs->mkdir("VtxDist3DPlots");
0457   VtxDist3DPlots.bookFromPSet(dirVtxDist3D, VtxDist3DConfiguration_);
0458 
0459   TFileDirectory dirVtxDistSig = fs->mkdir("VtxDistSigPlots");
0460   VtxDistSigPlots.bookFromPSet(dirVtxDistSig, VtxDistSigConfiguration_);
0461 
0462   TFileDirectory dirVtxDist3DSig = fs->mkdir("VtxDist3DSigPlots");
0463   VtxDist3DSigPlots.bookFromPSet(dirVtxDist3DSig, VtxDist3DSigConfiguration_);
0464 
0465   TFileDirectory dirInvariantMass = fs->mkdir("InvariantMassPlots");
0466   ZMassPlots.bookFromPSet(dirInvariantMass, DiMuMassConfiguration_);
0467 
0468   // cut flow 
0469 
0470   hCutFlow_ = fs->make<TH1F>("hCutFlow","cut flow;cut step;events left",6,-0.5,5.5);
0471   std::string names[6]={"Total","Mult.",">pT","<eta","hasVtx","VtxDist"};
0472   for(unsigned int i=0;i<6;i++){
0473     hCutFlow_->GetXaxis()->SetBinLabel(i+1,names[i].c_str());
0474   }
0475 
0476   myCounts.zeroAll();
0477 }
0478 
0479 // ------------ method called once each job just after ending the event loop  ------------
0480 void DiElectronVertexValidation::endJob() {
0481   myCounts.printCounts();
0482 
0483   hCutFlow_->SetBinContent(1,myCounts.eventsTotal);
0484   hCutFlow_->SetBinContent(2,myCounts.eventsAfterMult);
0485   hCutFlow_->SetBinContent(3,myCounts.eventsAfterPt);
0486   hCutFlow_->SetBinContent(4,myCounts.eventsAfterEta);
0487   hCutFlow_->SetBinContent(5,myCounts.eventsAfterVtx);
0488   hCutFlow_->SetBinContent(6,myCounts.eventsAfterDist);
0489 }
0490 
0491 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0492 void DiElectronVertexValidation::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0493   edm::ParameterSetDescription desc;
0494   desc.add<edm::InputTag>("gsfTracks",edm::InputTag("electronGsfTracks"));
0495   desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices"));
0496   desc.add<edm::InputTag>("electrons", edm::InputTag("gedGsfElectrons"));
0497   desc.add<std::vector<double>>("pTThresholds", {25., 15.});
0498   desc.add<double>("maxSVdist", 50.);
0499 
0500   {
0501     edm::ParameterSetDescription psDiMuMass;
0502     psDiMuMass.add<std::string>("name", "DiMuMass");
0503     psDiMuMass.add<std::string>("title", "M(#mu#mu)");
0504     psDiMuMass.add<std::string>("yUnits", "[GeV]");
0505     psDiMuMass.add<int>("NxBins", 24);
0506     psDiMuMass.add<int>("NyBins", 50);
0507     psDiMuMass.add<double>("ymin", 70.);
0508     psDiMuMass.add<double>("ymax", 120.);
0509     desc.add<edm::ParameterSetDescription>("DiMuMassConfig", psDiMuMass);
0510   }
0511   {
0512     edm::ParameterSetDescription psCosPhi;
0513     psCosPhi.add<std::string>("name", "CosPhi");
0514     psCosPhi.add<std::string>("title", "cos(#phi_{xy})");
0515     psCosPhi.add<std::string>("yUnits", "");
0516     psCosPhi.add<int>("NxBins", 50);
0517     psCosPhi.add<int>("NyBins", 50);
0518     psCosPhi.add<double>("ymin", -1.);
0519     psCosPhi.add<double>("ymax", 1.);
0520     desc.add<edm::ParameterSetDescription>("CosPhiConfig", psCosPhi);
0521   }
0522   {
0523     edm::ParameterSetDescription psCosPhi3D;
0524     psCosPhi3D.add<std::string>("name", "CosPhi3D");
0525     psCosPhi3D.add<std::string>("title", "cos(#phi_{3D})");
0526     psCosPhi3D.add<std::string>("yUnits", "");
0527     psCosPhi3D.add<int>("NxBins", 50);
0528     psCosPhi3D.add<int>("NyBins", 50);
0529     psCosPhi3D.add<double>("ymin", -1.);
0530     psCosPhi3D.add<double>("ymax", 1.);
0531     desc.add<edm::ParameterSetDescription>("CosPhi3DConfig", psCosPhi3D);
0532   }
0533   {
0534     edm::ParameterSetDescription psVtxProb;
0535     psVtxProb.add<std::string>("name", "VtxProb");
0536     psVtxProb.add<std::string>("title", "Prob(#chi^{2}_{SV})");
0537     psVtxProb.add<std::string>("yUnits", "");
0538     psVtxProb.add<int>("NxBins", 50);
0539     psVtxProb.add<int>("NyBins", 50);
0540     psVtxProb.add<double>("ymin", 0);
0541     psVtxProb.add<double>("ymax", 1.);
0542     desc.add<edm::ParameterSetDescription>("VtxProbConfig", psVtxProb);
0543   }
0544   {
0545     edm::ParameterSetDescription psVtxDist;
0546     psVtxDist.add<std::string>("name", "VtxDist");
0547     psVtxDist.add<std::string>("title", "d_{xy}(PV,SV)");
0548     psVtxDist.add<std::string>("yUnits", "[#mum]");
0549     psVtxDist.add<int>("NxBins", 50);
0550     psVtxDist.add<int>("NyBins", 100);
0551     psVtxDist.add<double>("ymin", 0);
0552     psVtxDist.add<double>("ymax", 300.);
0553     desc.add<edm::ParameterSetDescription>("VtxDistConfig", psVtxDist);
0554   }
0555   {
0556     edm::ParameterSetDescription psVtxDist3D;
0557     psVtxDist3D.add<std::string>("name", "VtxDist3D");
0558     psVtxDist3D.add<std::string>("title", "d_{3D}(PV,SV)");
0559     psVtxDist3D.add<std::string>("yUnits", "[#mum]");
0560     psVtxDist3D.add<int>("NxBins", 50);
0561     psVtxDist3D.add<int>("NyBins", 250);
0562     psVtxDist3D.add<double>("ymin", 0);
0563     psVtxDist3D.add<double>("ymax", 500.);
0564     desc.add<edm::ParameterSetDescription>("VtxDist3DConfig", psVtxDist3D);
0565   }
0566   {
0567     edm::ParameterSetDescription psVtxDistSig;
0568     psVtxDistSig.add<std::string>("name", "VtxDistSig");
0569     psVtxDistSig.add<std::string>("title", "d_{xy}(PV,SV)/#sigma_{dxy}(PV,SV)");
0570     psVtxDistSig.add<std::string>("yUnits", "");
0571     psVtxDistSig.add<int>("NxBins", 50);
0572     psVtxDistSig.add<int>("NyBins", 100);
0573     psVtxDistSig.add<double>("ymin", 0);
0574     psVtxDistSig.add<double>("ymax", 5.);
0575     desc.add<edm::ParameterSetDescription>("VtxDistSigConfig", psVtxDistSig);
0576   }
0577   {
0578     edm::ParameterSetDescription psVtxDist3DSig;
0579     psVtxDist3DSig.add<std::string>("name", "VtxDist3DSig");
0580     psVtxDist3DSig.add<std::string>("title", "d_{3D}(PV,SV)/#sigma_{d3D}(PV,SV)");
0581     psVtxDist3DSig.add<std::string>("yUnits", "");
0582     psVtxDist3DSig.add<int>("NxBins", 50);
0583     psVtxDist3DSig.add<int>("NyBins", 100);
0584     psVtxDist3DSig.add<double>("ymin", 0);
0585     psVtxDist3DSig.add<double>("ymax", 5.);
0586     desc.add<edm::ParameterSetDescription>("VtxDist3DSigConfig", psVtxDist3DSig);
0587   }
0588 
0589   descriptions.addWithDefaultLabel(desc);
0590 }
0591 
0592 //define this as a plug-in
0593 DEFINE_FWK_MODULE(DiElectronVertexValidation);