CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SeedGeneratorFromProtoTracksEDProducer.cc
Go to the documentation of this file.
2 
5 
10 
18 
23 
24 
26 #include <vector>
27 #include <cassert>
28 using namespace edm;
29 using namespace reco;
30 
31 template <class T> T sqr( T t) {return t*t;}
33 
34 struct HitLessByRadius { bool operator() (const Hit& h1, const Hit & h2) { return h1->globalPosition().perp2() < h2->globalPosition().perp2(); } };
35 
38  desc.add<InputTag>("InputCollection", InputTag("pixelTracks"));
39  desc.add<InputTag>("InputVertexCollection", InputTag(""));
40  desc.add<double>("originHalfLength", 1E9);
41  desc.add<double>("originRadius", 1E9);
42  desc.add<bool>("useProtoTrackKinematics", false);
43  desc.add<bool>("useEventsWithNoVertex", true);
44  desc.add<std::string>("TTRHBuilder", "TTRHBuilderWithoutAngle4PixelTriplets");
45  desc.add<bool>("usePV", false);
46  descriptions.add("SeedGeneratorFromProtoTracksEDProducer", desc);
47 }
48 
49 
51 
52 {
53  produces<TrajectorySeedCollection>();
54  theInputCollectionTag = consumes<reco::TrackCollection>(cfg.getParameter<InputTag>("InputCollection"));
55  theInputVertexCollectionTag = consumes<reco::VertexCollection>(cfg.getParameter<InputTag>("InputVertexCollection"));
56  originHalfLength = cfg.getParameter<double>("originHalfLength");
57  originRadius = cfg.getParameter<double>("originRadius");
58  useProtoTrackKinematics = cfg.getParameter<bool>("useProtoTrackKinematics");
59  useEventsWithNoVertex = cfg.getParameter<bool>("useEventsWithNoVertex");
60  builderName = cfg.getParameter<std::string>("TTRHBuilder");
61  usePV_ = cfg.getParameter<bool>( "usePV" );
62 }
63 
64 
66 {
67  std::auto_ptr<TrajectorySeedCollection> result(new TrajectorySeedCollection());
70 
71  const TrackCollection &protos = *(trks.product());
72 
74  bool foundVertices = ev.getByToken(theInputVertexCollectionTag, vertices);
75  //const reco::VertexCollection & vertices = *(h_vertices.product());
76 
80  for (TrackCollection::const_iterator it=protos.begin(); it!= protos.end(); ++it) {
81  const Track & proto = (*it);
82  GlobalPoint vtx(proto.vertex().x(), proto.vertex().y(), proto.vertex().z());
83 
84  // check the compatibility with a primary vertex
85  bool keepTrack = false;
86  if ( !foundVertices ) {
87  if (useEventsWithNoVertex) keepTrack = true;
88  }
89  else if (usePV_){
90 
91  GlobalPoint aPV(vertices->begin()->position().x(),vertices->begin()->position().y(),vertices->begin()->position().z());
92  double distR2 = sqr(vtx.x()-aPV.x()) +sqr(vtx.y()-aPV.y());
93  double distZ = fabs(vtx.z()-aPV.z());
94  if ( distR2 < sqr(originRadius) && distZ < originHalfLength ) {
95  keepTrack = true;
96  }
97  }
98 
99  else {
100  for (reco::VertexCollection::const_iterator iv=vertices->begin(); iv!= vertices->end(); ++iv) {
101  GlobalPoint aPV(iv->position().x(),iv->position().y(),iv->position().z());
102  double distR2 = sqr(vtx.x()-aPV.x()) +sqr(vtx.y()-aPV.y());
103  double distZ = fabs(vtx.z()-aPV.z());
104  if ( distR2 < sqr(originRadius) && distZ < originHalfLength ) {
105  keepTrack = true;
106  break;
107  }
108  }
109  }
110  if (!keepTrack) continue;
111 
112  if ( useProtoTrackKinematics ) {
113  SeedFromProtoTrack seedFromProtoTrack( proto, es);
114  if (seedFromProtoTrack.isValid()) (*result).push_back( seedFromProtoTrack.trajectorySeed() );
115  } else {
117  es.get<TransientRecHitRecord>().get(builderName,ttrhbESH);
118  std::vector<Hit> hits;
119  for (unsigned int iHit = 0, nHits = proto.recHitsSize(); iHit < nHits; ++iHit) {
120  TrackingRecHitRef refHit = proto.recHit(iHit);
121  if(refHit->isValid()) hits.push_back(ttrhbESH->build( &(*refHit) ));
122  }
123  sort(hits.begin(), hits.end(), HitLessByRadius());
124  assert(hits.size()<4);
125  if (hits.size() > 1) {
126  double mom_perp = sqrt(proto.momentum().x()*proto.momentum().x()+proto.momentum().y()*proto.momentum().y());
127  GlobalTrackingRegion region(mom_perp, vtx, 0.2, 0.2);
128  SeedFromConsecutiveHitsCreator seedCreator;
129  seedCreator.init(region, es, 0);
130  seedCreator.makeSeed(*result, SeedingHitSet(hits[0], hits[1], hits.size() >2 ? hits[2] : SeedingHitSet::nullPtr() ));
131  }
132  }
133  }
134 
135  ev.put(result);
136 }
137 
T getParameter(std::string const &) const
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:148
edm::EDGetTokenT< reco::VertexCollection > theInputVertexCollectionTag
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:68
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
edm::EDGetTokenT< reco::TrackCollection > theInputCollectionTag
TrajectorySeed trajectorySeed() const
std::vector< TrajectorySeed > TrajectorySeedCollection
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
T sqrt(T t)
Definition: SSEVec.h:48
virtual void init(const TrackingRegion &region, const edm::EventSetup &es, const SeedComparitor *filter) GCC11_FINAL
static ConstRecHitPointer nullPtr()
Definition: SeedingHitSet.h:11
tuple result
Definition: query.py:137
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:154
virtual void produce(edm::Event &ev, const edm::EventSetup &es) override
const T & get() const
Definition: EventSetup.h:55
void add(std::string const &label, ParameterSetDescription const &psetDescription)
T const * product() const
Definition: Handle.h:81
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Square< F >::type sqr(const F &f)
Definition: Square.h:13
TrackingRecHitRef recHit(size_t i) const
Get i-th hit on the track.
Definition: Track.h:66
virtual void makeSeed(TrajectorySeedCollection &seedCollection, const SeedingHitSet &hits) GCC11_FINAL
TransientTrackingRecHit::ConstRecHitPointer Hit
long double T