CMS 3D CMS Logo

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