CMS 3D CMS Logo

CandidateSeededTrackingRegionsProducer.h
Go to the documentation of this file.
1 #ifndef CandidateSeededTrackingRegionsProducer_h
2 #define CandidateSeededTrackingRegionsProducer_h
3 
7 
19 
46 {
47 public:
48 
50 
53  {
54  edm::ParameterSet regPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
55 
56  // operation mode
57  std::string modeString = regPSet.getParameter<std::string>("mode");
58  if (modeString == "BeamSpotFixed") m_mode = BEAM_SPOT_FIXED;
59  else if (modeString == "BeamSpotSigma") m_mode = BEAM_SPOT_SIGMA;
60  else if (modeString == "VerticesFixed") m_mode = VERTICES_FIXED;
61  else if (modeString == "VerticesSigma") m_mode = VERTICES_SIGMA;
62  else edm::LogError ("CandidateSeededTrackingRegionsProducer")<<"Unknown mode string: "<<modeString;
63 
64  // basic inputs
65  token_input = iC.consumes<reco::CandidateView>(regPSet.getParameter<edm::InputTag>("input"));
66  m_maxNRegions = regPSet.getParameter<int>("maxNRegions");
67  token_beamSpot = iC.consumes<reco::BeamSpot>(regPSet.getParameter<edm::InputTag>("beamSpot"));
68  m_maxNVertices = 1;
70  {
71  token_vertex = iC.consumes<reco::VertexCollection>(regPSet.getParameter<edm::InputTag>("vertexCollection"));
72  m_maxNVertices = regPSet.getParameter<int>("maxNVertices");
73  }
74 
75  // RectangularEtaPhiTrackingRegion parameters:
76  m_ptMin = regPSet.getParameter<double>("ptMin");
77  m_originRadius = regPSet.getParameter<double>("originRadius");
78  m_zErrorBeamSpot = regPSet.getParameter<double>("zErrorBeamSpot");
79  m_deltaEta = regPSet.getParameter<double>("deltaEta");
80  m_deltaPhi = regPSet.getParameter<double>("deltaPhi");
81  m_precise = regPSet.getParameter<bool>("precise");
83  if(m_whereToUseMeasurementTracker != RectangularEtaPhiTrackingRegion::UseMeasurementTracker::kNever) {
84  token_measurementTracker = iC.consumes<MeasurementTrackerEvent>(regPSet.getParameter<edm::InputTag>("measurementTrackerName"));
85  }
86  m_searchOpt = false;
87  if (regPSet.exists("searchOpt")) m_searchOpt = regPSet.getParameter<bool>("searchOpt");
88 
89  // mode-dependent z-halflength of tracking regions
90  if (m_mode == VERTICES_SIGMA) m_nSigmaZVertex = regPSet.getParameter<double>("nSigmaZVertex");
91  if (m_mode == VERTICES_FIXED) m_zErrorVetex = regPSet.getParameter<double>("zErrorVetex");
92  m_nSigmaZBeamSpot = -1.;
93  if (m_mode == BEAM_SPOT_SIGMA)
94  {
95  m_nSigmaZBeamSpot = regPSet.getParameter<double>("nSigmaZBeamSpot");
96  if (m_nSigmaZBeamSpot < 0.)
97  edm::LogError ("CandidateSeededTrackingRegionsProducer")<<"nSigmaZBeamSpot must be positive for BeamSpotSigma mode!";
98  }
99  }
100 
102 
105 
106  desc.add<std::string>("mode", "BeamSpotFixed");
107 
108  desc.add<edm::InputTag>("input", edm::InputTag(""));
109  desc.add<int>("maxNRegions", 10);
110  desc.add<edm::InputTag>("beamSpot", edm::InputTag("hltOnlineBeamSpot"));
111  desc.add<edm::InputTag>("vertexCollection", edm::InputTag("hltPixelVertices"));
112  desc.add<int>("maxNVertices", 1);
113 
114  desc.add<double>("ptMin", 0.9);
115  desc.add<double>("originRadius", 0.2);
116  desc.add<double>("zErrorBeamSpot", 24.2);
117  desc.add<double>("deltaEta", 0.5);
118  desc.add<double>("deltaPhi", 0.5);
119  desc.add<bool>("precise", true);
120 
121  desc.add<double>("nSigmaZVertex", 3.);
122  desc.add<double>("zErrorVetex", 0.2);
123  desc.add<double>("nSigmaZBeamSpot", 4.);
124 
125  desc.add<std::string>("whereToUseMeasurementTracker", "ForSiStrips");
126  desc.add<edm::InputTag>("measurementTrackerName", edm::InputTag(""));
127 
128  desc.add<bool>("searchOpt", false);
129 
130  // Only for backwards-compatibility
131  edm::ParameterSetDescription descRegion;
132  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
133 
134  descriptions.add("seededTrackingRegionsFromBeamSpotFixedZLength", descRegion);
135  }
136 
137 
138  std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& e, const edm::EventSetup& es) const override
139  {
140  std::vector<std::unique_ptr<TrackingRegion> > result;
141 
142  // pick up the candidate objects of interest
144  e.getByToken( token_input, objects );
145  size_t n_objects = objects->size();
146  if (n_objects == 0) return result;
147 
148  // always need the beam spot (as a fall back strategy for vertex modes)
150  e.getByToken( token_beamSpot, bs );
151  if( !bs.isValid() ) return result;
152 
153  // this is a default origin for all modes
154  GlobalPoint default_origin( bs->x0(), bs->y0(), bs->z0() );
155 
156  // vector of origin & halfLength pairs:
157  std::vector< std::pair< GlobalPoint, float > > origins;
158 
159  // fill the origins and halfLengths depending on the mode
161  {
162  origins.push_back( std::make_pair(
163  default_origin,
165  ));
166  }
167  else if (m_mode == VERTICES_FIXED || m_mode == VERTICES_SIGMA)
168  {
170  e.getByToken( token_vertex, vertices );
171  int n_vert = 0;
172  for (reco::VertexCollection::const_iterator v = vertices->begin(); v != vertices->end() && n_vert < m_maxNVertices; ++v)
173  {
174  if ( v->isFake() || !v->isValid() ) continue;
175 
176  origins.push_back( std::make_pair(
177  GlobalPoint( v->x(), v->y(), v->z() ),
179  ));
180  ++n_vert;
181  }
182  // no-vertex fall-back case:
183  if (origins.empty())
184  {
185  origins.push_back( std::make_pair(
186  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  // objects of interest (we expect that the collection was sorted in decreasing pt order)
201  int n_regions = 0;
202  for(size_t i = 0; i < n_objects && n_regions < m_maxNRegions; ++i )
203  {
204  const reco::Candidate & object = (*objects)[i];
205  GlobalVector direction( object.momentum().x(), object.momentum().y(), object.momentum().z() );
206 
207  for (size_t j=0; j<origins.size() && n_regions < m_maxNRegions; ++j)
208  {
209  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(
210  direction,
211  origins[j].first,
212  m_ptMin,
214  origins[j].second,
215  m_deltaEta,
216  m_deltaPhi,
218  m_precise,
219  measurementTracker,
221  ));
222  ++n_regions;
223  }
224  }
225  //std::cout<<"n_seeded_regions = "<<n_regions<<std::endl;
226  edm::LogInfo ("CandidateSeededTrackingRegionsProducer") << "produced "<<n_regions<<" regions";
227 
228  return result;
229  }
230 
231 private:
232 
233  Mode m_mode;
234 
240 
241  float m_ptMin;
244  float m_deltaEta;
245  float m_deltaPhi;
246  bool m_precise;
250 
254 };
255 
256 #endif
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:68
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
CandidateSeededTrackingRegionsProducer(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
bool exists(std::string const &parameterName) const
checks if a parameter exists
size_type size() const
edm::EDGetTokenT< reco::VertexCollection > token_vertex
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
U second(std::pair< T, U > const &p)
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override
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:74
double sigmaZ() const
sigma z
Definition: BeamSpot.h:80
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double y0() const
y coordinate
Definition: BeamSpot.h:66
bool isUninitialized() const
Definition: EDGetToken.h:70
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double x0() const
x coordinate
Definition: BeamSpot.h:64