CMS 3D CMS Logo

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