CMS 3D CMS Logo

ConversionTrackProducer.cc
Go to the documentation of this file.
1 //
2 // Package: RecoTracker/FinalTrackSelectors
3 // Class: ConversionTrackProducer
4 //
5 // Description: Trivial producer of ConversionTrack collection from an edm::View of a track collection
6 // (ConversionTrack is a simple wrappper class containing a TrackBaseRef and some additional flags)
7 //
8 // Original Author: J.Bendavid
9 //
10 //
11 
12 #include <memory>
13 #include <string>
14 #include <iostream>
15 #include <cmath>
16 #include <vector>
17 
19 
22 
25 
27 
29 
31 
34 
35 
37  conf_(conf),
38  trackProducer ( conf.getParameter<std::string>("TrackProducer") ),
39  useTrajectory ( conf.getParameter<bool>("useTrajectory") ),
40  setTrackerOnly ( conf.getParameter<bool>("setTrackerOnly") ),
41  setIsGsfTrackOpen ( conf.getParameter<bool>("setIsGsfTrackOpen") ),
42  setArbitratedEcalSeeded ( conf.getParameter<bool>("setArbitratedEcalSeeded") ),
43  setArbitratedMerged ( conf.getParameter<bool>("setArbitratedMerged") ),
44  setArbitratedMergedEcalGeneral ( conf.getParameter<bool>("setArbitratedMergedEcalGeneral") ),
45  beamSpotInputTag ( consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpotInputTag")) ),
46  filterOnConvTrackHyp( conf.getParameter<bool>("filterOnConvTrackHyp") ),
47  minConvRadius( conf.getParameter<double>("minConvRadius") )
48 {
51  consumes<edm::View<reco::Track> >(thetp);
53  consumes<TrajTrackAssociationCollection>(thetp);
55  consumes<TrajGsfTrackAssociationCollection>(thetp);
56  produces<reco::ConversionTrackCollection>();
57 
58 }
59 
60 
61  // Virtual destructor needed.
63 
64  // Functions that gets called by framework every event
66  {
67  //get input collection (through edm::View)
69  e.getByToken(genericTracks, hTrks);
70 
71  //get association maps between trajectories and tracks and build temporary maps
74 
75  std::map<reco::TrackRef,edm::Ref<std::vector<Trajectory> > > tracktrajmap;
76  std::map<reco::GsfTrackRef,edm::Ref<std::vector<Trajectory> > > gsftracktrajmap;
77 
78  if (useTrajectory) {
79  if (!hTrks->empty()) {
80  if (dynamic_cast<const reco::GsfTrack*>(&hTrks->at(0))) {
81  //fill map for gsf tracks
82  e.getByToken(gsfTrajectories, hTTAssGsf);
84  hTTAssGsf->begin();
85  iPair != hTTAssGsf->end(); ++iPair) {
86 
87  gsftracktrajmap[iPair->val] = iPair->key;
88 
89  }
90 
91  }
92  else {
93  //fill map for standard tracks
94  e.getByToken(kfTrajectories, hTTAss);
96  iPair != hTTAss->end();
97  ++iPair) {
98 
99  tracktrajmap[iPair->val] = iPair->key;
100 
101  }
102  }
103  }
104  }
105 
106  // Step B: create empty output collection
107  outputTrks = std::make_unique<reco::ConversionTrackCollection>();
108 
109  //--------------------------------------------------
110  //Added by D. Giordano
111  // 2011/08/05
112  // Reduction of the track sample based on geometric hypothesis for conversion tracks
113 
114  edm::Handle<reco::BeamSpot> beamSpotHandle;
115  e.getByToken(beamSpotInputTag,beamSpotHandle);
116 
117  edm::ESHandle<MagneticField> magFieldHandle;
118  es.get<IdealMagneticFieldRecord>().get( magFieldHandle );
119 
120 
121  if(filterOnConvTrackHyp && !beamSpotHandle.isValid()) {
122  edm::LogError("Invalid Collection")
123  << "invalid collection for the BeamSpot";
124  throw;
125  }
126 
127  ConvTrackPreSelector.setMagnField(magFieldHandle.product());
128 
129  //----------------------------------------------------------
130 
131 
132  // Simple conversion of tracks to conversion tracks, setting appropriate flags from configuration
133  for (size_t i = 0; i < hTrks->size(); ++i) {
134 
135  //--------------------------------------------------
136  //Added by D. Giordano
137  // 2011/08/05
138  // Reduction of the track sample based on geometric hypothesis for conversion tracks
139 
140  math::XYZVector beamSpot= math::XYZVector(beamSpotHandle->position());
141  edm::RefToBase<reco::Track> trackBaseRef = hTrks->refAt(i);
143  continue;
144  //--------------------------------------------------
145 
146  reco::ConversionTrack convTrack(trackBaseRef);
147  convTrack.setIsTrackerOnly(setTrackerOnly);
152 
153  //fill trajectory association if configured, using correct map depending on track type
154  if (useTrajectory) {
155  if (!gsftracktrajmap.empty()) {
156  convTrack.setTrajRef(gsftracktrajmap.find(trackBaseRef.castTo<reco::GsfTrackRef>())->second);
157  }
158  else {
159  convTrack.setTrajRef(tracktrajmap.find(trackBaseRef.castTo<reco::TrackRef>())->second);
160  }
161  }
162 
163  outputTrks->push_back(convTrack);
164  }
165 
167  return;
168 
169  }//end produce
170 
void setIsGsfTrackOpen(bool b)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
void setIsArbitratedMergedEcalGeneral(bool b)
const_iterator end() const
last iterator over the map (read only)
void setMagnField(const MagneticField *magnField)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< edm::View< reco::Track > > genericTracks
edm::EDGetTokenT< TrajTrackAssociationCollection > kfTrajectories
ConversionTrackProducer(const edm::ParameterSet &conf)
void setTrajRef(edm::Ref< std::vector< Trajectory > > tr)
void setIsArbitratedEcalSeeded(bool b)
std::unique_ptr< reco::ConversionTrackCollection > outputTrks
IdealHelixParameters ConvTrackPreSelector
void setIsTrackerOnly(bool b)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
edm::EDGetTokenT< reco::BeamSpot > beamSpotInputTag
fixed size matrix
HLT enums.
edm::EDGetTokenT< TrajGsfTrackAssociationCollection > gsfTrajectories
T get() const
Definition: EventSetup.h:71
bool isTangentPointDistanceLessThan(float rmax, const reco::Track *track, const math::XYZVector &refPoint)
void produce(edm::Event &e, const edm::EventSetup &c) override
const_iterator begin() const
first iterator over the map (read only)
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511
void setIsArbitratedMerged(bool b)