CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MtdTracksValidation.cc
Go to the documentation of this file.
1 #include <string>
2 
8 
11 
17 
23 
27 
35 
38 #include "HepMC/GenRanges.h"
39 #include "CLHEP/Units/PhysicalConstants.h"
40 
42 public:
43  explicit MtdTracksValidation(const edm::ParameterSet&);
44  ~MtdTracksValidation() override;
45 
46  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
47 
48 private:
49  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
50 
51  void analyze(const edm::Event&, const edm::EventSetup&) override;
52 
53  const bool mvaGenSel(const HepMC::GenParticle&, const float&);
54  const bool mvaRecSel(const reco::TrackBase&, const reco::Vertex&, const double&, const double&);
55  const bool mvaGenRecMatch(const HepMC::GenParticle&, const double&, const reco::TrackBase&);
56 
57  // ------------ member data ------------
58 
60  const float trackMinPt_;
61  const float trackMinEta_;
62  const float trackMaxEta_;
63 
64  static constexpr double etacutGEN_ = 4.; // |eta| < 4;
65  static constexpr double etacutREC_ = 3.; // |eta| < 3;
66  static constexpr double pTcut_ = 0.7; // PT > 0.7 GeV
67  static constexpr double deltaZcut_ = 0.1; // dz separation 1 mm
68  static constexpr double deltaPTcut_ = 0.05; // dPT < 5%
69  static constexpr double deltaDRcut_ = 0.03; // DeltaR separation
70 
74 
76 
79 
89 
92 
101 
110 
121 
131 };
132 
133 // ------------ constructor and destructor --------------
135  : folder_(iConfig.getParameter<std::string>("folder")),
136  trackMinPt_(iConfig.getParameter<double>("trackMinimumPt")),
137  trackMinEta_(iConfig.getParameter<double>("trackMinimumEta")),
138  trackMaxEta_(iConfig.getParameter<double>("trackMaximumEta")) {
139  GenRecTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTagG"));
140  RecTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTagT"));
141  RecVertexToken_ = consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("inputTagV"));
142  HepMCProductToken_ = consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("inputTagH"));
143  trackAssocToken_ = consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("trackAssocSrc"));
144  pathLengthToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("pathLengthSrc"));
145  tmtdToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("tmtd"));
146  SigmatmtdToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmatmtd"));
147  t0SrcToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0Src"));
148  Sigmat0SrcToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0Src"));
149  t0PidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0PID"));
150  Sigmat0PidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0PID"));
151  t0SafePidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0SafePID"));
152  Sigmat0SafePidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0SafePID"));
153  trackMVAQualToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("trackMVAQual"));
154  mtdtopoToken_ = esConsumes<MTDTopology, MTDTopologyRcd>();
155  particleTableToken_ = esConsumes<HepPDT::ParticleDataTable, edm::DefaultRecord>();
156 }
157 
159 
160 // ------------ method called for each event ------------
162  using namespace edm;
163  using namespace geant_units::operators;
164  using namespace std;
165 
166  auto topologyHandle = iSetup.getTransientHandle(mtdtopoToken_);
167  const MTDTopology* topology = topologyHandle.product();
168 
169  bool topo1Dis = false;
170  bool topo2Dis = false;
171  if (topology->getMTDTopologyMode() <= static_cast<int>(MTDTopologyMode::Mode::barphiflat)) {
172  topo1Dis = true;
173  }
174  if (topology->getMTDTopologyMode() > static_cast<int>(MTDTopologyMode::Mode::barphiflat)) {
175  topo2Dis = true;
176  }
177 
178  auto GenRecTrackHandle = makeValid(iEvent.getHandle(GenRecTrackToken_));
179  auto RecVertexHandle = makeValid(iEvent.getHandle(RecVertexToken_));
180 
181  const auto& tMtd = iEvent.get(tmtdToken_);
182  const auto& SigmatMtd = iEvent.get(SigmatmtdToken_);
183  const auto& t0Src = iEvent.get(t0SrcToken_);
184  const auto& Sigmat0Src = iEvent.get(Sigmat0SrcToken_);
185  const auto& t0Pid = iEvent.get(t0PidToken_);
186  const auto& Sigmat0Pid = iEvent.get(Sigmat0PidToken_);
187  const auto& t0Safe = iEvent.get(t0SafePidToken_);
188  const auto& Sigmat0Safe = iEvent.get(Sigmat0SafePidToken_);
189  const auto& mtdQualMVA = iEvent.get(trackMVAQualToken_);
190  const auto& trackAssoc = iEvent.get(trackAssocToken_);
191  const auto& pathLength = iEvent.get(pathLengthToken_);
192 
193  unsigned int index = 0;
194  // --- Loop over all RECO tracks ---
195  for (const auto& trackGen : *GenRecTrackHandle) {
196  const reco::TrackRef trackref(iEvent.getHandle(GenRecTrackToken_), index);
197  index++;
198 
199  if (trackAssoc[trackref] == -1) {
200  LogInfo("mtdTracks") << "Extended track not associated";
201  continue;
202  }
203 
204  const reco::TrackRef mtdTrackref = reco::TrackRef(iEvent.getHandle(RecTrackToken_), trackAssoc[trackref]);
205  const reco::Track track = *mtdTrackref;
206 
207  if (track.pt() < trackMinPt_)
208  continue;
209 
210  meTracktmtd_->Fill(tMtd[trackref]);
211  if (std::round(SigmatMtd[trackref] - Sigmat0Pid[trackref]) != 0) {
212  LogWarning("mtdTracks") << "TimeError associated to refitted track is different from TimeError stored in tofPID "
213  "sigmat0 ValueMap: this should not happen";
214  }
215 
216  meTrackt0Src_->Fill(t0Src[trackref]);
217  meTrackSigmat0Src_->Fill(Sigmat0Src[trackref]);
218 
219  meTrackt0Pid_->Fill(t0Pid[trackref]);
220  meTrackSigmat0Pid_->Fill(Sigmat0Pid[trackref]);
221  meTrackt0SafePid_->Fill(t0Safe[trackref]);
222  meTrackSigmat0SafePid_->Fill(Sigmat0Safe[trackref]);
223  meTrackMVAQual_->Fill(mtdQualMVA[trackref]);
224 
225  meTrackPathLenghtvsEta_->Fill(std::abs(track.eta()), pathLength[trackref]);
226 
227  if (std::abs(track.eta()) < trackMinEta_) {
228  // --- all BTL tracks (with and without hit in MTD) ---
229  meBTLTrackEffEtaTot_->Fill(track.eta());
230  meBTLTrackEffPhiTot_->Fill(track.phi());
231  meBTLTrackEffPtTot_->Fill(track.pt());
232 
233  bool MTDBtl = false;
234  int numMTDBtlvalidhits = 0;
235  for (const auto hit : track.recHits()) {
236  if (hit->isValid() == false)
237  continue;
238  MTDDetId Hit = hit->geographicalId();
239  if ((Hit.det() == 6) && (Hit.subdetId() == 1) && (Hit.mtdSubDetector() == 1)) {
240  MTDBtl = true;
241  numMTDBtlvalidhits++;
242  }
243  }
244  meTrackNumHits_->Fill(numMTDBtlvalidhits);
245 
246  // --- keeping only tracks with last hit in MTD ---
247  if (MTDBtl == true) {
248  meBTLTrackEffEtaMtd_->Fill(track.eta());
249  meBTLTrackEffPhiMtd_->Fill(track.phi());
250  meBTLTrackEffPtMtd_->Fill(track.pt());
251  meBTLTrackRPTime_->Fill(track.t0());
252  meBTLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
253  }
254  } //loop over (geometrical) BTL tracks
255 
256  else {
257  // --- all ETL tracks (with and without hit in MTD) ---
258  if ((track.eta() < -trackMinEta_) && (track.eta() > -trackMaxEta_)) {
259  meETLTrackEffEtaTot_[0]->Fill(track.eta());
260  meETLTrackEffPhiTot_[0]->Fill(track.phi());
261  meETLTrackEffPtTot_[0]->Fill(track.pt());
262  }
263 
264  if ((track.eta() > trackMinEta_) && (track.eta() < trackMaxEta_)) {
265  meETLTrackEffEtaTot_[1]->Fill(track.eta());
266  meETLTrackEffPhiTot_[1]->Fill(track.phi());
267  meETLTrackEffPtTot_[1]->Fill(track.pt());
268  }
269 
270  bool MTDEtlZnegD1 = false;
271  bool MTDEtlZnegD2 = false;
272  bool MTDEtlZposD1 = false;
273  bool MTDEtlZposD2 = false;
274  int numMTDEtlvalidhits = 0;
275  for (const auto hit : track.recHits()) {
276  if (hit->isValid() == false)
277  continue;
278  MTDDetId Hit = hit->geographicalId();
279  if ((Hit.det() == 6) && (Hit.subdetId() == 1) && (Hit.mtdSubDetector() == 2)) {
280  ETLDetId ETLHit = hit->geographicalId();
281 
282  if (topo2Dis) {
283  if ((ETLHit.zside() == -1) && (ETLHit.nDisc() == 1)) {
284  MTDEtlZnegD1 = true;
285  meETLTrackRPTime_->Fill(track.t0());
286  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
287  numMTDEtlvalidhits++;
288  }
289  if ((ETLHit.zside() == -1) && (ETLHit.nDisc() == 2)) {
290  MTDEtlZnegD2 = true;
291  meETLTrackRPTime_->Fill(track.t0());
292  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
293  numMTDEtlvalidhits++;
294  }
295  if ((ETLHit.zside() == 1) && (ETLHit.nDisc() == 1)) {
296  MTDEtlZposD1 = true;
297  meETLTrackRPTime_->Fill(track.t0());
298  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
299  numMTDEtlvalidhits++;
300  }
301  if ((ETLHit.zside() == 1) && (ETLHit.nDisc() == 2)) {
302  MTDEtlZposD2 = true;
303  meETLTrackRPTime_->Fill(track.t0());
304  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
305  numMTDEtlvalidhits++;
306  }
307  }
308 
309  if (topo1Dis) {
310  if (ETLHit.zside() == -1) {
311  MTDEtlZnegD1 = true;
312  meETLTrackRPTime_->Fill(track.t0());
313  numMTDEtlvalidhits++;
314  }
315  if (ETLHit.zside() == 1) {
316  MTDEtlZposD1 = true;
317  meETLTrackRPTime_->Fill(track.t0());
318  numMTDEtlvalidhits++;
319  }
320  }
321  }
322  }
323  meTrackNumHits_->Fill(-numMTDEtlvalidhits);
324 
325  // --- keeping only tracks with last hit in MTD ---
326  if ((track.eta() < -trackMinEta_) && (track.eta() > -trackMaxEta_)) {
327  if ((MTDEtlZnegD1 == true) || (MTDEtlZnegD2 == true)) {
328  meETLTrackEffEtaMtd_[0]->Fill(track.eta());
329  meETLTrackEffPhiMtd_[0]->Fill(track.phi());
330  meETLTrackEffPtMtd_[0]->Fill(track.pt());
331  }
332  }
333  if ((track.eta() > trackMinEta_) && (track.eta() < trackMaxEta_)) {
334  if ((MTDEtlZposD1 == true) || (MTDEtlZposD2 == true)) {
335  meETLTrackEffEtaMtd_[1]->Fill(track.eta());
336  meETLTrackEffPhiMtd_[1]->Fill(track.phi());
337  meETLTrackEffPtMtd_[1]->Fill(track.pt());
338  }
339  }
340  }
341  } //RECO tracks loop
342 
343  // reco-gen matching used for MVA quality flag
344  const auto& primRecoVtx = *(RecVertexHandle.product()->begin());
345 
346  auto GenEventHandle = makeValid(iEvent.getHandle(HepMCProductToken_));
347  const HepMC::GenEvent* mc = GenEventHandle->GetEvent();
348  double zsim = convertMmToCm((*(mc->vertices_begin()))->position().z());
349  double tsim = (*(mc->vertices_begin()))->position().t() * CLHEP::mm / CLHEP::c_light;
350 
351  auto pdt = iSetup.getHandle(particleTableToken_);
352  const HepPDT::ParticleDataTable* pdTable = pdt.product();
353 
354  // select events with reco vertex close to true simulated primary vertex
355  if (std::abs(primRecoVtx.z() - zsim) < deltaZcut_) {
356  index = 0;
357  for (const auto& trackGen : *GenRecTrackHandle) {
358  const reco::TrackRef trackref(iEvent.getHandle(GenRecTrackToken_), index);
359  index++;
360 
361  // select the reconstructed track
362 
363  if (trackAssoc[trackref] == -1) {
364  continue;
365  }
366 
367  if (mvaRecSel(trackGen, primRecoVtx, t0Safe[trackref], Sigmat0Safe[trackref])) {
368  meMVATrackEffPtTot_->Fill(trackGen.pt());
369  meMVATrackEffEtaTot_->Fill(std::abs(trackGen.eta()));
370 
371  double dZ = trackGen.vz() - zsim;
372  double dT(-9999.);
373  double pullT(-9999.);
374  if (Sigmat0Safe[trackref] != -1.) {
375  dT = t0Safe[trackref] - tsim;
376  pullT = dT / Sigmat0Safe[trackref];
377  }
378  for (const auto& genP : mc->particle_range()) {
379  // select status 1 genParticles and match them to the reconstructed track
380 
381  float charge = pdTable->particle(HepPDT::ParticleID(genP->pdg_id())) != nullptr
382  ? pdTable->particle(HepPDT::ParticleID(genP->pdg_id()))->charge()
383  : 0.f;
384  if (mvaGenSel(*genP, charge)) {
385  if (mvaGenRecMatch(*genP, zsim, trackGen)) {
387  meMVATrackMatchedEffPtTot_->Fill(trackGen.pt());
388  meMVATrackMatchedEffEtaTot_->Fill(std::abs(trackGen.eta()));
389  if (pullT > -9999.) {
390  meMVATrackResTot_->Fill(dT);
391  meMVATrackPullTot_->Fill(pullT);
392  meMVATrackMatchedEffPtMtd_->Fill(trackGen.pt());
393  meMVATrackMatchedEffEtaMtd_->Fill(std::abs(trackGen.eta()));
394  }
395  break;
396  }
397  }
398  }
399  }
400  }
401  }
402 }
403 
404 // ------------ method for histogram booking ------------
406  ibook.setCurrentFolder(folder_);
407 
408  // histogram booking
409  meBTLTrackRPTime_ = ibook.book1D("TrackBTLRPTime", "Track t0 with respect to R.P.;t0 [ns]", 100, -1, 3);
410  meBTLTrackEffEtaTot_ = ibook.book1D("TrackBTLEffEtaTot", "Track efficiency vs eta (Tot);#eta_{RECO}", 100, -1.6, 1.6);
412  ibook.book1D("TrackBTLEffPhiTot", "Track efficiency vs phi (Tot);#phi_{RECO} [rad]", 100, -3.2, 3.2);
413  meBTLTrackEffPtTot_ = ibook.book1D("TrackBTLEffPtTot", "Track efficiency vs pt (Tot);pt_{RECO} [GeV]", 50, 0, 10);
414  meBTLTrackEffEtaMtd_ = ibook.book1D("TrackBTLEffEtaMtd", "Track efficiency vs eta (Mtd);#eta_{RECO}", 100, -1.6, 1.6);
416  ibook.book1D("TrackBTLEffPhiMtd", "Track efficiency vs phi (Mtd);#phi_{RECO} [rad]", 100, -3.2, 3.2);
417  meBTLTrackEffPtMtd_ = ibook.book1D("TrackBTLEffPtMtd", "Track efficiency vs pt (Mtd);pt_{RECO} [GeV]", 50, 0, 10);
419  ibook.book1D("TrackBTLPtRes", "Track pT resolution ;pT_{Gentrack}-pT_{MTDtrack}/pT_{Gentrack} ", 100, -0.1, 0.1);
420  meETLTrackRPTime_ = ibook.book1D("TrackETLRPTime", "Track t0 with respect to R.P.;t0 [ns]", 100, -1, 3);
422  ibook.book1D("TrackETLEffEtaTotZneg", "Track efficiency vs eta (Tot) (-Z);#eta_{RECO}", 100, -3.2, -1.4);
424  ibook.book1D("TrackETLEffEtaTotZpos", "Track efficiency vs eta (Tot) (+Z);#eta_{RECO}", 100, 1.4, 3.2);
426  ibook.book1D("TrackETLEffPhiTotZneg", "Track efficiency vs phi (Tot) (-Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
428  ibook.book1D("TrackETLEffPhiTotZpos", "Track efficiency vs phi (Tot) (+Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
430  ibook.book1D("TrackETLEffPtTotZneg", "Track efficiency vs pt (Tot) (-Z);pt_{RECO} [GeV]", 50, 0, 10);
432  ibook.book1D("TrackETLEffPtTotZpos", "Track efficiency vs pt (Tot) (+Z);pt_{RECO} [GeV]", 50, 0, 10);
434  ibook.book1D("TrackETLEffEtaMtdZneg", "Track efficiency vs eta (Mtd) (-Z);#eta_{RECO}", 100, -3.2, -1.4);
436  ibook.book1D("TrackETLEffEtaMtdZpos", "Track efficiency vs eta (Mtd) (+Z);#eta_{RECO}", 100, 1.4, 3.2);
438  ibook.book1D("TrackETLEffPhiMtdZneg", "Track efficiency vs phi (Mtd) (-Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
440  ibook.book1D("TrackETLEffPhiMtdZpos", "Track efficiency vs phi (Mtd) (+Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
442  ibook.book1D("TrackETLEffPtMtdZneg", "Track efficiency vs pt (Mtd) (-Z);pt_{RECO} [GeV]", 50, 0, 10);
444  ibook.book1D("TrackETLEffPtMtdZpos", "Track efficiency vs pt (Mtd) (+Z);pt_{RECO} [GeV]", 50, 0, 10);
446  ibook.book1D("TrackETLPtRes", "Track pT resolution;pT_{Gentrack}-pT_{MTDtrack}/pT_{Gentrack} ", 100, -0.1, 0.1);
447 
448  meTracktmtd_ = ibook.book1D("Tracktmtd", "Track time from TrackExtenderWithMTD;tmtd [ns]", 150, 1, 16);
449  meTrackt0Src_ = ibook.book1D("Trackt0Src", "Track time from TrackExtenderWithMTD;t0Src [ns]", 100, -1.5, 1.5);
451  ibook.book1D("TrackSigmat0Src", "Time Error from TrackExtenderWithMTD; #sigma_{t0Src} [ns]", 100, 0, 0.1);
452 
453  meTrackt0Pid_ = ibook.book1D("Trackt0Pid", "Track t0 as stored in TofPid;t0 [ns]", 100, -1, 1);
454  meTrackSigmat0Pid_ = ibook.book1D("TrackSigmat0Pid", "Sigmat0 as stored in TofPid; #sigma_{t0} [ns]", 100, 0, 0.1);
455  meTrackt0SafePid_ = ibook.book1D("Trackt0SafePID", "Track t0 Safe as stored in TofPid;t0 [ns]", 100, -1, 1);
457  ibook.book1D("TrackSigmat0SafePID", "Sigmat0 Safe as stored in TofPid; #sigma_{t0} [ns]", 100, 0, 0.1);
458  meTrackNumHits_ = ibook.book1D("TrackNumHits", "Number of valid MTD hits per track ; Number of hits", 10, -5, 5);
459  meTrackMVAQual_ = ibook.book1D("TrackMVAQual", "Track MVA Quality as stored in Value Map ; MVAQual", 100, 0, 1);
461  "TrackPathLenghtvsEta", "MTD Track pathlength vs MTD track Eta;|#eta|;Pathlength", 100, 0, 3.2, 100.0, 400.0, "S");
462  meMVATrackEffPtTot_ = ibook.book1D("MVAEffPtTot", "Pt of tracks associated to LV; track pt [GeV] ", 110, 0., 11.);
464  ibook.book1D("MVAMatchedEffPtTot", "Pt of tracks associated to LV matched to GEN; track pt [GeV] ", 110, 0., 11.);
466  "MVAMatchedEffPtMtd", "Pt of tracks associated to LV matched to GEN with time; track pt [GeV] ", 110, 0., 11.);
467  meMVATrackEffEtaTot_ = ibook.book1D("MVAEffEtaTot", "Pt of tracks associated to LV; track eta ", 66, 0., 3.3);
469  ibook.book1D("MVAMatchedEffEtaTot", "Pt of tracks associated to LV matched to GEN; track eta ", 66, 0., 3.3);
471  "MVAMatchedEffEtaMtd", "Pt of tracks associated to LV matched to GEN with time; track eta ", 66, 0., 3.3);
472  meMVATrackResTot_ = ibook.book1D(
473  "MVATrackRes", "t_{rec} - t_{sim} for LV associated tracks; t_{rec} - t_{sim} [ns] ", 120, -0.15, 0.15);
475  ibook.book1D("MVATrackPull", "Pull for associated tracks; (t_{rec}-t_{sim})/#sigma_{t}", 50, -5., 5.);
477  "MVATrackZposResTot", "Z_{PCA} - Z_{sim} for associated tracks;Z_{PCA} - Z_{sim} [cm] ", 100, -0.1, 0.1);
478 }
479 
480 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
481 
484 
485  desc.add<std::string>("folder", "MTD/Tracks");
486  desc.add<edm::InputTag>("inputTagG", edm::InputTag("generalTracks"));
487  desc.add<edm::InputTag>("inputTagT", edm::InputTag("trackExtenderWithMTD"));
488  desc.add<edm::InputTag>("inputTagV", edm::InputTag("offlinePrimaryVertices4D"));
489  desc.add<edm::InputTag>("inputTagH", edm::InputTag("generatorSmeared"));
490  desc.add<edm::InputTag>("tmtd", edm::InputTag("trackExtenderWithMTD:generalTracktmtd"));
491  desc.add<edm::InputTag>("sigmatmtd", edm::InputTag("trackExtenderWithMTD:generalTracksigmatmtd"));
492  desc.add<edm::InputTag>("t0Src", edm::InputTag("trackExtenderWithMTD:generalTrackt0"));
493  desc.add<edm::InputTag>("sigmat0Src", edm::InputTag("trackExtenderWithMTD:generalTracksigmat0"));
494  desc.add<edm::InputTag>("trackAssocSrc", edm::InputTag("trackExtenderWithMTD:generalTrackassoc"))
495  ->setComment("Association between General and MTD Extended tracks");
496  desc.add<edm::InputTag>("pathLengthSrc", edm::InputTag("trackExtenderWithMTD:generalTrackPathLength"));
497  desc.add<edm::InputTag>("t0SafePID", edm::InputTag("tofPID:t0safe"));
498  desc.add<edm::InputTag>("sigmat0SafePID", edm::InputTag("tofPID:sigmat0safe"));
499  desc.add<edm::InputTag>("sigmat0PID", edm::InputTag("tofPID:sigmat0"));
500  desc.add<edm::InputTag>("t0PID", edm::InputTag("tofPID:t0"));
501  desc.add<edm::InputTag>("trackMVAQual", edm::InputTag("mtdTrackQualityMVA:mtdQualMVA"));
502  desc.add<double>("trackMinimumPt", 1.0); // [GeV]
503  desc.add<double>("trackMinimumEta", 1.5);
504  desc.add<double>("trackMaximumEta", 3.2);
505 
506  descriptions.add("mtdTracks", desc);
507 }
508 
509 const bool MtdTracksValidation::mvaGenSel(const HepMC::GenParticle& gp, const float& charge) {
510  bool match = false;
511  if (gp.status() != 1) {
512  return match;
513  }
514  match = charge != 0.f && gp.momentum().perp() > pTcut_ && std::abs(gp.momentum().eta()) < etacutGEN_;
515  return match;
516 }
517 
519  const reco::Vertex& vtx,
520  const double& t0,
521  const double& st0) {
522  bool match = false;
523  match = trk.pt() > pTcut_ && std::abs(trk.eta()) < etacutREC_ && std::abs(trk.vz() - vtx.z()) <= deltaZcut_;
524  if (st0 > 0.) {
525  match = match && std::abs(t0 - vtx.t()) < 3. * st0;
526  }
527  return match;
528 }
529 
531  const double& zsim,
532  const reco::TrackBase& trk) {
533  bool match = false;
534  double dR = reco::deltaR(genP.momentum(), trk.momentum());
535  double genPT = genP.momentum().perp();
536  match =
537  std::abs(genPT - trk.pt()) < trk.pt() * deltaPTcut_ && dR < deltaDRcut_ && std::abs(trk.vz() - zsim) < deltaZcut_;
538  return match;
539 }
540 
MonitorElement * meETLTrackEffEtaTot_[2]
edm::EDGetTokenT< edm::ValueMap< int > > trackAssocToken_
static constexpr double etacutGEN_
MonitorElement * meBTLTrackEffPhiMtd_
edm::EDGetTokenT< edm::ValueMap< float > > trackMVAQualToken_
edm::EDGetTokenT< edm::ValueMap< float > > pathLengthToken_
MonitorElement * meMVATrackResTot_
MonitorElement * meTrackNumHits_
MonitorElement * meTrackt0Pid_
static constexpr double pTcut_
MonitorElement * meETLTrackRPTime_
edm::EDGetTokenT< edm::ValueMap< float > > t0SafePidToken_
std::string folder_
edm::EDGetTokenT< edm::ValueMap< float > > t0SrcToken_
int nDisc() const
Definition: ETLDetId.h:122
edm::ESGetToken< MTDTopology, MTDTopologyRcd > mtdtopoToken_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
const std::string folder_
edm::EDGetTokenT< edm::HepMCProduct > HepMCProductToken_
HepPDT::ParticleDataTable ParticleDataTable
MonitorElement * meMVATrackMatchedEffEtaMtd_
const bool mvaRecSel(const reco::TrackBase &, const reco::Vertex &, const double &, const double &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static constexpr double etacutREC_
static constexpr double deltaPTcut_
MonitorElement * meBTLTrackPtRes_
edm::EDGetTokenT< std::vector< reco::Vertex > > RecVertexToken_
MonitorElement * meETLTrackEffPtTot_[2]
MonitorElement * meBTLTrackEffEtaMtd_
edm::EDGetTokenT< edm::ValueMap< float > > tmtdToken_
MonitorElement * meMVATrackMatchedEffPtTot_
MonitorElement * meMVATrackEffPtTot_
MonitorElement * meMVATrackEffEtaTot_
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:664
const bool mvaGenSel(const HepMC::GenParticle &, const float &)
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
edm::EDGetTokenT< reco::TrackCollection > GenRecTrackToken_
MonitorElement * meETLTrackEffPhiMtd_[2]
void Fill(long long x)
edm::EDGetTokenT< reco::TrackCollection > RecTrackToken_
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:563
MonitorElement * meBTLTrackEffPhiTot_
int iEvent
Definition: GenABIO.cc:224
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
MonitorElement * meBTLTrackRPTime_
MonitorElement * meETLTrackPtRes_
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int, double lowY, double highY, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:322
MonitorElement * meTrackt0SafePid_
MonitorElement * meBTLTrackEffEtaTot_
double pt() const
track transverse momentum
Definition: TrackBase.h:637
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
static constexpr double deltaZcut_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double z() const
z coordinate
Definition: Vertex.h:133
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
MonitorElement * meMVATrackMatchedEffEtaTot_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< edm::ValueMap< float > > t0PidToken_
MonitorElement * meTrackMVAQual_
MonitorElement * meETLTrackEffPhiTot_[2]
edm::EDGetTokenT< edm::ValueMap< float > > SigmatmtdToken_
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Log< level::Info, false > LogInfo
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:661
MonitorElement * meMVATrackPullTot_
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)
MonitorElement * meTrackSigmat0SafePid_
MonitorElement * meETLTrackEffEtaMtd_[2]
MonitorElement * meTrackSigmat0Pid_
MonitorElement * meBTLTrackEffPtMtd_
MonitorElement * meTrackt0Src_
edm::EDGetTokenT< edm::ValueMap< float > > Sigmat0PidToken_
static constexpr double deltaDRcut_
MonitorElement * meTracktmtd_
MonitorElement * meBTLTrackEffPtTot_
MonitorElement * meTrackPathLenghtvsEta_
Detector identifier class for the Endcap Timing Layer.
Definition: ETLDetId.h:15
ESTransientHandle< T > getTransientHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:168
int zside() const
Definition: MTDDetId.h:61
int mtdSubDetector() const
Definition: MTDDetId.h:56
constexpr NumType convertMmToCm(NumType millimeters)
Definition: GeantUnits.h:63
const bool mvaGenRecMatch(const HepMC::GenParticle &, const double &, const reco::TrackBase &)
edm::EDGetTokenT< edm::ValueMap< float > > Sigmat0SrcToken_
MonitorElement * meETLTrackEffPtMtd_[2]
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
MonitorElement * meMVATrackZposResTot_
Log< level::Warning, false > LogWarning
auto makeValid(const U &iOtherHandleType) noexcept(false)
Definition: ValidHandle.h:52
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 * meMVATrackMatchedEffPtMtd_
MonitorElement * meTrackSigmat0Src_
Definition: Run.h:45
double t() const
t coordinate
Definition: Vertex.h:135
edm::EDGetTokenT< edm::ValueMap< float > > Sigmat0SafePidToken_
MtdTracksValidation(const edm::ParameterSet &)
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > particleTableToken_
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46