ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/HiggsAnalysis/EarlyDataStudy/plugins/HedsEventInfoDumper.cc
Revision: 1.1
Committed: Thu May 6 17:56:01 2010 UTC (14 years, 11 months ago) by mangano
Content type: text/plain
Branch: MAIN
CVS Tags: V00-02-00, V00-01-03, V00-01-02, V00-01-01, V00-01-00, V00-00-01, V00-00-00, HEAD
Log Message:
add event dumper from Giovanni

File Contents

# User Rev Content
1 mangano 1.1 // -*- C++ -*-
2     //
3     // Package: HedsEventInfoDumper
4     // Class: HedsEventInfoDumper
5     //
6     /**\class HedsEventInfoDumper
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     // $Id: HedsEventInfoDumper.cc,v 1.5 2010/04/26 10:04:56 gpetrucc Exp $
17     //
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     #include "DataFormats/Common/interface/TriggerResults.h"
29     #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     #include "DataFormats/JetReco/interface/Jet.h"
37     #include "DataFormats/TrackReco/interface/Track.h"
38     #include "DataFormats/VertexReco/interface/Vertex.h"
39     #include "FWCore/Common/interface/TriggerNames.h"
40     #include "FWCore/Framework/interface/Event.h"
41     #include "FWCore/Framework/interface/MakerMacros.h"
42     #include "FWCore/ParameterSet/interface/ParameterSet.h"
43     #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
44     #include "Utilities/General/interface/ClassName.h"
45    
46     #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     #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
54     #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
55     #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
56     #include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
57     #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     #include <boost/regex.hpp>
76     #define foreach BOOST_FOREACH
77    
78    
79     //
80     // class decleration
81     //
82    
83     namespace edm { typedef std::vector<InputTag> VInputTag; }
84    
85     class HedsEventInfoDumper : public edm::EDAnalyzer {
86     public:
87     explicit HedsEventInfoDumper(const edm::ParameterSet&);
88     ~HedsEventInfoDumper();
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     void dumpIP(const reco::Track *tk) ;
96     int muonStations(const reco::TrackRef track, int subdet=0, bool validOnly=true) ;
97    
98     edm::InputTag primVtx_;
99     edm::VInputTag electrons_;
100     edm::VInputTag muons_;
101     edm::VInputTag mets_;
102     edm::VInputTag jets_;
103     StringCutObjectSelector<reco::Jet> jetCut_;
104     edm::VInputTag composites_;
105    
106     edm::VInputTag triggerResults_;
107     std::vector<boost::regex> triggerNames_;
108    
109     edm::ESHandle<TransientTrackBuilder> theTTBuilder;
110    
111     const reco::Vertex * primaryVertex;
112     };
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     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    
142     //
143     // constructors and destructor
144     //
145     HedsEventInfoDumper::HedsEventInfoDumper(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     jets_(iConfig.existsAs<edm::VInputTag>("jets") ? iConfig.getParameter<edm::VInputTag>("jets") : edm::VInputTag()),
151     jetCut_(jets_.empty() ? "" : iConfig.getParameter<std::string>("jetCut"), true),
152     composites_(iConfig.getParameter<edm::VInputTag>("composites"))
153     {
154     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     }
164    
165    
166     HedsEventInfoDumper::~HedsEventInfoDumper()
167     {
168     }
169    
170     //
171     // member functions
172     //
173    
174     // ------------ method called to for each event ------------
175     void
176     HedsEventInfoDumper::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
177     {
178     iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theTTBuilder);
179     KalmanVertexFitter vtxFitter;
180    
181     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     //std::cout << "========================================================================" << std::endl;
186     std::cout << "=========================== RECORD SEPARATOR =======================" << std::endl;
187     std::cout << "=== RUN " << iEvent.id().run() << " LS " << iEvent.id().luminosityBlock() << " EVENT " << iEvent.id().event() << " ( " << time << " GMT ) ===" << std::endl;
188     std::cout << "========================================================================" << std::endl;
189    
190     edm::Handle<std::vector<reco::Vertex> > vertices;
191     iEvent.getByLabel(primVtx_, vertices);
192     primaryVertex = &(*vertices)[0];
193     const reco::MET *metRef = 0;
194     if (!mets_.empty()) {
195     edm::Handle<edm::View<reco::MET> > meth;
196     iEvent.getByLabel(mets_.front(), meth);
197     metRef = &(*meth)[0];
198     }
199    
200     std::vector<const reco::Jet *> jetsRef;
201     if (!jets_.empty()) jetsRef = getSorted<reco::Jet>(iEvent, jets_.front(), jetCut_);
202    
203     foreach(const edm::InputTag tag, electrons_) {
204     std::vector<const reco::GsfElectron *> electrons(getSorted<reco::GsfElectron>(iEvent, tag));
205     if (electrons.empty()) continue;
206     std::cout << "Electrons from " << tag.encode() << std::endl;
207     foreach(const reco::GsfElectron * el, electrons) {
208     printf(" pt %.3f, eta %.3f, phi %.3f, charge %d\n", el->pt(), el->eta(), el->phi(), el->charge());
209     printf(" fbrem %.3f, eSuperClusterOverP %.3f, eSeedClusterOverP %.3f, hcalOverEcal %.3f, deltaPhiSuperClusterTrackAtVtx %.3f, deltaEtaSuperClusterTrackAtVtx %.3f, sigmaIetaIeta %.4f\n",
210     el->fbrem(), el->eSuperClusterOverP(), el->eSeedClusterOverP(), el->hcalOverEcal(),
211     el->deltaPhiSuperClusterTrackAtVtx(), el->deltaEtaSuperClusterTrackAtVtx(), el->sigmaIetaIeta());
212     printf(" isolation sums (R=0.3): tk %.2f, ecal %.2f, hcal %.2f, calo %.2f, all %.2f [GeV]\n",
213     el->dr03TkSumPt(), el->dr03EcalRecHitSumEt(), el->dr03HcalTowerSumEt(),
214     el->dr03EcalRecHitSumEt() + el->dr03HcalTowerSumEt(),
215     el->dr03TkSumPt()+el->dr03EcalRecHitSumEt()+el->dr03HcalTowerSumEt());
216     dumpIP(el->bestTrack());
217     }
218     std::cout << std::endl;
219     }
220    
221     foreach(const edm::InputTag tag, muons_) {
222     std::vector<const reco::Muon *> muons(getSorted<reco::Muon>(iEvent, tag));
223     if (muons.empty()) continue;
224     std::cout << "muons from " << tag.encode() << std::endl;
225     foreach(const reco::Muon * mu, muons) {
226     printf(" pt %.3f, eta %.3f, phi %.3f, charge %d\n", mu->pt(), mu->eta(), mu->phi(), mu->charge());
227     printf(" ");
228     if (mu->isGlobalMuon()) printf("global ");
229     if (mu->isTrackerMuon()) printf("tracker ");
230     if (mu->isStandAloneMuon()) printf("standalone ");
231     if (muon::isGoodMuon(*mu, muon::TMLastStationAngTight)) printf("TMLastStationAngTight ");
232     if (muon::isGoodMuon(*mu, muon::TMLastStationTight)) printf("TMLastStationTight ");
233     if (muon::isGoodMuon(*mu, muon::GlobalMuonPromptTight)) printf("GlobalMuonPromptTight ");
234     if (mu->isTrackerMuon()) {
235     printf(", matches %d (%d arbitrated), segment compat. %.2f (%.2f arbitrated)",
236     mu->numberOfMatches(reco::Muon::SegmentArbitration),
237     mu->numberOfMatches(reco::Muon::SegmentAndTrackArbitration),
238     muon::segmentCompatibility(*mu, reco::Muon::SegmentArbitration),
239     muon::segmentCompatibility(*mu, reco::Muon::SegmentAndTrackArbitration));
240     }
241     if (mu->isCaloCompatibilityValid()) {
242     printf(", calo compat. %.2f", mu->caloCompatibility());
243     }
244     if (mu->isStandAloneMuon()) {
245     printf(", stations %d (dt %d, csc %d, rpc %d)",
246     muonStations(mu->standAloneMuon()),
247     muonStations(mu->standAloneMuon(), MuonSubdetId::DT),
248     muonStations(mu->standAloneMuon(), MuonSubdetId::CSC),
249     muonStations(mu->standAloneMuon(), MuonSubdetId::RPC));
250     }
251     printf("\n");
252     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",
253     mu->isolationR03().sumPt , mu->isolationR03().emEt, mu->isolationR03().hadEt,
254     mu->isolationR03().emEt + mu->isolationR03().hadEt,
255     mu->isolationR03().sumPt + mu->isolationR03().emEt + mu->isolationR03().hadEt,
256     mu->isolationR03().emVetoEt, mu->isolationR03().hadVetoEt, mu->isolationR03().emVetoEt + mu->isolationR03().hadVetoEt );
257     reco::TrackRef track = mu->track();
258     if (track.isNonnull()) {
259     printf(" tracker track: hits %d, pixel layers %d, expected hits inner %d, outer %d; chi2/ndf %.2f, %s, algo %d\n",
260     track->numberOfValidHits(), track->hitPattern().pixelLayersWithMeasurement(),
261     track->trackerExpectedHitsInner().numberOfLostHits(), track->trackerExpectedHitsOuter().numberOfLostHits(),
262     track->normalizedChi2(), track->quality(reco::Track::highPurity) ? "highPurity" : "not highPurity", track->algo());
263     dumpIP(track.get());
264     }
265     reco::TrackRef gtrack = mu->globalTrack();
266     if (mu->isGlobalMuon()) {
267     printf(" global track: silicon hits %d, muon hits %d, chi2/ndf %.2f\n",
268     mu->globalTrack()->hitPattern().numberOfValidTrackerHits(),
269     mu->globalTrack()->hitPattern().numberOfValidMuonHits(),
270     mu->globalTrack()->normalizedChi2());
271     }
272     }
273     std::cout << std::endl;
274     }
275    
276     foreach(const edm::InputTag tag, jets_) {
277     std::vector<const reco::Jet *> jets(getSorted<reco::Jet>(iEvent, tag, jetCut_));
278     if (jets.empty()) continue;
279     std::cout << "jets from " << tag.encode() << std::endl;
280     foreach(const reco::Jet * jet, jets) {
281     printf(" pt %8.2f, eta % .3f, phi % .3f\n", jet->pt(), jet->eta(), jet->phi());
282     }
283     std::cout << std::endl;
284     }
285    
286    
287     foreach(const edm::InputTag tag, mets_) {
288     edm::Handle<edm::View<reco::MET> > meth;
289     iEvent.getByLabel(tag, meth);
290     const reco::MET *met = &(*meth)[0];
291     printf("missing et %s: %.1f (signif. %.1f), phi %.2f, sumEt %.1f\n", tag.encode().c_str(), met->pt(), met->mEtSig(), met->phi(), met->sumEt());
292     if (met != metRef) {
293     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()));
294     }
295     }
296     std::cout << std::endl;
297    
298     foreach(const edm::InputTag tag, composites_) {
299     std::vector<const reco::CompositeCandidate *> comps(getSorted<reco::CompositeCandidate>(iEvent, tag));
300     if (comps.empty()) continue;
301     std::cout << "composite objects from " << tag.encode() << std::endl;
302     foreach(const reco::CompositeCandidate *c, comps) {
303     std::vector<reco::TransientTrack> tracks;
304     dumpCompositeCandidate(c, 2, tracks, iEvent);
305     if (tracks.size() > 1) {
306     TransientVertex myVertex = vtxFitter.vertex(tracks);
307     if (myVertex.isValid()) {
308     float vChi2 = myVertex.totalChiSquared();
309     float vNDF = myVertex.degreesOfFreedom();
310     float vProb(TMath::Prob(vChi2,(int)vNDF));
311     reco::Vertex myVtx(myVertex);
312     float vtxCompXY = VertexDistanceXY().compatibility( *primaryVertex, myVtx);
313     Measurement1D vtxDistXY = VertexDistanceXY().signedDistance(*primaryVertex, myVtx, GlobalVector(c->px(), c->py(), c->pz()));
314     float vtxComp3D = VertexDistance3D().compatibility( *primaryVertex, myVtx);
315     Measurement1D vtxDist3D = VertexDistance3D().signedDistance(*primaryVertex, myVtx, GlobalVector(c->px(), c->py(), c->pz()));
316     printf(" common vertex: %d tracks, chi2/ndf = %.2f/%.1f = %.2f (p-val. %.3f)\n", tracks.size(), vChi2, vNDF, vChi2/vNDF, vProb);
317     printf(" signed xy distance from prim. vtx %.4f mm (signif. %.3f), chi2 %.2f\n", vtxDistXY.value(), vtxDistXY.significance(), vtxCompXY);
318     printf(" signed 3d distance from prim. vtx %.4f mm (signif. %.3f), chi2 %.2f\n", vtxDist3D.value(), vtxDist3D.significance(), vtxComp3D);
319     } else {
320     printf(" common vertex fit for %d tracks failed.\n",tracks.size());
321    
322     }
323     }
324     }
325     std::cout << std::endl;
326     }
327    
328     foreach(const edm::InputTag &tag, triggerResults_) {
329     edm::Handle<edm::TriggerResults> results;
330     if (!iEvent.getByLabel(tag, results)) {
331     std::cout << " oops, trigger results from " << tag.encode() << " not found" << std::endl;
332     continue;
333     }
334     std::cout << "passing triggers from " << tag.process() << ": ";
335     const edm::TriggerNames &names = iEvent.triggerNames(*results);
336     for (unsigned i = 0; i < names.size(); ++i) {
337     if (!results->accept(i)) continue;
338     foreach(const boost::regex &re, triggerNames_) {
339     if ( boost::regex_match(names.triggerName(i), re) ) {
340     std::cout << names.triggerName(i) << " ";
341     }
342     }
343     }
344     std::cout << std::endl;
345     }
346     std::cout << "=========================== RECORD SEPARATOR =======================" << std::endl;
347     std::cout << std::endl;
348    
349     }
350    
351     void HedsEventInfoDumper::dumpCompositeCandidate(const reco::CompositeCandidate *c, int indent, std::vector<reco::TransientTrack> &tracks, const edm::Event &iEvent) {
352     double scalarPtSum = 0;
353     for (size_t i = 0; i < c->numberOfDaughters(); ++i) { scalarPtSum += c->daughter(i)->masterClone()->pt(); }
354     double mt = std::sqrt(scalarPtSum*scalarPtSum - c->pt()*c->pt());
355     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());
356     if (c->numberOfDaughters() == 2) {
357     printf("%*sacoplanarity %.2f, delta eta %.2f\n", indent+2, "",
358     M_PI - std::abs(reco::deltaPhi(c->daughter(0)->phi(), c->daughter(1)->phi())),
359     std::abs(c->daughter(0)->eta() - c->daughter(1)->eta()));
360     }
361     for (size_t i = 0; i < c->numberOfDaughters(); ++i) {
362     reco::CandidateBaseRef r = c->daughter(i)->masterClone();
363     const edm::Provenance & prov = iEvent.getProvenance(r.id());
364     printf("%*sdaughter pt %6.3f, eta % 5.3f, phi % 5.3f, charge % 1d, type %s, collection %s:%s:%s\n", indent+2,"",
365     r->pt(), r->eta(), r->phi(), r->charge(),
366     className<reco::Candidate>(*r).c_str(), prov.moduleLabel().c_str(), prov.productInstanceName().c_str(), prov.processName().c_str());
367     const reco::RecoCandidate *rc = dynamic_cast<const reco::RecoCandidate *>(r.get());
368     if (rc) {
369     const reco::Track *tk = rc->bestTrack();
370     if (tk) tracks.push_back(theTTBuilder->build(tk));
371     }
372     if (typeid(*r) == typeid(reco::CompositeCandidate)) {
373     dumpCompositeCandidate(static_cast<const reco::CompositeCandidate *>(r.get()), indent+2, tracks, iEvent);
374     }
375     }
376    
377     }
378    
379     void HedsEventInfoDumper::dumpIP(const reco::Track *t) {
380     if (t == 0) return;
381     GlobalVector dir(t->px(), t->py(), t->pz());
382     reco::TransientTrack ttk = theTTBuilder->build(t);
383     Measurement1D tip = IPTools::signedTransverseImpactParameter(ttk, dir, *primaryVertex).second;
384     Measurement1D lip = Measurement1D(t->dz(primaryVertex->position()), hypot(t->dzError(), primaryVertex->zError()));
385     Measurement1D ip3d = IPTools::signedImpactParameter3D(ttk, dir, *primaryVertex).second;
386     printf(" primary vertex association: dz %.4f mm (signif. %.3f), dxy %.4f mm (signif. %.3f), 3d %.4f mm (signig. %.3f)\n",
387     10*lip.value(), std::abs(lip.significance()),
388     10*tip.value(), std::abs(tip.significance()),
389     10*ip3d.value(), std::abs(ip3d.significance()));
390     }
391    
392     int HedsEventInfoDumper::muonStations(const reco::TrackRef track, int subdet, bool validOnly) {
393     int stations[4] = { 0,0,0,0 };
394     if (track.isNull()) return 0;
395     for (trackingRecHit_iterator it = track->recHitsBegin(), ed = track->recHitsEnd(); it != ed; ++it) {
396     DetId id = (*it)->geographicalId();
397     if (id.det() != DetId::Muon) continue;
398     if (subdet != 0 && id.subdetId() != subdet) continue;
399     if (validOnly && !(*it)->isValid()) continue;
400     switch (id.subdetId()) {
401     case MuonSubdetId::DT:
402     stations[DTLayerId(id).station()-1] = 1;
403     break;
404     case MuonSubdetId::CSC:
405     stations[CSCDetId(id).station()-1] = 1;
406     break;
407     case MuonSubdetId::RPC:
408     stations[RPCDetId(id).station()-1] = 1;
409     break;
410     }
411     }
412     return stations[0]+stations[1]+stations[2]+stations[3];
413     }
414    
415     //define this as a plug-in
416     DEFINE_FWK_MODULE(HedsEventInfoDumper);