CMS 3D CMS Logo

SeedGeneratorFromProtoTracksEDProducer.cc
Go to the documentation of this file.
2 
5 
10 
16 
21 
22 #include <vector>
23 
24 using namespace edm;
25 using namespace reco;
26 
27 template <class T>
28 T sqr(T t) {
29  return t * t;
30 }
32 
34  bool operator()(const Hit& h1, const Hit& h2) { return h1->globalPosition().perp2() < h2->globalPosition().perp2(); }
35 };
36 
39  desc.add<InputTag>("InputCollection", InputTag("pixelTracks"));
40  desc.add<InputTag>("InputVertexCollection", InputTag(""));
41  desc.add<double>("originHalfLength", 1E9);
42  desc.add<double>("originRadius", 1E9);
43  desc.add<bool>("useProtoTrackKinematics", false);
44  desc.add<bool>("useEventsWithNoVertex", true);
45  desc.add<std::string>("TTRHBuilder", "TTRHBuilderWithoutAngle4PixelTriplets");
46  desc.add<bool>("usePV", false);
47  desc.add<bool>("includeFourthHit", false);
48  desc.add<bool>("produceComplement", false);
49 
51  psd0.add<std::string>("ComponentName", std::string("SeedFromConsecutiveHitsCreator"));
52  psd0.add<std::string>("propagator", std::string("PropagatorWithMaterial"));
53  psd0.add<double>("SeedMomentumForBOFF", 5.0);
54  psd0.add<double>("OriginTransverseErrorMultiplier", 1.0);
55  psd0.add<double>("MinOneOverPtError", 1.0);
56  psd0.add<std::string>("magneticField", std::string(""));
57  psd0.add<std::string>("TTRHBuilder", std::string("WithTrackAngle"));
58  psd0.add<bool>("forceKinematicWithRegionDirection", false);
59  desc.add<edm::ParameterSetDescription>("SeedCreatorPSet", psd0);
60 
61  descriptions.add("SeedGeneratorFromProtoTracksEDProducer", desc);
62 }
63 
65  : originHalfLength(cfg.getParameter<double>("originHalfLength")),
66  originRadius(cfg.getParameter<double>("originRadius")),
67  useProtoTrackKinematics(cfg.getParameter<bool>("useProtoTrackKinematics")),
68  useEventsWithNoVertex(cfg.getParameter<bool>("useEventsWithNoVertex")),
69  usePV_(cfg.getParameter<bool>("usePV")),
70  includeFourthHit_(cfg.getParameter<bool>("includeFourthHit")),
71  produceComplement_(cfg.getParameter<bool>("produceComplement")),
72  theInputCollectionTag(consumes<reco::TrackCollection>(cfg.getParameter<InputTag>("InputCollection"))),
73  theInputVertexCollectionTag(
74  consumes<reco::VertexCollection>(cfg.getParameter<InputTag>("InputVertexCollection"))),
75  seedCreator_(cfg.getParameter<edm::ParameterSet>("SeedCreatorPSet"), consumesCollector()),
76  config_(consumesCollector()) {
77  produces<TrajectorySeedCollection>();
78  if (produceComplement_) {
79  produces<reco::TrackCollection>();
80  }
81 }
82 
84  auto result = std::make_unique<TrajectorySeedCollection>();
85  auto leftTracks = std::make_unique<reco::TrackCollection>();
87  ev.getByToken(theInputCollectionTag, trks);
88 
89  const TrackCollection& protos = *(trks.product());
90 
92  bool foundVertices = ev.getByToken(theInputVertexCollectionTag, vertices);
93  //const reco::VertexCollection & vertices = *(h_vertices.product());
94 
98  for (TrackCollection::const_iterator it = protos.begin(); it != protos.end(); ++it) {
99  const Track& proto = (*it);
100  GlobalPoint vtx(proto.vertex().x(), proto.vertex().y(), proto.vertex().z());
101 
102  // check the compatibility with a primary vertex
103  bool keepTrack = false;
104  if ((!foundVertices) || vertices->empty()) {
106  keepTrack = true;
107  } else if (usePV_) {
108  GlobalPoint aPV(
109  vertices->begin()->position().x(), vertices->begin()->position().y(), vertices->begin()->position().z());
110  double distR2 = sqr(vtx.x() - aPV.x()) + sqr(vtx.y() - aPV.y());
111  double distZ = fabs(vtx.z() - aPV.z());
112  if (distR2 < sqr(originRadius) && distZ < originHalfLength) {
113  keepTrack = true;
114  }
115  } else {
116  for (reco::VertexCollection::const_iterator iv = vertices->begin(); iv != vertices->end(); ++iv) {
117  GlobalPoint aPV(iv->position().x(), iv->position().y(), iv->position().z());
118  double distR2 = sqr(vtx.x() - aPV.x()) + sqr(vtx.y() - aPV.y());
119  double distZ = fabs(vtx.z() - aPV.z());
120  if (distR2 < sqr(originRadius) && distZ < originHalfLength) {
121  keepTrack = true;
122  break;
123  }
124  }
125  }
126  if (produceComplement_ and !keepTrack)
127  (*leftTracks).push_back(proto);
128  if (!keepTrack)
129  continue;
130 
132  SeedFromProtoTrack seedFromProtoTrack(config_, proto, es);
133  if (seedFromProtoTrack.isValid())
134  (*result).push_back(seedFromProtoTrack.trajectorySeed());
135  } else {
136  std::vector<Hit> hits;
137  for (unsigned int iHit = 0, nHits = proto.recHitsSize(); iHit < nHits; ++iHit) {
138  TrackingRecHitRef refHit = proto.recHit(iHit);
139  if (refHit->isValid())
140  hits.push_back((Hit) & (*refHit));
141  }
142  sort(hits.begin(), hits.end(), HitLessByRadius());
143 
144  if (hits.size() > 1) {
145  double mom_perp =
146  sqrt(proto.momentum().x() * proto.momentum().x() + proto.momentum().y() * proto.momentum().y());
147  GlobalTrackingRegion region(mom_perp, vtx, 0.2, 0.2);
148 
149  seedCreator_.init(region, es, nullptr);
150  if (hits.size() > 3 and not includeFourthHit_)
151  seedCreator_.makeSeed(*result, {hits[0], hits[1], hits[2]});
152  else
153  seedCreator_.makeSeed(*result, hits);
154  }
155  }
156  }
157 
158  ev.put(std::move(result));
159  if (produceComplement_)
160  ev.put(std::move(leftTracks));
161 }
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:97
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
TrajectorySeed trajectorySeed() const
const edm::EDGetTokenT< reco::TrackCollection > theInputCollectionTag
bool operator()(const Hit &h1, const Hit &h2)
BaseTrackerRecHit const * ConstRecHitPointer
Definition: SeedingHitSet.h:14
T sqrt(T t)
Definition: SSEVec.h:23
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:676
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const edm::EDGetTokenT< reco::VertexCollection > theInputVertexCollectionTag
SeedingHitSet::ConstRecHitPointer Hit
void produce(edm::Event &ev, const edm::EventSetup &es) override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:664
TrackingRecHitRef recHit(size_t i) const
Get i-th hit on the track.
Definition: Track.h:94
fixed size matrix
HLT enums.
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
long double T
def move(src, dest)
Definition: eostools.py:511