CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 (regPSet.exists("measurementTrackerName"))
83  {
84  // FIXME: when next time altering the configuration of this
85  // class, please change the types of the following parameters:
86  // - whereToUseMeasurementTracker to at least int32 or to a string
87  // corresponding to the UseMeasurementTracker enumeration
88  // - measurementTrackerName to InputTag
89  if (regPSet.exists("whereToUseMeasurementTracker"))
90  m_whereToUseMeasurementTracker = RectangularEtaPhiTrackingRegion::doubleToUseMeasurementTracker(regPSet.getParameter<double>("whereToUseMeasurementTracker"));
91  if(m_whereToUseMeasurementTracker != RectangularEtaPhiTrackingRegion::UseMeasurementTracker::kNever)
92  token_measurementTracker = iC.consumes<MeasurementTrackerEvent>(regPSet.getParameter<std::string>("measurementTrackerName"));
93  }
94  m_searchOpt = false;
95  if (regPSet.exists("searchOpt")) m_searchOpt = regPSet.getParameter<bool>("searchOpt");
96 
97  // mode-dependent z-halflength of tracking regions
98  if (m_mode == VERTICES_SIGMA) m_nSigmaZVertex = regPSet.getParameter<double>("nSigmaZVertex");
99  if (m_mode == VERTICES_FIXED) m_zErrorVetex = regPSet.getParameter<double>("zErrorVetex");
100  m_nSigmaZBeamSpot = -1.;
101  if (m_mode == BEAM_SPOT_SIGMA)
102  {
103  m_nSigmaZBeamSpot = regPSet.getParameter<double>("nSigmaZBeamSpot");
104  if (m_nSigmaZBeamSpot < 0.)
105  edm::LogError ("CandidateSeededTrackingRegionsProducer")<<"nSigmaZBeamSpot must be positive for BeamSpotSigma mode!";
106  }
107  }
108 
110 
111 
112  virtual std::vector<TrackingRegion* > regions(const edm::Event& e, const edm::EventSetup& es) const
113  {
114  std::vector<TrackingRegion* > result;
115 
116  // pick up the candidate objects of interest
118  e.getByToken( token_input, objects );
119  size_t n_objects = objects->size();
120  if (n_objects == 0) return result;
121 
122  // always need the beam spot (as a fall back strategy for vertex modes)
124  e.getByToken( token_beamSpot, bs );
125  if( !bs.isValid() ) return result;
126 
127  // this is a default origin for all modes
128  GlobalPoint default_origin( bs->x0(), bs->y0(), bs->z0() );
129 
130  // vector of origin & halfLength pairs:
131  std::vector< std::pair< GlobalPoint, float > > origins;
132 
133  // fill the origins and halfLengths depending on the mode
135  {
136  origins.push_back( std::make_pair(
137  default_origin,
139  ));
140  }
141  else if (m_mode == VERTICES_FIXED || m_mode == VERTICES_SIGMA)
142  {
144  e.getByToken( token_vertex, vertices );
145  int n_vert = 0;
146  for (reco::VertexCollection::const_iterator v = vertices->begin(); v != vertices->end() && n_vert < m_maxNVertices; ++v)
147  {
148  if ( v->isFake() || !v->isValid() ) continue;
149 
150  origins.push_back( std::make_pair(
151  GlobalPoint( v->x(), v->y(), v->z() ),
153  ));
154  ++n_vert;
155  }
156  // no-vertex fall-back case:
157  if (origins.empty())
158  {
159  origins.push_back( std::make_pair(
160  default_origin,
161  (m_nSigmaZBeamSpot > 0.) ? m_nSigmaZBeamSpot*bs->z0Error() : m_zErrorBeamSpot
162  ));
163  }
164  }
165 
166  const MeasurementTrackerEvent *measurementTracker = nullptr;
170  measurementTracker = hmte.product();
171  }
172 
173  // create tracking regions (maximum MaxNRegions of them) in directions of the
174  // objects of interest (we expect that the collection was sorted in decreasing pt order)
175  int n_regions = 0;
176  for(size_t i = 0; i < n_objects && n_regions < m_maxNRegions; ++i )
177  {
178  const reco::Candidate & object = (*objects)[i];
179  GlobalVector direction( object.momentum().x(), object.momentum().y(), object.momentum().z() );
180 
181  for (size_t j=0; j<origins.size() && n_regions < m_maxNRegions; ++j)
182  {
183  result.push_back( new RectangularEtaPhiTrackingRegion(
184  direction,
185  origins[j].first,
186  m_ptMin,
188  origins[j].second,
189  m_deltaEta,
190  m_deltaPhi,
192  m_precise,
193  measurementTracker,
195  ));
196  ++n_regions;
197  }
198  }
199  //std::cout<<"n_seeded_regions = "<<n_regions<<std::endl;
200  edm::LogInfo ("CandidateSeededTrackingRegionsProducer") << "produced "<<n_regions<<" regions";
201 
202  return result;
203  }
204 
205 private:
206 
208 
214 
215  float m_ptMin;
218  float m_deltaEta;
219  float m_deltaPhi;
220  bool m_precise;
224 
228 };
229 
230 #endif
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
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
edm::EDGetTokenT< reco::VertexCollection > token_vertex
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
float float float z
virtual std::vector< TrackingRegion * > regions(const edm::Event &e, const edm::EventSetup &es) const
U second(std::pair< T, U > const &p)
tuple result
Definition: query.py:137
int j
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
bool first
Definition: L1TdeRCT.cc:75
bool isValid() const
Definition: HandleBase.h:76
tuple conf
Definition: dbtoconf.py:185
static UseMeasurementTracker doubleToUseMeasurementTracker(double value)
T const * product() const
Definition: Handle.h:81
bool isUninitialized() const
Definition: EDGetToken.h:71
Definition: DDAxes.h:10