ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/ShallowTools/plugins/ShallowSimTracksProducer.cc
Revision: 1.2
Committed: Fri Jul 10 16:57:49 2009 UTC (15 years, 10 months ago) by bbetchar
Content type: text/plain
Branch: MAIN
Changes since 1.1: +4 -0 lines
Log Message:
Add qoverp and qoverperr

File Contents

# Content
1 #include "UserCode/ShallowTools/interface/ShallowSimTracksProducer.h"
2 #include "UserCode/ShallowTools/interface/ShallowTools.h"
3
4 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
5 #include "SimTracker/Records/interface/TrackAssociatorRecord.h"
6 #include "SimTracker/TrackAssociation/interface/TrackAssociatorByChi2.h"
7 #include "SimTracker/TrackAssociation/interface/TrackAssociatorByHits.h"
8 #include "DataFormats/RecoCandidate/interface/TrackAssociation.h"
9
10 #include "FWCore/MessageLogger/interface/MessageLogger.h"
11
12 ShallowSimTracksProducer::ShallowSimTracksProducer(const edm::ParameterSet& conf)
13 : Prefix( conf.getParameter<std::string>("Prefix") ),
14 Suffix( conf.getParameter<std::string>("Suffix") ),
15 trackingParticles_tag( conf.getParameter<edm::InputTag>("TrackingParticles")),
16 associator_tag( conf.getParameter<edm::ESInputTag>("Associator")),
17 tracks_tag( conf.getParameter<edm::InputTag>("Tracks"))
18 {
19 produces <std::vector<unsigned> > ( Prefix + "multi" + Suffix );
20 produces <std::vector<int> > ( Prefix + "type" + Suffix );
21 produces <std::vector<float> > ( Prefix + "charge" + Suffix );
22 produces <std::vector<float> > ( Prefix + "momentum" + Suffix );
23 produces <std::vector<float> > ( Prefix + "pt" + Suffix );
24 produces <std::vector<double> > ( Prefix + "theta" + Suffix );
25 produces <std::vector<double> > ( Prefix + "phi" + Suffix );
26 produces <std::vector<double> > ( Prefix + "eta" + Suffix );
27 produces <std::vector<double> > ( Prefix + "qoverp" + Suffix );
28 produces <std::vector<double> > ( Prefix + "vx" + Suffix );
29 produces <std::vector<double> > ( Prefix + "vy" + Suffix );
30 produces <std::vector<double> > ( Prefix + "vz" + Suffix );
31 }
32
33
34 void ShallowSimTracksProducer::
35 produce(edm::Event& event, const edm::EventSetup& setup) {
36
37 edm::Handle<edm::View<reco::Track> > tracks ; event.getByLabel( tracks_tag, tracks);
38 edm::Handle<TrackingParticleCollection> trackingParticles ; event.getByLabel( trackingParticles_tag, trackingParticles );
39 edm::ESHandle<TrackAssociatorBase> associator ; setup.get<TrackAssociatorRecord>().get( associator_tag, associator);
40
41 unsigned size = tracks->size();
42 std::auto_ptr<std::vector<unsigned> > multi ( new std::vector<unsigned>(size, 0));
43 std::auto_ptr<std::vector<int> > type ( new std::vector<int> (size, 0));
44 std::auto_ptr<std::vector<float> > charge ( new std::vector<float> (size, 0));
45 std::auto_ptr<std::vector<float> > momentum ( new std::vector<float> (size, -1));
46 std::auto_ptr<std::vector<float> > pt ( new std::vector<float> (size, -1));
47 std::auto_ptr<std::vector<double> > theta ( new std::vector<double> (size,-1000));
48 std::auto_ptr<std::vector<double> > phi ( new std::vector<double> (size,-1000));
49 std::auto_ptr<std::vector<double> > eta ( new std::vector<double> (size,-1000));
50 std::auto_ptr<std::vector<double> > dxy ( new std::vector<double> (size,-1000));
51 std::auto_ptr<std::vector<double> > dsz ( new std::vector<double> (size,-1000));
52 std::auto_ptr<std::vector<double> > qoverp ( new std::vector<double> (size,-1000));
53 std::auto_ptr<std::vector<double> > vx ( new std::vector<double> (size,-1000));
54 std::auto_ptr<std::vector<double> > vy ( new std::vector<double> (size,-1000));
55 std::auto_ptr<std::vector<double> > vz ( new std::vector<double> (size,-1000));
56
57 reco::RecoToSimCollection associations = associator->associateRecoToSim( tracks, trackingParticles, &event );
58
59 for( reco::RecoToSimCollection::const_iterator association = associations.begin();
60 association != associations.end(); association++) {
61
62 const reco::Track* track = association->key.get();
63 const int matches = association->val.size();
64 if(matches>0) {
65 const TrackingParticle* tparticle = association->val[0].first.get();
66 unsigned i = shallow::findTrackIndex(tracks, track);
67
68 multi->at(i) = matches;
69 type->at(i) = tparticle->pdgId();
70 charge->at(i)= tparticle->charge();
71 momentum->at(i)=tparticle->p() ;
72 pt->at(i) = tparticle->pt() ;
73 theta->at(i) = tparticle->theta() ;
74 phi->at(i) = tparticle->phi() ;
75 eta->at(i) = tparticle->eta() ;
76 qoverp->at(i)= tparticle->charge()/tparticle->p();
77
78 const TrackingVertex* tvertex = tparticle->parentVertex().get();
79 vx->at(i) = tvertex->position().x();
80 vy->at(i) = tvertex->position().y();
81 vz->at(i) = tvertex->position().z();
82 }
83 }
84
85 event.put( multi ,Prefix + "multi" + Suffix );
86 event.put( type ,Prefix + "type" + Suffix );
87 event.put( charge ,Prefix + "charge" + Suffix );
88 event.put( momentum ,Prefix + "momentum" + Suffix );
89 event.put( pt ,Prefix + "pt" + Suffix );
90 event.put( theta ,Prefix + "theta" + Suffix );
91 event.put( phi ,Prefix + "phi" + Suffix );
92 event.put( eta ,Prefix + "eta" + Suffix );
93 event.put( qoverp ,Prefix + "qoverp" + Suffix );
94 event.put( vx ,Prefix + "vx" + Suffix );
95 event.put( vy ,Prefix + "vy" + Suffix );
96 event.put( vz ,Prefix + "vz" + Suffix );
97
98 }