CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
HitPairGeneratorFromLayerPairForPhotonConversion Class Reference

#include <HitPairGeneratorFromLayerPairForPhotonConversion.h>

Public Types

typedef SeedingLayerSetsHits::SeedingLayer Layer
 
typedef LayerHitMapCache LayerCacheType
 
typedef SeedingLayerSetsHits::SeedingLayerSet Layers
 

Public Member Functions

bool checkBoundaries (const DetLayer &layer, const ConversionRegion &convRegion, float maxSearchR, float maxSearchZ)
 
bool checkRZCompatibilityWithSeedTrack (const RecHitsSortedInPhi::Hit &hit, const DetLayer &layer, const ConversionRegion &convRegion)
 
float getLayerRadius (const DetLayer &layer)
 
float getLayerZ (const DetLayer &layer)
 
bool getPhiRange (float &Phimin, float &Phimax, const DetLayer &layer, const ConversionRegion &convRegion, const MagneticField &field)
 
bool getPhiRange (float &Phimin, float &Phimax, const float &layerR, const ConversionRegion &convRegion, const MagneticField &field)
 
 HitPairGeneratorFromLayerPairForPhotonConversion (edm::ConsumesCollector iC, unsigned int inner, unsigned int outer, LayerCacheType *layerCache, unsigned int nSize=30000, unsigned int max=0)
 
void hitPairs (const ConversionRegion &convRegion, const TrackingRegion &reg, OrderedHitPairs &prs, const Layers &layers, const edm::Event &ev, const edm::EventSetup &es)
 

Private Member Functions

double getCot (double dz, double dr)
 

Private Attributes

std::stringstream ss
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordtheFieldToken
 
const unsigned int theInnerLayer
 
LayerCacheTypetheLayerCache
 
const unsigned int theMaxElement
 
const unsigned int theOuterLayer
 

Detailed Description

Definition at line 17 of file HitPairGeneratorFromLayerPairForPhotonConversion.h.

Member Typedef Documentation

◆ Layer

◆ LayerCacheType

◆ Layers

Constructor & Destructor Documentation

◆ HitPairGeneratorFromLayerPairForPhotonConversion()

HitPairGeneratorFromLayerPairForPhotonConversion::HitPairGeneratorFromLayerPairForPhotonConversion ( edm::ConsumesCollector  iC,
unsigned int  inner,
unsigned int  outer,
LayerCacheType layerCache,
unsigned int  nSize = 30000,
unsigned int  max = 0 
)

Member Function Documentation

◆ checkBoundaries()

bool HitPairGeneratorFromLayerPairForPhotonConversion::checkBoundaries ( const DetLayer layer,
const ConversionRegion convRegion,
float  maxSearchR,
float  maxSearchZ 
)

Definition at line 231 of file HitPairGeneratorFromLayerPairForPhotonConversion.cc.

References GeomDetEnumerators::barrel, ConversionRegion::convPoint(), gather_cfg::cout, GeomDetEnumerators::endcap, getLayerRadius(), getLayerZ(), phase1PixelTopology::layer, PV3DBase< T, PVType, FrameType >::perp(), dttmaxenums::R, ss, BeamSpotPI::Z, and PV3DBase< T, PVType, FrameType >::z().

Referenced by hitPairs().

