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");
83  if (regPSet.exists("measurementTrackerName"))
84  {
85  m_measurementTrackerName = regPSet.getParameter<std::string>("measurementTrackerName");
86  if (regPSet.exists("whereToUseMeasurementTracker"))
87  m_whereToUseMeasurementTracker = regPSet.getParameter<double>("whereToUseMeasurementTracker");
88  }
89  m_searchOpt = false;
90  if (regPSet.exists("searchOpt")) m_searchOpt = regPSet.getParameter<bool>("searchOpt");
91 
92  // mode-dependent z-halflength of tracking regions
93  if (m_mode == VERTICES_SIGMA) m_nSigmaZVertex = regPSet.getParameter<double>("nSigmaZVertex");
94  if (m_mode == VERTICES_FIXED) m_zErrorVetex = regPSet.getParameter<double>("zErrorVetex");
95  m_nSigmaZBeamSpot = -1.;
96  if (m_mode == BEAM_SPOT_SIGMA)
97  {
98  m_nSigmaZBeamSpot = regPSet.getParameter<double>("nSigmaZBeamSpot");
99  if (m_nSigmaZBeamSpot < 0.)
100  edm::LogError ("CandidateSeededTrackingRegionsProducer")<<"nSigmaZBeamSpot must be positive for BeamSpotSigma mode!";
101  }
102  }
103 
105 
106 
107  virtual std::vector<TrackingRegion* > regions(const edm::Event& e, const edm::EventSetup& es) const
108  {
109  std::vector<TrackingRegion* > result;
110 
111  // pick up the candidate objects of interest
113  e.getByToken( token_input, objects );
114  size_t n_objects = objects->size();
115  if (n_objects == 0) return result;
116 
117  // always need the beam spot (as a fall back strategy for vertex modes)
119  e.getByToken( token_beamSpot, bs );
120  if( !bs.isValid() ) return result;
121 
122  // this is a default origin for all modes
123  GlobalPoint default_origin( bs->x0(), bs->y0(), bs->z0() );
124 
125  // vector of origin & halfLength pairs:
126  std::vector< std::pair< GlobalPoint, float > > origins;
127 
128  // fill the origins and halfLengths depending on the mode
130  {
131  origins.push_back( std::make_pair(
132  default_origin,
134  ));
135  }
136  else if (m_mode == VERTICES_FIXED || m_mode == VERTICES_SIGMA)
137  {
139  e.getByToken( token_vertex, vertices );
140  int n_vert = 0;
141  for (reco::VertexCollection::const_iterator v = vertices->begin(); v != vertices->end() && n_vert < m_maxNVertices; ++v)
142  {
143  if ( v->isFake() || !v->isValid() ) continue;
144 
145  origins.push_back( std::make_pair(
146  GlobalPoint( v->x(), v->y(), v->z() ),
148  ));
149  ++n_vert;
150  }
151  // no-vertex fall-back case:
152  if (origins.empty())
153  {
154  origins.push_back( std::make_pair(
155  default_origin,
156  (m_nSigmaZBeamSpot > 0.) ? m_nSigmaZBeamSpot*bs->z0Error() : m_zErrorBeamSpot
157  ));
158  }
159  }
160 
161  // create tracking regions (maximum MaxNRegions of them) in directions of the
162  // objects of interest (we expect that the collection was sorted in decreasing pt order)
163  int n_regions = 0;
164  for(size_t i = 0; i < n_objects && n_regions < m_maxNRegions; ++i )
165  {
166  const reco::Candidate & object = (*objects)[i];
167  GlobalVector direction( object.momentum().x(), object.momentum().y(), object.momentum().z() );
168 
169  for (size_t j=0; j<origins.size() && n_regions < m_maxNRegions; ++j)
170  {
171  result.push_back( new RectangularEtaPhiTrackingRegion(
172  direction,
173  origins[j].first,
174  m_ptMin,
176  origins[j].second,
177  m_deltaEta,
178  m_deltaPhi,
180  m_precise,
183  ));
184  ++n_regions;
185  }
186  }
187  //std::cout<<"n_seeded_regions = "<<n_regions<<std::endl;
188  edm::LogInfo ("CandidateSeededTrackingRegionsProducer") << "produced "<<n_regions<<" regions";
189 
190  return result;
191  }
192 
193 private:
194 
196 
202 
203  float m_ptMin;
206  float m_deltaEta;
207  float m_deltaPhi;
208  bool m_precise;
212 
216 };
217 
218 #endif
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
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
bool first
Definition: L1TdeRCT.cc:79
bool isValid() const
Definition: HandleBase.h:76
tuple conf
Definition: dbtoconf.py:185
Definition: DDAxes.h:10