CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Attributes | Private Attributes
BDHadronTrackMonitoringAnalyzer Class Reference

#include <BDHadronTrackMonitoringAnalyzer.h>

Inheritance diagram for BDHadronTrackMonitoringAnalyzer:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Types

enum  HistoryClasses {
  BCWeakDecay = 0, BWeakDecay = 1, CWeakDecay = 2, PU = 3,
  Other = 4, Fake = 5
}
 

Public Member Functions

void analyze (const edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
 BDHadronTrackMonitoringAnalyzer (const edm::ParameterSet &pSet)
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
 ~BDHadronTrackMonitoringAnalyzer () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Static Public Attributes

static const std::vector< std::string > TrkHistCat
 

Private Attributes

TrackClassifier classifier_
 
edm::InputTag ClusterTPMapSrc_
 
edm::EDGetTokenT< ClusterTPAssociationclusterTPMapToken_
 
double decayLength_
 
double distJetAxis_
 
std::string ipTagInfos_
 
double maxJetEta_
 
double minJetPt_
 
MonitorElementnTrk_bjet [6]
 
MonitorElementnTrk_cjet [6]
 
MonitorElementnTrk_dusgjet [6]
 
MonitorElementnTrkAll_bjet
 
MonitorElementnTrkAll_cjet
 
MonitorElementnTrkAll_dusgjet
 
edm::EDGetTokenT< pat::JetCollectionPatJetCollectionTag_
 
edm::InputTag PatJetSrc_
 
edm::EDGetTokenT< reco::VertexCollectionPrimaryVertexColl_
 
edm::InputTag PVSrc_
 
edm::EDGetTokenT< reco::TrackCollectionTrackCollectionTag_
 
edm::InputTag TrackSrc_
 
MonitorElementTrkDxy_alljets [6]
 
MonitorElementTrkDz_alljets [6]
 
MonitorElementTrkEta_alljets [6]
 
MonitorElementTrkHitAll_alljets [6]
 
MonitorElementTrkHitPixel_alljets [6]
 
MonitorElementTrkHitStrip_alljets [6]
 
MonitorElementTrkPhi_alljets [6]
 
MonitorElementTrkPt_alljets [6]
 
MonitorElementTrkTruthDxy_alljets [5]
 
MonitorElementTrkTruthDz_alljets [5]
 
MonitorElementTrkTruthEta_alljets [5]
 
MonitorElementTrkTruthHitAll_alljets [5]
 
MonitorElementTrkTruthHitPixel_alljets [5]
 
MonitorElementTrkTruthHitStrip_alljets [5]
 
MonitorElementTrkTruthPhi_alljets [5]
 
MonitorElementTrkTruthPt_alljets [5]
 

Detailed Description

Top level steering routine for B + D hadron track monitoring tool from RECODEBUG samples.

Definition at line 51 of file BDHadronTrackMonitoringAnalyzer.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

BDHadronTrackMonitoringAnalyzer::BDHadronTrackMonitoringAnalyzer ( const edm::ParameterSet pSet)
explicit

Definition at line 20 of file BDHadronTrackMonitoringAnalyzer.cc.

References ClusterTPMapSrc_, clusterTPMapToken_, PatJetCollectionTag_, PatJetSrc_, PrimaryVertexColl_, PVSrc_, TrackCollectionTag_, and TrackSrc_.

21  : distJetAxis_(pSet.getParameter<double>("distJetAxisCut")),
22  decayLength_(pSet.getParameter<double>("decayLengthCut")),
23  minJetPt_(pSet.getParameter<double>("minJetPt")),
24  maxJetEta_(pSet.getParameter<double>("maxJetEta")),
25  ipTagInfos_(pSet.getParameter<std::string>("ipTagInfos")),
26  PatJetSrc_(pSet.getParameter<InputTag>("PatJetSource")),
27  TrackSrc_(pSet.getParameter<InputTag>("TrackSource")),
28  PVSrc_(pSet.getParameter<InputTag>("PrimaryVertexSource")),
29  ClusterTPMapSrc_(pSet.getParameter<InputTag>("clusterTPMap")),
30  classifier_(pSet, consumesCollector())
31 
32 {
33  PatJetCollectionTag_ = consumes<pat::JetCollection>(PatJetSrc_);
34  TrackCollectionTag_ = consumes<reco::TrackCollection>(TrackSrc_);
35  PrimaryVertexColl_ = consumes<reco::VertexCollection>(PVSrc_);
36  clusterTPMapToken_ = consumes<ClusterTPAssociation>(ClusterTPMapSrc_);
37  // TrkHistCat = {"BCWeakDecay", "BWeakDecay", "CWeakDecay", "PU", "Other",
38  // "Fake"};
39 }
T getParameter(std::string const &) const
edm::EDGetTokenT< ClusterTPAssociation > clusterTPMapToken_
edm::EDGetTokenT< reco::TrackCollection > TrackCollectionTag_
edm::EDGetTokenT< reco::VertexCollection > PrimaryVertexColl_
edm::EDGetTokenT< pat::JetCollection > PatJetCollectionTag_
BDHadronTrackMonitoringAnalyzer::~BDHadronTrackMonitoringAnalyzer ( )
override

Definition at line 193 of file BDHadronTrackMonitoringAnalyzer.cc.

193 {}

Member Function Documentation

void BDHadronTrackMonitoringAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 201 of file BDHadronTrackMonitoringAnalyzer.cc.

References funct::abs(), BCWeakDecay, TransientTrackBuilder::build(), TrackCategories::BWeakDecay, BWeakDecay, classifier_, IPTools::closestApproachToJet(), clusterTPMapToken_, compare(), RecoVertex::convertPos(), TrackCategories::CWeakDecay, CWeakDecay, listHistos::decayLength, decayLength_, DEFINE_FWK_MODULE, distJetAxis_, reco::TrackBase::dxy(), reco::TrackBase::dz(), MillePedeFileConverter_cfg::e, reco::TrackBase::eta(), TrackClassifier::evaluate(), TrackCategories::Fake, Fake, reco::TransientTrack::field(), MonitorElement::Fill(), TrackCategories::flags(), edm::EventSetup::get(), edm::Event::getByToken(), TrackHistory::getMatchedTrackingParticle(), TrajectoryStateOnSurface::globalPosition(), TrackClassifier::history(), reco::TrackBase::hitPattern(), mps_fire::i, reco::TransientTrack::impactPointState(), ipTagInfos_, OmniClusterRef::isPixel(), OmniClusterRef::isStrip(), TrajectoryStateOnSurface::isValid(), OmniClusterRef::isValid(), metsig::jet, IPTools::jetTrackDistance(), ClusterTPAssociation::map(), maxJetEta_, minJetPt_, TrackClassifier::newEvent(), nTrk_bjet, nTrk_cjet, nTrk_dusgjet, nTrkAll_bjet, nTrkAll_cjet, nTrkAll_dusgjet, reco::TrackBase::numberOfValidHits(), reco::HitPattern::numberOfValidPixelHits(), reco::HitPattern::numberOfValidStripHits(), Other, AlCaHLTBitMon_ParallelJobs::p, PatJetCollectionTag_, reco::TrackBase::phi(), reco::Vertex::position(), impactParameterTagInfos_cfi::primaryVertex, PrimaryVertexColl_, reco::TrackBase::pt(), PU, MetAnalyzer::pv(), trackingPlots::selectedTracks, TrackCategories::SignalEvent, mathSSE::sqrt(), reco::btag::toTrack(), toTrackRef(), TrackCollectionTag_, TrkDxy_alljets, TrkDz_alljets, TrkEta_alljets, TrkHistCat, TrkHitAll_alljets, TrkHitPixel_alljets, TrkHitStrip_alljets, TrkPhi_alljets, TrkPt_alljets, TrkTruthDxy_alljets, TrkTruthDz_alljets, TrkTruthEta_alljets, TrkTruthHitAll_alljets, TrkTruthHitPixel_alljets, TrkTruthHitStrip_alljets, TrkTruthPhi_alljets, and TrkTruthPt_alljets.

201  {
203  iEvent.getByToken(PatJetCollectionTag_, patJetsColl);
204 
206  iEvent.getByToken(TrackCollectionTag_, tracksHandle);
207 
208  edm::Handle<ClusterTPAssociation> pCluster2TPListH;
209  iEvent.getByToken(clusterTPMapToken_, pCluster2TPListH);
210  const ClusterTPAssociation &clusterToTPMap = *pCluster2TPListH;
211 
213  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", trackBuilder);
214 
215  classifier_.newEvent(iEvent, iSetup);
216 
217  // -----Primary Vertex-----
218  const reco::Vertex *pv;
219 
221  iEvent.getByToken(PrimaryVertexColl_, primaryVertex);
222 
223  bool pvFound = (!primaryVertex->empty());
224  if (pvFound) {
225  pv = &(*primaryVertex->begin());
226  } else {
228  e(0, 0) = 0.0015 * 0.0015;
229  e(1, 1) = 0.0015 * 0.0015;
230  e(2, 2) = 15. * 15.;
231  reco::Vertex::Point p(0, 0, 0);
232  pv = new reco::Vertex(p, e, 1, 1, 1);
233  }
234  // -----------------------
235 
236  // -------- Loop Over Jets ----------
237  for (pat::JetCollection::const_iterator jet = patJetsColl->begin(); jet != patJetsColl->end(); ++jet) {
238  if (jet->pt() < minJetPt_ || std::fabs(jet->eta()) > maxJetEta_)
239  continue;
240 
241  unsigned int flav = abs(jet->hadronFlavour());
242 
243  // std::cout << "patJet collection has pfImpactParameterTagInfo?: " <<
244  // jet->hasTagInfo("pfImpactParameter") << std::endl;
245  const CandIPTagInfo *trackIpTagInfo = jet->tagInfoCandIP(ipTagInfos_);
246  const std::vector<edm::Ptr<reco::Candidate>> &selectedTracks(trackIpTagInfo->selectedTracks());
247 
248  unsigned int nseltracks = 0;
249  std::vector<int> nseltracksCat(TrkHistCat.size(),
250  0); // following the order of TrkHistCat
251 
252  unsigned int nTrackSize = selectedTracks.size(); // number of tracks from IPInfos to loop over
253  // -------- Loop Over (selected) Tracks ----------
254  for (unsigned int itt = 0; itt < nTrackSize; ++itt) {
255  const TrackBaseRef ptrackRef = toTrackRef(selectedTracks[itt]);
256  const reco::Track *ptrackPtr = reco::btag::toTrack(ptrackRef);
257  const reco::Track &ptrack = *ptrackPtr;
258 
259  reco::TransientTrack transientTrack = trackBuilder->build(ptrackPtr);
260  GlobalVector direction(jet->px(), jet->py(), jet->pz());
261 
262  Double_t distJetAxis = IPTools::jetTrackDistance(transientTrack, direction, *pv).second.value();
263 
264  Double_t decayLength = 999;
265  TrajectoryStateOnSurface closest =
266  IPTools::closestApproachToJet(transientTrack.impactPointState(), *pv, direction, transientTrack.field());
267  if (closest.isValid())
268  decayLength = (closest.globalPosition() - RecoVertex::convertPos(pv->position())).mag();
269  else
270  decayLength = 999;
271 
272  // extra cut ons the tracks
273  if (std::fabs(distJetAxis) > distJetAxis_ || decayLength > decayLength_) {
274  continue;
275  }
276  nseltracks += 1; // if it passed these cuts, nselectedtracks +1
277 
279 
280  double TrkPt = ptrack.pt();
281  double TrkEta = ptrack.eta();
282  double TrkPhi = ptrack.phi();
283  double TrkDxy = ptrack.dxy(pv->position());
284  double TrkDz = ptrack.dz(pv->position());
285  int TrknHitAll = ptrack.numberOfValidHits();
286  int TrknHitPixel = ptrack.hitPattern().numberOfValidPixelHits();
287  int TrknHitStrip = ptrack.hitPattern().numberOfValidStripHits();
288 
289  double TrkTruthPt = -99;
290  double TrkTruthEta = -99;
291  double TrkTruthPhi = -99;
292  double TrkTruthDxy = -1;
293  double TrkTruthDz = -1;
294  int TrkTruthnHitAll = -1;
295  int TrkTruthnHitPixel = -1;
296  int TrkTruthnHitStrip = -1;
297 
298  // Get corresponding Trackingparticle
299  std::pair<TrackingParticleRef, double> res = classifier_.history().getMatchedTrackingParticle();
300  TrackingParticleRef tpr = res.first;
301  double quality_tpr = res.second;
302 
303  // Match TP to hit-cluster (re-ordering according to TP rather than
304  // clusters and look for equal_range of a given tpr)
305  auto clusterTPmap = clusterToTPMap.map();
306  std::sort(clusterTPmap.begin(), clusterTPmap.end(), compare);
307  auto clusterRange =
308  std::equal_range(clusterTPmap.begin(), clusterTPmap.end(), std::make_pair(OmniClusterRef(), tpr), compare);
309  if (quality_tpr != 0) {
310  TrkTruthPt = tpr->pt();
311  TrkTruthEta = tpr->eta();
312  TrkTruthPhi = tpr->phi();
313 
314  const TrackingParticle::Point &vertex_pv = pv->position();
315  TrackingParticle::Point vertex_tpr = tpr->vertex();
316  TrackingParticle::Vector momentum_tpr = tpr->momentum();
317  TrkTruthDxy = (-(vertex_tpr.x() - vertex_pv.x()) * momentum_tpr.y() +
318  (vertex_tpr.y() - vertex_pv.y()) * momentum_tpr.x()) /
319  tpr->pt();
320  TrkTruthDz = (vertex_tpr.z() - vertex_pv.z()) - ((vertex_tpr.x() - vertex_pv.x()) * momentum_tpr.x() +
321  (vertex_tpr.y() - vertex_pv.y()) * momentum_tpr.y()) /
322  sqrt(momentum_tpr.perp2()) * momentum_tpr.z() /
323  sqrt(momentum_tpr.perp2());
324 
325  TrkTruthnHitAll = 0;
326  TrkTruthnHitPixel = 0;
327  TrkTruthnHitStrip = 0;
328  if (clusterRange.first != clusterRange.second) {
329  for (auto ip = clusterRange.first; ip != clusterRange.second; ++ip) {
330  const OmniClusterRef &cluster = ip->first;
331  if (cluster.isPixel() && cluster.isValid()) {
332  TrkTruthnHitPixel += 1;
333  }
334  if (cluster.isStrip() && cluster.isValid()) {
335  TrkTruthnHitStrip += 1;
336  }
337  }
338  }
339  TrkTruthnHitAll = TrkTruthnHitPixel + TrkTruthnHitStrip;
340  }
341 
342  // ----------- Filling the correct histograms based on jet flavour and
343  // Track history Category --------
344 
345  // BCWeakDecay
347  theFlag[TrackCategories::CWeakDecay]) {
357  if (quality_tpr != 0) {
366  }
367  }
368  // BWeakDecay
369  else if (theFlag[TrackCategories::SignalEvent] && theFlag[TrackCategories::BWeakDecay] &&
370  !theFlag[TrackCategories::CWeakDecay]) {
380  if (quality_tpr != 0) {
389  }
390  }
391  // CWeakDecay
392  else if (theFlag[TrackCategories::SignalEvent] && !theFlag[TrackCategories::BWeakDecay] &&
393  theFlag[TrackCategories::CWeakDecay]) {
403  if (quality_tpr != 0) {
412  }
413  }
414  // PU
415  else if (!theFlag[TrackCategories::SignalEvent] && !theFlag[TrackCategories::Fake]) {
416  nseltracksCat[BDHadronTrackMonitoringAnalyzer::PU] += 1;
425  if (quality_tpr != 0) {
434  }
435  }
436  // Other
437  else if (theFlag[TrackCategories::SignalEvent] && !theFlag[TrackCategories::BWeakDecay] &&
438  !theFlag[TrackCategories::CWeakDecay]) {
439  nseltracksCat[BDHadronTrackMonitoringAnalyzer::Other] += 1;
448  if (quality_tpr != 0) {
457  }
458  }
459  // Fake
460  else if (!theFlag[TrackCategories::SignalEvent] && theFlag[TrackCategories::Fake]) {
461  nseltracksCat[BDHadronTrackMonitoringAnalyzer::Fake] += 1;
470  // NO TRUTH FOR FAKES!!!
471  }
472  }
473  // -------- END Loop Over (selected) Tracks ----------
474  // Still have to fill some jet-flavour specific variables
475  if (flav == 5) {
476  nTrkAll_bjet->Fill(nseltracks);
477  for (unsigned int i = 0; i < TrkHistCat.size(); i++) {
478  nTrk_bjet[i]->Fill(nseltracksCat[i]);
479  }
480  } else if (flav == 4) {
481  nTrkAll_cjet->Fill(nseltracks);
482  for (unsigned int i = 0; i < TrkHistCat.size(); i++) {
483  nTrk_cjet[i]->Fill(nseltracksCat[i]);
484  }
485  } else {
486  nTrkAll_dusgjet->Fill(nseltracks);
487  for (unsigned int i = 0; i < TrkHistCat.size(); i++) {
488  nTrk_dusgjet[i]->Fill(nseltracksCat[i]);
489  }
490  }
491  }
492  // -------- END Loop Over Jets ----------
493 
494  if (!pvFound) {
495  delete pv;
496  }
497 }
reco::Vertex::Point convertPos(const GlobalPoint &p)
bool isStrip() const
bool compare(const P &i, const P &j)
edm::EDGetTokenT< ClusterTPAssociation > clusterTPMapToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< reco::TrackCollection > TrackCollectionTag_
edm::EDGetTokenT< reco::VertexCollection > PrimaryVertexColl_
reco::TransientTrack build(const reco::Track *p) const
int numberOfValidStripHits() const
Definition: HitPattern.h:931
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
GlobalPoint globalPosition() const
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:678
TrajectoryStateOnSurface closestApproachToJet(const TrajectoryStateOnSurface &state, const reco::Vertex &vertex, const GlobalVector &aJetDirection, const MagneticField *field)
Definition: IPTools.cc:177
const MagneticField * field() const
TrackHistory const & history() const
Returns a reference to the track history used in the classification.
const reco::Track * toTrack(const reco::TrackBaseRef &t)
Definition: IPTagInfo.h:24
const Point & position() const
position
Definition: Vertex.h:109
Definition: Electron.h:6
std::pair< double, Measurement1D > jetTrackDistance(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:200
void Fill(long long x)
edm::EDGetTokenT< pat::JetCollection > PatJetCollectionTag_
math::XYZPointD Point
point in the space
bool isValid() const
TrackClassifier const & evaluate(reco::TrackBaseRef const &)
Classify the RecoTrack in categories.
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:684
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:654
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:889
static const std::vector< std::string > TrkHistCat
void newEvent(edm::Event const &, edm::EventSetup const &)
Pre-process event information (for accessing reconstraction information)
bool isPixel() const
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:642
const Flags & flags() const
Returns flags with the category descriptions.
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:479
const reco::TrackBaseRef toTrackRef(const edm::Ptr< reco::Candidate > &cnd)
T get() const
Definition: EventSetup.h:71
const map_type & map() const
int numberOfValidPixelHits() const
Definition: HitPattern.h:916
std::vector< bool > Flags
Main types associated to the class.
TrajectoryStateOnSurface impactPointState() const
const std::pair< TrackingParticleRef, double > getMatchedTrackingParticle() const
Definition: TrackHistory.h:58
math::XYZVectorD Vector
point in the space
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:624
void BDHadronTrackMonitoringAnalyzer::bookHistograms ( DQMStore::IBooker ibook,
edm::Run const &  run,
edm::EventSetup const &  es 
)
override

Definition at line 43 of file BDHadronTrackMonitoringAnalyzer.cc.

References DQMStore::IBooker::book1D(), mps_fire::i, nTrk_bjet, nTrk_cjet, nTrk_dusgjet, nTrkAll_bjet, nTrkAll_cjet, nTrkAll_dusgjet, DQMStore::IBooker::setCurrentFolder(), RecoBTag::setTDRStyle(), TrkDxy_alljets, TrkDz_alljets, TrkEta_alljets, TrkHistCat, TrkHitAll_alljets, TrkHitPixel_alljets, TrkHitStrip_alljets, TrkPhi_alljets, TrkPt_alljets, TrkTruthDxy_alljets, TrkTruthDz_alljets, TrkTruthEta_alljets, TrkTruthHitAll_alljets, TrkTruthHitPixel_alljets, TrkTruthHitStrip_alljets, TrkTruthPhi_alljets, and TrkTruthPt_alljets.

45  {
46  ibook.setCurrentFolder("BDHadronTracks/JetContent");
47  //
48  // Book all histograms.
49  //
51 
52  nTrkAll_bjet = ibook.book1D(
53  "nTrkAll_bjet", "Number of selected tracks in b jets;number of selected tracks;jets", 16, -0.5, 15.5);
54 
55  nTrkAll_cjet = ibook.book1D(
56  "nTrkAll_cjet", "Number of selected tracks in c jets;number of selected tracks;jets", 16, -0.5, 15.5);
57 
58  nTrkAll_dusgjet = ibook.book1D(
59  "nTrkAll_dusgjet", "Number of selected tracks in dusg jets;number of selected tracks;jets", 16, -0.5, 15.5);
60 
61  // Loop over different Track History Categories
62  for (unsigned int i = 0; i < TrkHistCat.size(); i++) {
63  ibook.setCurrentFolder("BDHadronTracks/JetContent");
64  // b jets
65  nTrk_bjet[i] = ibook.book1D("nTrk_bjet_" + TrkHistCat[i],
66  "Number of selected tracks in b jets (" + TrkHistCat[i] +
67  ");number of selected tracks (" + TrkHistCat[i] + ");jets",
68  16,
69  -0.5,
70  15.5);
71 
72  // c jets
73  nTrk_cjet[i] = ibook.book1D("nTrk_cjet_" + TrkHistCat[i],
74  "Number of selected tracks in c jets (" + TrkHistCat[i] +
75  ");number of selected tracks (" + TrkHistCat[i] + ");jets",
76  16,
77  -0.5,
78  15.5);
79 
80  // dusg jets
81  nTrk_dusgjet[i] = ibook.book1D("nTrk_dusgjet_" + TrkHistCat[i],
82  "Number of selected tracks in dusg jets (" + TrkHistCat[i] +
83  ");number of selected tracks (" + TrkHistCat[i] + ");jets",
84  16,
85  -0.5,
86  15.5);
87 
88  ibook.setCurrentFolder("BDHadronTracks/TrackInfo");
89  // track properties for all flavours combined
90  TrkPt_alljets[i] = ibook.book1D("TrkPt_" + TrkHistCat[i],
91  "Track pT (" + TrkHistCat[i] + ");track p_{T} (" + TrkHistCat[i] + ");tracks",
92  30,
93  0,
94  100);
95  TrkEta_alljets[i] = ibook.book1D("TrkEta_" + TrkHistCat[i],
96  "Track #eta (" + TrkHistCat[i] + ");track #eta (" + TrkHistCat[i] + ");tracks",
97  30,
98  -2.5,
99  2.5);
100  TrkPhi_alljets[i] = ibook.book1D("TrkPhi_" + TrkHistCat[i],
101  "Track #phi (" + TrkHistCat[i] + ");track #phi (" + TrkHistCat[i] + ");tracks",
102  30,
103  -3.15,
104  3.15);
105  TrkDxy_alljets[i] = ibook.book1D("TrkDxy_" + TrkHistCat[i],
106  "Track dxy (" + TrkHistCat[i] + ");track dxy (" + TrkHistCat[i] + ");tracks",
107  30,
108  -0.1,
109  0.1);
110  TrkDz_alljets[i] = ibook.book1D("TrkDz_" + TrkHistCat[i],
111  "Track dz (" + TrkHistCat[i] + ");track dz (" + TrkHistCat[i] + ");tracks",
112  30,
113  -0.1,
114  0.1);
115  TrkHitAll_alljets[i] = ibook.book1D(
116  "TrkHitAll_" + TrkHistCat[i],
117  "Number of tracker hits (" + TrkHistCat[i] + ");track number of all hits (" + TrkHistCat[i] + ");tracks",
118  31,
119  -0.5,
120  30.5);
121  TrkHitStrip_alljets[i] = ibook.book1D(
122  "TrkHitStrip_" + TrkHistCat[i],
123  "Number of strip hits (" + TrkHistCat[i] + ");track number of strip hits (" + TrkHistCat[i] + ");tracks",
124  31,
125  -0.5,
126  30.5);
127  TrkHitPixel_alljets[i] = ibook.book1D(
128  "TrkHitPixel_" + TrkHistCat[i],
129  "Number of pixel hits (" + TrkHistCat[i] + ");track number of pixel hits (" + TrkHistCat[i] + ");tracks",
130  9,
131  -0.5,
132  8.5);
133 
134  ibook.setCurrentFolder("BDHadronTracks/TrackTruthInfo");
135  if (i < 5) { // Fakes (i == 5) have no truth by definition!
137  ibook.book1D("TrkTruthPt_" + TrkHistCat[i],
138  "Track pT (" + TrkHistCat[i] + " Truth);track p_{T} (" + TrkHistCat[i] + " Truth);tracks",
139  30,
140  0,
141  100);
143  ibook.book1D("TrkTruthEta_" + TrkHistCat[i],
144  "Track #eta (" + TrkHistCat[i] + " Truth);track #eta (" + TrkHistCat[i] + " Truth);tracks",
145  30,
146  -2.5,
147  2.5);
149  ibook.book1D("TrkTruthPhi_" + TrkHistCat[i],
150  "Track #phi (" + TrkHistCat[i] + " Truth);track #phi (" + TrkHistCat[i] + " Truth);tracks",
151  30,
152  -3.15,
153  3.15);
155  ibook.book1D("TrkTruthDxy_" + TrkHistCat[i],
156  "Track dxy (" + TrkHistCat[i] + " Truth);track dxy (" + TrkHistCat[i] + " Truth);tracks",
157  30,
158  -0.1,
159  0.1);
161  ibook.book1D("TrkTruthDz_" + TrkHistCat[i],
162  "Track dz (" + TrkHistCat[i] + " Truth);track dz (" + TrkHistCat[i] + " Truth);tracks",
163  30,
164  -0.1,
165  0.1);
167  ibook.book1D("TrkTruthHitAll_" + TrkHistCat[i],
168  "Number of tracker hits (" + TrkHistCat[i] + " Truth);track number of all hits (" +
169  TrkHistCat[i] + " Truth);tracks",
170  31,
171  -0.5,
172  30.5);
174  ibook.book1D("TrkTruthHitStrip_" + TrkHistCat[i],
175  "Number of strip hits (" + TrkHistCat[i] + " Truth);track number of strip hits (" +
176  TrkHistCat[i] + " Truth);tracks",
177  31,
178  -0.5,
179  30.5);
181  ibook.book1D("TrkTruthHitPixel_" + TrkHistCat[i],
182  "Number of pixel hits (" + TrkHistCat[i] + " Truth);track number of pixel hits (" +
183  TrkHistCat[i] + " Truth);tracks",
184  9,
185  -0.5,
186  8.5);
187  }
188  }
189 }
TStyle * setTDRStyle()
Definition: Tools.cc:346
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
static const std::vector< std::string > TrkHistCat

Member Data Documentation

TrackClassifier BDHadronTrackMonitoringAnalyzer::classifier_
private

Definition at line 86 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze().

edm::InputTag BDHadronTrackMonitoringAnalyzer::ClusterTPMapSrc_
private

Definition at line 77 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by BDHadronTrackMonitoringAnalyzer().

edm::EDGetTokenT<ClusterTPAssociation> BDHadronTrackMonitoringAnalyzer::clusterTPMapToken_
private

Definition at line 83 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and BDHadronTrackMonitoringAnalyzer().

double BDHadronTrackMonitoringAnalyzer::decayLength_
private

Definition at line 66 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze().

double BDHadronTrackMonitoringAnalyzer::distJetAxis_
private

Definition at line 65 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze().

std::string BDHadronTrackMonitoringAnalyzer::ipTagInfos_
private

Definition at line 71 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze().

double BDHadronTrackMonitoringAnalyzer::maxJetEta_
private

Definition at line 68 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze().

double BDHadronTrackMonitoringAnalyzer::minJetPt_
private

Definition at line 67 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze().

MonitorElement* BDHadronTrackMonitoringAnalyzer::nTrk_bjet[6]
private

Definition at line 91 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BDHadronTrackMonitoringAnalyzer::nTrk_cjet[6]
private

Definition at line 95 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BDHadronTrackMonitoringAnalyzer::nTrk_dusgjet[6]
private

Definition at line 99 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BDHadronTrackMonitoringAnalyzer::nTrkAll_bjet
private

Definition at line 90 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BDHadronTrackMonitoringAnalyzer::nTrkAll_cjet
private

Definition at line 94 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BDHadronTrackMonitoringAnalyzer::nTrkAll_dusgjet
private

Definition at line 98 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

edm::EDGetTokenT<pat::JetCollection> BDHadronTrackMonitoringAnalyzer::PatJetCollectionTag_
private

Definition at line 80 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and BDHadronTrackMonitoringAnalyzer().

edm::InputTag BDHadronTrackMonitoringAnalyzer::PatJetSrc_
private

Definition at line 74 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by BDHadronTrackMonitoringAnalyzer().

edm::EDGetTokenT<reco::VertexCollection> BDHadronTrackMonitoringAnalyzer::PrimaryVertexColl_
private

Definition at line 82 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and BDHadronTrackMonitoringAnalyzer().

edm::InputTag BDHadronTrackMonitoringAnalyzer::PVSrc_
private

Definition at line 76 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by BDHadronTrackMonitoringAnalyzer().

edm::EDGetTokenT<reco::TrackCollection> BDHadronTrackMonitoringAnalyzer::TrackCollectionTag_
private

Definition at line 81 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and BDHadronTrackMonitoringAnalyzer().

edm::InputTag BDHadronTrackMonitoringAnalyzer::TrackSrc_
private

Definition at line 75 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by BDHadronTrackMonitoringAnalyzer().

MonitorElement* BDHadronTrackMonitoringAnalyzer::TrkDxy_alljets[6]
private

Definition at line 109 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BDHadronTrackMonitoringAnalyzer::TrkDz_alljets[6]
private

Definition at line 112 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BDHadronTrackMonitoringAnalyzer::TrkEta_alljets[6]
private

Definition at line 105 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

const std::vector< std::string > BDHadronTrackMonitoringAnalyzer::TrkHistCat
static
Initial value:
= {
"BCWeakDecay", "BWeakDecay", "CWeakDecay", "PU", "Other", "Fake"}

Definition at line 61 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), bookHistograms(), BDHadronTrackMonitoringHarvester::dqmEndJob(), and toTrackRef().

