Imperial Analysis
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ICVertexProducer.cc
Go to the documentation of this file.
2 #include <memory>
3 #include <string>
4 #include <vector>
5 #include "FWCore/Framework/interface/Event.h"
6 #include "FWCore/Framework/interface/EventSetup.h"
7 #include "FWCore/Framework/interface/MakerMacros.h"
8 #include "FWCore/ParameterSet/interface/ParameterSet.h"
9 #include "FWCore/Utilities/interface/InputTag.h"
10 #include "DataFormats/Common/interface/Handle.h"
11 #include "DataFormats/Common/interface/View.h"
12 #include "DataFormats/VertexReco/interface/Vertex.h"
13 #include "DataFormats/VertexReco/interface/VertexFwd.h"
14 #include "DataFormats/TrackReco/interface/Track.h"
15 #include "DataFormats/TrackReco/interface/TrackFwd.h"
21 
22 ICVertexProducer::ICVertexProducer(const edm::ParameterSet& config)
23  : input_(config.getParameter<edm::InputTag>("input")),
24  branch_(config.getParameter<std::string>("branch")),
25  first_only_(config.getParameter<bool>("firstVertexOnly")),
26  track_pt_threshold_(config.getParameter<double>("trackPtThreshold")),
27  request_trks_(config.getParameter<bool>("requestTracks")) {
28  consumes<edm::View<reco::Vertex>>(input_);
29  vertices_ = new std::vector<ic::Vertex>();
30  if (request_trks_) {
31  produces<reco::TrackRefVector>("requestedTracks");
32  }
33  PrintHeaderWithProduces(config, input_, branch_);
34  PrintOptional(1, first_only_, "firstVertexOnly");
35  PrintOptional(1, request_trks_, "requestTracks");
36 }
37 
38 ICVertexProducer::~ICVertexProducer() { delete vertices_; }
39 
40 void ICVertexProducer::produce(edm::Event& event,
41  const edm::EventSetup& setup) {
42  edm::Handle<edm::View<reco::Vertex> > vtxs_handle;
43  event.getByLabel(input_, vtxs_handle);
44 
45  std::auto_ptr<reco::TrackRefVector> trk_requests(new reco::TrackRefVector());
46 
47  vertices_->clear();
48  if (!first_only_) {
49  vertices_->resize(vtxs_handle->size());
50  } else {
51  vertices_->resize(vtxs_handle->size() > 0 ? 1 : 0);
52  }
53 
54  for (unsigned i = 0; i < vtxs_handle->size(); ++i) {
55  reco::Vertex const& src = vtxs_handle->at(i);
56  ic::Vertex & dest = vertices_->at(i);
57  dest.set_id(vertex_hasher_(&src));
58  dest.set_vx(src.x());
59  dest.set_vy(src.y());
60  dest.set_vz(src.z());
61  dest.set_chi2(src.chi2());
62  dest.set_ndof(src.ndof());
63  if (request_trks_) {
64  for (reco::Vertex::trackRef_iterator trk_iter = src.tracks_begin();
65  trk_iter != src.tracks_end(); ++trk_iter) {
66  float weight = src.trackWeight(*trk_iter);
67  reco::TrackRef trk_ref = trk_iter->castTo<reco::TrackRef>();
68  if (trk_ref->pt() < track_pt_threshold_) continue;
69  dest.AddTrack(track_hasher_(&(*trk_ref)), weight);
70  trk_requests->push_back(trk_ref);
71  }
72  }
73  if (first_only_) break;
74  }
75  if (request_trks_) event.put(trk_requests, "requestedTracks");
76 }
77 
78 void ICVertexProducer::beginJob() {
79  ic::StaticTree::tree_->Branch(branch_.c_str(), &vertices_);
80 }
81 
82 void ICVertexProducer::endJob() {}
83 
84 // define this as a plug-in
void set_vz(float const &z)
Vertex z-coordinate.
Definition: Vertex.hh:71
void set_vx(float const &x)
Vertex x-coordinate.
Definition: Vertex.hh:65
ICVertexProducer(const edm::ParameterSet &)
void set_ndof(float const &ndof)
The number-of-degrees-of-freedom in the vertex fit.
Definition: Vertex.hh:77
static TTree * tree_
Definition: StaticTree.hh:13
See documentation here.
void set_id(std::size_t const &id)
Unique identifier.
Definition: Vertex.hh:62
void set_chi2(float const &chi2)
The of the vertex fit.
Definition: Vertex.hh:74
void AddTrack(std::size_t id, float weight)
Add an ic::Track::id with a vertex fit weight.
Definition: Vertex.hh:86
void PrintOptional(unsigned depth, bool value, std::string text)
Definition: Consumes.h:19
void set_vy(float const &y)
Vertex y-coordinate.
Definition: Vertex.hh:68
DEFINE_FWK_MODULE(ICVertexProducer)
void PrintHeaderWithProduces(edm::ParameterSet const &config, edm::InputTag const &in, std::string branch)
Stores information about the position of an event vertex and the quality of the track fit...
Definition: Vertex.hh:15