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 
34 
35 #include <string>
36 #include <vector>
37 
41 
42 public:
43  explicit ConversionTrackProducer(const edm::ParameterSet& conf);
44 
45  ~ConversionTrackProducer() override;
46 
47  void produce(edm::Event& e, const edm::EventSetup& c) override;
48 
49 private:
59 
60  //--------------------------------------------------
61  //Added by D. Giordano
62  // 2011/08/05
63  // Reduction of the track sample based on geometric hypothesis for conversion tracks
64 
68  double minConvRadius;
70 };
71 
74 
76  : useTrajectory(conf.getParameter<bool>("useTrajectory")),
77  setTrackerOnly(conf.getParameter<bool>("setTrackerOnly")),
78  setIsGsfTrackOpen(conf.getParameter<bool>("setIsGsfTrackOpen")),
79  setArbitratedEcalSeeded(conf.getParameter<bool>("setArbitratedEcalSeeded")),
80  setArbitratedMerged(conf.getParameter<bool>("setArbitratedMerged")),
81  setArbitratedMergedEcalGeneral(conf.getParameter<bool>("setArbitratedMergedEcalGeneral")),
82  beamSpotInputTag(consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpotInputTag"))),
83  filterOnConvTrackHyp(conf.getParameter<bool>("filterOnConvTrackHyp")),
84  minConvRadius(conf.getParameter<double>("minConvRadius")) {
85  edm::InputTag thetp(conf.getParameter<std::string>("TrackProducer"));
86  genericTracks = consumes<edm::View<reco::Track> >(thetp);
87  if (useTrajectory) {
88  kfTrajectories = consumes<TrajTrackAssociationCollection>(thetp);
89  gsfTrajectories = consumes<TrajGsfTrackAssociationCollection>(thetp);
90  }
92  produces<reco::ConversionTrackCollection>();
93 }
94 
95 // Virtual destructor needed.
97 
98 // Functions that gets called by framework every event
100  //get input collection (through edm::View)
102 
103  //get association maps between trajectories and tracks and build temporary maps
104  std::map<reco::TrackRef, edm::Ref<std::vector<Trajectory> > > tracktrajmap;
105  std::map<reco::GsfTrackRef, edm::Ref<std::vector<Trajectory> > > gsftracktrajmap;
106 
107  if (useTrajectory) {
108  if (!trks.empty()) {
109  if (dynamic_cast<const reco::GsfTrack*>(&trks.at(0))) {
110  //fill map for gsf tracks
111  for (auto const& pair : e.get(gsfTrajectories)) {
112  gsftracktrajmap[pair.val] = pair.key;
113  }
114  } else {
115  //fill map for standard tracks
116  for (auto const& pair : e.get(kfTrajectories)) {
117  tracktrajmap[pair.val] = pair.key;
118  }
119  }
120  }
121  }
122 
123  // Step B: create empty output collection
124  auto outputTrks = std::make_unique<reco::ConversionTrackCollection>();
125 
126  //--------------------------------------------------
127  //Added by D. Giordano
128  // 2011/08/05
129  // Reduction of the track sample based on geometric hypothesis for conversion tracks
130 
131  math::XYZVector beamSpot{e.get(beamSpotInputTag).position()};
132 
134 
135  //----------------------------------------------------------
136 
137  // Simple conversion of tracks to conversion tracks, setting appropriate flags from configuration
138  for (size_t i = 0; i < trks.size(); ++i) {
139  //--------------------------------------------------
140  //Added by D. Giordano
141  // 2011/08/05
142  // Reduction of the track sample based on geometric hypothesis for conversion tracks
143 
144  edm::RefToBase<reco::Track> trackBaseRef = trks.refAt(i);
145  if (filterOnConvTrackHyp &&
147  continue;
148  //--------------------------------------------------
149 
150  reco::ConversionTrack convTrack(trackBaseRef);
151  convTrack.setIsTrackerOnly(setTrackerOnly);
156 
157  //fill trajectory association if configured, using correct map depending on track type
158  if (useTrajectory) {
159  if (!gsftracktrajmap.empty()) {
160  convTrack.setTrajRef(gsftracktrajmap.find(trackBaseRef.castTo<reco::GsfTrackRef>())->second);
161  } else {
162  convTrack.setTrajRef(tracktrajmap.find(trackBaseRef.castTo<reco::TrackRef>())->second);
163  }
164  }
165 
166  outputTrks->push_back(convTrack);
167  }
168 
169  e.put(std::move(outputTrks));
170  return;
171 
172 } //end produce
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void setIsGsfTrackOpen(bool b)
void setIsArbitratedMergedEcalGeneral(bool b)
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
void setMagnField(const MagneticField *magnField)
edm::EDGetTokenT< edm::View< reco::Track > > genericTracks
edm::EDGetTokenT< TrajTrackAssociationCollection > kfTrajectories
ConversionTrackProducer(const edm::ParameterSet &conf)
REF castTo() const
Definition: RefToBase.h:243
edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::GsfTrackCollection, unsigned short > > TrajGsfTrackAssociationCollection
void setTrajRef(edm::Ref< std::vector< Trajectory > > tr)
void setIsArbitratedEcalSeeded(bool b)
std::vector< GsfTrack > GsfTrackCollection
collection of GsfTracks
Definition: GsfTrackFwd.h:9
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
IdealHelixParameters ConvTrackPreSelector
void setIsTrackerOnly(bool b)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
edm::EDGetTokenT< reco::BeamSpot > beamSpotInputTag
fixed size matrix
HLT enums.
edm::EDGetTokenT< TrajGsfTrackAssociationCollection > gsfTrajectories
bool isTangentPointDistanceLessThan(float rmax, const reco::Track *track, const math::XYZVector &refPoint)
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magFieldToken
void produce(edm::Event &e, const edm::EventSetup &c) override
value_type const * get() const
Definition: RefToBase.h:216
def move(src, dest)
Definition: eostools.py:511
void setIsArbitratedMerged(bool b)