CMS 3D CMS Logo

PhysicsObjectsMonitor.cc
Go to the documentation of this file.
1 
10 
11 // Collaborating Class Header
15 
17 
21 
23 
27 
28 using namespace std;
29 using namespace edm;
30 
33  theSTAMuonLabel = pset.getUntrackedParameter<string>("StandAloneTrackCollectionLabel");
34  theSeedCollectionLabel = pset.getUntrackedParameter<string>("MuonSeedCollectionLabel");
35  theDataType = pset.getUntrackedParameter<string>("DataType");
36 
37  if (theDataType != "RealData" && theDataType != "SimData")
38  edm::LogInfo("PhysicsObjectsMonitor") << "Error in Data Type!!" << endl;
39 
40  if (theDataType == "SimData") {
41  edm::LogInfo("PhysicsObjectsMonitor") << "Sorry! Running this package on simulation is no longer supported! ";
42  }
43 
44  // set Token(-s)
45  theSTAMuonToken_ =
46  consumes<reco::TrackCollection>(pset.getUntrackedParameter<string>("StandAloneTrackCollectionLabel"));
47 }
48 
51 
53  iBooker.setCurrentFolder("PhysicsObjects/MuonReconstruction");
54 
55  hPres = iBooker.book1D("pTRes", "pT Resolution", 100, -2, 2);
56  h1_Pres = iBooker.book1D("invPTRes", "1/pT Resolution", 100, -2, 2);
57 
58  charge = iBooker.book1D("charge", "track charge", 5, -2.5, 2.5);
59  ptot = iBooker.book1D("ptot", "track momentum", 50, 0, 50);
60  pt = iBooker.book1D("pt", "track pT", 100, 0, 50);
61  px = iBooker.book1D("px ", "track px", 100, -50, 50);
62  py = iBooker.book1D("py", "track py", 100, -50, 50);
63  pz = iBooker.book1D("pz", "track pz", 100, -50, 50);
64  Nmuon = iBooker.book1D("Nmuon", "Number of muon tracks", 11, -.5, 10.5);
65  Nrechits = iBooker.book1D("Nrechits", "Number of RecHits/Segments on track", 21, -.5, 21.5);
66  NDThits = iBooker.book1D("NDThits", "Number of DT Hits/Segments on track", 31, -.5, 31.5);
67  NCSChits = iBooker.book1D("NCSChits", "Number of CSC Hits/Segments on track", 31, -.5, 31.5);
68  NRPChits = iBooker.book1D("NRPChits", "Number of RPC hits on track", 11, -.5, 11.5);
69 
70  DTvsCSC = iBooker.book2D("DTvsCSC", "Number of DT vs CSC hits on track", 29, -.5, 28.5, 29, -.5, 28.5);
71  TH2F *root_ob = DTvsCSC->getTH2F();
72  root_ob->SetXTitle("Number of DT hits");
73  root_ob->SetYTitle("Number of CSC hits");
74 }
75 
76 void PhysicsObjectsMonitor::analyze(const Event &event, const EventSetup &eventSetup) {
77  edm::LogInfo("PhysicsObjectsMonitor") << "Run: " << event.id().run() << " Event: " << event.id().event();
79 
80  // Get the RecTrack collection from the event
82  event.getByToken(theSTAMuonToken_, staTracks);
83 
84  ESHandle<MagneticField> theMGField;
85  eventSetup.get<IdealMagneticFieldRecord>().get(theMGField);
86 
87  double recPt = 0.;
88  double simPt = 0.;
89 
90  reco::TrackCollection::const_iterator staTrack;
91 
92  edm::LogInfo("PhysicsObjectsMonitor") << "Reconstructed tracks: " << staTracks->size() << endl;
93  Nmuon->Fill(staTracks->size());
94  for (staTrack = staTracks->begin(); staTrack != staTracks->end(); ++staTrack) {
95  reco::TransientTrack track(*staTrack, &*theMGField);
96 
97  int nrechits = 0;
98  int nDThits = 0;
99  int nCSChits = 0;
100  int nRPChits = 0;
101 
102  for (trackingRecHit_iterator it = track.recHitsBegin(); it != track.recHitsEnd(); it++) {
103  if ((*it)->isValid()) {
104  edm::LogInfo("PhysicsObjectsMonitor") << "Analyzer: Aha this looks like a Rechit!" << std::endl;
105  if ((*it)->geographicalId().subdetId() == MuonSubdetId::DT) {
106  nDThits++;
107  } else if ((*it)->geographicalId().subdetId() == MuonSubdetId::CSC) {
108  nCSChits++;
109  } else if ((*it)->geographicalId().subdetId() == MuonSubdetId::RPC) {
110  nRPChits++;
111  } else {
112  edm::LogInfo("PhysicsObjectsMonitor") << "This is an UNKNOWN hit !! " << std::endl;
113  }
114  nrechits++;
115  }
116  }
117 
118  Nrechits->Fill(nrechits);
119  NDThits->Fill(nDThits);
120  NCSChits->Fill(nCSChits);
121  DTvsCSC->Fill(nDThits, nCSChits);
122  NRPChits->Fill(nRPChits);
123 
124  debug.dumpFTS(track.impactPointTSCP().theState());
125 
126  recPt = track.impactPointTSCP().momentum().perp();
127  edm::LogInfo("PhysicsObjectsMonitor")
128  << " p: " << track.impactPointTSCP().momentum().mag() << " pT: " << recPt << endl;
129  pt->Fill(recPt);
130  ptot->Fill(track.impactPointTSCP().momentum().mag());
131  charge->Fill(track.impactPointTSCP().charge());
132  px->Fill(track.impactPointTSCP().momentum().x());
133  py->Fill(track.impactPointTSCP().momentum().y());
134  pz->Fill(track.impactPointTSCP().momentum().z());
135  }
136  edm::LogInfo("PhysicsObjectsMonitor") << "---" << endl;
137  if (recPt && theDataType == "SimData") {
138  hPres->Fill((recPt - simPt) / simPt);
139  h1_Pres->Fill((1 / recPt - 1 / simPt) / (1 / simPt));
140  }
141 }
142 
T getUntrackedParameter(std::string const &, T const &) const
T perp() const
Definition: PV3DBase.h:72
TrajectoryStateClosestToPoint impactPointTSCP() const
const FreeTrajectoryState & theState() const
T y() const
Definition: PV3DBase.h:63
PhysicsObjectsMonitor(const edm::ParameterSet &pset)
Constructor.
std::string dumpFTS(const FreeTrajectoryState &fts) const
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T mag() const
Definition: PV3DBase.h:67
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
T z() const
Definition: PV3DBase.h:64
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
TH2F * getTH2F() const
trackingRecHit_iterator recHitsEnd() const
last iterator to RecHits
#define debug
Definition: HDRShower.cc:19
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
static constexpr int RPC
Definition: MuonSubdetId.h:14
~PhysicsObjectsMonitor() override
Destructor.
HLT enums.
T get() const
Definition: EventSetup.h:71
double recPt
static constexpr int DT
Definition: MuonSubdetId.h:12
static constexpr int CSC
Definition: MuonSubdetId.h:13
T x() const
Definition: PV3DBase.h:62
trackingRecHit_iterator recHitsBegin() const
first iterator to RecHits
Definition: event.py:1
Definition: Run.h:45