CMS 3D CMS Logo

PatBJetVertexAnalyzer.cc
Go to the documentation of this file.
1 #include <algorithm>
2 #include <iostream>
3 #include <cmath>
4 #include <vector>
5 #include <string>
6 
7 #include <TH1.h>
8 #include <TProfile.h>
9 
10 #include <Math/VectorUtil.h>
11 #include <Math/GenVector/PxPyPzE4D.h>
12 #include <Math/GenVector/PxPyPzM4D.h>
13 
21 
23 
32 
33 class PatBJetVertexAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
34 public:
37  ~PatBJetVertexAnalyzer() override;
38 
39  // virtual methods called from base class EDAnalyzer
40  void beginJob() override;
41  void analyze(const edm::Event &event, const edm::EventSetup &es) override;
42 
43 private:
44  // configuration parameters
46 
47  double jetPtCut_; // minimum (uncorrected) jet energy
48  double jetEtaCut_; // maximum |eta| for jet
49  double maxDeltaR_; // angle between jet and tracks
50 
51  // jet flavour constants
52 
54 
55  TH1 *flavours_;
56 
57  // one group of plots per jet flavour;
58  struct Plots {
59  TH1 *nVertices;
60  TH1 *deltaR, *mass, *dist, *distErr, *distSig;
61  TH1 *nTracks, *chi2;
63 };
64 
66  : jetsToken_(consumes<pat::JetCollection>(params.getParameter<edm::InputTag>("jets"))),
67  jetPtCut_(params.getParameter<double>("jetPtCut")),
68  jetEtaCut_(params.getParameter<double>("jetEtaCut")) {
69  usesResource(TFileService::kSharedResource);
70 }
71 
73 
75  // retrieve handle to auxiliary service
76  // used for storing histograms into ROOT file
78 
79  flavours_ = fs->make<TH1F>("flavours", "jet flavours", 5, 0, 5);
80 
81  // book histograms for all jet flavours
82  for (unsigned int i = 0; i < N_JET_TYPES; i++) {
83  Plots &plots = plots_[i];
84  const char *flavour, *name;
85 
86  switch ((Flavour)i) {
87  case ALL_JETS:
88  flavour = "all jets";
89  name = "all";
90  break;
91  case UDSG_JETS:
92  flavour = "light flavour jets";
93  name = "udsg";
94  break;
95  case C_JETS:
96  flavour = "charm jets";
97  name = "c";
98  break;
99  case B_JETS:
100  flavour = "bottom jets";
101  name = "b";
102  break;
103  default:
104  flavour = "unidentified jets";
105  name = "ni";
106  break;
107  }
108 
109  plots.nVertices =
110  fs->make<TH1F>(Form("nVertices_%s", name), Form("number of secondary vertices in %s", flavour), 5, 0, 5);
111  plots.deltaR = fs->make<TH1F>(Form("deltaR_%s", name),
112  Form("\\DeltaR between vertex direction and jet direction in %s", flavour),
113  100,
114  0,
115  0.5);
116  plots.mass = fs->make<TH1F>(Form("mass_%s", name), Form("vertex mass in %s", flavour), 100, 0, 10);
117  plots.dist =
118  fs->make<TH1F>(Form("dist_%s", name), Form("Transverse distance between PV and SV in %s", flavour), 100, 0, 2);
119  plots.distErr = fs->make<TH1F>(
120  Form("distErr_%s", name), Form("Transverse distance error between PV and SV in %s", flavour), 100, 0, 0.5);
121  plots.distSig = fs->make<TH1F>(
122  Form("distSig_%s", name), Form("Transverse distance significance between PV and SV in %s", flavour), 100, 0, 50);
123  plots.nTracks = fs->make<TH1F>(
124  Form("nTracks_%s", name), Form("number of tracks at secondary vertex in %s", flavour), 20, 0, 20);
125  plots.chi2 =
126  fs->make<TH1F>(Form("chi2_%s", name), Form("secondary vertex fit \\chi^{2} in %s", flavour), 100, 0, 50);
127  }
128 }
129 
130 // helper function to sort the tracks by impact parameter significance
131 
133  // handle to the jets collection
135  event.getByToken(jetsToken_, jetsHandle);
136 
137  // now go through all jets
138  for (pat::JetCollection::const_iterator jet = jetsHandle->begin(); jet != jetsHandle->end(); ++jet) {
139  // only look at jets that pass the pt and eta cut
140  if (jet->pt() < jetPtCut_ || std::abs(jet->eta()) > jetEtaCut_)
141  continue;
142 
144  // find out the jet flavour (differs between quark and anti-quark)
145  switch (std::abs(jet->partonFlavour())) {
146  case 1:
147  case 2:
148  case 3:
149  case 21:
150  flavour = UDSG_JETS;
151  break;
152  case 4:
153  flavour = C_JETS;
154  break;
155  case 5:
156  flavour = B_JETS;
157  break;
158  default:
160  }
161 
162  // simply count the number of accepted jets
163  flavours_->Fill(ALL_JETS);
164  flavours_->Fill(flavour);
165 
166  // retrieve the "secondary vertex tag infos"
167  // this is the output of the b-tagging reconstruction code
168  // and contains secondary vertices in the jets
169  const reco::SecondaryVertexTagInfo &svTagInfo = *jet->tagInfoSecondaryVertex();
170 
171  // count the number of secondary vertices
172  plots_[ALL_JETS].nVertices->Fill(svTagInfo.nVertices());
173  plots_[flavour].nVertices->Fill(svTagInfo.nVertices());
174 
175  // ignore jets without SV from now on
176  if (svTagInfo.nVertices() < 1)
177  continue;
178 
179  // pick the first secondary vertex (the "best" one)
180  const reco::Vertex &sv = svTagInfo.secondaryVertex(0);
181 
182  // and plot number of tracks and chi^2
183  plots_[ALL_JETS].nTracks->Fill(sv.tracksSize());
184  plots_[flavour].nTracks->Fill(sv.tracksSize());
185 
186  plots_[ALL_JETS].chi2->Fill(sv.chi2());
187  plots_[flavour].chi2->Fill(sv.chi2());
188 
189  // the precomputed transverse distance to the primary vertex
190  Measurement1D distance = svTagInfo.flightDistance(0, true);
191 
192  plots_[ALL_JETS].dist->Fill(distance.value());
193  plots_[flavour].dist->Fill(distance.value());
194 
195  plots_[ALL_JETS].distErr->Fill(distance.error());
196  plots_[flavour].distErr->Fill(distance.error());
197 
198  plots_[ALL_JETS].distSig->Fill(distance.significance());
199  plots_[flavour].distSig->Fill(distance.significance());
200 
201  // the precomputed direction with respect to the primary vertex
202  const GlobalVector &dir = svTagInfo.flightDirection(0);
203 
204  // unfortunately CMSSW hsa all kinds of vectors,
205  // and sometimes we need to convert them *sigh*
206  math::XYZVector dir2(dir.x(), dir.y(), dir.z());
207 
208  // compute a few variables that we are plotting
209  double deltaR = ROOT::Math::VectorUtil::DeltaR(jet->momentum(), dir2);
210 
211  plots_[ALL_JETS].deltaR->Fill(deltaR);
212  plots_[flavour].deltaR->Fill(deltaR);
213 
214  // compute the invariant mass from a four-vector sum
215  math::XYZTLorentzVector trackFourVectorSum;
216 
217  // loop over all tracks in the vertex
218  for (reco::Vertex::trackRef_iterator track = sv.tracks_begin(); track != sv.tracks_end(); ++track) {
219  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec;
220  vec.SetPx((*track)->px());
221  vec.SetPy((*track)->py());
222  vec.SetPz((*track)->pz());
223  vec.SetM(0.13957); // pion mass
224  trackFourVectorSum += vec;
225  }
226 
227  // get the invariant mass: sqrt(E² - px² - py² - pz²)
228  double vertexMass = trackFourVectorSum.M();
229 
230  plots_[ALL_JETS].mass->Fill(vertexMass);
231  plots_[flavour].mass->Fill(vertexMass);
232  }
233 }
234 
236 
static const std::string kSharedResource
Definition: TFileService.h:76
std::vector< Jet > JetCollection
Definition: Jet.h:53
edm::EDGetTokenT< pat::JetCollection > jetsToken_
void analyze(const edm::Event &event, const edm::EventSetup &es) override
Definition: HeavyIon.h:7
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
struct PatBJetVertexAnalyzer::Plots plots_[N_JET_TYPES]
PatBJetVertexAnalyzer(const edm::ParameterSet &params)
constructor and destructor
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
HLT enums.
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:38
Definition: event.py:1