CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LightPFTrackProducer.cc
Go to the documentation of this file.
1 #include <memory>
11 
12 using namespace std;
13 using namespace edm;
15  pfTransformer_(0)
16 {
17  produces<reco::PFRecTrackCollection>();
18 
19 
20 
21  std::vector<InputTag> tags =
22  iConfig.getParameter< vector < InputTag > >("TkColList");
23 
24  for (unsigned int i=0;i<tags.size();++i)
25  tracksContainers_.push_back(consumes<reco::TrackCollection>(tags[i]));
26 
27  useQuality_ = iConfig.getParameter<bool>("UseQuality");
29 
30 }
31 
33 {
34  delete pfTransformer_;
35 }
36 
37 void
39 {
40 
41  //create the empty collections
42  auto_ptr< reco::PFRecTrackCollection >
43  PfTrColl (new reco::PFRecTrackCollection);
44 
45  for (unsigned int istr=0; istr<tracksContainers_.size();istr++){
46 
47  //Track collection
48  Handle<reco::TrackCollection> tkRefCollection;
49  iEvent.getByToken(tracksContainers_[istr], tkRefCollection);
50  reco::TrackCollection Tk=*(tkRefCollection.product());
51  for(unsigned int i=0;i<Tk.size();i++){
52  if (useQuality_ &&
53  (!(Tk[i].quality(trackQuality_)))) continue;
54  reco::TrackRef trackRef(tkRefCollection, i);
55  reco::PFRecTrack pftrack( trackRef->charge(),
57  i, trackRef );
58  Trajectory FakeTraj;
59  bool mymsgwarning = false;
60  bool valid = pfTransformer_->addPoints( pftrack, *trackRef, FakeTraj, mymsgwarning);
61  if(valid)
62  PfTrColl->push_back(pftrack);
63 
64  }
65  }
66  iEvent.put(PfTrColl);
67 }
68 
69 // ------------ method called once each job just before starting event loop ------------
70 void
72  const EventSetup& iSetup)
73 {
75  iSetup.get<IdealMagneticFieldRecord>().get(magneticField);
76  pfTransformer_= new PFTrackTransformer(math::XYZVector(magneticField->inTesla(GlobalPoint(0,0,0))));
78 }
79 
80 // ------------ method called once each job just after ending the event loop ------------
81 void
83  const EventSetup& iSetup) {
84  delete pfTransformer_;
85  pfTransformer_=nullptr;
86 }
T getParameter(std::string const &) const
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:22
int i
Definition: DBlmapReader.cc:9
LightPFTrackProducer(const edm::ParameterSet &)
Constructor.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
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.
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
~LightPFTrackProducer()
Destructor.
tuple magneticField
virtual void produce(edm::Event &, const edm::EventSetup &) override
Produce the PFRecTrack collection.
PFTrackTransformer * pfTransformer_
PFTrackTransformer.
int iEvent
Definition: GenABIO.cc:230
bool useQuality_
TRACK QUALITY.
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
tuple tags
Definition: o2o.py:248
virtual void endRun(const edm::Run &, const edm::EventSetup &) override
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
T const * product() const
Definition: Handle.h:81
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
const T & get() const
Definition: EventSetup.h:56
virtual void beginRun(const edm::Run &, const edm::EventSetup &) override
std::vector< edm::EDGetTokenT< reco::TrackCollection > > tracksContainers_
reco::TrackBase::TrackQuality trackQuality_
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
Definition: Run.h:43