CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GenericTripletGenerator.cc
Go to the documentation of this file.
4 
5 
8 #include <map>
9 using namespace ctfseeding;
10 
11 
13  theSeedingLayerToken(iC.consumes<SeedingLayerSetsHits>(conf.getParameter<edm::InputTag>("LayerSrc"))) {
14  edm::LogInfo("CtfSpecialSeedGenerator|GenericTripletGenerator") << "Constructing GenericTripletGenerator";
15 }
16 
17 
19  const edm::Event& e,
20  const edm::EventSetup& es){
21  hitTriplets.clear();
22  hitTriplets.reserve(0);
24  e.getByToken(theSeedingLayerToken, hlayers);
25  const SeedingLayerSetsHits& layers = *hlayers;
26  if(layers.numberOfLayersInSet() != 3)
27  throw cms::Exception("CtfSpecialSeedGenerator") << "You are using " << layers.numberOfLayersInSet() <<" layers in set instead of 3 ";
28  std::map<float, OrderedHitTriplet> radius_triplet_map;
30  auto innerHits = region.hits(e, es, ls[0]);
31  //std::cout << "innerHits.size()=" << innerHits.size() << std::endl;
32  auto middleHits = region.hits(e, es, ls[1]);
33  //std::cout << "middleHits.size()=" << middleHits.size() << std::endl;
34  auto outerHits = region.hits(e, es, ls[2]);
35  //std::cout << "outerHits.size()=" << outerHits.size() << std::endl;
36  //std::cout << "trying " << innerHits.size()*middleHits.size()*outerHits.size() << " combinations "<<std::endl;
37  for (auto iOuterHit = outerHits.begin(); iOuterHit != outerHits.end(); iOuterHit++){
38  for (auto iMiddleHit = middleHits.begin(); iMiddleHit != middleHits.end(); iMiddleHit++){
39  for (auto iInnerHit = innerHits.begin(); iInnerHit != innerHits.end(); iInnerHit++){
40  //GlobalPoint innerpos = ls[0].hitBuilder()->build(&(**iInnerHit))->globalPosition();
41  //GlobalPoint middlepos = ls[1].hitBuilder()->build(&(**iMiddleHit))->globalPosition();
42  //GlobalPoint outerpos = ls[2].hitBuilder()->build(&(**iOuterHit))->globalPosition();
43  //FastCircle circle(innerpos,
44  // middlepos,
45  // outerpos);
46  //do a first check on the radius of curvature to reduce the combinatorics
47  OrderedHitTriplet oht(&(**iInnerHit),&(**iMiddleHit),&(**iOuterHit));
48  std::pair<bool,float> val_radius = qualityFilter(oht,radius_triplet_map,ls);
49  if (val_radius.first){
50  //if (circle.rho() > 200 || circle.rho() == 0) { //0 radius means straight line
51  //hitTriplets.push_back(OrderedHitTriplet(*iInnerHit,
52  // *iMiddleHit,
53  // *iOuterHit));
54  hitTriplets.push_back(oht);
55  radius_triplet_map.insert(std::make_pair(val_radius.second,oht));
56  }
57  }
58  }
59  }
60  }
61  //std::cout << "ending with " << hitTriplets.size() << " triplets" << std::endl;
62  return hitTriplets;
63 }
64 
66  const std::map<float, OrderedHitTriplet>& map,
68  //first check the radius
69  GlobalPoint innerpos = oht.inner()->globalPosition();
70  GlobalPoint middlepos = oht.middle()->globalPosition();
71  GlobalPoint outerpos = oht.outer()->globalPosition();
72  std::vector<const TrackingRecHit*> ohttrh;
73  ohttrh.push_back(&(*(oht.inner()))); ohttrh.push_back(&(*(oht.middle()))); ohttrh.push_back(&(*(oht.outer())));
74  std::vector<const TrackingRecHit*>::const_iterator ioht;
75  //now chech that the change in phi is reasonable. the limiting angular distance is the one in case
76  //one of the two points is a tangent point.
77  float limit_phi_distance1 = sqrt((middlepos.x()-outerpos.x())*(middlepos.x()-outerpos.x()) +
78  (middlepos.y()-outerpos.y())*(middlepos.y()-outerpos.y()))/middlepos.mag();//actually this is the tangent of the limiting angle
79  float limit_phi_distance2 = sqrt((middlepos.x()-innerpos.x())*(middlepos.x()-innerpos.x()) +
80  (middlepos.y()-innerpos.y())*(middlepos.y()-innerpos.y()))/innerpos.mag();
81  //if (fabs(tan(outerpos.phi()-middlepos.phi()))>limit_phi_distance1 ||
82  // fabs(tan(innerpos.phi()-middlepos.phi()))>limit_phi_distance2) {
83  if (fabs(outerpos.phi()-middlepos.phi())>fabs(atan(limit_phi_distance1)) ||
84  fabs(innerpos.phi()-middlepos.phi())>fabs(atan(limit_phi_distance2)) ) {
85  //std::cout << "rejected because phi" << std::endl;
86  return std::make_pair(false, 0.);
87  }
88  //should we add a control on the r-z plane too?
89  /*
90  //now check that there is no big change in the r-z projection
91  float dz1 = outerpos.z()-middlepos.z();
92  float dr1 = sqrt(outerpos.x()*outerpos.x()+outerpos.y()*outerpos.y())-
93  sqrt(middlepos.x()*middlepos.x()+middlepos.y()*middlepos.y());
94  float dz2 = middlepos.z()-innerpos.z();
95  float dr2 = sqrt(middlepos.x()*middlepos.x()+middlepos.y()*middlepos.y())-
96  sqrt(innerpos.x()*innerpos.x()+innerpos.y()*innerpos.y());
97  float tan1 = dz1/dr1;
98  float tan2 = dz2/dr2;
99  //how much should we allow? should we make it configurable?
100  if (fabs(tan1-tan2)/tan1>0.5){
101  //std::cout << "rejected because z" << std::endl;
102  return std::make_pair(false, 0.);
103  }
104  */
105  //now check the radius is not too small
106  FastCircle circle(innerpos, middlepos, outerpos);
107  if (circle.rho() < 200 && circle.rho() != 0) return std::make_pair(false, circle.rho()); //to small radius
108  //now check if at least 2 hits are shared with an existing triplet
109  //look for similar radii in the map
110  std::map<float, OrderedHitTriplet>::const_iterator lower_bound = map.lower_bound((1-0.01)*circle.rho());
111  std::map<float, OrderedHitTriplet>::const_iterator upper_bound = map.upper_bound((1+0.01)*circle.rho());
112  std::map<float, OrderedHitTriplet>::const_iterator iter;
113  for (iter = lower_bound; iter != upper_bound && iter->first <= upper_bound->first; iter++){
114  int shared=0;
115  std::vector<const TrackingRecHit*> curtrh;
116  curtrh.push_back(&*(iter->second.inner()));curtrh.push_back(&*(iter->second.middle()));curtrh.push_back(&*(iter->second.outer()));
117  std::vector<const TrackingRecHit*>::const_iterator curiter;
118  for (curiter = curtrh.begin(); curiter != curtrh.end(); curiter++){
119  for (ioht = ohttrh.begin(); ioht != ohttrh.end(); ioht++){
120  if ((*ioht)->geographicalId()==(*curiter)->geographicalId() &&
121  ((*ioht)->localPosition()-(*curiter)->localPosition()).mag()<1e-5) shared++;
122  }
123  }
124  if (shared>1) return std::make_pair(false, circle.rho());
125  }
126 
127  return std::make_pair(true,circle.rho());
128 }
unsigned short numberOfLayersInSet() const
Get number of layers in each SeedingLayerSets.
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
MiddleRecHit middle() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
virtual Hits hits(const edm::Event &ev, const edm::EventSetup &es, const SeedingLayerSetsHits::SeedingLayer &layer) const =0
get hits from layer compatible with region constraints
double rho() const
Definition: FastCircle.h:54
edm::EDGetTokenT< SeedingLayerSetsHits > theSeedingLayerToken
InnerRecHit inner() const
BaseTrackerRecHit const * ConstRecHitPointer
Definition: SeedingHitSet.h:11
T mag() const
Definition: PV3DBase.h:67
T sqrt(T t)
Definition: SSEVec.h:48
std::pair< bool, float > qualityFilter(const OrderedHitTriplet &oht, const std::map< float, OrderedHitTriplet > &map, const SeedingLayerSetsHits::SeedingLayerSet &ls) const
tuple conf
Definition: dbtoconf.py:185
BaseTrackerRecHit const * SeedingHit
GenericTripletGenerator(const edm::ParameterSet &conf, edm::ConsumesCollector &iC)
virtual const OrderedSeedingHits & run(const TrackingRegion &region, const edm::Event &ev, const edm::EventSetup &es)
OuterRecHit outer() const
T x() const
Definition: PV3DBase.h:62