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