CMS 3D CMS Logo

MiniAODSVAnalyzer.cc
Go to the documentation of this file.
1 
9 
17 public:
18  explicit MiniAODSVAnalyzer(const edm::ParameterSet& pSet);
19  ~MiniAODSVAnalyzer() override = default;
20 
21  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
22 
23 private:
24  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
25 
28  const double ptMin_;
29  const double etaMax_;
30 
32 
38 
39  // relation to jet
43 
48 };
49 
51  : jetToken_(consumes<std::vector<pat::Jet>>(pSet.getParameter<edm::InputTag>("JetTag"))),
52  svTagInfo_(pSet.getParameter<std::string>("svTagInfo")),
53  ptMin_(pSet.getParameter<double>("ptMin")),
54  etaMax_(pSet.getParameter<double>("etaMax")) {}
55 
57  ibook.setCurrentFolder("Btag/SV");
58 
59  n_sv_ = ibook.book1D("n_sv", "number of SV in jet", 5, 0, 5);
60  n_sv_->setAxisTitle("number of SV in jet");
61 
62  sv_mass_ = ibook.book1D("sv_mass", "SV mass", 30, 0., 6.);
63  sv_mass_->setAxisTitle("SV mass");
64 
65  sv_pt_ = ibook.book1D("sv_pt", "SV transverse momentum", 40, 0., 120.);
66  sv_pt_->setAxisTitle("SV pt");
67 
68  sv_ntracks_ = ibook.book1D("sv_ntracks", "SV number of daugthers", 10, 0, 10);
69  sv_ntracks_->setAxisTitle("number of tracks at SV");
70 
71  sv_chi2norm_ = ibook.book1D("sv_chi2norm", "normalized Chi2 of vertex", 30, 0, 15);
72  sv_chi2norm_->setAxisTitle("normalized Chi2 of SV");
73 
74  sv_chi2prob_ = ibook.book1D("sv_chi2prob", "Chi2 probability of vertex", 20, 0., 1.);
75  sv_chi2prob_->setAxisTitle("Chi2 probability of SV");
76 
77  sv_ptrel_ = ibook.book1D("sv_ptrel", "SV jet transverse momentum ratio", 25, 0., 1.);
78  sv_ptrel_->setAxisTitle("pt(SV)/pt(jet)");
79 
80  sv_energyratio_ = ibook.book1D("sv_energyratio", "SV jet energy ratio", 25, 0., 1.);
81  sv_energyratio_->setAxisTitle("E(SV)/E(jet)");
82 
83  sv_deltaR_ = ibook.book1D("sv_deltaR", "SV jet deltaR", 40, 0., 0.4);
84  sv_deltaR_->setAxisTitle("deltaR(jet, SV)");
85 
86  sv_dxy_ = ibook.book1D("sv_dxy", "2D flight distance", 40, 0., 8.);
87  sv_dxy_->setAxisTitle("dxy");
88 
89  sv_dxysig_ = ibook.book1D("sv_dxysig", "2D flight distance significance", 25, 0., 250.);
90  sv_dxysig_->setAxisTitle("dxy significance");
91 
92  sv_d3d_ = ibook.book1D("sv_d3d", "3D flight distance", 40, 0., 8.);
93  sv_d3d_->setAxisTitle("d3d");
94 
95  sv_d3dsig_ = ibook.book1D("sv_d3dsig", "3D flight distance significance", 25, 0., 250.);
96  sv_d3dsig_->setAxisTitle("d3d significance");
97 }
98 
101  iEvent.getByToken(jetToken_, jetCollection);
102 
103  // Loop over the pat::Jets
104  for (std::vector<pat::Jet>::const_iterator jet = jetCollection->begin(); jet != jetCollection->end(); ++jet) {
105  // jet selection
106  if (jet->hasTagInfo(svTagInfo_) && jet->pt() > ptMin_ && std::abs(jet->eta()) < etaMax_) {
107  const reco::CandSecondaryVertexTagInfo* taginfo =
108  static_cast<const reco::CandSecondaryVertexTagInfo*>(jet->tagInfo(svTagInfo_));
109  n_sv_->Fill(taginfo->nVertices());
110 
111  // loop secondary vertices
112  for (unsigned int i = 0; i < taginfo->nVertices(); i++) {
114 
115  sv_mass_->Fill(sv.mass());
116  sv_pt_->Fill(sv.pt());
117  sv_ntracks_->Fill(sv.numberOfDaughters());
118  sv_chi2norm_->Fill(sv.vertexNormalizedChi2());
119  sv_chi2prob_->Fill(ChiSquaredProbability(sv.vertexChi2(), sv.vertexNdof()));
120 
121  sv_ptrel_->Fill(sv.pt() / jet->pt());
122  sv_energyratio_->Fill(sv.energy() / jet->energy());
123  sv_deltaR_->Fill(reco::deltaR(sv, jet->momentum()));
124 
125  sv_dxy_->Fill(taginfo->flightDistance(i, 2).value());
126  sv_dxysig_->Fill(taginfo->flightDistance(i, 2).significance());
127  sv_d3d_->Fill(taginfo->flightDistance(i, 3).value());
128  sv_d3dsig_->Fill(taginfo->flightDistance(i, 3).significance());
129  }
130  }
131  }
132 }
133 
134 //define this as a plug-in
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
MonitorElement * sv_d3d_
MonitorElement * sv_dxy_
MiniAODSVAnalyzer(const edm::ParameterSet &pSet)
~MiniAODSVAnalyzer() override=default
Definition: HeavyIon.h:7
void Fill(long long x)
MonitorElement * sv_mass_
MonitorElement * sv_dxysig_
MonitorElement * sv_chi2norm_
MonitorElement * sv_ptrel_
int iEvent
Definition: GenABIO.cc:224
MonitorElement * sv_energyratio_
const edm::EDGetTokenT< std::vector< pat::Jet > > jetToken_
Definition: Jet.py:1
MonitorElement * sv_pt_
const std::string svTagInfo_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float ChiSquaredProbability(double chiSquared, double nrDOF)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
Measurement1D flightDistance(unsigned int index, int dim=0) const
MonitorElement * sv_deltaR_
MonitorElement * sv_chi2prob_
const VTX & secondaryVertex(unsigned int index) const
double value() const
Definition: Measurement1D.h:25
double significance() const
Definition: Measurement1D.h:29
HLT enums.
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * n_sv_
MonitorElement * sv_ntracks_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * sv_d3dsig_
Definition: Run.h:45
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)