ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/GPetrucc/plugins/EventInfoDumper.cc
Revision: 1.5
Committed: Mon Apr 26 10:04:56 2010 UTC (15 years ago) by gpetrucc
Content type: text/plain
Branch: MAIN
CVS Tags: HEAD
Changes since 1.4: +105 -40 lines
Log Message:
updates

File Contents

# User Rev Content
1 gpetrucc 1.1 // -*- C++ -*-
2     //
3     // Package: EventInfoDumper
4     // Class: EventInfoDumper
5     //
6     /**\class EventInfoDumper EventInfoDumper.cc Gio/EventInfoDumper/src/EventInfoDumper.cc
7    
8     Description: <one line class summary>
9    
10     Implementation:
11     <Notes on implementation>
12     */
13     //
14     // Original Author: Sep 15 09:45
15     // Created: Mon Sep 15 09:49:08 CEST 2008
16 gpetrucc 1.5 // $Id: EventInfoDumper.cc,v 1.4 2010/04/19 07:57:45 gpetrucc Exp $
17 gpetrucc 1.1 //
18     //
19    
20    
21     // system include files
22     #include <memory>
23    
24     // user include files
25     #include "FWCore/Framework/interface/Frameworkfwd.h"
26     #include "FWCore/Framework/interface/EDAnalyzer.h"
27    
28 gpetrucc 1.2 #include "DataFormats/Common/interface/TriggerResults.h"
29 gpetrucc 1.1 #include "DataFormats/Common/interface/View.h"
30     #include "DataFormats/Math/interface/deltaR.h"
31     #include "DataFormats/MuonReco/interface/Muon.h"
32     #include "DataFormats/MuonReco/interface/MuonSelectors.h"
33     #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
34     #include "DataFormats/METReco/interface/MET.h"
35     #include "DataFormats/Candidate/interface/CompositeCandidate.h"
36 gpetrucc 1.3 #include "DataFormats/JetReco/interface/Jet.h"
37 gpetrucc 1.1 #include "DataFormats/TrackReco/interface/Track.h"
38     #include "DataFormats/VertexReco/interface/Vertex.h"
39 gpetrucc 1.2 #include "FWCore/Common/interface/TriggerNames.h"
40 gpetrucc 1.1 #include "FWCore/Framework/interface/Event.h"
41     #include "FWCore/Framework/interface/MakerMacros.h"
42     #include "FWCore/ParameterSet/interface/ParameterSet.h"
43 gpetrucc 1.3 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
44 gpetrucc 1.1 #include "Utilities/General/interface/ClassName.h"
45    
46 gpetrucc 1.5 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
47     #include "DataFormats/MuonDetId/interface/CSCDetId.h"
48     #include "DataFormats/MuonDetId/interface/DTLayerId.h"
49     #include "DataFormats/MuonDetId/interface/MuonSubdetId.h"
50     #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
51     #include "TrackingTools/IPTools/interface/IPTools.h"
52     #include "TrackingTools/Records/interface/TransientTrackRecord.h"
53 gpetrucc 1.1 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
54     #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
55 gpetrucc 1.5 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
56     #include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
57 gpetrucc 1.1 #if 0
58     #include "Geometry/CommonDetUnit/interface/GeomDet.h"
59     #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
60     #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
61     #include "MagneticField/Engine/interface/MagneticField.h"
62     #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
63     #include "TrackingTools/GeomPropagators/interface/Propagator.h"
64     #include "TrackingTools/IPTools/interface/IPTools.h"
65     #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
66     #include "TrackingTools/Records/interface/TransientTrackRecord.h"
67     #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
68     #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
69     #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
70     #include <TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h>
71     #include <TrackingTools/PatternTools/interface/TwoTrackMinimumDistance.h>
72     #endif
73    
74     #include <boost/foreach.hpp>
75 gpetrucc 1.2 #include <boost/regex.hpp>
76 gpetrucc 1.1 #define foreach BOOST_FOREACH
77    
78    
79     //
80     // class decleration
81     //
82    
83     namespace edm { typedef std::vector<InputTag> VInputTag; }
84    
85     class EventInfoDumper : public edm::EDAnalyzer {
86     public:
87     explicit EventInfoDumper(const edm::ParameterSet&);
88     ~EventInfoDumper();
89    
90    
91     private:
92     virtual void analyze(const edm::Event&, const edm::EventSetup&);
93    
94     void dumpCompositeCandidate(const reco::CompositeCandidate *c, int indent, std::vector<reco::TransientTrack> &tracks, const edm::Event &iEvent) ;
95 gpetrucc 1.5 void dumpIP(const reco::Track *tk) ;
96     int muonStations(const reco::TrackRef track, int subdet=0, bool validOnly=true) ;
97 gpetrucc 1.1
98     edm::InputTag primVtx_;
99     edm::VInputTag electrons_;
100     edm::VInputTag muons_;
101     edm::VInputTag mets_;
102 gpetrucc 1.3 edm::VInputTag jets_;
103     StringCutObjectSelector<reco::Jet> jetCut_;
104 gpetrucc 1.1 edm::VInputTag composites_;
105    
106 gpetrucc 1.2 edm::VInputTag triggerResults_;
107     std::vector<boost::regex> triggerNames_;
108    
109 gpetrucc 1.1 edm::ESHandle<TransientTrackBuilder> theTTBuilder;
110 gpetrucc 1.5
111     const reco::Vertex * primaryVertex;
112 gpetrucc 1.1 };
113    
114     struct PtrPtSort {
115     bool operator()(const reco::Candidate *c1, const reco::Candidate *c2) const { return c1->pt() > c2->pt(); }
116     };
117     template<typename T>
118     const std::vector<const T*> getSorted(const edm::Event &iev, const edm::InputTag tag) {
119     edm::Handle<edm::View<T> > view;
120     std::vector<const T *> pointers;
121     if (iev.getByLabel(tag, view)) {
122     foreach(const T &t, *view) { pointers.push_back(&t); }
123     std::sort(pointers.begin(), pointers.end(), PtrPtSort());
124     } else {
125     std::cout << " oops, " << tag << " is missing." << std::endl;
126     }
127     return pointers;
128     }
129 gpetrucc 1.3 template<typename T, typename S>
130     const std::vector<const T*> getSorted(const edm::Event &iev, const edm::InputTag tag, const S &selector) {
131     edm::Handle<edm::View<T> > view;
132     std::vector<const T *> pointers;
133     if (iev.getByLabel(tag, view)) {
134     foreach(const T &t, *view) { if (selector(t)) pointers.push_back(&t); }
135     std::sort(pointers.begin(), pointers.end(), PtrPtSort());
136     } else {
137     std::cout << " oops, " << tag << " is missing." << std::endl;
138     }
139     return pointers;
140     }
141 gpetrucc 1.1
142     //
143     // constructors and destructor
144     //
145     EventInfoDumper::EventInfoDumper(const edm::ParameterSet& iConfig) :
146     primVtx_(iConfig.getParameter<edm::InputTag>("primaryVertices")),
147     electrons_(iConfig.getParameter<edm::VInputTag>("electrons")),
148     muons_(iConfig.getParameter<edm::VInputTag>("muons")),
149     mets_(iConfig.getParameter<edm::VInputTag>("mets")),
150 gpetrucc 1.3 jets_(iConfig.existsAs<edm::VInputTag>("jets") ? iConfig.getParameter<edm::VInputTag>("jets") : edm::VInputTag()),
151     jetCut_(jets_.empty() ? "" : iConfig.getParameter<std::string>("jetCut"), true),
152 gpetrucc 1.1 composites_(iConfig.getParameter<edm::VInputTag>("composites"))
153     {
154 gpetrucc 1.2 if (iConfig.existsAs<edm::VInputTag>("triggerResults")) {
155     triggerResults_ = iConfig.getParameter<edm::VInputTag>("triggerResults");
156     }
157     if (iConfig.existsAs<std::vector<std::string> >("triggerNames")) {
158     std::vector<std::string> names = iConfig.getParameter<std::vector<std::string> >("triggerNames");
159     foreach(const std::string &n, names) {
160     triggerNames_.push_back(boost::regex(n));
161     }
162     }
163 gpetrucc 1.1 }
164    
165    
166     EventInfoDumper::~EventInfoDumper()
167     {
168     }
169    
170     //
171     // member functions
172     //
173    
174     // ------------ method called to for each event ------------
175     void
176     EventInfoDumper::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
177     {
178     iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theTTBuilder);
179     KalmanVertexFitter vtxFitter;
180    
181 gpetrucc 1.2 time_t t(iEvent.time().value() >> 32);
182     std::string time( asctime( gmtime(&t) ) );
183     size_t pos = time.find('\n');
184     if ( pos != std::string::npos ) time = time.substr(0,pos);
185 gpetrucc 1.1 std::cout << "========================================================================" << std::endl;
186 gpetrucc 1.2 std::cout << "=== RUN " << iEvent.id().run() << " LS " << iEvent.id().luminosityBlock() << " EVENT " << iEvent.id().event() << " ( " << time << " GMT ) ===" << std::endl;
187 gpetrucc 1.1 std::cout << "========================================================================" << std::endl;
188    
189     edm::Handle<std::vector<reco::Vertex> > vertices;
190     iEvent.getByLabel(primVtx_, vertices);
191 gpetrucc 1.5 primaryVertex = &(*vertices)[0];
192 gpetrucc 1.1 const reco::MET *metRef = 0;
193     if (!mets_.empty()) {
194     edm::Handle<edm::View<reco::MET> > meth;
195     iEvent.getByLabel(mets_.front(), meth);
196     metRef = &(*meth)[0];
197     }
198    
199 gpetrucc 1.3 std::vector<const reco::Jet *> jetsRef;
200     if (!jets_.empty()) jetsRef = getSorted<reco::Jet>(iEvent, jets_.front(), jetCut_);
201    
202 gpetrucc 1.1 foreach(const edm::InputTag tag, electrons_) {
203     std::vector<const reco::GsfElectron *> electrons(getSorted<reco::GsfElectron>(iEvent, tag));
204     if (electrons.empty()) continue;
205     std::cout << "Electrons from " << tag.encode() << std::endl;
206     foreach(const reco::GsfElectron * el, electrons) {
207     printf(" pt %.3f, eta %.3f, phi %.3f, charge %d\n", el->pt(), el->eta(), el->phi(), el->charge());
208 gpetrucc 1.5 printf(" fbrem %.3f, eSuperClusterOverP %.3f, eSeedClusterOverP %.3f, hcalOverEcal %.3f, deltaPhiSuperClusterTrackAtVtx %.3f, deltaEtaSuperClusterTrackAtVtx %.3f, sigmaIetaIeta %.4f\n",
209     el->fbrem(), el->eSuperClusterOverP(), el->eSeedClusterOverP(), el->hcalOverEcal(),
210     el->deltaPhiSuperClusterTrackAtVtx(), el->deltaEtaSuperClusterTrackAtVtx(), el->sigmaIetaIeta());
211 gpetrucc 1.1 printf(" isolation sums (R=0.3): tk %.2f, ecal %.2f, hcal %.2f, calo %.2f, all %.2f [GeV]\n",
212     el->dr03TkSumPt(), el->dr03EcalRecHitSumEt(), el->dr03HcalTowerSumEt(),
213     el->dr03EcalRecHitSumEt() + el->dr03HcalTowerSumEt(),
214     el->dr03TkSumPt()+el->dr03EcalRecHitSumEt()+el->dr03HcalTowerSumEt());
215 gpetrucc 1.5 dumpIP(el->bestTrack());
216 gpetrucc 1.1 }
217     std::cout << std::endl;
218     }
219    
220     foreach(const edm::InputTag tag, muons_) {
221     std::vector<const reco::Muon *> muons(getSorted<reco::Muon>(iEvent, tag));
222     if (muons.empty()) continue;
223     std::cout << "muons from " << tag.encode() << std::endl;
224     foreach(const reco::Muon * mu, muons) {
225     printf(" pt %.3f, eta %.3f, phi %.3f, charge %d\n", mu->pt(), mu->eta(), mu->phi(), mu->charge());
226     printf(" ");
227     if (mu->isGlobalMuon()) printf("global ");
228     if (mu->isTrackerMuon()) printf("tracker ");
229     if (mu->isStandAloneMuon()) printf("standalone ");
230     if (muon::isGoodMuon(*mu, muon::TMLastStationAngTight)) printf("TMLastStationAngTight ");
231 gpetrucc 1.5 if (muon::isGoodMuon(*mu, muon::TMLastStationTight)) printf("TMLastStationTight ");
232     if (muon::isGoodMuon(*mu, muon::GlobalMuonPromptTight)) printf("GlobalMuonPromptTight ");
233     if (mu->isTrackerMuon()) {
234     printf(", matches %d (%d arbitrated), segment compat. %.2f (%.2f arbitrated)",
235 gpetrucc 1.3 mu->numberOfMatches(reco::Muon::SegmentArbitration),
236     mu->numberOfMatches(reco::Muon::SegmentAndTrackArbitration),
237     muon::segmentCompatibility(*mu, reco::Muon::SegmentArbitration),
238 gpetrucc 1.5 muon::segmentCompatibility(*mu, reco::Muon::SegmentAndTrackArbitration));
239     }
240     if (mu->isCaloCompatibilityValid()) {
241     printf(", calo compat. %.2f", mu->caloCompatibility());
242     }
243     if (mu->isStandAloneMuon()) {
244     printf(", stations %d (dt %d, csc %d, rpc %d)",
245     muonStations(mu->standAloneMuon()),
246     muonStations(mu->standAloneMuon(), MuonSubdetId::DT),
247     muonStations(mu->standAloneMuon(), MuonSubdetId::CSC),
248     muonStations(mu->standAloneMuon(), MuonSubdetId::RPC));
249     }
250     printf("\n");
251 gpetrucc 1.4 printf(" isolation sums (R=0.3): tk %.2f, ecal %.2f, hcal %.2f, calo %.2f, all %.2f; in veto cone ecal %.2f, hcal %.2f, calo %.2f [GeV]\n",
252 gpetrucc 1.1 mu->isolationR03().sumPt , mu->isolationR03().emEt, mu->isolationR03().hadEt,
253     mu->isolationR03().emEt + mu->isolationR03().hadEt,
254     mu->isolationR03().sumPt + mu->isolationR03().emEt + mu->isolationR03().hadEt,
255     mu->isolationR03().emVetoEt, mu->isolationR03().hadVetoEt, mu->isolationR03().emVetoEt + mu->isolationR03().hadVetoEt );
256     reco::TrackRef track = mu->track();
257     if (track.isNonnull()) {
258     printf(" tracker track: hits %d, pixel layers %d, expected hits inner %d, outer %d; chi2/ndf %.2f, %s, algo %d\n",
259     track->numberOfValidHits(), track->hitPattern().pixelLayersWithMeasurement(),
260     track->trackerExpectedHitsInner().numberOfLostHits(), track->trackerExpectedHitsOuter().numberOfLostHits(),
261     track->normalizedChi2(), track->quality(reco::Track::highPurity) ? "highPurity" : "not highPurity", track->algo());
262 gpetrucc 1.5 dumpIP(track.get());
263 gpetrucc 1.1 }
264 gpetrucc 1.2 reco::TrackRef gtrack = mu->globalTrack();
265     if (mu->isGlobalMuon()) {
266     printf(" global track: silicon hits %d, muon hits %d, chi2/ndf %.2f\n",
267     mu->globalTrack()->hitPattern().numberOfValidTrackerHits(),
268     mu->globalTrack()->hitPattern().numberOfValidMuonHits(),
269     mu->globalTrack()->normalizedChi2());
270     }
271 gpetrucc 1.1 }
272     std::cout << std::endl;
273     }
274    
275 gpetrucc 1.3 foreach(const edm::InputTag tag, jets_) {
276     std::vector<const reco::Jet *> jets(getSorted<reco::Jet>(iEvent, tag, jetCut_));
277     if (jets.empty()) continue;
278     std::cout << "jets from " << tag.encode() << std::endl;
279     foreach(const reco::Jet * jet, jets) {
280     printf(" pt %8.2f, eta % .3f, phi % .3f\n", jet->pt(), jet->eta(), jet->phi());
281     }
282     std::cout << std::endl;
283     }
284    
285    
286 gpetrucc 1.1 foreach(const edm::InputTag tag, mets_) {
287     edm::Handle<edm::View<reco::MET> > meth;
288     iEvent.getByLabel(tag, meth);
289     const reco::MET *met = &(*meth)[0];
290     printf("missing et %s: %.1f (signif. %.1f), phi %.2f, sumEt %.1f\n", tag.encode().c_str(), met->pt(), met->mEtSig(), met->phi(), met->sumEt());
291     if (met != metRef) {
292     printf(" with respect to reference %s: diff/sum %.2f, deltaPhi %.2f\n", mets_.front().encode().c_str(), (met->pt()-metRef->pt())/(met->pt()+metRef->pt()), deltaPhi(met->phi(), metRef->phi()));
293     }
294     }
295     std::cout << std::endl;
296    
297     foreach(const edm::InputTag tag, composites_) {
298     std::vector<const reco::CompositeCandidate *> comps(getSorted<reco::CompositeCandidate>(iEvent, tag));
299     if (comps.empty()) continue;
300     std::cout << "composite objects from " << tag.encode() << std::endl;
301     foreach(const reco::CompositeCandidate *c, comps) {
302     std::vector<reco::TransientTrack> tracks;
303     dumpCompositeCandidate(c, 2, tracks, iEvent);
304     if (tracks.size() > 1) {
305     TransientVertex myVertex = vtxFitter.vertex(tracks);
306     if (myVertex.isValid()) {
307     float vChi2 = myVertex.totalChiSquared();
308     float vNDF = myVertex.degreesOfFreedom();
309     float vProb(TMath::Prob(vChi2,(int)vNDF));
310 gpetrucc 1.5 reco::Vertex myVtx(myVertex);
311     float vtxCompXY = VertexDistanceXY().compatibility( *primaryVertex, myVtx);
312     Measurement1D vtxDistXY = VertexDistanceXY().signedDistance(*primaryVertex, myVtx, GlobalVector(c->px(), c->py(), c->pz()));
313     float vtxComp3D = VertexDistance3D().compatibility( *primaryVertex, myVtx);
314     Measurement1D vtxDist3D = VertexDistance3D().signedDistance(*primaryVertex, myVtx, GlobalVector(c->px(), c->py(), c->pz()));
315 gpetrucc 1.1 printf(" common vertex: %d tracks, chi2/ndf = %.2f/%.1f = %.2f (p-val. %.3f)\n", tracks.size(), vChi2, vNDF, vChi2/vNDF, vProb);
316 gpetrucc 1.5 printf(" signed xy distance from prim. vtx %.4f mm (signif. %.3f), chi2 %.2f\n", vtxDistXY.value(), vtxDistXY.significance(), vtxCompXY);
317     printf(" signed 3d distance from prim. vtx %.4f mm (signif. %.3f), chi2 %.2f\n", vtxDist3D.value(), vtxDist3D.significance(), vtxComp3D);
318 gpetrucc 1.1 } else {
319     printf(" common vertex fit for %d tracks failed.\n",tracks.size());
320    
321     }
322     }
323     }
324     std::cout << std::endl;
325     }
326    
327 gpetrucc 1.2 foreach(const edm::InputTag &tag, triggerResults_) {
328     edm::Handle<edm::TriggerResults> results;
329     if (!iEvent.getByLabel(tag, results)) {
330     std::cout << " oops, trigger results from " << tag.encode() << " not found" << std::endl;
331     continue;
332     }
333     std::cout << "passing triggers from " << tag.process() << ": ";
334     const edm::TriggerNames &names = iEvent.triggerNames(*results);
335     for (unsigned i = 0; i < names.size(); ++i) {
336     if (!results->accept(i)) continue;
337     foreach(const boost::regex &re, triggerNames_) {
338     if ( boost::regex_match(names.triggerName(i), re) ) {
339     std::cout << names.triggerName(i) << " ";
340     }
341     }
342     }
343     std::cout << std::endl;
344     }
345 gpetrucc 1.1 std::cout << std::endl;
346    
347     }
348    
349     void EventInfoDumper::dumpCompositeCandidate(const reco::CompositeCandidate *c, int indent, std::vector<reco::TransientTrack> &tracks, const edm::Event &iEvent) {
350 gpetrucc 1.5 double scalarPtSum = 0;
351     for (size_t i = 0; i < c->numberOfDaughters(); ++i) { scalarPtSum += c->daughter(i)->masterClone()->pt(); }
352     double mt = std::sqrt(scalarPtSum*scalarPtSum - c->pt()*c->pt());
353     printf("%*smass %.3f, mt %.2f, pt %.3f, eta %.3f, phi %.3f, charge %d\n", indent,"",c->mass(), mt, c->pt(), c->eta(), c->phi(), c->charge());
354     if (c->numberOfDaughters() == 2) {
355     printf("%*sacoplanarity %.2f, delta eta %.2f\n", indent+2, "",
356     M_PI - std::abs(reco::deltaPhi(c->daughter(0)->phi(), c->daughter(1)->phi())),
357     std::abs(c->daughter(0)->eta() - c->daughter(1)->eta()));
358     }
359     for (size_t i = 0; i < c->numberOfDaughters(); ++i) {
360     reco::CandidateBaseRef r = c->daughter(i)->masterClone();
361     const edm::Provenance & prov = iEvent.getProvenance(r.id());
362     printf("%*sdaughter pt %6.3f, eta % 5.3f, phi % 5.3f, charge % 1d, type %s, collection %s:%s:%s\n", indent+2,"",
363     r->pt(), r->eta(), r->phi(), r->charge(),
364     className<reco::Candidate>(*r).c_str(), prov.moduleLabel().c_str(), prov.productInstanceName().c_str(), prov.processName().c_str());
365     const reco::RecoCandidate *rc = dynamic_cast<const reco::RecoCandidate *>(r.get());
366     if (rc) {
367     const reco::Track *tk = rc->bestTrack();
368     if (tk) tracks.push_back(theTTBuilder->build(tk));
369     }
370     if (typeid(*r) == typeid(reco::CompositeCandidate)) {
371     dumpCompositeCandidate(static_cast<const reco::CompositeCandidate *>(r.get()), indent+2, tracks, iEvent);
372     }
373     }
374    
375     }
376    
377     void EventInfoDumper::dumpIP(const reco::Track *t) {
378     if (t == 0) return;
379     GlobalVector dir(t->px(), t->py(), t->pz());
380     reco::TransientTrack ttk = theTTBuilder->build(t);
381     Measurement1D tip = IPTools::signedTransverseImpactParameter(ttk, dir, *primaryVertex).second;
382     Measurement1D lip = Measurement1D(t->dz(primaryVertex->position()), hypot(t->dzError(), primaryVertex->zError()));
383     Measurement1D ip3d = IPTools::signedImpactParameter3D(ttk, dir, *primaryVertex).second;
384     printf(" primary vertex association: dz %.4f mm (signif. %.3f), dxy %.4f mm (signif. %.3f), 3d %.4f mm (signig. %.3f)\n",
385     10*lip.value(), std::abs(lip.significance()),
386     10*tip.value(), std::abs(tip.significance()),
387     10*ip3d.value(), std::abs(ip3d.significance()));
388     }
389 gpetrucc 1.1
390 gpetrucc 1.5 int EventInfoDumper::muonStations(const reco::TrackRef track, int subdet, bool validOnly) {
391     int stations[4] = { 0,0,0,0 };
392     if (track.isNull()) return 0;
393     for (trackingRecHit_iterator it = track->recHitsBegin(), ed = track->recHitsEnd(); it != ed; ++it) {
394     DetId id = (*it)->geographicalId();
395     if (id.det() != DetId::Muon) continue;
396     if (subdet != 0 && id.subdetId() != subdet) continue;
397     if (validOnly && !(*it)->isValid()) continue;
398     switch (id.subdetId()) {
399     case MuonSubdetId::DT:
400     stations[DTLayerId(id).station()-1] = 1;
401     break;
402     case MuonSubdetId::CSC:
403     stations[CSCDetId(id).station()-1] = 1;
404     break;
405     case MuonSubdetId::RPC:
406     stations[RPCDetId(id).station()-1] = 1;
407     break;
408     }
409     }
410     return stations[0]+stations[1]+stations[2]+stations[3];
411 gpetrucc 1.1 }
412 gpetrucc 1.5
413 gpetrucc 1.1 //define this as a plug-in
414     DEFINE_FWK_MODULE(EventInfoDumper);