CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ShallowSimTracksProducer.cc
Go to the documentation of this file.
3 
7 
9 
11  : Prefix(conf.getParameter<std::string>("Prefix")),
12  Suffix(conf.getParameter<std::string>("Suffix")),
13  trackingParticles_token_(
14  consumes<TrackingParticleCollection>(conf.getParameter<edm::InputTag>("TrackingParticles"))),
15  associator_token_(
16  consumes<reco::TrackToTrackingParticleAssociator>(conf.getParameter<edm::InputTag>("Associator"))),
17  tracks_token_(consumes<edm::View<reco::Track>>(conf.getParameter<edm::InputTag>("Tracks"))) {
18  produces<std::vector<unsigned>>(Prefix + "multi" + Suffix);
19  produces<std::vector<int>>(Prefix + "type" + Suffix);
20  produces<std::vector<float>>(Prefix + "charge" + Suffix);
21  produces<std::vector<float>>(Prefix + "momentum" + Suffix);
22  produces<std::vector<float>>(Prefix + "pt" + Suffix);
23  produces<std::vector<double>>(Prefix + "theta" + Suffix);
24  produces<std::vector<double>>(Prefix + "phi" + Suffix);
25  produces<std::vector<double>>(Prefix + "eta" + Suffix);
26  produces<std::vector<double>>(Prefix + "qoverp" + Suffix);
27  produces<std::vector<double>>(Prefix + "vx" + Suffix);
28  produces<std::vector<double>>(Prefix + "vy" + Suffix);
29  produces<std::vector<double>>(Prefix + "vz" + Suffix);
30 }
31 
34  event.getByToken(tracks_token_, tracks);
36  event.getByToken(trackingParticles_token_, trackingParticles);
38  event.getByToken(associator_token_, associator);
39 
40  unsigned size = tracks->size();
41  auto multi = std::make_unique<std::vector<unsigned>>(size, 0);
42  auto type = std::make_unique<std::vector<int>>(size, 0);
43  auto charge = std::make_unique<std::vector<float>>(size, 0);
44  auto momentum = std::make_unique<std::vector<float>>(size, -1);
45  auto pt = std::make_unique<std::vector<float>>(size, -1);
46  auto theta = std::make_unique<std::vector<double>>(size, -1000);
47  auto phi = std::make_unique<std::vector<double>>(size, -1000);
48  auto eta = std::make_unique<std::vector<double>>(size, -1000);
49  auto dxy = std::make_unique<std::vector<double>>(size, -1000);
50  auto dsz = std::make_unique<std::vector<double>>(size, -1000);
51  auto qoverp = std::make_unique<std::vector<double>>(size, -1000);
52  auto vx = std::make_unique<std::vector<double>>(size, -1000);
53  auto vy = std::make_unique<std::vector<double>>(size, -1000);
54  auto vz = std::make_unique<std::vector<double>>(size, -1000);
55 
56  reco::RecoToSimCollection associations = associator->associateRecoToSim(tracks, trackingParticles);
57 
58  for (reco::RecoToSimCollection::const_iterator association = associations.begin(); association != associations.end();
59  association++) {
60  const reco::Track* track = association->key.get();
61  const int matches = association->val.size();
62  if (matches > 0) {
63  const TrackingParticle* tparticle = association->val[0].first.get();
64  unsigned i = shallow::findTrackIndex(tracks, track);
65 
66  multi->at(i) = matches;
67  type->at(i) = tparticle->pdgId();
68  charge->at(i) = tparticle->charge();
69  momentum->at(i) = tparticle->p();
70  pt->at(i) = tparticle->pt();
71  theta->at(i) = tparticle->theta();
72  phi->at(i) = tparticle->phi();
73  eta->at(i) = tparticle->eta();
74  qoverp->at(i) = tparticle->charge() / tparticle->p();
75 
76  const TrackingVertex* tvertex = tparticle->parentVertex().get();
77  vx->at(i) = tvertex->position().x();
78  vy->at(i) = tvertex->position().y();
79  vz->at(i) = tvertex->position().z();
80  }
81  }
82 
83  event.put(std::move(multi), Prefix + "multi" + Suffix);
84  event.put(std::move(type), Prefix + "type" + Suffix);
85  event.put(std::move(charge), Prefix + "charge" + Suffix);
86  event.put(std::move(momentum), Prefix + "momentum" + Suffix);
87  event.put(std::move(pt), Prefix + "pt" + Suffix);
88  event.put(std::move(theta), Prefix + "theta" + Suffix);
89  event.put(std::move(phi), Prefix + "phi" + Suffix);
90  event.put(std::move(eta), Prefix + "eta" + Suffix);
91  event.put(std::move(qoverp), Prefix + "qoverp" + Suffix);
92  event.put(std::move(vx), Prefix + "vx" + Suffix);
93  event.put(std::move(vy), Prefix + "vy" + Suffix);
94  event.put(std::move(vz), Prefix + "vz" + Suffix);
95 }
const edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > associator_token_
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
int findTrackIndex(const edm::Handle< edm::View< reco::Track > > &h, const reco::Track *t)
Definition: ShallowTools.cc:25
int pdgId() const
PDG ID.
ShallowSimTracksProducer(const edm::ParameterSet &)
Geom::Theta< T > theta() const
double pt() const
Transverse momentum. Note this is taken from the first SimTrack only.
auto const & tracks
cannot be loose
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
std::tuple< layerClusterToCaloParticle, caloParticleToLayerCluster > association
def move
Definition: eostools.py:511
double p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
const TrackingVertexRef & parentVertex() const
double eta() const
Momentum pseudorapidity. Note this is taken from the first SimTrack only.
Monte Carlo truth information used for tracking validation.
std::vector< TrackingParticle > TrackingParticleCollection
const edm::EDGetTokenT< edm::View< reco::Track > > tracks_token_
const edm::EDGetTokenT< TrackingParticleCollection > trackingParticles_token_
tuple size
Write out results.
double phi() const
Momentum azimuthal angle. Note this is taken from the first SimTrack only.
double theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.
const LorentzVector & position() const