CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFDisplacedTrackerVertexProducer.cc
Go to the documentation of this file.
1 #include <memory>
8 using namespace std;
9 using namespace edm;
11  pfTransformer_(0)
12 {
13  produces<reco::PFRecTrackCollection>();
14  produces<reco::PFDisplacedTrackerVertexCollection>();
15 
16  pfDisplacedVertexContainer_ = consumes<reco::PFDisplacedVertexCollection>( iConfig.getParameter< InputTag >("displacedTrackerVertexColl"));
17 
18 
19  pfTrackContainer_ =consumes<reco::TrackCollection>(
20  iConfig.getParameter< InputTag >("trackColl"));
21 
22 }
23 
25 {
26  delete pfTransformer_;
27 }
28 
29 void
31 {
32 
33  //create the empty collections
34  auto_ptr< reco::PFDisplacedTrackerVertexCollection >
35  pfDisplacedTrackerVertexColl (new reco::PFDisplacedTrackerVertexCollection);
36  auto_ptr< reco::PFRecTrackCollection >
37  pfRecTrackColl (new reco::PFRecTrackCollection);
38 
40 
41 
42 
44  iEvent.getByToken(pfDisplacedVertexContainer_, nuclCollH);
45  const reco::PFDisplacedVertexCollection& nuclColl = *(nuclCollH.product());
46 
48  iEvent.getByToken(pfTrackContainer_, trackColl);
49 
50  int idx = 0;
51 
52  // cout << "Size of Displaced Vertices "
53  // << nuclColl.size() << endl;
54 
55  // loop on all NuclearInteraction
56  for( unsigned int icoll=0; icoll < nuclColl.size(); icoll++) {
57 
58  reco::PFRecTrackRefVector pfRecTkcoll;
59 
60  std::vector<reco::Track> refittedTracks = nuclColl[icoll].refittedTracks();
61 
62  // convert the secondary tracks
63  for(unsigned it = 0; it < refittedTracks.size(); it++){
64 
65  reco::TrackBaseRef trackBaseRef = nuclColl[icoll].originalTrack(refittedTracks[it]);
66 
67  // cout << "track base pt = " << trackBaseRef->pt() << endl;
68 
69  reco::TrackRef trackRef(trackColl, trackBaseRef.key());
70 
71  // cout << "track pt = " << trackRef->pt() << endl;
72 
73 
74  reco::PFRecTrack pfRecTrack( trackBaseRef->charge(),
76  trackBaseRef.key(),
77  trackRef );
78 
79  // cout << pfRecTrack << endl;
80 
81  Trajectory FakeTraj;
82  bool valid = pfTransformer_->addPoints( pfRecTrack, *trackBaseRef, FakeTraj);
83  if(valid) {
84  pfRecTkcoll.push_back(reco::PFRecTrackRef( pfTrackRefProd, idx++));
85  pfRecTrackColl->push_back(pfRecTrack);
86  // cout << "after "<< pfRecTrack << endl;
87 
88  }
89  }
90  reco::PFDisplacedVertexRef niRef(nuclCollH, icoll);
91  pfDisplacedTrackerVertexColl->push_back( reco::PFDisplacedTrackerVertex( niRef, pfRecTkcoll ));
92  }
93 
94  iEvent.put(pfRecTrackColl);
95  iEvent.put(pfDisplacedTrackerVertexColl);
96 }
97 
98 // ------------ method called once each job just before starting event loop ------------
99 void
101  const EventSetup& iSetup)
102 {
104  iSetup.get<IdealMagneticFieldRecord>().get(magneticField);
105  pfTransformer_= new PFTrackTransformer(math::XYZVector(magneticField->inTesla(GlobalPoint(0,0,0))));
107 }
108 
109 // ------------ method called once each job just after ending the event loop ------------
110 void
112  const EventSetup& iSetup) {
113  delete pfTransformer_;
114  pfTransformer_=nullptr;
115 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::TrackCollection > pfTrackContainer_
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:22
std::vector< PFDisplacedVertex > PFDisplacedVertexCollection
collection of PFDisplacedVertex objects
virtual void produce(edm::Event &, const edm::EventSetup &) override
Produce the PFRecTrack collection.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
PFDisplacedTrackerVertexProducer(const edm::ParameterSet &)
Constructor.
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.
edm::EDGetTokenT< reco::PFDisplacedVertexCollection > pfDisplacedVertexContainer_
virtual void beginRun(const edm::Run &, const edm::EventSetup &) override
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
size_t key() const
Definition: RefToBase.h:241
virtual void endRun(const edm::Run &, const edm::EventSetup &) override
PFTrackTransformer * pfTransformer_
PFTrackTransformer.
RefProd< PROD > getRefBeforePut()
Definition: Event.h:140
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
const T & get() const
Definition: EventSetup.h:56
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:62
int charge() const
track electric charge
Definition: TrackBase.h:554
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
std::vector< PFDisplacedTrackerVertex > PFDisplacedTrackerVertexCollection
collection of DisplacedTrackerVertexs
Definition: Run.h:43
int icoll
Definition: AMPTWrapper.h:136