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 
34  public:
38 
39  // virtual methods called from base class EDAnalyzer
40  virtual void beginJob() override;
41  virtual 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 
53  enum Flavour {
54  ALL_JETS = 0,
60  };
61 
62  TH1 * flavours_;
63 
64  // one group of plots per jet flavour;
65  struct Plots {
66  TH1 *nVertices;
67  TH1 *deltaR, *mass, *dist, *distErr, *distSig;
68  TH1 *nTracks, *chi2;
70 };
71 
73  jetsToken_(consumes<pat::JetCollection>(params.getParameter<edm::InputTag>("jets"))),
74  jetPtCut_(params.getParameter<double>("jetPtCut")),
75  jetEtaCut_(params.getParameter<double>("jetEtaCut"))
76 {
77 }
78 
80 {
81 }
82 
84 {
85  // retrieve handle to auxiliary service
86  // used for storing histograms into ROOT file
88 
89  flavours_ = fs->make<TH1F>("flavours", "jet flavours", 5, 0, 5);
90 
91  // book histograms for all jet flavours
92  for(unsigned int i = 0; i < N_JET_TYPES; i++) {
93  Plots &plots = plots_[i];
94  const char *flavour, *name;
95 
96  switch((Flavour)i) {
97  case ALL_JETS:
98  flavour = "all jets";
99  name = "all";
100  break;
101  case UDSG_JETS:
102  flavour = "light flavour jets";
103  name = "udsg";
104  break;
105  case C_JETS:
106  flavour = "charm jets";
107  name = "c";
108  break;
109  case B_JETS:
110  flavour = "bottom jets";
111  name = "b";
112  break;
113  default:
114  flavour = "unidentified jets";
115  name = "ni";
116  break;
117  }
118 
119  plots.nVertices = fs->make<TH1F>(Form("nVertices_%s", name),
120  Form("number of secondary vertices in %s", flavour),
121  5, 0, 5);
122  plots.deltaR = fs->make<TH1F>(Form("deltaR_%s", name),
123  Form("\\DeltaR between vertex direction and jet direction in %s", flavour),
124  100, 0, 0.5);
125  plots.mass = fs->make<TH1F>(Form("mass_%s", name),
126  Form("vertex mass in %s", flavour),
127  100, 0, 10);
128  plots.dist = fs->make<TH1F>(Form("dist_%s", name),
129  Form("Transverse distance between PV and SV in %s", flavour),
130  100, 0, 2);
131  plots.distErr = fs->make<TH1F>(Form("distErr_%s", name),
132  Form("Transverse distance error between PV and SV in %s", flavour),
133  100, 0, 0.5);
134  plots.distSig = fs->make<TH1F>(Form("distSig_%s", name),
135  Form("Transverse distance significance between PV and SV in %s", flavour),
136  100, 0, 50);
137  plots.nTracks = fs->make<TH1F>(Form("nTracks_%s", name),
138  Form("number of tracks at secondary vertex in %s", flavour),
139  20, 0, 20);
140  plots.chi2 = fs->make<TH1F>(Form("chi2_%s", name),
141  Form("secondary vertex fit \\chi^{2} in %s", flavour),
142  100, 0, 50);
143  }
144 }
145 
146 // helper function to sort the tracks by impact parameter significance
147 
149 {
150  // handle to the jets collection
152  event.getByToken(jetsToken_, jetsHandle);
153 
154  // now go through all jets
155  for(pat::JetCollection::const_iterator jet = jetsHandle->begin();
156  jet != jetsHandle->end(); ++jet) {
157 
158  // only look at jets that pass the pt and eta cut
159  if (jet->pt() < jetPtCut_ ||
160  std::abs(jet->eta()) > jetEtaCut_)
161  continue;
162 
164  // find out the jet flavour (differs between quark and anti-quark)
165  switch(std::abs(jet->partonFlavour())) {
166  case 1:
167  case 2:
168  case 3:
169  case 21:
170  flavour = UDSG_JETS;
171  break;
172  case 4:
173  flavour = C_JETS;
174  break;
175  case 5:
176  flavour = B_JETS;
177  break;
178  default:
179  flavour = NONID_JETS;
180  }
181 
182  // simply count the number of accepted jets
183  flavours_->Fill(ALL_JETS);
184  flavours_->Fill(flavour);
185 
186  // retrieve the "secondary vertex tag infos"
187  // this is the output of the b-tagging reconstruction code
188  // and contains secondary vertices in the jets
189  const reco::SecondaryVertexTagInfo &svTagInfo =
190  *jet->tagInfoSecondaryVertex();
191 
192  // count the number of secondary vertices
193  plots_[ALL_JETS].nVertices->Fill(svTagInfo.nVertices());
194  plots_[flavour].nVertices->Fill(svTagInfo.nVertices());
195 
196  // ignore jets without SV from now on
197  if (svTagInfo.nVertices() < 1)
198  continue;
199 
200  // pick the first secondary vertex (the "best" one)
201  const reco::Vertex &sv = svTagInfo.secondaryVertex(0);
202 
203  // and plot number of tracks and chi^2
204  plots_[ALL_JETS].nTracks->Fill(sv.tracksSize());
205  plots_[flavour].nTracks->Fill(sv.tracksSize());
206 
207  plots_[ALL_JETS].chi2->Fill(sv.chi2());
208  plots_[flavour].chi2->Fill(sv.chi2());
209 
210  // the precomputed transverse distance to the primary vertex
211  Measurement1D distance = svTagInfo.flightDistance(0, true);
212 
213  plots_[ALL_JETS].dist->Fill(distance.value());
214  plots_[flavour].dist->Fill(distance.value());
215 
216  plots_[ALL_JETS].distErr->Fill(distance.error());
217  plots_[flavour].distErr->Fill(distance.error());
218 
219  plots_[ALL_JETS].distSig->Fill(distance.significance());
220  plots_[flavour].distSig->Fill(distance.significance());
221 
222 
223  // the precomputed direction with respect to the primary vertex
224  GlobalVector dir = svTagInfo.flightDirection(0);
225 
226  // unfortunately CMSSW hsa all kinds of vectors,
227  // and sometimes we need to convert them *sigh*
228  math::XYZVector dir2(dir.x(), dir.y(), dir.z());
229 
230  // compute a few variables that we are plotting
232  jet->momentum(), dir2);
233 
234  plots_[ALL_JETS].deltaR->Fill(deltaR);
235  plots_[flavour].deltaR->Fill(deltaR);
236 
237  // compute the invariant mass from a four-vector sum
238  math::XYZTLorentzVector trackFourVectorSum;
239 
240  // loop over all tracks in the vertex
242  track != sv.tracks_end(); ++track) {
243  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec;
244  vec.SetPx((*track)->px());
245  vec.SetPy((*track)->py());
246  vec.SetPz((*track)->pz());
247  vec.SetM(0.13957); // pion mass
248  trackFourVectorSum += vec;
249  }
250 
251  // get the invariant mass: sqrt(E² - px² - py² - pz²)
252  double vertexMass = trackFourVectorSum.M();
253 
254  plots_[ALL_JETS].mass->Fill(vertexMass);
255  plots_[flavour].mass->Fill(vertexMass);
256  }
257 }
258 
260 
std::vector< Jet > JetCollection
Definition: Jet.h:53
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:81
const VTX & secondaryVertex(unsigned int index) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
T y() const
Definition: PV3DBase.h:63
double error() const
Definition: Measurement1D.h:30
edm::EDGetTokenT< pat::JetCollection > jetsToken_
virtual void analyze(const edm::Event &event, const edm::EventSetup &es) override
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
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
T z() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double chi2() const
chi-squares
Definition: Vertex.h:98
const GlobalVector & flightDirection(unsigned int index) const
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
double significance() const
Definition: Measurement1D.h:32
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
double value() const
Definition: Measurement1D.h:28
HLT enums.
Measurement1D flightDistance(unsigned int index, int dim=0) const
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
dbl *** dir
Definition: mlp_gen.cc:35
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76
T x() const
Definition: PV3DBase.h:62
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:71
Definition: event.py:1
virtual void beginJob() override