CMS 3D CMS Logo

LightPFTrackProducer.cc
Go to the documentation of this file.
11 
13 public:
15  explicit LightPFTrackProducer(const edm::ParameterSet&);
16 
18  ~LightPFTrackProducer() override;
19 
20 private:
21  void beginRun(const edm::Run&, const edm::EventSetup&) override;
22  void endRun(const edm::Run&, const edm::EventSetup&) override;
23 
25  void produce(edm::Event&, const edm::EventSetup&) override;
26 
29  std::vector<edm::EDGetTokenT<reco::TrackCollection> > tracksContainers_;
30 
35 };
36 
39 
40 using namespace std;
41 using namespace edm;
43  : pfTransformer_(nullptr), magneticFieldToken_(esConsumes<edm::Transition::BeginRun>()) {
44  produces<reco::PFRecTrackCollection>();
45 
46  std::vector<InputTag> tags = iConfig.getParameter<vector<InputTag> >("TkColList");
47 
48  for (unsigned int i = 0; i < tags.size(); ++i)
49  tracksContainers_.push_back(consumes<reco::TrackCollection>(tags[i]));
50 
51  useQuality_ = iConfig.getParameter<bool>("UseQuality");
53 }
54 
56 
58  //create the empty collections
59  auto PfTrColl = std::make_unique<reco::PFRecTrackCollection>();
60 
61  for (unsigned int istr = 0; istr < tracksContainers_.size(); istr++) {
62  //Track collection
63  Handle<reco::TrackCollection> tkRefCollection;
64  iEvent.getByToken(tracksContainers_[istr], tkRefCollection);
65  reco::TrackCollection Tk = *(tkRefCollection.product());
66  for (unsigned int i = 0; i < Tk.size(); i++) {
67  if (useQuality_ && (!(Tk[i].quality(trackQuality_))))
68  continue;
69  reco::TrackRef trackRef(tkRefCollection, i);
70  reco::PFRecTrack pftrack(trackRef->charge(), reco::PFRecTrack::KF, i, trackRef);
71  Trajectory FakeTraj;
72  bool mymsgwarning = false;
73  bool valid = pfTransformer_->addPoints(pftrack, *trackRef, FakeTraj, mymsgwarning);
74  if (valid)
75  PfTrColl->push_back(pftrack);
76  }
77  }
78  iEvent.put(std::move(PfTrColl));
79 }
80 
81 // ------------ method called once each job just before starting event loop ------------
83  auto const& magneticField = &iSetup.getData(magneticFieldToken_);
86 }
87 
88 // ------------ method called once each job just after ending the event loop ------------
89 void LightPFTrackProducer::endRun(const edm::Run& run, const EventSetup& iSetup) {
90  delete pfTransformer_;
91  pfTransformer_ = nullptr;
92 }
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:20
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
LightPFTrackProducer(const edm::ParameterSet &)
Constructor.
TrackQuality
track quality
Definition: TrackBase.h:150
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T const * product() const
Definition: Handle.h:70
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
~LightPFTrackProducer() override
Destructor.
void produce(edm::Event &, const edm::EventSetup &) override
Produce the PFRecTrack collection.
PFTrackTransformer * pfTransformer_
PFTrackTransformer.
int iEvent
Definition: GenABIO.cc:224
bool useQuality_
TRACK QUALITY.
Transition
Definition: Transition.h:12
bool getData(T &iHolder) const
Definition: EventSetup.h:122
void endRun(const edm::Run &, const edm::EventSetup &) override
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:126
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
bool addPoints(reco::PFRecTrack &pftrack, const reco::Track &track, const Trajectory &traj, bool msgwarning=true) const
Add points to a PFTrack. return false if a TSOS is invalid.
void beginRun(const edm::Run &, const edm::EventSetup &) override
HLT enums.
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
std::vector< edm::EDGetTokenT< reco::TrackCollection > > tracksContainers_
string quality
reco::TrackBase::TrackQuality trackQuality_
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45