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