234  {
235  if (layer.location() == GeomDetEnumerators::barrel) {
236  float minZEndCap = 130;
237  if (fabs(convRegion.convPoint().z()) > minZEndCap) {
238 #ifdef mydebug_Seed
239  ss << "\tthe conversion seems to be in the endcap. Zconv " << convRegion.convPoint().z() << std::endl;
240  std::cout << ss.str();
241 #endif
242  return false;
243  }
244 
245  float R = getLayerRadius(layer);
246 
247  if (convRegion.convPoint().perp() > R) {
248 #ifdef mydebug_Seed
249  ss << "\tthis layer is before the conversion : R layer " << R << " [ Rconv " << convRegion.convPoint().perp()
250  << " Zconv " << convRegion.convPoint().z() << std::endl;
251  std::cout << ss.str();
252 #endif
253  return false;
254  }
255 
256  if (R - convRegion.convPoint().perp() > maxSearchR) {
257 #ifdef mydebug_Seed
258  ss << "\tthis layer is far from the conversion more than cut " << maxSearchR << " cm. R layer " << R
259  << " [ Rconv " << convRegion.convPoint().perp() << " Zconv " << convRegion.convPoint().z() << std::endl;
260  std::cout << ss.str();
261 #endif
262  return false;
263  }
264 
265  } else if (layer.location() == GeomDetEnumerators::endcap) {
266  float Z = getLayerZ(layer);
267  if ((convRegion.convPoint().z() > 0 && convRegion.convPoint().z() > Z) ||
268  (convRegion.convPoint().z() < 0 && convRegion.convPoint().z() < Z)) {
269 #ifdef mydebug_Seed
270  ss << "\tthis layer is before the conversion : Z layer " << Z << " [ Rconv " << convRegion.convPoint().perp()
271  << " Zconv " << convRegion.convPoint().z() << std::endl;
272  std::cout << ss.str();
273 #endif
274  return false;
275  }
276 
277  if (fabs(Z - convRegion.convPoint().z()) > maxSearchZ) {
278 #ifdef mydebug_Seed
279  ss << "\tthis layer is far from the conversion more than cut " << maxSearchZ << " cm. Z layer " << Z
280  << " [ Rconv " << convRegion.convPoint().perp() << " Zconv " << convRegion.convPoint().z() << std::endl;
281  std::cout << ss.str();
282 #endif
283  return false;
284  }
285  }
286  return true;
287 }
T perp() const
Definition: PV3DBase.h:69
T z() const
Definition: PV3DBase.h:61
constexpr std::array< uint8_t, layerIndexSize > layer
GlobalPoint convPoint() const

◆ checkRZCompatibilityWithSeedTrack()

bool HitPairGeneratorFromLayerPairForPhotonConversion::checkRZCompatibilityWithSeedTrack ( const RecHitsSortedInPhi::Hit hit,
const DetLayer layer,
const ConversionRegion convRegion 
)

Definition at line 336 of file HitPairGeneratorFromLayerPairForPhotonConversion.cc.

References GeomDetEnumerators::barrel, ConversionRegion::convPoint(), ConversionRegion::cotTheta(), PVValHelper::dz, ConversionRegion::errTheta(), f, getCot(), PixelRecoRange< T >::intersection(), phase1PixelTopology::layer, nSigmaRZ, ConversionRegion::pvtxPoint(), funct::sqr(), mathSSE::sqrt(), ss, PV3DBase< T, PVType, FrameType >::x(), hit::x, PV3DBase< T, PVType, FrameType >::y(), hit::y, PV3DBase< T, PVType, FrameType >::z(), and hit::z.

Referenced by hitPairs().

337  {
338  static const float nSigmaRZ = std::sqrt(12.f);
339  Range hitCotTheta;
340 
341  double sigmaCotTheta = convRegion.errTheta() *
342  (1 + convRegion.cotTheta() * convRegion.cotTheta()); //Error Propagation from sigma theta.
343  Range allowedCotTheta(convRegion.cotTheta() - nSigmaRZ * sigmaCotTheta,
344  convRegion.cotTheta() + nSigmaRZ * sigmaCotTheta);
345 
346  double dz = hit->globalPosition().z() - convRegion.pvtxPoint().z();
347  double r_reduced = std::sqrt(sqr(hit->globalPosition().x() - convRegion.pvtxPoint().x()) +
348  sqr(hit->globalPosition().y() - convRegion.pvtxPoint().y()));
349 
350  if (layer.location() == GeomDetEnumerators::barrel) {
351  float zErr = nSigmaRZ * hit->errorGlobalZ();
352  hitCotTheta = Range(getCot(dz - zErr, r_reduced), getCot(dz + zErr, r_reduced));
353  } else {
354  float rErr = nSigmaRZ * hit->errorGlobalR();
355  if (dz > 0)
356  hitCotTheta = Range(getCot(dz, r_reduced + rErr), getCot(dz, r_reduced - rErr));
357  else
358  hitCotTheta = Range(getCot(dz, r_reduced - rErr), getCot(dz, r_reduced + rErr));
359  }
360 
361  Range crossRange = allowedCotTheta.intersection(hitCotTheta);
362 
363 #ifdef mydebug_Seed
364  ss << "\n\t\t cotTheta allowed Range " << allowedCotTheta.min() << " \t, " << allowedCotTheta.max()
365  << "\n\t\t hitCotTheta Range " << hitCotTheta.min() << " \t, " << hitCotTheta.max() << "\n\t\t Cross Range "
366  << crossRange.min() << " \t, " << crossRange.max() << "\n\t\t the seed track has origin " << convRegion.convPoint()
367  << " \t cotTheta " << convRegion.cotTheta() << std::endl;
368 #endif
369 
370  return crossRange.empty();
371 }
GlobalPoint pvtxPoint() const
PixelRecoRange< float > Range
T z() const
Definition: PV3DBase.h:61
constexpr std::array< uint8_t, layerIndexSize > layer
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
float cotTheta() const
constexpr double nSigmaRZ
T sqrt(T t)
Definition: SSEVec.h:19
double f[11][100]
PixelRecoRange< T > intersection(const PixelRecoRange< T > &r) const
Square< F >::type sqr(const F &f)
Definition: Square.h:14
GlobalPoint convPoint() const
double errTheta() const

