CMS 3D CMS Logo

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