Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 
0002 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0003 #include "DataFormats/PatCandidates/interface/Jet.h"
0004 #include "DataFormats/BTauReco/interface/SecondaryVertexTagInfo.h"
0005 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 
0010 /** \class MiniAODSVAnalyzer
0011  *
0012  *  Secondary Vertex Analyzer to run on MiniAOD
0013  *
0014  */
0015 
0016 class MiniAODSVAnalyzer : public DQMEDAnalyzer {
0017 public:
0018   explicit MiniAODSVAnalyzer(const edm::ParameterSet& pSet);
0019   ~MiniAODSVAnalyzer() override = default;
0020 
0021   void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0022 
0023 private:
0024   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0025 
0026   const edm::EDGetTokenT<std::vector<pat::Jet>> jetToken_;
0027   const std::string svTagInfo_;
0028   const double ptMin_;
0029   const double etaMax_;
0030 
0031   MonitorElement* n_sv_;
0032 
0033   MonitorElement* sv_mass_;
0034   MonitorElement* sv_pt_;
0035   MonitorElement* sv_ntracks_;
0036   MonitorElement* sv_chi2norm_;
0037   MonitorElement* sv_chi2prob_;
0038 
0039   // relation to jet
0040   MonitorElement* sv_ptrel_;
0041   MonitorElement* sv_energyratio_;
0042   MonitorElement* sv_deltaR_;
0043 
0044   MonitorElement* sv_dxy_;
0045   MonitorElement* sv_dxysig_;
0046   MonitorElement* sv_d3d_;
0047   MonitorElement* sv_d3dsig_;
0048 };
0049 
0050 MiniAODSVAnalyzer::MiniAODSVAnalyzer(const edm::ParameterSet& pSet)
0051     : jetToken_(consumes<std::vector<pat::Jet>>(pSet.getParameter<edm::InputTag>("JetTag"))),
0052       svTagInfo_(pSet.getParameter<std::string>("svTagInfo")),
0053       ptMin_(pSet.getParameter<double>("ptMin")),
0054       etaMax_(pSet.getParameter<double>("etaMax")) {}
0055 
0056 void MiniAODSVAnalyzer::bookHistograms(DQMStore::IBooker& ibook, edm::Run const& run, edm::EventSetup const& es) {
0057   ibook.setCurrentFolder("Btag/SV");
0058 
0059   n_sv_ = ibook.book1D("n_sv", "number of SV in jet", 5, 0, 5);
0060   n_sv_->setAxisTitle("number of SV in jet");
0061 
0062   sv_mass_ = ibook.book1D("sv_mass", "SV mass", 30, 0., 6.);
0063   sv_mass_->setAxisTitle("SV mass");
0064 
0065   sv_pt_ = ibook.book1D("sv_pt", "SV transverse momentum", 40, 0., 120.);
0066   sv_pt_->setAxisTitle("SV pt");
0067 
0068   sv_ntracks_ = ibook.book1D("sv_ntracks", "SV number of daugthers", 10, 0, 10);
0069   sv_ntracks_->setAxisTitle("number of tracks at SV");
0070 
0071   sv_chi2norm_ = ibook.book1D("sv_chi2norm", "normalized Chi2 of vertex", 30, 0, 15);
0072   sv_chi2norm_->setAxisTitle("normalized Chi2 of SV");
0073 
0074   sv_chi2prob_ = ibook.book1D("sv_chi2prob", "Chi2 probability of vertex", 20, 0., 1.);
0075   sv_chi2prob_->setAxisTitle("Chi2 probability of SV");
0076 
0077   sv_ptrel_ = ibook.book1D("sv_ptrel", "SV jet transverse momentum ratio", 25, 0., 1.);
0078   sv_ptrel_->setAxisTitle("pt(SV)/pt(jet)");
0079 
0080   sv_energyratio_ = ibook.book1D("sv_energyratio", "SV jet energy ratio", 25, 0., 1.);
0081   sv_energyratio_->setAxisTitle("E(SV)/E(jet)");
0082 
0083   sv_deltaR_ = ibook.book1D("sv_deltaR", "SV jet deltaR", 40, 0., 0.4);
0084   sv_deltaR_->setAxisTitle("deltaR(jet, SV)");
0085 
0086   sv_dxy_ = ibook.book1D("sv_dxy", "2D flight distance", 40, 0., 8.);
0087   sv_dxy_->setAxisTitle("dxy");
0088 
0089   sv_dxysig_ = ibook.book1D("sv_dxysig", "2D flight distance significance", 25, 0., 250.);
0090   sv_dxysig_->setAxisTitle("dxy significance");
0091 
0092   sv_d3d_ = ibook.book1D("sv_d3d", "3D flight distance", 40, 0., 8.);
0093   sv_d3d_->setAxisTitle("d3d");
0094 
0095   sv_d3dsig_ = ibook.book1D("sv_d3dsig", "3D flight distance significance", 25, 0., 250.);
0096   sv_d3dsig_->setAxisTitle("d3d significance");
0097 }
0098 
0099 void MiniAODSVAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0100   edm::Handle<std::vector<pat::Jet>> jetCollection;
0101   iEvent.getByToken(jetToken_, jetCollection);
0102 
0103   // Loop over the pat::Jets
0104   for (std::vector<pat::Jet>::const_iterator jet = jetCollection->begin(); jet != jetCollection->end(); ++jet) {
0105     // jet selection
0106     if (jet->hasTagInfo(svTagInfo_) && jet->pt() > ptMin_ && std::abs(jet->eta()) < etaMax_) {
0107       const reco::CandSecondaryVertexTagInfo* taginfo =
0108           static_cast<const reco::CandSecondaryVertexTagInfo*>(jet->tagInfo(svTagInfo_));
0109       n_sv_->Fill(taginfo->nVertices());
0110 
0111       // loop secondary vertices
0112       for (unsigned int i = 0; i < taginfo->nVertices(); i++) {
0113         const reco::VertexCompositePtrCandidate& sv = taginfo->secondaryVertex(i);
0114 
0115         sv_mass_->Fill(sv.mass());
0116         sv_pt_->Fill(sv.pt());
0117         sv_ntracks_->Fill(sv.numberOfDaughters());
0118         sv_chi2norm_->Fill(sv.vertexNormalizedChi2());
0119         sv_chi2prob_->Fill(ChiSquaredProbability(sv.vertexChi2(), sv.vertexNdof()));
0120 
0121         sv_ptrel_->Fill(sv.pt() / jet->pt());
0122         sv_energyratio_->Fill(sv.energy() / jet->energy());
0123         sv_deltaR_->Fill(reco::deltaR(sv, jet->momentum()));
0124 
0125         sv_dxy_->Fill(taginfo->flightDistance(i, 2).value());
0126         sv_dxysig_->Fill(taginfo->flightDistance(i, 2).significance());
0127         sv_d3d_->Fill(taginfo->flightDistance(i, 3).value());
0128         sv_d3dsig_->Fill(taginfo->flightDistance(i, 3).significance());
0129       }
0130     }
0131   }
0132 }
0133 
0134 //define this as a plug-in
0135 DEFINE_FWK_MODULE(MiniAODSVAnalyzer);