◆ getCot()

double HitPairGeneratorFromLayerPairForPhotonConversion::getCot ( double  dz,
double  dr 
)
private

Definition at line 373 of file HitPairGeneratorFromLayerPairForPhotonConversion.cc.

References funct::abs(), l1ctLayer1_cff::dr, PVValHelper::dz, MillePedeFileConverter_cfg::e, and f.

Referenced by checkRZCompatibilityWithSeedTrack().

373  {
374  if (std::abs(dr) > 1.e-4f)
375  return dz / dr;
376  else if (dz > 0)
377  return 99999.f;
378  else
379  return -99999.f;
380 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]

◆ getLayerRadius()

float HitPairGeneratorFromLayerPairForPhotonConversion::getLayerRadius ( const DetLayer layer)

Definition at line 206 of file HitPairGeneratorFromLayerPairForPhotonConversion.cc.

References GeomDetEnumerators::barrel, phase1PixelTopology::layer, and BarrelDetLayer::specificSurface().

Referenced by checkBoundaries(), and getPhiRange().

206  {
207  if (layer.location() == GeomDetEnumerators::barrel) {
208  const BarrelDetLayer& bl = static_cast<const BarrelDetLayer&>(layer);
209  float rLayer = bl.specificSurface().radius();
210 
211  // the maximal delta phi will be for the innermost hits
212  float theThickness = layer.surface().bounds().thickness();
213  return rLayer + 0.5f * theThickness;
214  }
215 
216  //Fixme
217  return 0;
218 }
constexpr std::array< uint8_t, layerIndexSize > layer
virtual const BoundCylinder & specificSurface() const final
Extension of the interface.

◆ getLayerZ()

float HitPairGeneratorFromLayerPairForPhotonConversion::getLayerZ ( const DetLayer layer)

Definition at line 220 of file HitPairGeneratorFromLayerPairForPhotonConversion.cc.

References GeomDetEnumerators::endcap, and phase1PixelTopology::layer.

Referenced by checkBoundaries(), and getPhiRange().

220  {
221  if (layer.location() == GeomDetEnumerators::endcap) {
222  float layerZ = layer.position().z();
223  float theThickness = layer.surface().bounds().thickness();
224  float layerZmax = layerZ > 0 ? layerZ + 0.5f * theThickness : layerZ - 0.5f * theThickness;
225  return layerZmax;
226  } else {
227  //Fixme
228  return 0;
229  }
230 }
constexpr std::array< uint8_t, layerIndexSize > layer

◆ getPhiRange() [1/2]

bool HitPairGeneratorFromLayerPairForPhotonConversion::getPhiRange ( float &  Phimin,
float &  Phimax,
const DetLayer layer,
const ConversionRegion convRegion,
const MagneticField field 
)

Definition at line 289 of file HitPairGeneratorFromLayerPairForPhotonConversion.cc.

References GeomDetEnumerators::barrel, ConversionRegion::cotTheta(), GeomDetEnumerators::endcap, getLayerRadius(), getLayerZ(), phase1PixelTopology::layer, ElectronMcFakeValidator_cfi::Phimax, ElectronMcFakeValidator_cfi::Phimin, dttmaxenums::R, and BeamSpotPI::Z.

Referenced by hitPairs().

293  {
294  if (layer.location() == GeomDetEnumerators::barrel) {
295  return getPhiRange(Phimin, Phimax, getLayerRadius(layer), convRegion, field);
296  } else if (layer.location() == GeomDetEnumerators::endcap) {
297  float Z = getLayerZ(layer);
298  float R = Z / convRegion.cotTheta();
299  return getPhiRange(Phimin, Phimax, R, convRegion, field); //FIXME
300  }
301  return false;
302 }
constexpr std::array< uint8_t, layerIndexSize > layer
float cotTheta() const
bool getPhiRange(float &Phimin, float &Phimax, const DetLayer &layer, const ConversionRegion &convRegion, const MagneticField &field)

