CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
AreaSeededTrackingRegionsBuilder.h
Go to the documentation of this file.
1 #ifndef RecoTracker_TkTrackingRegions_AreaSeededTrackingRegionsBuilder_h
2 #define RecoTracker_TkTrackingRegions_AreaSeededTrackingRegionsBuilder_h
3 
5 
14 
17 
20 
21 #include "TrackingSeedCandidates.h"
23 public:
24  using Origin = std::pair<GlobalPoint, float>; // (origin, half-length in z)
25  using Origins = std::vector<Origin>;
26 
27  class Area {
28  public:
29  Area() {}
30 
31  // phimin and phimax, and hence xmin+xmax and ymin+ymax are
32  // ordered by which way one goes around the unit circle, so it may
33  // happen that actually phimax < phimin
34  Area(float rmin, float rmax, float phimin, float phimax, float zmin, float zmax) : m_zmin(zmin), m_zmax(zmax) {
35  auto cosphimin = std::cos(phimin);
36  auto sinphimin = std::sin(phimin);
37  auto cosphimax = std::cos(phimax);
38  auto sinphimax = std::sin(phimax);
39 
40  m_x_rmin_phimin = rmin * cosphimin;
41  m_x_rmin_phimax = rmin * cosphimax;
42  m_x_rmax_phimin = rmax * cosphimin;
43  m_x_rmax_phimax = rmax * cosphimax;
44 
45  m_y_rmin_phimin = rmin * sinphimin;
46  m_y_rmin_phimax = rmin * sinphimax;
47  m_y_rmax_phimin = rmax * sinphimin;
48  m_y_rmax_phimax = rmax * sinphimax;
49  }
50 
51  float x_rmin_phimin() const { return m_x_rmin_phimin; }
52  float x_rmin_phimax() const { return m_x_rmin_phimax; }
53  float x_rmax_phimin() const { return m_x_rmax_phimin; }
54  float x_rmax_phimax() const { return m_x_rmax_phimax; }
55  float y_rmin_phimin() const { return m_y_rmin_phimin; }
56  float y_rmin_phimax() const { return m_y_rmin_phimax; }
57  float y_rmax_phimin() const { return m_y_rmax_phimin; }
58  float y_rmax_phimax() const { return m_y_rmax_phimax; }
59 
60  float zmin() const { return m_zmin; }
61  float zmax() const { return m_zmax; }
62 
63  private:
64  // all of these are in global coordinates
65  float m_x_rmin_phimin = 0;
66  float m_x_rmin_phimax = 0;
67  float m_x_rmax_phimin = 0;
68  float m_x_rmax_phimax = 0;
69 
70  float m_y_rmin_phimin = 0;
71  float m_y_rmin_phimax = 0;
72  float m_y_rmax_phimin = 0;
73  float m_y_rmax_phimax = 0;
74 
75  float m_zmin = 0;
76  float m_zmax = 0;
77  };
78 
79  class Builder {
80  public:
82  const MagneticField* field,
84  : m_conf(conf), m_field(field), m_msmaker(msmaker) {}
85  ~Builder() = default;
86 
89 
90  std::vector<std::unique_ptr<TrackingRegion> > regions(const Origins& origins, const std::vector<Area>& areas) const;
91  std::unique_ptr<TrackingRegion> region(const Origin& origin, const std::vector<Area>& areas) const;
92  std::unique_ptr<TrackingRegion> region(const Origin& origin, const edm::VecArray<Area, 2>& areas) const;
93 
94  private:
95  template <typename T>
96  std::unique_ptr<TrackingRegion> regionImpl(const Origin& origin, const T& areas) const;
97 
100  const MagneticField* m_field = nullptr;
103  };
104 
106  : AreaSeededTrackingRegionsBuilder(regPSet, iC) {}
109 
111 
112  Builder beginEvent(const edm::Event& e, const edm::EventSetup& es) const;
113 
114 private:
115  std::vector<Area> m_areas;
117  float m_extraPhi;
118  float m_extraEta;
119  float m_ptMin;
121  bool m_precise;
127 };
128 
129 #endif
static void fillDescriptions(edm::ParameterSetDescription &desc)
const AreaSeededTrackingRegionsBuilder * m_conf
std::pair< const reco::CandidateView *, std::pair< float, float > > Objects
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void setCandidates(const TrackingSeedCandidates::Objects cands)
Builder beginEvent(const edm::Event &e, const edm::EventSetup &es) const
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
Builder(const AreaSeededTrackingRegionsBuilder *conf, const MagneticField *field, const MultipleScatteringParametrisationMaker *msmaker)
std::unique_ptr< TrackingRegion > region(const Origin &origin, const std::vector< Area > &areas) const
Area(float rmin, float rmax, float phimin, float phimax, float zmin, float zmax)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > token_field
std::vector< std::unique_ptr< TrackingRegion > > regions(const Origins &origins, const std::vector< Area > &areas) const
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > token_msmaker
std::unique_ptr< TrackingRegion > regionImpl(const Origin &origin, const T &areas) const
void setMeasurementTracker(const MeasurementTrackerEvent *mte)
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
long double T
AreaSeededTrackingRegionsBuilder(const edm::ParameterSet &regPSet, edm::ConsumesCollector &&iC)
const MultipleScatteringParametrisationMaker * m_msmaker