CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HitPairGeneratorFromLayerPairForPhotonConversion.cc
Go to the documentation of this file.
3 
7 
13 
16 
17 using namespace GeomDetEnumerators;
18 using namespace std;
20 
21 // #define mydebug_Seed
22 
24 
25 namespace {
26  template<class T> inline T sqr( T t) {return t*t;}
27 }
28 
29 
36 
38  unsigned int inner,
39  unsigned int outer,
40  LayerCacheType* layerCache,
41  unsigned int nSize,
42  unsigned int max)
43  :
44  theLayerCache(*layerCache), theOuterLayer(outer), theInnerLayer(inner), theMaxElement(max) {
45 
46 }
47 
50  const Layers& layers,
51  const edm::Event& event, const edm::EventSetup& es)
52 {
53  auto oldSize = result.size();
54  auto maxNum = result.size()+theMaxElement;
55 
56 #ifdef mydebug_Seed
57  ss.str("");
58 #endif
59 
60  typedef OrderedHitPair::InnerRecHit InnerHit;
61  typedef OrderedHitPair::OuterRecHit OuterHit;
63 
64  Layer innerLayerObj = layers[theInnerLayer];
65  Layer outerLayerObj = layers[theOuterLayer];
66 
67 #ifdef mydebug_Seed
68  ss << "In " << innerLayerObj.name() << " Out " << outerLayerObj.name() << std::endl;
69 #endif
70 
71  if(!checkBoundaries(*innerLayerObj.detLayer(),convRegion,40.,60.)) return; //FIXME, the maxSearchR(Z) are not optimized
72  if(!checkBoundaries(*outerLayerObj.detLayer(),convRegion,50.,60.)) return; //FIXME, the maxSearchR(Z) are not optimized
73 
74  /*get hit sorted in phi for each layer: NB: doesn't apply any region cut*/
75  const RecHitsSortedInPhi & innerHitsMap = theLayerCache(innerLayerObj, region, es);
76  if (innerHitsMap.empty()) return;
77 
78  const RecHitsSortedInPhi& outerHitsMap = theLayerCache(outerLayerObj, region, es);
79  if (outerHitsMap.empty()) return;
80  /*----------------*/
81 
82  /*This object will check the compatibility of the his in phi among the two layers. */
83  //InnerDeltaPhi deltaPhi(innerLayerObj.detLayer(), region, es);
84 
85  static const float nSigmaRZ = std::sqrt(12.f);
86  // static const float nSigmaPhi = 3.f;
87  vector<RecHitsSortedInPhi::Hit> innerHits, outerHits;
88  float outerPhimin, outerPhimax;
89  float innerPhimin, innerPhimax;
90 
91  /*Getting only the Hits in the outer layer that are compatible with the conversion region*/
92  if(!getPhiRange(outerPhimin,outerPhimax,*outerLayerObj.detLayer(),convRegion,es)) return;
93  outerHitsMap.hits( outerPhimin, outerPhimax, outerHits);
94 
95 #ifdef mydebug_Seed
96  ss << "\tophimin, ophimax " << outerPhimin << " " << outerPhimax << std::endl;
97 #endif
98 
99  /* loop on outer hits*/
100  for (vector<RecHitsSortedInPhi::Hit>::const_iterator oh = outerHits.begin(); oh!= outerHits.end(); ++oh) {
101  RecHitsSortedInPhi::Hit ohit = (*oh);
102 #ifdef mydebug_Seed
103  GlobalPoint oPos = ohit->globalPosition();
104 
105  ss << "\toPos " << oPos << " r " << oPos.perp() << " phi " << oPos.phi() << " cotTheta " << oPos.z()/oPos.perp() << std::endl;
106 #endif
107 
108  /*Check the compatibility of the ohit with the eta of the seeding track*/
109  if(checkRZCompatibilityWithSeedTrack(ohit,*outerLayerObj.detLayer(),convRegion)) continue;
110 
111  /*
112  //Do I need this? it uses a compatibility that probably I wouldn't
113  //Removing for the time being
114 
115  PixelRecoRange<float> phiRange = deltaPhi( oPos.perp(), oPos.phi(), oPos.z(), nSigmaPhi*(ohit->errorGlobalRPhi()));
116  if (phiRange.empty()) continue;
117  */
118 
119  const HitRZCompatibility *checkRZ = region.checkRZ(innerLayerObj.detLayer(), ohit, es);
120  if(!checkRZ) {
121 #ifdef mydebug_Seed
122  ss << "*******\nNo valid checkRZ\n*******" << std::endl;
123 #endif
124  continue;
125  }
126 
127  /*Get only the inner hits compatible with the conversion region*/
128  innerHits.clear();
129  if(!getPhiRange(innerPhimin,innerPhimax,*innerLayerObj.detLayer(),convRegion,es)) continue;
130  innerHitsMap.hits(innerPhimin, innerPhimax, innerHits);
131 
132 #ifdef mydebug_Seed
133  ss << "\tiphimin, iphimax " << innerPhimin << " " << innerPhimax << std::endl;
134 #endif
135 
136  /*Loop on inner hits*/
137  for ( vector<RecHitsSortedInPhi::Hit>::const_iterator ih=innerHits.begin(), ieh = innerHits.end(); ih < ieh; ++ih) {
138  GlobalPoint innPos = (*ih)->globalPosition();
139 
140 
141 #ifdef mydebug_Seed
142  ss << "\tinnPos " << innPos << " r " << innPos.perp() << " phi " << innPos.phi() << " cotTheta " << innPos.z()/innPos.perp() << std::endl;
143 #endif
144 
145  /*Check the compatibility of the ohit with the eta of the seeding track*/
146  if(checkRZCompatibilityWithSeedTrack(*ih,*innerLayerObj.detLayer(),convRegion)) continue;
147 
148  float r_reduced = std::sqrt( sqr(innPos.x()-region.origin().x())+sqr(innPos.y()-region.origin().y()));
149  Range allowed;
150  Range hitRZ;
151  if (innerLayerObj.detLayer()->location() == barrel) {
152  allowed = checkRZ->range(r_reduced);
153  float zErr = nSigmaRZ * (*ih)->errorGlobalZ();
154  hitRZ = Range(innPos.z()-zErr, innPos.z()+zErr);
155  } else {
156  allowed = checkRZ->range(innPos.z());
157  float rErr = nSigmaRZ * (*ih)->errorGlobalR();
158  hitRZ = Range(r_reduced-rErr, r_reduced+rErr);
159  }
160  Range crossRange = allowed.intersection(hitRZ);
161 
162 #ifdef mydebug_Seed
163  ss
164  << "\n\t\t allowed Range " << allowed.min() << " \t, " << allowed.max()
165  << "\n\t\t hitRz Range " << hitRZ.min() << " \t, " << hitRZ.max()
166  << "\n\t\t Cross Range " << crossRange.min() << " \t, " << crossRange.max()
167  << "\n\t\t the seed track has origin " << convRegion.convPoint() << " \t cotTheta " << convRegion.cotTheta()
168  << std::endl;
169 #endif
170 
171  if (! crossRange.empty() ) {
172 #ifdef mydebug_Seed
173  ss
174  << "\n\t\t !!!!ACCEPTED!!! \n\n";
175 #endif
176  if (theMaxElement!=0 && result.size() >= maxNum ){
177  result.resize(oldSize);
178 #ifdef mydebug_Seed
179  edm::LogError("TooManySeeds")<<"number of pairs exceed maximum, no pairs produced";
180 #endif
181  delete checkRZ;
182 
183 #ifdef mydebug_Seed
184  std::cout << ss.str();
185 #endif
186  return;
187  }
188  result.push_back( OrderedHitPair( *ih, ohit) );
189  }
190  }
191  delete checkRZ;
192  }
193 
194 #ifdef mydebug_Seed
195  std::cout << ss.str();
196 #endif
197 }
198 
199 
201  if (layer.location() == GeomDetEnumerators::barrel){
202  const BarrelDetLayer& bl = static_cast<const BarrelDetLayer&>(layer);
203  float rLayer = bl.specificSurface().radius();
204 
205  // the maximal delta phi will be for the innermost hits
206  float theThickness = layer.surface().bounds().thickness();
207  return rLayer + 0.5f*theThickness;
208  }
209 
210  //Fixme
211  return 0;
212 }
213 
215  if (layer.location() == GeomDetEnumerators::endcap){
216  float layerZ = layer.position().z();
217  float theThickness = layer.surface().bounds().thickness();
218  float layerZmax = layerZ > 0 ? layerZ+0.5f*theThickness: layerZ-0.5f*theThickness;
219  return layerZmax;
220  }else{
221  //Fixme
222  return 0;
223  }
224 }
225 bool HitPairGeneratorFromLayerPairForPhotonConversion::checkBoundaries(const DetLayer& layer,const ConversionRegion& convRegion,float maxSearchR, float maxSearchZ){
226 
227  if(layer.location() == GeomDetEnumerators::barrel){
228 
229  float minZEndCap=130;
230  if(fabs(convRegion.convPoint().z()) > minZEndCap){
231 #ifdef mydebug_Seed
232  ss << "\tthe conversion seems to be in the endcap. Zconv " << convRegion.convPoint().z() << std::endl;
233  std::cout << ss.str();
234 #endif
235  return false;
236  }
237 
238  float R=getLayerRadius(layer);
239 
240  if(convRegion.convPoint().perp()>R){
241 #ifdef mydebug_Seed
242  ss << "\tthis layer is before the conversion : R layer " << R << " [ Rconv " << convRegion.convPoint().perp() << " Zconv " << convRegion.convPoint().z()<< std::endl;
243  std::cout << ss.str();
244 #endif
245  return false;
246  }
247 
248  if(R - convRegion.convPoint().perp() > maxSearchR ){
249 #ifdef mydebug_Seed
250  ss << "\tthis layer is far from the conversion more than cut " << maxSearchR << " cm. R layer " << R << " [ Rconv " << convRegion.convPoint().perp() << " Zconv " << convRegion.convPoint().z()<< std::endl;
251  std::cout << ss.str();
252 #endif
253  return false;
254  }
255 
256  }else if (layer.location() == GeomDetEnumerators::endcap){
257 
258  float Z=getLayerZ(layer);
259  if(
260  (convRegion.convPoint().z()>0 && convRegion.convPoint().z()>Z)
261  ||
262  (convRegion.convPoint().z()<0 && convRegion.convPoint().z()<Z)
263  ) {
264 #ifdef mydebug_Seed
265  ss << "\tthis layer is before the conversion : Z layer " << Z << " [ Rconv " << convRegion.convPoint().perp()<< " Zconv " << convRegion.convPoint().z() << std::endl;
266  std::cout << ss.str();
267 #endif
268  return false;
269  }
270 
271 
272  if(fabs(Z - convRegion.convPoint().z()) > maxSearchZ ){
273 #ifdef mydebug_Seed
274  ss << "\tthis layer is far from the conversion more than cut " << maxSearchZ << " cm. Z layer " << Z << " [ Rconv " << convRegion.convPoint().perp()<< " Zconv " << convRegion.convPoint().z() << std::endl;
275  std::cout << ss.str();
276 #endif
277  return false;
278  }
279 
280  }
281  return true;
282 }
283 
284 bool HitPairGeneratorFromLayerPairForPhotonConversion::getPhiRange(float& Phimin, float& Phimax, const DetLayer& layer, const ConversionRegion &convRegion, const edm::EventSetup& es){
285  if(layer.location() == GeomDetEnumerators::barrel){
286  return getPhiRange(Phimin,Phimax,getLayerRadius(layer),convRegion,es);
287  }else if (layer.location() == GeomDetEnumerators::endcap){
288  float Z=getLayerZ(layer);
289  float R=Z/convRegion.cotTheta();
290  return getPhiRange(Phimin,Phimax,R,convRegion,es); //FIXME
291  }
292  return false;
293 }
294 
295 bool HitPairGeneratorFromLayerPairForPhotonConversion::getPhiRange(float& Phimin, float& Phimax, const float& layerR, const ConversionRegion &convRegion, const edm::EventSetup& es){
296  Phimin = reco::deltaPhi(convRegion.convPoint().phi(),0.);
297 
298  float dphi;
299  float ptmin=0.1;
300  float DeltaL=layerR-convRegion.convPoint().perp();
301 
302  if(DeltaL<0){
303  Phimin=0;
304  Phimax=0;
305  return false;
306  }
307 
308  float theRCurvatureMin = PixelRecoUtilities::bendingRadius(ptmin,es);
309 
310  if(theRCurvatureMin<DeltaL)
311  dphi = atan(DeltaL/layerR);
312  else
313  dphi = atan(theRCurvatureMin/layerR * ( 1 - sqrt(1-sqr(DeltaL/theRCurvatureMin)) ) );
314 
315  if(convRegion.charge()>0){
316  Phimax=Phimin;
317  Phimin=Phimax-dphi;
318  }else{
319  Phimax=Phimin+dphi;
320  }
321 
322  //std::cout << dphi << " " << Phimin << " " << Phimax << " " << layerR << " " << DeltaL << " " << convRegion.convPoint().phi() << " " << convRegion.convPoint().perp()<< std::endl;
323  return true;
324 }
325 
328  static const float nSigmaRZ = std::sqrt(12.f);
329  Range hitCotTheta;
330 
331  double sigmaCotTheta = convRegion.errTheta() * (1+convRegion.cotTheta()*convRegion.cotTheta()); //Error Propagation from sigma theta.
332  Range allowedCotTheta(convRegion.cotTheta()-nSigmaRZ*sigmaCotTheta,convRegion.cotTheta()+nSigmaRZ*sigmaCotTheta);
333 
334  double dz = hit->globalPosition().z()-convRegion.pvtxPoint().z();
335  double r_reduced = std::sqrt( sqr(hit->globalPosition().x()-convRegion.pvtxPoint().x())+sqr(hit->globalPosition().y()-convRegion.pvtxPoint().y()));
336 
337  if (layer.location() == GeomDetEnumerators::barrel){
338  float zErr = nSigmaRZ * hit->errorGlobalZ();
339  hitCotTheta = Range(getCot(dz-zErr,r_reduced),getCot(dz+zErr,r_reduced));
340  }else{
341  float rErr = nSigmaRZ * hit->errorGlobalR();
342  if(dz>0)
343  hitCotTheta = Range(getCot(dz,r_reduced+rErr),getCot(dz,r_reduced-rErr));
344  else
345  hitCotTheta = Range(getCot(dz,r_reduced-rErr), getCot(dz,r_reduced+rErr));
346  }
347 
348  Range crossRange = allowedCotTheta.intersection(hitCotTheta);
349 
350 #ifdef mydebug_Seed
351  ss
352  << "\n\t\t cotTheta allowed Range " << allowedCotTheta.min() << " \t, " << allowedCotTheta.max()
353  << "\n\t\t hitCotTheta Range " << hitCotTheta.min() << " \t, " << hitCotTheta.max()
354  << "\n\t\t Cross Range " << crossRange.min() << " \t, " << crossRange.max()
355  << "\n\t\t the seed track has origin " << convRegion.convPoint() << " \t cotTheta " << convRegion.cotTheta()
356  << std::endl;
357 #endif
358 
359  return crossRange.empty();
360 }
361 
362 
364 getCot(double dz, double dr){
365  if ( std::abs(dr) > 1.e-4f ) return dz/dr;
366  else
367  if(dz>0) return 99999.f;
368  else return -99999.f;
369 }
370 
const double Z[kNumberCalorimeter]
std::vector< Hit > hits(float phiMin, float phiMax) const
float cotTheta() const
virtual const BoundSurface & surface() const =0
The surface of the GeometricSearchDet.
bool checkRZCompatibilityWithSeedTrack(const RecHitsSortedInPhi::Hit &hit, const DetLayer &layer, const ConversionRegion &convRegion)
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
T perp() const
Definition: PV3DBase.h:72
GlobalPoint const & origin() const
virtual Location location() const =0
Which part of the detector (barrel, endcap)
virtual Range range(const float &rORz) const =0
GlobalPoint pvtxPoint() const
int charge() const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
SeedingHitSet::ConstRecHitPointer InnerRecHit
GlobalPoint convPoint() const
const Bounds & bounds() const
Definition: Surface.h:120
virtual unsigned int size() const
void hitPairs(const ConversionRegion &convRegion, const TrackingRegion &reg, OrderedHitPairs &prs, const Layers &layers, const edm::Event &ev, const edm::EventSetup &es)
tuple result
Definition: mps_fire.py:84
virtual float thickness() const =0
PixelRecoRange< float > Range
virtual HitRZCompatibility * checkRZ(const DetLayer *layer, const Hit &outerHit, const edm::EventSetup &iSetup, const DetLayer *outerlayer=0, float lr=0, float gz=0, float dr=0, float dz=0) const =0
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
bool getPhiRange(float &Phimin, float &Phimax, const DetLayer &layer, const ConversionRegion &convRegion, const edm::EventSetup &es)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
const std::string & name() const
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
SeedingHitSet::ConstRecHitPointer Hit
BaseTrackerRecHit const * Hit
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
bool checkBoundaries(const DetLayer &layer, const ConversionRegion &convRegion, float maxSearchR, float maxSearchZ)
virtual const Surface::PositionType & position() const
Returns position of the surface.
virtual const BoundCylinder & specificSurface() const final
Extension of the interface.
const DetLayer * detLayer() const
double errTheta() const
T bendingRadius(T pt, const edm::EventSetup &iSetup)
double ptmin
Definition: HydjetWrapper.h:90
PixelRecoRange< T > intersection(const PixelRecoRange< T > &r) const
Square< F >::type sqr(const F &f)
Definition: Square.h:13
HitPairGeneratorFromLayerPairForPhotonConversion(unsigned int inner, unsigned int outer, LayerCacheType *layerCache, unsigned int nSize=30000, unsigned int max=0)
tuple cout
Definition: gather_cfg.py:145
long double T
T x() const
Definition: PV3DBase.h:62
SeedingHitSet::ConstRecHitPointer OuterRecHit