◆ getPhiRange() [2/2]

bool HitPairGeneratorFromLayerPairForPhotonConversion::getPhiRange ( float &  Phimin,
float &  Phimax,
const float &  layerR,
const ConversionRegion convRegion,
const MagneticField field 
)

Definition at line 304 of file HitPairGeneratorFromLayerPairForPhotonConversion.cc.

References PixelRecoUtilities::bendingRadius(), ConversionRegion::charge(), ConversionRegion::convPoint(), reco::deltaPhi(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), ElectronMcFakeValidator_cfi::Phimax, ElectronMcFakeValidator_cfi::Phimin, ptmin, funct::sqr(), and mathSSE::sqrt().

305  {
306  Phimin = reco::deltaPhi(convRegion.convPoint().phi(), 0.);
307 
308  float dphi;
309  float ptmin = 0.1;
310  float DeltaL = layerR - convRegion.convPoint().perp();
311 
312  if (DeltaL < 0) {
313  Phimin = 0;
314  Phimax = 0;
315  return false;
316  }
317 
318  float theRCurvatureMin = PixelRecoUtilities::bendingRadius(ptmin, field);
319 
320  if (theRCurvatureMin < DeltaL)
321  dphi = atan(DeltaL / layerR);
322  else
323  dphi = atan(theRCurvatureMin / layerR * (1 - sqrt(1 - sqr(DeltaL / theRCurvatureMin))));
324 
325  if (convRegion.charge() > 0) {
326  Phimax = Phimin;
327  Phimin = Phimax - dphi;
328  } else {
329  Phimax = Phimin + dphi;
330  }
331 
332  //std::cout << dphi << " " << Phimin << " " << Phimax << " " << layerR << " " << DeltaL << " " << convRegion.convPoint().phi() << " " << convRegion.convPoint().perp()<< std::endl;
333  return true;
334 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
T perp() const
Definition: PV3DBase.h:69
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
int charge() const
T sqrt(T t)
Definition: SSEVec.h:19
double ptmin
Definition: HydjetWrapper.h:84
Square< F >::type sqr(const F &f)
Definition: Square.h:14
T bendingRadius(T pt, const MagneticField &field)
GlobalPoint convPoint() const

◆ hitPairs()

void HitPairGeneratorFromLayerPairForPhotonConversion::hitPairs ( const ConversionRegion convRegion,
const TrackingRegion reg,
OrderedHitPairs prs,
const Layers layers,
const edm::Event ev,
const edm::EventSetup es 
)

Definition at line 51 of file HitPairGeneratorFromLayerPairForPhotonConversion.cc.

References Reference_intrackfit_cff::barrel, checkBoundaries(), checkRZCompatibilityWithSeedTrack(), ConversionRegion::convPoint(), ConversionRegion::cotTheta(), gather_cfg::cout, SeedingLayerSetsHits::SeedingLayer::detLayer(), RecHitsSortedInPhi::empty(), f, edm::EventSetup::getData(), getPhiRange(), RecHitsSortedInPhi::hits(), PixelRecoRange< T >::intersection(), hgcalTopologyTester_cfi::layers, DetLayer::location(), SeedingLayerSetsHits::SeedingLayer::name(), nSigmaRZ, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), HLT_2022v15_cff::region, mps_fire::result, funct::sqr(), mathSSE::sqrt(), ss, theFieldToken, theInnerLayer, theLayerCache, theMaxElement, theOuterLayer, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

