CMS 3D CMS Logo

CAHitTripletGenerator.h
Go to the documentation of this file.
1 #ifndef RecoTracker_PixelSeeding_CAHitTripletGenerator_h
2 #define RecoTracker_PixelSeeding_CAHitTripletGenerator_h
3 
7 
12 
17 
19 
22 class TrackingRegion;
24 
25 namespace edm {
26  class Event;
27  class EventSetup;
29 } // namespace edm
30 
32 public:
34 
35  static constexpr unsigned int minLayers = 3;
37 
38 public:
41 
42  ~CAHitTripletGenerator() = default;
43 
45  static const char* fillDescriptionsLabel() { return "caHitTriplet"; }
46 
47  void initEvent(const edm::Event& ev, const edm::EventSetup& es);
48 
49  void hitNtuplets(const IntermediateHitDoublets& regionDoublets,
50  std::vector<OrderedHitSeeds>& result,
52 
53 private:
55 
56  std::unique_ptr<SeedComparitor> theComparitor;
57 
59  const MagneticField* theField = nullptr;
60 
62  public:
63  QuantityDependsPtEval(float v1, float v2, float c1, float c2)
64  : value1_(v1), value2_(v2), curvature1_(c1), curvature2_(c2) {}
65 
66  float value(float curvature) const {
67  if (value1_ == value2_) // not enabled
68  return value1_;
69 
70  if (curvature1_ < curvature)
71  return value1_;
74  return value2_;
75  }
76 
77  private:
78  const float value1_;
79  const float value2_;
80  const float curvature1_;
81  const float curvature2_;
82  };
83 
84  // Linear interpolation (in curvature) between value1 at pt1 and
85  // value2 at pt2. If disabled, value2 is given (the point is to
86  // allow larger/smaller values of the quantity at low pt, so it
87  // makes more sense to have the high-pt value as the default).
88 
90  public:
92  : value1_(pset.getParameter<double>("value1")),
93  value2_(pset.getParameter<double>("value2")),
94  pt1_(pset.getParameter<double>("pt1")),
95  pt2_(pset.getParameter<double>("pt2")),
96  enabled_(pset.getParameter<bool>("enabled")) {
97  if (enabled_ && pt1_ >= pt2_)
98  throw cms::Exception("Configuration") << "CAHitTripletGenerator::QuantityDependsPt: pt1 (" << pt1_
99  << ") needs to be smaller than pt2 (" << pt2_ << ")";
100  if (pt1_ <= 0)
101  throw cms::Exception("Configuration")
102  << "CAHitTripletGenerator::QuantityDependsPt: pt1 needs to be > 0; is " << pt1_;
103  if (pt2_ <= 0)
104  throw cms::Exception("Configuration")
105  << "CAHitTripletGenerator::QuantityDependsPt: pt2 needs to be > 0; is " << pt2_;
106  }
107 
109  if (enabled_) {
111  value2_,
113  PixelRecoUtilities::curvature(1.f / pt2_, field));
114  }
115  return QuantityDependsPtEval(value2_, value2_, 0.f, 0.f);
116  }
117 
118  private:
119  const float value1_;
120  const float value2_;
121  const float pt1_;
122  const float pt2_;
123  const bool enabled_;
124  };
125 
127 
130 
133  const float caHardPtCut = 0.f;
134 };
135 
136 #endif
QuantityDependsPt(const edm::ParameterSet &pset)
const QuantityDependsPt maxChi2
static const char * fillDescriptionsLabel()
CAHitTripletGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
T curvature(T InversePt, const MagneticField &field)
Definition: CACut.h:16
const MagneticField * theField
QuantityDependsPtEval evaluator(const MagneticField &field) const
LayerHitMapCache LayerCacheType
double f[11][100]
void initEvent(const edm::Event &ev, const edm::EventSetup &es)
QuantityDependsPtEval(float v1, float v2, float c1, float c2)
static constexpr unsigned int minLayers
~CAHitTripletGenerator()=default
HLT enums.
static void fillDescriptions(edm::ParameterSetDescription &desc)
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > theFieldToken
std::unique_ptr< SeedComparitor > theComparitor
void hitNtuplets(const IntermediateHitDoublets &regionDoublets, std::vector< OrderedHitSeeds > &result, const SeedingLayerSetsHits &layers)