CMS 3D CMS Logo

PointSeededTrackingRegionsProducer.h
Go to the documentation of this file.
1 #ifndef PointSeededTrackingRegionsProducer_h
2 #define PointSeededTrackingRegionsProducer_h
3 
4 
7 
20 
21 #include <TMath.h>
22 #include <TLorentzVector.h>
23 
49 {
50 public:
51 
53 
56  {
57  edm::ParameterSet regPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
58 
59  // operation mode
60  std::string modeString = regPSet.getParameter<std::string>("mode");
61  if (modeString == "BeamSpotFixed") m_mode = BEAM_SPOT_FIXED;
62  else if (modeString == "BeamSpotSigma") m_mode = BEAM_SPOT_SIGMA;
63  else if (modeString == "VerticesFixed") m_mode = VERTICES_FIXED;
64  else if (modeString == "VerticesSigma") m_mode = VERTICES_SIGMA;
65  else edm::LogError ("PointSeededTrackingRegionsProducer")<<"Unknown mode string: "<<modeString;
66 
67  // basic inputsi
69  etaPoints = points.getParameter<std::vector<double>>("eta");
70  phiPoints = points.getParameter<std::vector<double>>("phi");
71  if (!(etaPoints.size() == phiPoints.size())) throw edm::Exception(edm::errors::Configuration) << "The parameters 'eta' and 'phi' must have the same size";;
72  m_maxNRegions = regPSet.getParameter<int>("maxNRegions");
73  token_beamSpot = iC.consumes<reco::BeamSpot>(regPSet.getParameter<edm::InputTag>("beamSpot"));
74  m_maxNVertices = 1;
76  {
77  token_vertex = iC.consumes<reco::VertexCollection>(regPSet.getParameter<edm::InputTag>("vertexCollection"));
78  m_maxNVertices = regPSet.getParameter<int>("maxNVertices");
79  }
80 
81  // RectangularEtaPhiTrackingRegion parameters:
82  m_ptMin = regPSet.getParameter<double>("ptMin");
83  m_originRadius = regPSet.getParameter<double>("originRadius");
84  m_zErrorBeamSpot = regPSet.getParameter<double>("zErrorBeamSpot");
85  m_deltaEta = regPSet.getParameter<double>("deltaEta");
86  m_deltaPhi = regPSet.getParameter<double>("deltaPhi");
87  m_precise = regPSet.getParameter<bool>("precise");
89  if(m_whereToUseMeasurementTracker != RectangularEtaPhiTrackingRegion::UseMeasurementTracker::kNever) {
90  token_measurementTracker = iC.consumes<MeasurementTrackerEvent>(regPSet.getParameter<edm::InputTag>("measurementTrackerName"));
91  }
92  m_searchOpt = false;
93  if (regPSet.exists("searchOpt")) m_searchOpt = regPSet.getParameter<bool>("searchOpt");
94 
95  // mode-dependent z-halflength of tracking regions
96  if (m_mode == VERTICES_SIGMA) m_nSigmaZVertex = regPSet.getParameter<double>("nSigmaZVertex");
97  if (m_mode == VERTICES_FIXED) m_zErrorVetex = regPSet.getParameter<double>("zErrorVetex");
98  m_nSigmaZBeamSpot = -1.;
99  if (m_mode == BEAM_SPOT_SIGMA)
100  {
101  m_nSigmaZBeamSpot = regPSet.getParameter<double>("nSigmaZBeamSpot");
102  if (m_nSigmaZBeamSpot < 0.)
103  edm::LogError ("PointSeededTrackingRegionsProducer")<<"nSigmaZBeamSpot must be positive for BeamSpotSigma mode!";
104  }
105  }
106 
108 
111 
112  edm::ParameterSetDescription descPoints;
113  descPoints.add<std::vector<double>> ("eta", {0.} );
114  descPoints.add<std::vector<double>> ("phi", {0.} );
115  desc.add<edm::ParameterSetDescription>("points", descPoints);
116 
117  desc.add<std::string>("mode", "BeamSpotFixed");
118  desc.add<int>("maxNRegions", 10);
119  desc.add<edm::InputTag>("beamSpot", edm::InputTag("hltOnlineBeamSpot"));
120  desc.add<edm::InputTag>("vertexCollection", edm::InputTag("hltPixelVertices"));
121  desc.add<int>("maxNVertices", 1);
122 
123  desc.add<double>("ptMin", 0.9);
124  desc.add<double>("originRadius", 0.2);
125  desc.add<double>("zErrorBeamSpot", 24.2);
126  desc.add<double>("deltaEta", 0.5);
127  desc.add<double>("deltaPhi", 0.5);
128  desc.add<bool>("precise", true);
129 
130  desc.add<double>("nSigmaZVertex", 3.);
131  desc.add<double>("zErrorVetex", 0.2);
132  desc.add<double>("nSigmaZBeamSpot", 4.);
133 
134  desc.add<std::string>("whereToUseMeasurementTracker", "ForSiStrips");
135  desc.add<edm::InputTag>("measurementTrackerName", edm::InputTag(""));
136 
137  desc.add<bool>("searchOpt", false);
138 
139  // Only for backwards-compatibility
140  edm::ParameterSetDescription descRegion;
141  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
142  //edm::ParameterSetDescription descPoint;
143  //descPoint.add<edm::ParameterSetDescription>("point_input", desc);
144 
145 
146  descriptions.add("pointSeededTrackingRegion", descRegion);
147  }
148 
149 
150 
151  std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& e, const edm::EventSetup& es) const override
152  {
153  std::vector<std::unique_ptr<TrackingRegion> > result;
154 
155  // always need the beam spot (as a fall back strategy for vertex modes)
157  e.getByToken( token_beamSpot, bs );
158  if( !bs.isValid() ) return result;
159 
160  // this is a default origin for all modes
161  GlobalPoint default_origin( bs->x0(), bs->y0(), bs->z0() );
162 
163  // vector of origin & halfLength pairs:
164  std::vector< std::pair< GlobalPoint, float > > origins;
165 
166  // fill the origins and halfLengths depending on the mode
168  origins.push_back( std::make_pair( default_origin,
170  ));
171  } else if (m_mode == VERTICES_FIXED || m_mode == VERTICES_SIGMA) {
173  e.getByToken( token_vertex, vertices );
174  int n_vert = 0;
175  for (reco::VertexCollection::const_iterator iv = vertices->begin(), ev = vertices->end();
176  iv != ev && n_vert < m_maxNVertices; ++iv) {
177  if ( iv->isFake() || !iv->isValid() ) continue;
178 
179  origins.push_back( std::make_pair( GlobalPoint( iv->x(), iv->y(), iv->z() ),
180  (m_mode == VERTICES_FIXED) ? m_zErrorVetex : m_nSigmaZVertex*iv->zError()
181  ));
182  ++n_vert;
183  }
184  // no-vertex fall-back case:
185  if ( origins.empty() ) {
186  origins.push_back( std::make_pair( default_origin,
188  ));
189  }
190  }
191 
192  const MeasurementTrackerEvent *measurementTracker = nullptr;
196  measurementTracker = hmte.product();
197  }
198 
199  // create tracking regions (maximum MaxNRegions of them) in directions of the
200  // points of interest
201  size_t n_points = etaPoints.size();
202  int n_regions = 0;
203  for(size_t i = 0; i < n_points && n_regions < m_maxNRegions; ++i ) {
204 
205  double x = std::cos(phiPoints[i]);
206  double y = std::sin(phiPoints[i]);
207  double theta = 2*std::atan(std::exp(-etaPoints[i]));
208  double z = 1./std::tan(theta);
209 
210  GlobalVector direction( x,y,z );
211 
212  for (size_t j=0; j<origins.size() && n_regions < m_maxNRegions; ++j) {
213 
214  result.push_back( std::make_unique<RectangularEtaPhiTrackingRegion>(
215  direction, // GlobalVector
216  origins[j].first, // GlobalPoint
217  m_ptMin,
219  origins[j].second,
220  m_deltaEta,
221  m_deltaPhi,
223  m_precise,
224  measurementTracker,
226  ));
227  ++n_regions;
228  }
229  }
230  edm::LogInfo ("PointSeededTrackingRegionsProducer") << "produced "<<n_regions<<" regions";
231 
232  return result;
233  }
234 
235 private:
236 
237  Mode m_mode;
238 
243 
244  std::vector<double> etaPoints;
245  std::vector<double> phiPoints;
246 
247  float m_ptMin;
250  float m_deltaEta;
251  float m_deltaPhi;
252  bool m_precise;
256 
260 };
261 
262 #endif
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:68
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
PointSeededTrackingRegionsProducer(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Geom::Theta< T > theta() const
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool ev
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
U second(std::pair< T, U > const &p)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::VertexCollection > token_vertex
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
double z0Error() const
error on z
Definition: BeamSpot.h:94
T const * product() const
Definition: Handle.h:81
double sigmaZ() const
sigma z
Definition: BeamSpot.h:80
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double y0() const
y coordinate
Definition: BeamSpot.h:66
bool isUninitialized() const
Definition: EDGetToken.h:73
double x0() const
x coordinate
Definition: BeamSpot.h:64