56  {
57  auto oldSize = result.size();
58  auto maxNum = result.size() + theMaxElement;
59 
60 #ifdef mydebug_Seed
61  ss.str("");
62 #endif
63 
64  typedef OrderedHitPair::InnerRecHit InnerHit;
65  typedef OrderedHitPair::OuterRecHit OuterHit;
67 
68  Layer innerLayerObj = layers[theInnerLayer];
69  Layer outerLayerObj = layers[theOuterLayer];
70 
71 #ifdef mydebug_Seed
72  ss << "In " << innerLayerObj.name() << " Out " << outerLayerObj.name() << std::endl;
73 #endif
74 
75  if (!checkBoundaries(*innerLayerObj.detLayer(), convRegion, 40., 60.))
76  return; //FIXME, the maxSearchR(Z) are not optimized
77  if (!checkBoundaries(*outerLayerObj.detLayer(), convRegion, 50., 60.))
78  return; //FIXME, the maxSearchR(Z) are not optimized
79 
80  /*get hit sorted in phi for each layer: NB: doesn't apply any region cut*/
81  const RecHitsSortedInPhi& innerHitsMap = theLayerCache(innerLayerObj, region);
82  if (innerHitsMap.empty())
83  return;
84 
85  const RecHitsSortedInPhi& outerHitsMap = theLayerCache(outerLayerObj, region);
86  if (outerHitsMap.empty())
87  return;
88  /*----------------*/
89 
90  /*This object will check the compatibility of the his in phi among the two layers. */
91  //InnerDeltaPhi deltaPhi(innerLayerObj.detLayer(), region, es);
92 
93  static const float nSigmaRZ = std::sqrt(12.f);
94  // static const float nSigmaPhi = 3.f;
95  vector<RecHitsSortedInPhi::Hit> innerHits, outerHits;
96  float outerPhimin, outerPhimax;
97  float innerPhimin, innerPhimax;
98 
99  /*Getting only the Hits in the outer layer that are compatible with the conversion region*/
100  const auto& field = es.getData(theFieldToken);
101  if (!getPhiRange(outerPhimin, outerPhimax, *outerLayerObj.detLayer(), convRegion, field))
102  return;
103  outerHitsMap.hits(outerPhimin, outerPhimax, outerHits);
104 
105 #ifdef mydebug_Seed
106  ss << "\tophimin, ophimax " << outerPhimin << " " << outerPhimax << std::endl;
107 #endif
108 
109  /* loop on outer hits*/
110  for (vector<RecHitsSortedInPhi::Hit>::const_iterator oh = outerHits.begin(); oh != outerHits.end(); ++oh) {
111  RecHitsSortedInPhi::Hit ohit = (*oh);
112 #ifdef mydebug_Seed
113  GlobalPoint oPos = ohit->globalPosition();
114 
115  ss << "\toPos " << oPos << " r " << oPos.perp() << " phi " << oPos.phi() << " cotTheta " << oPos.z() / oPos.perp()
116  << std::endl;
117 #endif
118 
119  /*Check the compatibility of the ohit with the eta of the seeding track*/
120  if (checkRZCompatibilityWithSeedTrack(ohit, *outerLayerObj.detLayer(), convRegion))
121  continue;
122 
123  /*
124  //Do I need this? it uses a compatibility that probably I wouldn't
125  //Removing for the time being
126 
127  PixelRecoRange<float> phiRange = deltaPhi( oPos.perp(), oPos.phi(), oPos.z(), nSigmaPhi*(ohit->errorGlobalRPhi()));
128  if (phiRange.empty()) continue;
129  */
130 
131  std::unique_ptr<const HitRZCompatibility> checkRZ = region.checkRZ(innerLayerObj.detLayer(), ohit);
132  if (!checkRZ) {
133 #ifdef mydebug_Seed
134  ss << "*******\nNo valid checkRZ\n*******" << std::endl;
135 #endif
136  continue;
137  }
138 
139  /*Get only the inner hits compatible with the conversion region*/
140  innerHits.clear();
141  if (!getPhiRange(innerPhimin, innerPhimax, *innerLayerObj.detLayer(), convRegion, field))
142  continue;
143  innerHitsMap.hits(innerPhimin, innerPhimax, innerHits);
144 
145 #ifdef mydebug_Seed
146  ss << "\tiphimin, iphimax " << innerPhimin << " " << innerPhimax << std::endl;
147 #endif
148 
149  /*Loop on inner hits*/
150  for (vector<RecHitsSortedInPhi::Hit>::const_iterator ih = innerHits.begin(), ieh = innerHits.end(); ih < ieh;
151  ++ih) {
152  GlobalPoint innPos = (*ih)->globalPosition();
153 
154 #ifdef mydebug_Seed
155  ss << "\tinnPos " << innPos << " r " << innPos.perp() << " phi " << innPos.phi() << " cotTheta "
156  << innPos.z() / innPos.perp() << std::endl;
157 #endif
158 
159  /*Check the compatibility of the ohit with the eta of the seeding track*/
160  if (checkRZCompatibilityWithSeedTrack(*ih, *innerLayerObj.detLayer(), convRegion))
161  continue;
162 
163  float r_reduced = std::sqrt(sqr(innPos.x() - region.origin().x()) + sqr(innPos.y() - region.origin().y()));
164  Range allowed;
165  Range hitRZ;
166  if (innerLayerObj.detLayer()->location() == barrel) {
167  allowed = checkRZ->range(r_reduced);
168  float zErr = nSigmaRZ * (*ih)->errorGlobalZ();
169  hitRZ = Range(innPos.z() - zErr, innPos.z() + zErr);
170  } else {
171  allowed = checkRZ->range(innPos.z());
172  float rErr = nSigmaRZ * (*ih)->errorGlobalR();
173  hitRZ = Range(r_reduced - rErr, r_reduced + rErr);
174  }
175  Range crossRange = allowed.intersection(hitRZ);
176 
177 #ifdef mydebug_Seed
178  ss << "\n\t\t allowed Range " << allowed.min() << " \t, " << allowed.max() << "\n\t\t hitRz Range "
179  << hitRZ.min() << " \t, " << hitRZ.max() << "\n\t\t Cross Range " << crossRange.min() << " \t, "
180  << crossRange.max() << "\n\t\t the seed track has origin " << convRegion.convPoint() << " \t cotTheta "
181  << convRegion.cotTheta() << std::endl;
182 #endif
183 
184  if (!crossRange.empty()) {
185 #ifdef mydebug_Seed
186  ss << "\n\t\t !!!!ACCEPTED!!! \n\n";
187 #endif
188  if (theMaxElement != 0 && result.size() >= maxNum) {
189  result.resize(oldSize);
190 #ifdef mydebug_Seed
191  edm::LogError("TooManySeeds") << "number of pairs exceed maximum, no pairs produced";
192  std::cout << ss.str();
193 #endif
194  return;
195  }
196  result.push_back(OrderedHitPair(*ih, ohit));
197  }
198  }
199  }
200 
201 #ifdef mydebug_Seed
202  std::cout << ss.str();
203 #endif
204 }
bool checkRZCompatibilityWithSeedTrack(const RecHitsSortedInPhi::Hit &hit, const DetLayer &layer, const ConversionRegion &convRegion)
T perp() const
Definition: PV3DBase.h:69
PixelRecoRange< float > Range
T z() const
Definition: PV3DBase.h:61
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
SeedingHitSet::ConstRecHitPointer InnerRecHit
Definition: OrderedHitPair.h:9
Log< level::Error, false > LogError
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
float cotTheta() const
constexpr double nSigmaRZ
T sqrt(T t)
Definition: SSEVec.h:19
double f[11][100]
bool getData(T &iHolder) const
Definition: EventSetup.h:122
std::vector< Hit > hits(float phiMin, float phiMax) const
SeedingHitSet::ConstRecHitPointer Hit
BaseTrackerRecHit const * Hit
bool checkBoundaries(const DetLayer &layer, const ConversionRegion &convRegion, float maxSearchR, float maxSearchZ)
PixelRecoRange< T > intersection(const PixelRecoRange< T > &r) const
Square< F >::type sqr(const F &f)
Definition: Square.h:14
bool getPhiRange(float &Phimin, float &Phimax, const DetLayer &layer, const ConversionRegion &convRegion, const MagneticField &field)
GlobalPoint convPoint() const
SeedingHitSet::ConstRecHitPointer OuterRecHit
Definition: OrderedHitPair.h:8
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > theFieldToken