MonitorElement* BDHadronTrackMonitoringAnalyzer::TrkHitAll_alljets[6]
private

Definition at line 115 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BDHadronTrackMonitoringAnalyzer::TrkHitPixel_alljets[6]
private

Definition at line 121 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BDHadronTrackMonitoringAnalyzer::TrkHitStrip_alljets[6]
private

Definition at line 118 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BDHadronTrackMonitoringAnalyzer::TrkPhi_alljets[6]
private

Definition at line 107 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BDHadronTrackMonitoringAnalyzer::TrkPt_alljets[6]
private

Definition at line 103 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * BDHadronTrackMonitoringAnalyzer::TrkTruthDxy_alljets[5]
private

Definition at line 109 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * BDHadronTrackMonitoringAnalyzer::TrkTruthDz_alljets[5]
private

Definition at line 112 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * BDHadronTrackMonitoringAnalyzer::TrkTruthEta_alljets[5]
private

Definition at line 105 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * BDHadronTrackMonitoringAnalyzer::TrkTruthHitAll_alljets[5]
private

Definition at line 115 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * BDHadronTrackMonitoringAnalyzer::TrkTruthHitPixel_alljets[5]
private

Definition at line 121 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * BDHadronTrackMonitoringAnalyzer::TrkTruthHitStrip_alljets[5]
private

Definition at line 118 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * BDHadronTrackMonitoringAnalyzer::TrkTruthPhi_alljets[5]
private

Definition at line 107 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().

MonitorElement * BDHadronTrackMonitoringAnalyzer::TrkTruthPt_alljets[5]
private

Definition at line 103 of file BDHadronTrackMonitoringAnalyzer.h.

Referenced by analyze(), and bookHistograms().