CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PrimaryVertexAnalyzer4PUSlimmed.cc
Go to the documentation of this file.
2 
7 
8 // reco track and vertex
11 
12 // TrackingParticle
15 
16 // associator
18 
19 // DQM
21 
22 //
23 // constructors and destructor
24 //
26  const edm::ParameterSet& iConfig)
27  : verbose_(iConfig.getUntrackedParameter<bool>("verbose", false)),
28  use_only_charged_tracks_(iConfig.getUntrackedParameter<bool>(
29  "use_only_charged_tracks", true)),
30  use_TP_associator_(
31  iConfig.getUntrackedParameter<bool>("use_TP_associator", false)),
32  sigma_z_match_(
33  iConfig.getUntrackedParameter<double>("sigma_z_match", 3.0)),
34  abs_z_match_(
35  iConfig.getUntrackedParameter<double>("abs_z_match", 0.1)),
36  root_folder_(
37  iConfig.getUntrackedParameter<std::string>("root_folder",
38  "Validation/Vertices")),
39  vecPileupSummaryInfoToken_(consumes<std::vector<PileupSummaryInfo> >(
40  edm::InputTag(std::string("addPileupInfo")))),
41  recoTrackCollectionToken_(consumes<reco::TrackCollection>(edm::InputTag(
42  iConfig.getUntrackedParameter<std::string>("recoTrackProducer")))),
43  edmView_recoTrack_Token_(consumes<edm::View<reco::Track> >(edm::InputTag(
44  iConfig.getUntrackedParameter<std::string>("recoTrackProducer")))),
45  trackingParticleCollectionToken_(consumes<TrackingParticleCollection>(
46  edm::InputTag(std::string("mix"), std::string("MergedTrackTruth")))),
47  trackingVertexCollectionToken_(consumes<TrackingVertexCollection>(
48  edm::InputTag(std::string("mix"), std::string("MergedTrackTruth")))) {
49  reco_vertex_collections_ = iConfig.getParameter<std::vector<edm::InputTag> >(
50  "vertexRecoCollections");
51  for (auto const& l : reco_vertex_collections_) {
54  consumes<reco::VertexCollection>(l)));
55  }
56  if(use_TP_associator_) {
57  recoTrackToTrackingParticleAssociatorToken_ = consumes<reco::TrackToTrackingParticleAssociator>(edm::InputTag("trackAssociatorByHits"));
58  }
59 }
60 
62 
63 //
64 // member functions
65 //
67  DQMStore::IBooker& i, edm::Run const& iRun, edm::EventSetup const& iSetup) {
68  // TODO(rovere) make this booking method shorter and smarter,
69  // factorizing similar histograms with different prefix in a single
70  // method call.
71  float log_bins[31] = {
72  0.0, 0.0002, 0.0004, 0.0006, 0.0008, 0.001, 0.002,
73  0.004, 0.006, 0.008, 0.01, 0.02,
74  0.04, 0.06, 0.08, 0.1, 0.2,
75  0.4, 0.6, 0.8, 1.0, 2.0,
76  4.0, 6.0, 8.0, 10.0, 20.0,
77  40.0, 60.0, 80.0, 100.0
78  };
79  float log_mergez_bins[18] = {
80  0.0, 0.0025, 0.005, 0.0075, 0.01, 0.025, 0.05, 0.075, 0.1,
81  0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0
82  };
83  float log_pt2_bins[16] = {
84  0.0, 0.1, 0.5,
85  1.0, 2.0, 5.0,
86  10.0, 20.0, 50.0,
87  100.0, 200.0, 500.0,
88  1000.0, 2000.0, 5000.0,10000.0
89  };
90  float log_ntrk_bins[25] = {
91  0., 2.0, 4.0, 6.0, 8.0, 10.,
92  12.0, 14.0, 16.0, 18.0, 22.0,
93  26.0, 30.0, 35.0, 40.0,
94  45.0, 50.0, 55.0, 60.0, 70.0,
95  80.0, 90.0, 100.0, 150.0, 200.0
96  };
97  // TODO(rovere) Possibly change or add the main DQMStore booking
98  // interface to allow booking a TProfile with variable bin-width
99  // using an array of floats, as done for the TH1F case, not of
100  // doubles.
101  double log_pt2_bins_double[16] = {
102  0.0, 0.1, 0.5,
103  1.0, 2.0, 5.0,
104  10.0, 20.0, 50.0,
105  100.0, 200.0, 500.0,
106  1000.0, 2000.0, 5000.0,10000.0
107  };
108 
110  mes_["root_folder"]["GenVtx_vs_BX"] =
111  i.book2D("GenVtx_vs_BX", "GenVtx_vs_BX", 16, -12.5, 3.5, 200, 0., 200.);
112  // Generated Primary Vertex Plots
113  mes_["root_folder"]["GenPV_X"] =
114  i.book1D("GenPV_X", "GeneratedPV_X", 120, -0.6, 0.6);
115  mes_["root_folder"]["GenPV_Y"] =
116  i.book1D("GenPV_Y", "GeneratedPV_Y", 120, -0.6, 0.6);
117  mes_["root_folder"]["GenPV_Z"] =
118  i.book1D("GenPV_Z", "GeneratedPV_Z", 120, -60., 60.);
119  mes_["root_folder"]["GenPV_R"] =
120  i.book1D("GenPV_R", "GeneratedPV_R", 120, 0, 0.6);
121  mes_["root_folder"]["GenPV_Pt2"] =
122  i.book1D("GenPV_Pt2", "GeneratedPV_Sum-pt2", 15, &log_pt2_bins[0]);
123  mes_["root_folder"]["GenPV_NumTracks"] =
124  i.book1D("GenPV_NumTracks", "GeneratedPV_NumTracks", 24, &log_ntrk_bins[0]);
125  mes_["root_folder"]["GenPV_ClosestDistanceZ"] =
126  i.book1D("GenPV_ClosestDistanceZ", "GeneratedPV_ClosestDistanceZ", 30,
127  &log_bins[0]);
128 
129  // All Generated Vertices, used for efficiency plots
130  mes_["root_folder"]["GenAllV_NumVertices"] = i.book1D(
131  "GenAllV_NumVertices", "GeneratedAllV_NumVertices", 100, 0., 200.);
132  mes_["root_folder"]["GenAllV_X"] =
133  i.book1D("GenAllV_X", "GeneratedAllV_X", 120, -0.6, 0.6);
134  mes_["root_folder"]["GenAllV_Y"] =
135  i.book1D("GenAllV_Y", "GeneratedAllV_Y", 120, -0.6, 0.6);
136  mes_["root_folder"]["GenAllV_Z"] =
137  i.book1D("GenAllV_Z", "GeneratedAllV_Z", 120, -60, 60);
138  mes_["root_folder"]["GenAllV_R"] =
139  i.book1D("GenAllV_R", "GeneratedAllV_R", 120, 0, 0.6);
140  mes_["root_folder"]["GenAllV_Pt2"] =
141  i.book1D("GenAllV_Pt2", "GeneratedAllV_Sum-pt2", 15, &log_pt2_bins[0]);
142  mes_["root_folder"]["GenAllV_NumTracks"] =
143  i.book1D("GenAllV_NumTracks", "GeneratedAllV_NumTracks", 24, &log_ntrk_bins[0]);
144  mes_["root_folder"]["GenAllV_ClosestDistanceZ"] =
145  i.book1D("GenAllV_ClosestDistanceZ", "GeneratedAllV_ClosestDistanceZ", 30,
146  &log_bins[0]);
147  mes_["root_folder"]["GenAllV_PairDistanceZ"] =
148  i.book1D("GenAllV_PairDistanceZ", "GeneratedAllV_PairDistanceZ",
149  1000, 0, 20);
150  mes_["root_folder"]["SignalIsHighestPt2"] =
151  i.book1D("SignalIsHighestPt2", "SignalIsHighestPt2", 2, -0.5, 1.5);
152 
153  for (auto const& l : reco_vertex_collections_) {
154  std::string label = l.label();
155  std::string current_folder = root_folder_ + "/" + label;
156  i.setCurrentFolder(current_folder);
157 
158  mes_[label]["RecoVtx_vs_GenVtx"] = i.bookProfile(
159  "RecoVtx_vs_GenVtx", "RecoVtx_vs_GenVtx", 125, 0., 250., 250, 0., 250.);
160  mes_[label]["MatchedRecoVtx_vs_GenVtx"] =
161  i.bookProfile("MatchedRecoVtx_vs_GenVtx", "MatchedRecoVtx_vs_GenVtx",
162  125, 0., 250., 250, 0., 250.);
163  mes_[label]["KindOfSignalPV"] =
164  i.book1D("KindOfSignalPV", "KindOfSignalPV", 9, -0.5, 8.5);
165  mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(1, "!Highest!Assoc2Any");
166  mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(2, "Highest!Assoc2Any");
167  mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(3, "!HighestAssoc2First");
168  mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(4, "HighestAssoc2First");
169  mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(5, "!HighestAssoc2!First");
170  mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(6, "HighestAssoc2!First");
171  mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(7, "!HighestAssoc2First");
172  mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(8, "HighestAssoc2First");
173  mes_[label]["MisTagRate"] =
174  i.book1D("MisTagRate", "MisTagRate", 2, -0.5, 1.5);
175  mes_[label]["MisTagRate_vs_PU"] =
176  i.bookProfile("MisTagRate_vs_PU", "MisTagRate_vs_PU", 125, 0., 250.,
177  2, 0., 1.);
178  mes_[label]["MisTagRate_vs_sum-pt2"] =
179  i.bookProfile("MisTagRate_vs_sum-pt2", "MisTagRate_vs_sum-pt2",
180  15, &log_pt2_bins_double[0], 2, 0., 1.);
181  mes_[label]["MisTagRate_vs_Z"] =
182  i.bookProfile("MisTagRate_vs_Z", "MisTagRate_vs_Z",
183  120, -60., 60., 2, 0., 1.);
184  mes_[label]["MisTagRate_vs_R"] =
185  i.bookProfile("MisTagRate_vs_R", "MisTagRate_vs_R",
186  120, 0., 0.6, 2, 0., 1.);
187  mes_[label]["MisTagRate_vs_NumTracks"] =
188  i.bookProfile("MisTagRate_vs_NumTracks", "MisTagRate_vs_NumTracks",
189  100, 0., 200, 2, 0., 1.);
190  mes_[label]["MisTagRateSignalIsHighest"] =
191  i.book1D("MisTagRateSignalIsHighest",
192  "MisTagRateSignalIsHighest", 2, -0.5, 1.5);
193  mes_[label]["MisTagRateSignalIsHighest_vs_PU"] =
194  i.bookProfile("MisTagRateSignalIsHighest_vs_PU",
195  "MisTagRateSignalIsHighest_vs_PU", 125, 0., 250.,
196  2, 0., 1.);
197  mes_[label]["MisTagRateSignalIsHighest_vs_sum-pt2"] =
198  i.bookProfile("MisTagRateSignalIsHighest_vs_sum-pt2",
199  "MisTagRateSignalIsHighest_vs_sum-pt2",
200  15, &log_pt2_bins_double[0], 2, 0., 1.);
201  mes_[label]["MisTagRateSignalIsHighest_vs_Z"] =
202  i.bookProfile("MisTagRateSignalIsHighest_vs_Z",
203  "MisTagRateSignalIsHighest_vs_Z",
204  120, -60., 60., 2, 0., 1.);
205  mes_[label]["MisTagRateSignalIsHighest_vs_R"] =
206  i.bookProfile("MisTagRateSignalIsHighest_vs_R",
207  "MisTagRateSignalIsHighest_vs_R",
208  120, 0., 0.6, 2, 0., 1.);
209  mes_[label]["MisTagRateSignalIsHighest_vs_NumTracks"] =
210  i.bookProfile("MisTagRateSignalIsHighest_vs_NumTracks",
211  "MisTagRateSignalIsHighest_vs_NumTracks",
212  100, 0., 200, 2, 0., 1.);
213  mes_[label]["MisTagRateSignalIsNotHighest"] =
214  i.book1D("MisTagRateSignalIsNotHighest",
215  "MisTagRateSignalIsNotHighest", 2, -0.5, 1.5);
216  mes_[label]["MisTagRateSignalIsNotHighest_vs_PU"] =
217  i.bookProfile("MisTagRateSignalIsNotHighest_vs_PU",
218  "MisTagRateSignalIsNotHighest_vs_PU", 125, 0., 250.,
219  2, 0., 1.);
220  mes_[label]["MisTagRateSignalIsNotHighest_vs_sum-pt2"] =
221  i.bookProfile("MisTagRateSignalIsNotHighest_vs_sum-pt2",
222  "MisTagRateSignalIsNotHighest_vs_sum-pt2",
223  15, &log_pt2_bins_double[0], 2, 0., 1.);
224  mes_[label]["MisTagRateSignalIsNotHighest_vs_Z"] =
225  i.bookProfile("MisTagRateSignalIsNotHighest_vs_Z",
226  "MisTagRateSignalIsNotHighest_vs_Z",
227  120, -60., 60., 2, 0., 1.);
228  mes_[label]["MisTagRateSignalIsNotHighest_vs_R"] =
229  i.bookProfile("MisTagRateSignalIsNotHighest_vs_R",
230  "MisTagRateSignalIsNotHighest_vs_R",
231  120, 0., 0.6, 2, 0., 1.);
232  mes_[label]["MisTagRateSignalIsNotHighest_vs_NumTracks"] =
233  i.bookProfile("MisTagRateSignalIsNotHighest_vs_NumTracks",
234  "MisTagRateSignalIsNotHighest_vs_NumTracks",
235  100, 0., 200, 2, 0., 1.);
236  mes_[label]["TruePVLocationIndex"] =
237  i.book1D("TruePVLocationIndex",
238  "TruePVLocationIndexInRecoVertexCollection", 12, -1.5, 10.5);
239  mes_[label]["TruePVLocationIndexCumulative"] =
240  i.book1D("TruePVLocationIndexCumulative",
241  "TruePVLocationIndexInRecoVertexCollectionCumulative",
242  3, -1.5, 1.5);
243  mes_[label]["TruePVLocationIndexSignalIsHighest"] =
244  i.book1D("TruePVLocationIndexSignalIsHighest",
245  "TruePVLocationIndexSignalIsHighestInRecoVertexCollection",
246  12, -1.5, 10.5);
247  mes_[label]["TruePVLocationIndexSignalIsNotHighest"] =
248  i.book1D("TruePVLocationIndexSignalIsNotHighest",
249  "TruePVLocationIndexSignalIsNotHighestInRecoVertexCollection",
250  12, -1.5, 10.5);
251  // All Generated Vertices. Used for Efficiency plots We kind of
252  // duplicate plots here in case we want to perform more detailed
253  // studies on a selection of generated vertices, not on all of them.
254  mes_[label]["GenAllAssoc2Reco_NumVertices"] =
255  i.book1D("GenAllAssoc2Reco_NumVertices",
256  "GeneratedAllAssoc2Reco_NumVertices", 100, 0., 200.);
257  mes_[label]["GenAllAssoc2Reco_X"] = i.book1D(
258  "GenAllAssoc2Reco_X", "GeneratedAllAssoc2Reco_X", 120, -0.6, 0.6);
259  mes_[label]["GenAllAssoc2Reco_Y"] = i.book1D(
260  "GenAllAssoc2Reco_Y", "GeneratedAllAssoc2Reco_Y", 120, -0.6, 0.6);
261  mes_[label]["GenAllAssoc2Reco_Z"] = i.book1D(
262  "GenAllAssoc2Reco_Z", "GeneratedAllAssoc2Reco_Z", 120, -60, 60);
263  mes_[label]["GenAllAssoc2Reco_R"] =
264  i.book1D("GenAllAssoc2Reco_R", "GeneratedAllAssoc2Reco_R", 120, 0, 0.6);
265  mes_[label]["GenAllAssoc2Reco_Pt2"] =
266  i.book1D("GenAllAssoc2Reco_Pt2", "GeneratedAllAssoc2Reco_Sum-pt2", 15,
267  &log_pt2_bins[0]);
268  mes_[label]["GenAllAssoc2Reco_NumTracks"] =
269  i.book1D("GenAllAssoc2Reco_NumTracks",
270  "GeneratedAllAssoc2Reco_NumTracks", 24, &log_ntrk_bins[0]);
271  mes_[label]["GenAllAssoc2Reco_ClosestDistanceZ"] =
272  i.book1D("GenAllAssoc2Reco_ClosestDistanceZ",
273  "GeneratedAllAssoc2Reco_ClosestDistanceZ", 30, &log_bins[0]);
274 
275  // All Generated Vertices Matched to a Reconstructed vertex. Used
276  // for Efficiency plots
277  mes_[label]["GenAllAssoc2RecoMatched_NumVertices"] =
278  i.book1D("GenAllAssoc2RecoMatched_NumVertices",
279  "GeneratedAllAssoc2RecoMatched_NumVertices", 100, 0., 200.);
280  mes_[label]["GenAllAssoc2RecoMatched_X"] =
281  i.book1D("GenAllAssoc2RecoMatched_X", "GeneratedAllAssoc2RecoMatched_X",
282  120, -0.6, 0.6);
283  mes_[label]["GenAllAssoc2RecoMatched_Y"] =
284  i.book1D("GenAllAssoc2RecoMatched_Y", "GeneratedAllAssoc2RecoMatched_Y",
285  120, -0.6, 0.6);
286  mes_[label]["GenAllAssoc2RecoMatched_Z"] =
287  i.book1D("GenAllAssoc2RecoMatched_Z", "GeneratedAllAssoc2RecoMatched_Z",
288  120, -60, 60);
289  mes_[label]["GenAllAssoc2RecoMatched_R"] =
290  i.book1D("GenAllAssoc2RecoMatched_R", "GeneratedAllAssoc2RecoMatched_R",
291  120, 0, 0.6);
292  mes_[label]["GenAllAssoc2RecoMatched_Pt2"] =
293  i.book1D("GenAllAssoc2RecoMatched_Pt2",
294  "GeneratedAllAssoc2RecoMatched_Sum-pt2", 15, &log_pt2_bins[0]);
295  mes_[label]["GenAllAssoc2RecoMatched_NumTracks"] =
296  i.book1D("GenAllAssoc2RecoMatched_NumTracks",
297  "GeneratedAllAssoc2RecoMatched_NumTracks", 24, &log_ntrk_bins[0]);
298  mes_[label]["GenAllAssoc2RecoMatched_ClosestDistanceZ"] = i.book1D(
299  "GenAllAssoc2RecoMatched_ClosestDistanceZ",
300  "GeneratedAllAssoc2RecoMatched_ClosestDistanceZ", 30, &log_bins[0]);
301 
302  // All Generated Vertices Multi-Matched to a Reconstructed vertex. Used
303  // for Duplicate rate plots
304  mes_[label]["GenAllAssoc2RecoMultiMatched_NumVertices"] = i.book1D(
305  "GenAllAssoc2RecoMultiMatched_NumVertices",
306  "GeneratedAllAssoc2RecoMultiMatched_NumVertices", 100, 0., 200.);
307  mes_[label]["GenAllAssoc2RecoMultiMatched_X"] =
308  i.book1D("GenAllAssoc2RecoMultiMatched_X",
309  "GeneratedAllAssoc2RecoMultiMatched_X", 120, -0.6, 0.6);
310  mes_[label]["GenAllAssoc2RecoMultiMatched_Y"] =
311  i.book1D("GenAllAssoc2RecoMultiMatched_Y",
312  "GeneratedAllAssoc2RecoMultiMatched_Y", 120, -0.6, 0.6);
313  mes_[label]["GenAllAssoc2RecoMultiMatched_Z"] =
314  i.book1D("GenAllAssoc2RecoMultiMatched_Z",
315  "GeneratedAllAssoc2RecoMultiMatched_Z", 120, -60, 60);
316  mes_[label]["GenAllAssoc2RecoMultiMatched_R"] =
317  i.book1D("GenAllAssoc2RecoMultiMatched_R",
318  "GeneratedAllAssoc2RecoMultiMatched_R", 120, 0, 0.6);
319  mes_[label]["GenAllAssoc2RecoMultiMatched_Pt2"] =
320  i.book1D("GenAllAssoc2RecoMultiMatched_Pt2",
321  "GeneratedAllAssoc2RecoMultiMatched_Sum-pt2",
322  15, &log_pt2_bins[0]);
323  mes_[label]["GenAllAssoc2RecoMultiMatched_NumTracks"] =
324  i.book1D("GenAllAssoc2RecoMultiMatched_NumTracks",
325  "GeneratedAllAssoc2RecoMultiMatched_NumTracks", 24, &log_ntrk_bins[0]);
326  mes_[label]["GenAllAssoc2RecoMultiMatched_ClosestDistanceZ"] = i.book1D(
327  "GenAllAssoc2RecoMultiMatched_ClosestDistanceZ",
328  "GeneratedAllAssoc2RecoMultiMatched_ClosestDistanceZ",
329  30, &log_bins[0]);
330 
331  // All Reco Vertices. Used for {Fake,Duplicate}-Rate plots
332  mes_[label]["RecoAllAssoc2Gen_NumVertices"] =
333  i.book1D("RecoAllAssoc2Gen_NumVertices",
334  "ReconstructedAllAssoc2Gen_NumVertices", 100, 0., 200.);
335  mes_[label]["RecoAllAssoc2Gen_X"] = i.book1D(
336  "RecoAllAssoc2Gen_X", "ReconstructedAllAssoc2Gen_X", 120, -0.6, 0.6);
337  mes_[label]["RecoAllAssoc2Gen_Y"] = i.book1D(
338  "RecoAllAssoc2Gen_Y", "ReconstructedAllAssoc2Gen_Y", 120, -0.6, 0.6);
339  mes_[label]["RecoAllAssoc2Gen_Z"] = i.book1D(
340  "RecoAllAssoc2Gen_Z", "ReconstructedAllAssoc2Gen_Z", 120, -60, 60);
341  mes_[label]["RecoAllAssoc2Gen_R"] = i.book1D(
342  "RecoAllAssoc2Gen_R", "ReconstructedAllAssoc2Gen_R", 120, 0, 0.6);
343  mes_[label]["RecoAllAssoc2Gen_Pt2"] =
344  i.book1D("RecoAllAssoc2Gen_Pt2", "ReconstructedAllAssoc2Gen_Sum-pt2",
345  15, &log_pt2_bins[0]);
346  mes_[label]["RecoAllAssoc2Gen_Ndof"] =
347  i.book1D("RecoAllAssoc2Gen_Ndof",
348  "ReconstructedAllAssoc2Gen_Ndof", 120, 0., 240.);
349  mes_[label]["RecoAllAssoc2Gen_NumTracks"] =
350  i.book1D("RecoAllAssoc2Gen_NumTracks",
351  "ReconstructedAllAssoc2Gen_NumTracks", 24, &log_ntrk_bins[0]);
352  mes_[label]["RecoAllAssoc2Gen_PU"] =
353  i.book1D("RecoAllAssoc2Gen_PU",
354  "ReconstructedAllAssoc2Gen_PU", 125, 0., 250.);
355  mes_[label]["RecoAllAssoc2Gen_ClosestDistanceZ"] =
356  i.book1D("RecoAllAssoc2Gen_ClosestDistanceZ",
357  "ReconstructedAllAssoc2Gen_ClosestDistanceZ",
358  30, &log_bins[0]);
359  mes_[label]["RecoAllAssoc2GenProperties"] =
360  i.book1D("RecoAllAssoc2GenProperties",
361  "ReconstructedAllAssoc2Gen_Properties", 8, -0.5, 7.5);
362  mes_[label]["RecoAllAssoc2Gen_PairDistanceZ"] =
363  i.book1D("RecoAllAssoc2Gen_PairDistanceZ",
364  "RecoAllAssoc2Gen_PairDistanceZ", 1000, 0, 20);
365 
366  // All Reconstructed Vertices Matched to a Generated vertex. Used
367  // for Fake-Rate plots
368  mes_[label]["RecoAllAssoc2GenMatched_NumVertices"] =
369  i.book1D("RecoAllAssoc2GenMatched_NumVertices",
370  "ReconstructedAllAssoc2GenMatched_NumVertices", 100, 0., 200.);
371  mes_[label]["RecoAllAssoc2GenMatched_X"] =
372  i.book1D("RecoAllAssoc2GenMatched_X",
373  "ReconstructedAllAssoc2GenMatched_X", 120, -0.6, 0.6);
374  mes_[label]["RecoAllAssoc2GenMatched_Y"] =
375  i.book1D("RecoAllAssoc2GenMatched_Y",
376  "ReconstructedAllAssoc2GenMatched_Y", 120, -0.6, 0.6);
377  mes_[label]["RecoAllAssoc2GenMatched_Z"] =
378  i.book1D("RecoAllAssoc2GenMatched_Z",
379  "ReconstructedAllAssoc2GenMatched_Z", 120, -60, 60);
380  mes_[label]["RecoAllAssoc2GenMatched_R"] =
381  i.book1D("RecoAllAssoc2GenMatched_R",
382  "ReconstructedAllAssoc2GenMatched_R", 120, 0, 0.6);
383  mes_[label]["RecoAllAssoc2GenMatched_Pt2"] =
384  i.book1D("RecoAllAssoc2GenMatched_Pt2",
385  "ReconstructedAllAssoc2GenMatched_Sum-pt2",
386  15, &log_pt2_bins[0]);
387  mes_[label]["RecoAllAssoc2GenMatched_Ndof"] =
388  i.book1D("RecoAllAssoc2GenMatched_Ndof",
389  "ReconstructedAllAssoc2GenMatched_Ndof", 120, 0., 240.);
390  mes_[label]["RecoAllAssoc2GenMatched_NumTracks"] =
391  i.book1D("RecoAllAssoc2GenMatched_NumTracks",
392  "ReconstructedAllAssoc2GenMatched_NumTracks", 24, &log_ntrk_bins[0]);
393  mes_[label]["RecoAllAssoc2GenMatched_PU"] =
394  i.book1D("RecoAllAssoc2GenMatched_PU",
395  "ReconstructedAllAssoc2GenMatched_PU", 125, 0., 250.);
396  mes_[label]["RecoAllAssoc2GenMatched_ClosestDistanceZ"] = i.book1D(
397  "RecoAllAssoc2GenMatched_ClosestDistanceZ",
398  "ReconstructedAllAssoc2GenMatched_ClosestDistanceZ", 30, &log_bins[0]);
399 
400  // All Reconstructed Vertices Multi-Matched to a Generated vertex. Used
401  // for Merge-Rate plots
402  mes_[label]["RecoAllAssoc2GenMultiMatched_NumVertices"] = i.book1D(
403  "RecoAllAssoc2GenMultiMatched_NumVertices",
404  "ReconstructedAllAssoc2GenMultiMatched_NumVertices", 100, 0., 200.);
405  mes_[label]["RecoAllAssoc2GenMultiMatched_X"] =
406  i.book1D("RecoAllAssoc2GenMultiMatched_X",
407  "ReconstructedAllAssoc2GenMultiMatched_X", 120, -0.6, 0.6);
408  mes_[label]["RecoAllAssoc2GenMultiMatched_Y"] =
409  i.book1D("RecoAllAssoc2GenMultiMatched_Y",
410  "ReconstructedAllAssoc2GenMultiMatched_Y", 120, -0.6, 0.6);
411  mes_[label]["RecoAllAssoc2GenMultiMatched_Z"] =
412  i.book1D("RecoAllAssoc2GenMultiMatched_Z",
413  "ReconstructedAllAssoc2GenMultiMatched_Z", 120, -60, 60);
414  mes_[label]["RecoAllAssoc2GenMultiMatched_R"] =
415  i.book1D("RecoAllAssoc2GenMultiMatched_R",
416  "ReconstructedAllAssoc2GenMultiMatched_R", 120, 0, 0.6);
417  mes_[label]["RecoAllAssoc2GenMultiMatched_Pt2"] = i.book1D(
418  "RecoAllAssoc2GenMultiMatched_Pt2",
419  "ReconstructedAllAssoc2GenMultiMatched_Sum-pt2", 15, &log_pt2_bins[0]);
420  mes_[label]["RecoAllAssoc2GenMultiMatched_NumTracks"] = i.book1D(
421  "RecoAllAssoc2GenMultiMatched_NumTracks",
422  "ReconstructedAllAssoc2GenMultiMatched_NumTracks", 24, &log_ntrk_bins[0]);
423  mes_[label]["RecoAllAssoc2GenMultiMatched_PU"] =
424  i.book1D("RecoAllAssoc2GenMultiMatched_PU",
425  "ReconstructedAllAssoc2GenMultiMatched_PU", 125, 0., 250.);
426  mes_[label]["RecoAllAssoc2GenMultiMatched_ClosestDistanceZ"] =
427  i.book1D("RecoAllAssoc2GenMultiMatched_ClosestDistanceZ",
428  "ReconstructedAllAssoc2GenMultiMatched_ClosestDistanceZ",
429  17, &log_mergez_bins[0]);
430 
431  // All Reconstructed Vertices Matched to a Multi-Matched Gen
432  // Vertex. Used for Duplicate rate plots done w.r.t. Reco
433  // Quantities. We basically want to ask how many times a RecoVTX
434  // has been reconstructed and associated to a SimulatedVTX that
435  // has been linked to at least another RecoVTX. In this sense this
436  // RecoVTX is a duplicate of the same, real GenVTX.
437  mes_[label]["RecoAllAssoc2MultiMatchedGen_NumVertices"] = i.book1D(
438  "RecoAllAssoc2MultiMatchedGen_NumVertices",
439  "RecoAllAssoc2MultiMatchedGen_NumVertices", 100, 0., 200.);
440  mes_[label]["RecoAllAssoc2MultiMatchedGen_X"] =
441  i.book1D("RecoAllAssoc2MultiMatchedGen_X",
442  "RecoAllAssoc2MultiMatchedGen_X", 120, -0.6, 0.6);
443  mes_[label]["RecoAllAssoc2MultiMatchedGen_Y"] =
444  i.book1D("RecoAllAssoc2MultiMatchedGen_Y",
445  "RecoAllAssoc2MultiMatchedGen_Y", 120, -0.6, 0.6);
446  mes_[label]["RecoAllAssoc2MultiMatchedGen_Z"] =
447  i.book1D("RecoAllAssoc2MultiMatchedGen_Z",
448  "RecoAllAssoc2MultiMatchedGen_Z", 120, -60, 60);
449  mes_[label]["RecoAllAssoc2MultiMatchedGen_R"] =
450  i.book1D("RecoAllAssoc2MultiMatchedGen_R",
451  "RecoAllAssoc2MultiMatchedGen_R", 120, 0, 0.6);
452  mes_[label]["RecoAllAssoc2MultiMatchedGen_Pt2"] =
453  i.book1D("RecoAllAssoc2MultiMatchedGen_Pt2",
454  "RecoAllAssoc2MultiMatchedGen_Sum-pt2", 15, &log_pt2_bins[0]);
455  mes_[label]["RecoAllAssoc2MultiMatchedGen_NumTracks"] =
456  i.book1D("RecoAllAssoc2MultiMatchedGen_NumTracks",
457  "RecoAllAssoc2MultiMatchedGen_NumTracks", 24, &log_ntrk_bins[0]);
458  mes_[label]["RecoAllAssoc2MultiMatchedGen_PU"] =
459  i.book1D("RecoAllAssoc2MultiMatchedGen_PU",
460  "RecoAllAssoc2MultiMatchedGen_PU", 125, 0., 250.);
461  mes_[label]["RecoAllAssoc2MultiMatchedGen_ClosestDistanceZ"] = i.book1D(
462  "RecoAllAssoc2MultiMatchedGen_ClosestDistanceZ",
463  "RecoAllAssoc2MultiMatchedGen_ClosestDistanceZ", 30, &log_bins[0]);
464  mes_[label]["RecoAllAssoc2GenSimForMerge_ClosestDistanceZ"] = i.book1D(
465  "RecoAllAssoc2GenSimForMerge_ClosestDistanceZ",
466  "RecoAllAssoc2GenSimForMerge_ClosestDistanceZ",
467  17, &log_mergez_bins[0]);
468  }
469 }
470 
472  const simPrimaryVertex& v) {
473  if (v.eventId.event() == 0) {
474  mes_["root_folder"]["GenPV_X"]->Fill(v.x);
475  mes_["root_folder"]["GenPV_Y"]->Fill(v.y);
476  mes_["root_folder"]["GenPV_Z"]->Fill(v.z);
477  mes_["root_folder"]["GenPV_R"]->Fill(v.r);
478  mes_["root_folder"]["GenPV_Pt2"]->Fill(v.ptsq);
479  mes_["root_folder"]["GenPV_NumTracks"]->Fill(v.nGenTrk);
480  if (v.closest_vertex_distance_z > 0.)
481  mes_["root_folder"]["GenPV_ClosestDistanceZ"]
482  ->Fill(v.closest_vertex_distance_z);
483  }
484  mes_["root_folder"]["GenAllV_X"]->Fill(v.x);
485  mes_["root_folder"]["GenAllV_Y"]->Fill(v.y);
486  mes_["root_folder"]["GenAllV_Z"]->Fill(v.z);
487  mes_["root_folder"]["GenAllV_R"]->Fill(v.r);
488  mes_["root_folder"]["GenAllV_Pt2"]->Fill(v.ptsq);
489  mes_["root_folder"]["GenAllV_NumTracks"]->Fill(v.nGenTrk);
490  if (v.closest_vertex_distance_z > 0.)
491  mes_["root_folder"]["GenAllV_ClosestDistanceZ"]
492  ->Fill(v.closest_vertex_distance_z);
493 }
494 
496  const std::string& label,
498  mes_[label]["GenAllAssoc2Reco_X"]->Fill(v.x);
499  mes_[label]["GenAllAssoc2Reco_Y"]->Fill(v.y);
500  mes_[label]["GenAllAssoc2Reco_Z"]->Fill(v.z);
501  mes_[label]["GenAllAssoc2Reco_R"]->Fill(v.r);
502  mes_[label]["GenAllAssoc2Reco_Pt2"]->Fill(v.ptsq);
503  mes_[label]["GenAllAssoc2Reco_NumTracks"]->Fill(v.nGenTrk);
504  if (v.closest_vertex_distance_z > 0.)
505  mes_[label]["GenAllAssoc2Reco_ClosestDistanceZ"]
506  ->Fill(v.closest_vertex_distance_z);
507  if (v.rec_vertices.size()) {
508  mes_[label]["GenAllAssoc2RecoMatched_X"]->Fill(v.x);
509  mes_[label]["GenAllAssoc2RecoMatched_Y"]->Fill(v.y);
510  mes_[label]["GenAllAssoc2RecoMatched_Z"]->Fill(v.z);
511  mes_[label]["GenAllAssoc2RecoMatched_R"]->Fill(v.r);
512  mes_[label]["GenAllAssoc2RecoMatched_Pt2"]->Fill(v.ptsq);
513  mes_[label]["GenAllAssoc2RecoMatched_NumTracks"]->Fill(v.nGenTrk);
514  if (v.closest_vertex_distance_z > 0.)
515  mes_[label]["GenAllAssoc2RecoMatched_ClosestDistanceZ"]
516  ->Fill(v.closest_vertex_distance_z);
517  }
518  if (v.rec_vertices.size() > 1) {
519  mes_[label]["GenAllAssoc2RecoMultiMatched_X"]->Fill(v.x);
520  mes_[label]["GenAllAssoc2RecoMultiMatched_Y"]->Fill(v.y);
521  mes_[label]["GenAllAssoc2RecoMultiMatched_Z"]->Fill(v.z);
522  mes_[label]["GenAllAssoc2RecoMultiMatched_R"]->Fill(v.r);
523  mes_[label]["GenAllAssoc2RecoMultiMatched_Pt2"]->Fill(v.ptsq);
524  mes_[label]["GenAllAssoc2RecoMultiMatched_NumTracks"]->Fill(v.nGenTrk);
525  if (v.closest_vertex_distance_z > 0.)
526  mes_[label]["GenAllAssoc2RecoMultiMatched_ClosestDistanceZ"]
527  ->Fill(v.closest_vertex_distance_z);
528  }
529 }
530 
532  const std::string& label,
533  int num_pileup_vertices,
535  mes_[label]["RecoAllAssoc2Gen_X"]->Fill(v.x);
536  mes_[label]["RecoAllAssoc2Gen_Y"]->Fill(v.y);
537  mes_[label]["RecoAllAssoc2Gen_Z"]->Fill(v.z);
538  mes_[label]["RecoAllAssoc2Gen_R"]->Fill(v.r);
539  mes_[label]["RecoAllAssoc2Gen_Pt2"]->Fill(v.ptsq);
540  mes_[label]["RecoAllAssoc2Gen_Ndof"]->Fill(v.recVtx->ndof());
541  mes_[label]["RecoAllAssoc2Gen_NumTracks"]->Fill(v.nRecoTrk);
542  mes_[label]["RecoAllAssoc2Gen_PU"]->Fill(num_pileup_vertices);
543  if (v.closest_vertex_distance_z > 0.)
544  mes_[label]["RecoAllAssoc2Gen_ClosestDistanceZ"]
545  ->Fill(v.closest_vertex_distance_z);
546  if (v.sim_vertices.size()) {
548  mes_[label]["RecoAllAssoc2GenMatched_X"]->Fill(v.x);
549  mes_[label]["RecoAllAssoc2GenMatched_Y"]->Fill(v.y);
550  mes_[label]["RecoAllAssoc2GenMatched_Z"]->Fill(v.z);
551  mes_[label]["RecoAllAssoc2GenMatched_R"]->Fill(v.r);
552  mes_[label]["RecoAllAssoc2GenMatched_Pt2"]->Fill(v.ptsq);
553  mes_[label]["RecoAllAssoc2GenMatched_Ndof"]->Fill(v.recVtx->ndof());
554  mes_[label]["RecoAllAssoc2GenMatched_NumTracks"]->Fill(v.nRecoTrk);
555  mes_[label]["RecoAllAssoc2GenMatched_PU"]->Fill(num_pileup_vertices);
556  if (v.closest_vertex_distance_z > 0.)
557  mes_[label]["RecoAllAssoc2GenMatched_ClosestDistanceZ"]
558  ->Fill(v.closest_vertex_distance_z);
559  // Now keep track of all RecoVTX associated to a SimVTX that
560  // itself is associated to more than one RecoVTX, for
561  // duplicate-rate plots on reco quantities.
562  if (v.sim_vertices_internal[0]->rec_vertices.size() > 1) {
564  mes_[label]["RecoAllAssoc2MultiMatchedGen_X"]->Fill(v.x);
565  mes_[label]["RecoAllAssoc2MultiMatchedGen_Y"]->Fill(v.y);
566  mes_[label]["RecoAllAssoc2MultiMatchedGen_Z"]->Fill(v.z);
567  mes_[label]["RecoAllAssoc2MultiMatchedGen_R"]->Fill(v.r);
568  mes_[label]["RecoAllAssoc2MultiMatchedGen_Pt2"]->Fill(v.ptsq);
569  mes_[label]["RecoAllAssoc2MultiMatchedGen_NumTracks"]->Fill(v.nRecoTrk);
570  mes_[label]["RecoAllAssoc2MultiMatchedGen_PU"]->Fill(num_pileup_vertices);
571  if (v.closest_vertex_distance_z > 0.)
572  mes_[label]["RecoAllAssoc2MultiMatchedGen_ClosestDistanceZ"]
573  ->Fill(v.closest_vertex_distance_z);
574  }
575  // This is meant to be used as "denominator" for the merge-rate
576  // plots produced starting from reco quantities. We enter here
577  // only if the reco vertex has been associated, since we need info
578  // from the SimVTX associated to it. In this regard, the final
579  // merge-rate plot coming from reco is not to be intended as a
580  // pure efficiency-like plot, since the normalization is biased.
581  if (v.sim_vertices_internal[0]->closest_vertex_distance_z > 0.)
582  mes_[label]["RecoAllAssoc2GenSimForMerge_ClosestDistanceZ"]
583  ->Fill(v.sim_vertices_internal[0]->closest_vertex_distance_z);
584  }
585  // this plots are meant to be used to compute the merge rate
586  if (v.sim_vertices.size() > 1) {
588  mes_[label]["RecoAllAssoc2GenMultiMatched_X"]->Fill(v.x);
589  mes_[label]["RecoAllAssoc2GenMultiMatched_Y"]->Fill(v.y);
590  mes_[label]["RecoAllAssoc2GenMultiMatched_Z"]->Fill(v.z);
591  mes_[label]["RecoAllAssoc2GenMultiMatched_R"]->Fill(v.r);
592  mes_[label]["RecoAllAssoc2GenMultiMatched_Pt2"]->Fill(v.ptsq);
593  mes_[label]["RecoAllAssoc2GenMultiMatched_NumTracks"]->Fill(v.nRecoTrk);
594  mes_[label]["RecoAllAssoc2GenMultiMatched_PU"]->Fill(num_pileup_vertices);
595  if (v.sim_vertices_internal[0]->closest_vertex_distance_z > 0.)
596  mes_[label]["RecoAllAssoc2GenMultiMatched_ClosestDistanceZ"]
597  ->Fill(v.sim_vertices_internal[0]->closest_vertex_distance_z);
598  }
599  mes_[label]["RecoAllAssoc2GenProperties"]->Fill(v.kind_of_vertex);
600 }
601 
602 /* Extract information form TrackingParticles/TrackingVertex and fill
603  * the helper class simPrimaryVertex with proper generation-level
604  * information */
605 std::vector<PrimaryVertexAnalyzer4PUSlimmed::simPrimaryVertex>
608  std::vector<PrimaryVertexAnalyzer4PUSlimmed::simPrimaryVertex> simpv;
609  int current_event = -1;
610 
611  if (verbose_) {
612  std::cout << "getSimPVs TrackingVertexCollection " << std::endl;
613  }
614 
615  for (TrackingVertexCollection::const_iterator v = tVC->begin();
616  v != tVC->end(); ++v) {
617  if (verbose_) {
618  std::cout << "BunchX.EventId: " << v->eventId().bunchCrossing() << "."
619  << (v->eventId()).event() << " Position: " << v->position()
620  << " G4/HepMC Vertices: " << v->g4Vertices().size() << "/"
621  << v->genVertices().size()
622  << " t = " << v->position().t() * 1.e12
623  << " == 0:" << (v->position().t() > 0) << std::endl;
624  for (TrackingVertex::g4v_iterator gv = v->g4Vertices_begin();
625  gv != v->g4Vertices_end(); gv++) {
626  std::cout << *gv << std::endl;
627  }
628  std::cout << "----------" << std::endl;
629  } // end of verbose_ session
630 
631  // I'd rather change this and select only vertices that come from
632  // BX=0. We should keep only the first vertex from all the events
633  // at BX=0.
634  if (v->eventId().bunchCrossing() != 0) continue;
635  if (v->eventId().event() != current_event) {
636  current_event = v->eventId().event();
637  } else {
638  continue;
639  }
640  // TODO(rovere) is this really necessary?
641  if (fabs(v->position().z()) > 1000) continue; // skip funny junk vertices
642 
643  // could be a new vertex, check all primaries found so far to avoid
644  // multiple entries
645  simPrimaryVertex sv(v->position().x(), v->position().y(),
646  v->position().z());
647  sv.eventId = v->eventId();
648  sv.sim_vertex = &(*v);
649 
650  for (TrackingParticleRefVector::iterator iTrack = v->daughterTracks_begin();
651  iTrack != v->daughterTracks_end(); ++iTrack) {
652  // TODO(rovere) isn't it always the case? Is it really worth
653  // checking this out?
654  // sv.eventId = (**iTrack).eventId();
655  assert((**iTrack).eventId().bunchCrossing() == 0);
656  }
657  // TODO(rovere) maybe get rid of this old logic completely ... ?
658  simPrimaryVertex* vp = NULL; // will become non-NULL if a vertex
659  // is found and then point to it
660  for (std::vector<simPrimaryVertex>::iterator v0 = simpv.begin();
661  v0 != simpv.end(); v0++) {
662  if ((sv.eventId == v0->eventId) && (fabs(sv.x - v0->x) < 1e-5) &&
663  (fabs(sv.y - v0->y) < 1e-5) && (fabs(sv.z - v0->z) < 1e-5)) {
664  vp = &(*v0);
665  break;
666  }
667  }
668  if (!vp) {
669  // this is a new vertex, add it to the list of sim-vertices
670  simpv.push_back(sv);
671  vp = &simpv.back();
672  if (verbose_) {
673  std::cout << "this is a new vertex " << sv.eventId.event() << " "
674  << sv.x << " " << sv.y << " " << sv.z << std::endl;
675  }
676  } else {
677  if (verbose_) {
678  std::cout << "this is not a new vertex" << sv.x << " " << sv.y << " "
679  << sv.z << std::endl;
680  }
681  }
682 
683  // Loop over daughter track(s) as Tracking Particles
684  for (TrackingVertex::tp_iterator iTP = v->daughterTracks_begin();
685  iTP != v->daughterTracks_end(); ++iTP) {
686  auto momentum = (*(*iTP)).momentum();
687  const reco::Track* matched_best_reco_track = nullptr;
688  double match_quality = -1;
689  if (use_only_charged_tracks_ && (**iTP).charge() == 0)
690  continue;
691  if (s2r_.find(*iTP) != s2r_.end()) {
692  matched_best_reco_track = s2r_[*iTP][0].first.get();
693  match_quality = s2r_[*iTP][0].second;
694  }
695  if (verbose_) {
696  std::cout << " Daughter momentum: " << momentum;
697  std::cout << " Daughter type " << (*(*iTP)).pdgId();
698  std::cout << " matched: " << (matched_best_reco_track != nullptr);
699  std::cout << " match-quality: " << match_quality;
700  std::cout << std::endl;
701  }
702  vp->ptot.setPx(vp->ptot.x() + momentum.x());
703  vp->ptot.setPy(vp->ptot.y() + momentum.y());
704  vp->ptot.setPz(vp->ptot.z() + momentum.z());
705  vp->ptot.setE(vp->ptot.e() + (**iTP).energy());
706  vp->ptsq += ((**iTP).pt() * (**iTP).pt());
707  // TODO(rovere) only select charged sim-particles? If so, maybe
708  // put it as a configuration parameter?
709  if (matched_best_reco_track) {
711  vp->average_match_quality += match_quality;
712  }
713  // TODO(rovere) get rid of cuts on sim-tracks
714  // TODO(rovere) be consistent between simulated tracks and
715  // reconstructed tracks selection
716  // count relevant particles
717  if (((**iTP).pt() > 0.2) && (fabs((**iTP).eta()) < 2.5) &&
718  (**iTP).charge() != 0) {
719  vp->nGenTrk++;
720  }
721  } // End of for loop on daughters sim-particles
722  if (vp->num_matched_reco_tracks)
723  vp->average_match_quality /=
724  static_cast<float>(vp->num_matched_reco_tracks);
725  if (verbose_) {
726  std::cout << "average number of associated tracks: "
727  << vp->num_matched_reco_tracks / static_cast<float>(vp->nGenTrk)
728  << " with average quality: " << vp->average_match_quality
729  << std::endl;
730  }
731  } // End of for loop on tracking vertices
732 
733  if (verbose_) {
734  std::cout << "------- PrimaryVertexAnalyzer4PUSlimmed simPVs from "
735  "TrackingVertices "
736  "-------" << std::endl;
737  for (std::vector<simPrimaryVertex>::iterator v0 = simpv.begin();
738  v0 != simpv.end(); v0++) {
739  std::cout << "z=" << v0->z << " event=" << v0->eventId.event()
740  << std::endl;
741  }
742  std::cout << "-----------------------------------------------" << std::endl;
743  } // End of for summary on discovered simulated primary vertices.
744 
745  // Now compute the closest distance in z between all simulated vertex
746  for (std::vector<simPrimaryVertex>::iterator vsim = simpv.begin();
747  vsim != simpv.end(); vsim++) {
748  std::vector<simPrimaryVertex>::iterator vsim2 = vsim;
749  vsim2++;
750  for (; vsim2 != simpv.end(); vsim2++) {
751  double distance_z = fabs(vsim->z - vsim2->z);
752  // Initialize with the next-sibling in the vector: minimization
753  // is performed by the next if.
754  if (vsim->closest_vertex_distance_z < 0) {
755  vsim->closest_vertex_distance_z = distance_z;
756  continue;
757  }
758  if (distance_z < vsim->closest_vertex_distance_z)
759  vsim->closest_vertex_distance_z = distance_z;
760  }
761  }
762  return simpv;
763 }
764 
765 /* Extract information form recoVertex and fill the helper class
766  * recoPrimaryVertex with proper reco-level information */
767 std::vector<PrimaryVertexAnalyzer4PUSlimmed::recoPrimaryVertex>
770  std::vector<PrimaryVertexAnalyzer4PUSlimmed::recoPrimaryVertex> recopv;
771 
772  if (verbose_) {
773  std::cout << "getRecoPVs TrackingVertexCollection " << std::endl;
774  }
775 
776  for (std::vector<reco::Vertex>::const_iterator v = tVC->begin();
777  v != tVC->end(); ++v) {
778  if (verbose_) {
779  std::cout << " Position: " << v->position() << std::endl;
780  }
781 
782  // Skip junk vertices
783  if (fabs(v->z()) > 1000) continue;
784  if (v->isFake() || !v->isValid()) continue;
785 
786  recoPrimaryVertex sv(v->position().x(), v->position().y(),
787  v->position().z());
788  sv.recVtx = &(*v);
789  // this is a new vertex, add it to the list of sim-vertices
790  recopv.push_back(sv);
792 
793  // Loop over daughter track(s)
794  for (auto iTrack = v->tracks_begin(); iTrack != v->tracks_end(); ++iTrack) {
795  auto momentum = (*(*iTrack)).innerMomentum();
796  // TODO(rovere) better handle the pixelVerticies, whose tracks
797  // do not have the innerMomentum defined. This is a temporary
798  // hack to overcome this problem.
799  if (momentum.mag2() == 0)
800  momentum = (*(*iTrack)).momentum();
801  if (verbose_) {
802  std::cout << " Daughter momentum: " << momentum;
803  std::cout << std::endl;
804  }
805  vp->ptsq += (momentum.perp2());
806  vp->nRecoTrk++;
807  } // End of for loop on daughters reconstructed tracks
808  } // End of for loop on tracking vertices
809 
810  if (verbose_) {
811  std::cout << "------- PrimaryVertexAnalyzer4PUSlimmed recoPVs from "
812  "VertexCollection "
813  "-------" << std::endl;
814  for (std::vector<recoPrimaryVertex>::iterator v0 = recopv.begin();
815  v0 != recopv.end(); v0++) {
816  std::cout << "z=" << v0->z << std::endl;
817  }
818  std::cout << "-----------------------------------------------" << std::endl;
819  } // End of for summary on reconstructed primary vertices.
820 
821  // Now compute the closest distance in z between all reconstructed vertex
822  for (std::vector<recoPrimaryVertex>::iterator vreco = recopv.begin();
823  vreco != recopv.end(); vreco++) {
824  std::vector<recoPrimaryVertex>::iterator vreco2 = vreco;
825  vreco2++;
826  for (; vreco2 != recopv.end(); vreco2++) {
827  double distance_z = fabs(vreco->z - vreco2->z);
828  // Initialize with the next-sibling in the vector: minimization
829  // is performed by the next if.
830  if (vreco->closest_vertex_distance_z < 0) {
831  vreco->closest_vertex_distance_z = distance_z;
832  continue;
833  }
834  if (distance_z < vreco->closest_vertex_distance_z)
835  vreco->closest_vertex_distance_z = distance_z;
836  }
837  }
838  return recopv;
839 }
840 
842  std::vector<simPrimaryVertex> & simpv) {
843  for (auto & v : simpv) {
844  v.num_matched_reco_tracks = 0;
845  v.average_match_quality = 0;
846  v.rec_vertices.clear();
847  }
848 }
849 
850 // ------------ method called to produce the data ------------
852  std::vector<simPrimaryVertex>& simpv,
853  const reco::VertexCollection& reco_vertices) {
854  if (verbose_) {
855  std::cout << "PrimaryVertexAnalyzer4PUSlimmed::matchSim2RecoVertices "
856  << reco_vertices.size() << std::endl;
857  }
858  for (std::vector<simPrimaryVertex>::iterator vsim = simpv.begin();
859  vsim != simpv.end(); vsim++) {
860  for (std::vector<reco::Vertex>::const_iterator vrec = reco_vertices.begin();
861  vrec != reco_vertices.end(); vrec++) {
862  if (vrec->isFake() || vrec->ndof() < 0.) {
863  continue; // skip fake vertices
864  }
865  if (verbose_) {
866  std::cout << "Considering reconstructed vertex at Z:" << vrec->z()
867  << std::endl;
868  }
869  if (((fabs(vrec->z() - vsim->z) / vrec->zError()) < sigma_z_match_)
870  && (fabs(vrec->z() - vsim->z) < abs_z_match_)) {
871  vsim->rec_vertices.push_back(&(*vrec));
872  if (verbose_) {
873  std::cout << "Trying a matching vertex for " << vsim->z << " at "
874  << vrec->z() << " with sign: "
875  << fabs(vrec->z() - vsim->z) / vrec->zError() << std::endl;
876  }
877  }
878  } // end of loop on reconstructed vertices
879  if (verbose_) {
880  if (vsim->rec_vertices.size()) {
881  for (auto const& v : vsim->rec_vertices) {
882  std::cout << "Found a matching vertex for genVtx "
883  << vsim->z << " at " << v->z()
884  << " with sign: " << fabs(v->z() - vsim->z) / v->zError()
885  << std::endl;
886  }
887  } else {
888  std::cout << "No matching vertex for " << vsim->z << std::endl;
889  }
890  }
891  } // end for loop on simulated vertices
892  if (verbose_) {
893  std::cout << "Done with matching sim vertices" << std::endl;
894  }
895 }
896 
898  std::vector<recoPrimaryVertex>& recopv,
899  const TrackingVertexCollection& gen_vertices,
900  const std::vector<simPrimaryVertex>& simpv) {
901  for (std::vector<recoPrimaryVertex>::iterator vrec = recopv.begin();
902  vrec != recopv.end(); vrec++) {
903  int current_event = -1;
904  for (std::vector<TrackingVertex>::const_iterator vsim =
905  gen_vertices.begin();
906  vsim != gen_vertices.end(); vsim++) {
907  // Keep only signal events
908  if (vsim->eventId().bunchCrossing() != 0) continue;
909 
910  // Keep only the primary vertex for each PU event
911  if (vsim->eventId().event() != current_event) {
912  current_event = vsim->eventId().event();
913  } else {
914  continue;
915  }
916 
917  // if the matching criteria are fulfilled, accept all the
918  // gen-vertices that are close in z, in unit of sigma_z of the
919  // reconstructed vertex, at least of sigma_z_match_. Require
920  // also a maximum absolute distance between the 2 vertices of at
921  // most abs_z_match_ along the Z axis(in cm).
922  if (((fabs(vrec->z - vsim->position().z()) / vrec->recVtx->zError()) <
924  && (fabs(vrec->z - vsim->position().z()) < abs_z_match_)) {
925  vrec->sim_vertices.push_back(&(*vsim));
926  for (std::vector<simPrimaryVertex>::const_iterator vv = simpv.begin();
927  vv != simpv.end(); vv++) {
928  if (vv->sim_vertex == &(*vsim)) {
929  vrec->sim_vertices_internal.push_back(&(*vv));
930  }
931  }
932 
933  if (verbose_) {
934  std::cout << "Matching Reco vertex for " << vrec->z
935  << " at Sim Vertex " << vsim->position().z()
936  << " with sign: " << fabs(vsim->position().z() - vrec->z) /
937  vrec->recVtx->zError()
938  << std::endl;
939  }
940  }
941  } // end of loop on simulated vertices
942  if (verbose_) {
943  for (auto v : vrec->sim_vertices) {
944  std::cout << "Found a matching vertex for reco: " << vrec->z
945  << " at gen:" << v->position().z() << " with sign: "
946  << fabs(vrec->z - v->position().z()) / vrec->recVtx->zError()
947  << std::endl;
948  }
949  }
950  } // end for loop on reconstructed vertices
951 }
952 
954  const edm::EventSetup& iSetup) {
955  using std::vector;
956  using std::cout;
957  using std::endl;
958  using edm::Handle;
959  using edm::View;
960  using edm::LogInfo;
961  using namespace reco;
962 
963  std::vector<float> pileUpInfo_z;
964 
965  // get the pileup information
967  if (iEvent.getByToken(vecPileupSummaryInfoToken_, puinfoH)) {
968  for (auto const& pu_info : *puinfoH.product()) {
969  mes_["root_folder"]["GenVtx_vs_BX"]
970  ->Fill(pu_info.getBunchCrossing(), pu_info.getPU_NumInteractions());
971  if (pu_info.getBunchCrossing() == 0) {
972  pileUpInfo_z = pu_info.getPU_zpositions();
973  if (verbose_) {
974  for (auto const& p : pileUpInfo_z) {
975  std::cout << "PileUpInfo on Z vertex: " << p << std::endl;
976  }
977  }
978  break;
979  }
980  }
981  }
982 
983  Handle<reco::TrackCollection> recTrks;
984  iEvent.getByToken(recoTrackCollectionToken_, recTrks);
985 
986  // for the associator
987  Handle<View<Track> > trackCollectionH;
988  iEvent.getByToken(edmView_recoTrack_Token_, trackCollectionH);
989 
992  bool gotTP =
993  iEvent.getByToken(trackingParticleCollectionToken_, TPCollectionH);
994  bool gotTV = iEvent.getByToken(trackingVertexCollectionToken_, TVCollectionH);
995 
996  // TODO(rovere) the idea is to put in case a track-selector in front
997  // of this module and then use its label to get the selected tracks
998  // out of the event instead of making an hard-coded selection in the
999  // code.
1000 
1001  if (gotTP) {
1002  // TODO(rovere) fetch an already existing collection from the
1003  // event instead of making another association on the fly???
1004  if (use_TP_associator_) {
1007  theHitsAssociator);
1008  associatorByHits_ = theHitsAssociator.product();
1010  trackCollectionH, TPCollectionH);
1012  trackCollectionH, TPCollectionH);
1013  }
1014  }
1015 
1016  std::vector<simPrimaryVertex> simpv; // a list of simulated primary
1017  // MC vertices
1018  // TODO(rovere) use move semantic?
1019  simpv = getSimPVs(TVCollectionH);
1020  // TODO(rovere) 1 vertex is not, by definition, pileup, and should
1021  // probably be subtracted?
1022  int kind_of_signal_vertex = 0;
1023  int num_pileup_vertices = simpv.size();
1024  mes_["root_folder"]["GenAllV_NumVertices"]->Fill(simpv.size());
1025  bool signal_is_highest_pt = std::max_element(simpv.begin(), simpv.end(),
1026  [](const simPrimaryVertex& lhs,
1027  const simPrimaryVertex& rhs) {
1028  return lhs.ptsq < rhs.ptsq;
1029  }) == simpv.begin();
1030  kind_of_signal_vertex = (kind_of_signal_vertex & ~(1<<HIGHEST_PT)) |
1031  (signal_is_highest_pt << HIGHEST_PT);
1032  mes_["root_folder"]["SignalIsHighestPt2"]->Fill(
1033  signal_is_highest_pt ? 1. : 0.);
1034  computePairDistance(simpv,
1035  mes_["root_folder"]["GenAllV_PairDistanceZ"]);
1036 
1037  int label_index = -1;
1038  for (auto const& vertex_token : reco_vertex_collection_tokens_) {
1039  std::vector<recoPrimaryVertex> recopv; // a list of reconstructed
1040  // primary MC vertices
1041  std::string label = reco_vertex_collections_[++label_index].label();
1042  Handle<reco::VertexCollection> recVtxs;
1043  if (!iEvent.getByToken(vertex_token, recVtxs)) {
1044  LogInfo("PrimaryVertexAnalyzer4PUSlimmed")
1045  << "Skipping vertex collection: " << label << " since it is missing."
1046  << std::endl;
1047  continue;
1048  }
1049  if (gotTV) {
1050  resetSimPVAssociation(simpv);
1051  matchSim2RecoVertices(simpv, *recVtxs.product());
1052  recopv = getRecoPVs(recVtxs);
1053  computePairDistance(recopv,
1054  mes_[label]["RecoAllAssoc2Gen_PairDistanceZ"]);
1055  matchReco2SimVertices(recopv, *TVCollectionH.product(), simpv);
1056  }
1057 
1058  int num_total_gen_vertices_assoc2reco = 0;
1059  int num_total_reco_vertices_assoc2gen = 0;
1060  int num_total_gen_vertices_multiassoc2reco = 0;
1061  int num_total_reco_vertices_multiassoc2gen = 0;
1062  int num_total_reco_vertices_duplicate = 0;
1063  for (auto const& v : simpv) {
1064  float mistag = 1.;
1065  // TODO(rovere) put selectors here in front of fill* methods.
1066  if (v.eventId.event() == 0) {
1067  if (std::find(v.rec_vertices.begin(), v.rec_vertices.end(),
1068  &((*recVtxs.product())[0])) != v.rec_vertices.end()) {
1069  mistag = 0.;
1070  kind_of_signal_vertex =
1071  (kind_of_signal_vertex & ~(1<<IS_ASSOC2FIRST_RECO)) |
1072  (signal_is_highest_pt << IS_ASSOC2FIRST_RECO);
1073  } else {
1074  if (v.rec_vertices.size()) {
1075  kind_of_signal_vertex =
1076  (kind_of_signal_vertex & ~(1<<IS_ASSOC2ANY_RECO)) |
1077  (signal_is_highest_pt << IS_ASSOC2ANY_RECO);
1078  }
1079  }
1080  mes_[label]["KindOfSignalPV"]->Fill(kind_of_signal_vertex);
1081  mes_[label]["MisTagRate"]->Fill(mistag);
1082  mes_[label]["MisTagRate_vs_PU"]->Fill(simpv.size(), mistag);
1083  mes_[label]["MisTagRate_vs_sum-pt2"]->Fill(v.ptsq, mistag);
1084  mes_[label]["MisTagRate_vs_Z"]->Fill(v.z, mistag);
1085  mes_[label]["MisTagRate_vs_R"]->Fill(v.r, mistag);
1086  mes_[label]["MisTagRate_vs_NumTracks"]->Fill(v.nGenTrk, mistag);
1087  if (signal_is_highest_pt) {
1088  mes_[label]["MisTagRateSignalIsHighest"]->Fill(mistag);
1089  mes_[label]["MisTagRateSignalIsHighest_vs_PU"]->Fill(simpv.size(),
1090  mistag);
1091  mes_[label]["MisTagRateSignalIsHighest_vs_sum-pt2"]->Fill(v.ptsq,
1092  mistag);
1093  mes_[label]["MisTagRateSignalIsHighest_vs_Z"]->Fill(v.z, mistag);
1094  mes_[label]["MisTagRateSignalIsHighest_vs_R"]->Fill(v.r, mistag);
1095  mes_[label]["MisTagRateSignalIsHighest_vs_NumTracks"]->Fill(v.nGenTrk,
1096  mistag);
1097  } else {
1098  mes_[label]["MisTagRateSignalIsNotHighest"]->Fill(mistag);
1099  mes_[label]["MisTagRateSignalIsNotHighest_vs_PU"]->Fill(simpv.size(),
1100  mistag);
1101  mes_[label]["MisTagRateSignalIsNotHighest_vs_sum-pt2"]->Fill(v.ptsq,
1102  mistag);
1103  mes_[label]["MisTagRateSignalIsNotHighest_vs_Z"]->Fill(v.z, mistag);
1104  mes_[label]["MisTagRateSignalIsNotHighest_vs_R"]->Fill(v.r, mistag);
1105  mes_[label]["MisTagRateSignalIsNotHighest_vs_NumTracks"]->
1106  Fill(v.nGenTrk, mistag);
1107  }
1108  // Now check at which location the Simulated PV has been
1109  // reconstructed in the primary vertex collection
1110  // at-hand. Mark it with fake index -1 if it was not
1111  // reconstructed at all.
1112 
1113  auto iv = (*recVtxs.product()).begin();
1114  for (int pv_position_in_reco_collection = 0;
1115  iv != (*recVtxs.product()).end();
1116  ++pv_position_in_reco_collection, ++iv) {
1117  if (std::find(v.rec_vertices.begin(), v.rec_vertices.end(),
1118  &(*iv)) != v.rec_vertices.end()) {
1119  mes_[label]["TruePVLocationIndex"]
1120  ->Fill(pv_position_in_reco_collection);
1121  mes_[label]["TruePVLocationIndexCumulative"]
1122  ->Fill(pv_position_in_reco_collection > 0 ? 1 : 0);
1123  if (signal_is_highest_pt) {
1124  mes_[label]["TruePVLocationIndexSignalIsHighest"]
1125  ->Fill(pv_position_in_reco_collection);
1126  } else {
1127  mes_[label]["TruePVLocationIndexSignalIsNotHighest"]
1128  ->Fill(pv_position_in_reco_collection);
1129  }
1130  break;
1131  }
1132  }
1133 
1134  // If we reached the end, it means that the Simulated PV has not
1135  // been associated to any reconstructed vertex: mark it as
1136  // missing in the reconstructed vertex collection using the fake
1137  // index -1.
1138  if (iv == (*recVtxs.product()).end()) {
1139  mes_[label]["TruePVLocationIndex"]->Fill(-1.);
1140  mes_[label]["TruePVLocationIndexCumulative"]->Fill(-1.);
1141  if (signal_is_highest_pt)
1142  mes_[label]["TruePVLocationIndexSignalIsHighest"]->Fill(-1.);
1143  else
1144  mes_[label]["TruePVLocationIndexSignalIsNotHighest"]->Fill(-1.);
1145  }
1146  }
1147 
1148  if (v.rec_vertices.size()) num_total_gen_vertices_assoc2reco++;
1149  if (v.rec_vertices.size() > 1) num_total_gen_vertices_multiassoc2reco++;
1150  // No need to N-tplicate the Gen-related cumulative histograms:
1151  // fill them only at the first iteration
1152  if (label_index == 0) fillGenericGenVertexHistograms(v);
1154  }
1155  mes_[label]["GenAllAssoc2Reco_NumVertices"]
1156  ->Fill(simpv.size(), simpv.size());
1157  mes_[label]["GenAllAssoc2RecoMatched_NumVertices"]
1158  ->Fill(simpv.size(), num_total_gen_vertices_assoc2reco);
1159  mes_[label]["GenAllAssoc2RecoMultiMatched_NumVertices"]
1160  ->Fill(simpv.size(), num_total_gen_vertices_multiassoc2reco);
1161  for (auto & v : recopv) {
1162  fillGenAssociatedRecoVertexHistograms(label, num_pileup_vertices, v);
1163  if (v.sim_vertices.size()) {
1164  num_total_reco_vertices_assoc2gen++;
1165  if (v.sim_vertices_internal[0]->rec_vertices.size() > 1) {
1166  num_total_reco_vertices_duplicate++;
1167  }
1168  }
1169  if (v.sim_vertices.size() > 1) num_total_reco_vertices_multiassoc2gen++;
1170  }
1171  mes_[label]["RecoAllAssoc2Gen_NumVertices"]
1172  ->Fill(recopv.size(), recopv.size());
1173  mes_[label]["RecoAllAssoc2GenMatched_NumVertices"]
1174  ->Fill(recopv.size(), num_total_reco_vertices_assoc2gen);
1175  mes_[label]["RecoAllAssoc2GenMultiMatched_NumVertices"]
1176  ->Fill(recopv.size(), num_total_reco_vertices_multiassoc2gen);
1177  mes_[label]["RecoAllAssoc2MultiMatchedGen_NumVertices"]
1178  ->Fill(recopv.size(), num_total_reco_vertices_duplicate);
1179  mes_[label]["RecoVtx_vs_GenVtx"]->Fill(simpv.size(), recopv.size());
1180  mes_[label]["MatchedRecoVtx_vs_GenVtx"]
1181  ->Fill(simpv.size(), num_total_reco_vertices_assoc2gen);
1182  }
1183 } // end of analyze
1184 
1185 template<class T>
1187  MonitorElement *me) {
1188  for (unsigned int i = 0; i < collection.size(); ++i) {
1189  for (unsigned int j = i+1; j < collection.size(); ++j) {
1190  me->Fill(
1191  std::abs(collection[i].z-collection[j].z));
1192  }
1193  }
1194 }
void matchSim2RecoVertices(std::vector< simPrimaryVertex > &, const reco::VertexCollection &)
T getParameter(std::string const &) const
void fillGenAssociatedRecoVertexHistograms(const std::string &, int, recoPrimaryVertex &v)
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< edm::View< reco::Track > > edmView_recoTrack_Token_
int event() const
get the contents of the subdetector field (should be protected?)
void resetSimPVAssociation(std::vector< simPrimaryVertex > &)
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > vecPileupSummaryInfoToken_
std::vector< TrackingParticle > TrackingParticleCollection
void computePairDistance(const T &collection, MonitorElement *me)
const_iterator end() const
last iterator over the map (read only)
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:157
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
const reco::TrackToTrackingParticleAssociator * associatorByHits_
std::vector< SimVertex >::const_iterator g4v_iterator
assert(m_qm.get())
const_iterator find(const key_type &k) const
find element with specified reference key
virtual void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
PrimaryVertexAnalyzer4PUSlimmed(const edm::ParameterSet &)
#define NULL
Definition: scimark2.h:8
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::TrackCollection > recoTrackCollectionToken_
edm::EDGetTokenT< TrackingVertexCollection > trackingVertexCollectionToken_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
float float float z
std::vector< PrimaryVertexAnalyzer4PUSlimmed::simPrimaryVertex > getSimPVs(const edm::Handle< TrackingVertexCollection >)
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:230
std::vector< PrimaryVertexAnalyzer4PUSlimmed::recoPrimaryVertex > getRecoPVs(const edm::Handle< reco::VertexCollection >)
std::vector< edm::EDGetTokenT< reco::VertexCollection > > reco_vertex_collection_tokens_
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
#define end
Definition: vmac.h:37
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
double ndof() const
Definition: Vertex.h:102
void fillGenericGenVertexHistograms(const simPrimaryVertex &v)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
void matchReco2SimVertices(std::vector< recoPrimaryVertex > &, const TrackingVertexCollection &, const std::vector< simPrimaryVertex > &)
T const * product() const
Definition: Handle.h:81
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
reco::RecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
std::vector< TrackingVertex > TrackingVertexCollection
void fillRecoAssociatedGenVertexHistograms(const std::string &, const simPrimaryVertex &v)
std::vector< edm::InputTag > reco_vertex_collections_
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > recoTrackToTrackingParticleAssociatorToken_
#define begin
Definition: vmac.h:30
std::map< std::string, std::map< std::string, MonitorElement * > > mes_
tuple cout
Definition: gather_cfg.py:121
volatile std::atomic< bool > shutdown_flag false
reco::SimToRecoCollection associateSimToReco(const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
long double T
Definition: Run.h:41
virtual void analyze(const edm::Event &, const edm::EventSetup &)
edm::EDGetTokenT< TrackingParticleCollection > trackingParticleCollectionToken_