Member Data Documentation

◆ ss

std::stringstream HitPairGeneratorFromLayerPairForPhotonConversion::ss
private

◆ theFieldToken

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> HitPairGeneratorFromLayerPairForPhotonConversion::theFieldToken
private

Definition at line 65 of file HitPairGeneratorFromLayerPairForPhotonConversion.h.

Referenced by hitPairs().

◆ theInnerLayer

const unsigned int HitPairGeneratorFromLayerPairForPhotonConversion::theInnerLayer
private

Definition at line 63 of file HitPairGeneratorFromLayerPairForPhotonConversion.h.

Referenced by hitPairs().

◆ theLayerCache

LayerCacheType& HitPairGeneratorFromLayerPairForPhotonConversion::theLayerCache
private

Definition at line 61 of file HitPairGeneratorFromLayerPairForPhotonConversion.h.

Referenced by hitPairs().

◆ theMaxElement

const unsigned int HitPairGeneratorFromLayerPairForPhotonConversion::theMaxElement
private

Definition at line 64 of file HitPairGeneratorFromLayerPairForPhotonConversion.h.

Referenced by hitPairs().

◆ theOuterLayer

const unsigned int HitPairGeneratorFromLayerPairForPhotonConversion::theOuterLayer
private

Definition at line 62 of file HitPairGeneratorFromLayerPairForPhotonConversion.h.

Referenced by hitPairs().