CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PointSeededTrackingRegionsProducer.h
Go to the documentation of this file.
1 #ifndef PointSeededTrackingRegionsProducer_h
2 #define PointSeededTrackingRegionsProducer_h
3 
4 
7 
18 
19 #include <TMath.h>
20 #include <TLorentzVector.h>
21 
47 
49 {
50 public:
51 
53 
56  {
57  edm::ParameterSet regPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
58 
59  // operation mode
60  std::string modeString = regPSet.getParameter<std::string>("mode");
61  if (modeString == "BeamSpotFixed") m_mode = BEAM_SPOT_FIXED;
62  else if (modeString == "BeamSpotSigma") m_mode = BEAM_SPOT_SIGMA;
63  else if (modeString == "VerticesFixed") m_mode = VERTICES_FIXED;
64  else if (modeString == "VerticesSigma") m_mode = VERTICES_SIGMA;
65  else edm::LogError ("PointSeededTrackingRegionsProducer")<<"Unknown mode string: "<<modeString;
66 
67  // basic inputs
68  edm::ParameterSet point_input = regPSet.getParameter<edm::ParameterSet>("point_input");
69  eta_input = point_input.getParameter<double>("eta");
70  phi_input = point_input.getParameter<double>("phi");
71  m_maxNRegions = regPSet.getParameter<int>("maxNRegions");
72  token_beamSpot = iC.consumes<reco::BeamSpot>(regPSet.getParameter<edm::InputTag>("beamSpot"));
73  m_maxNVertices = 1;
75  {
76  token_vertex = iC.consumes<reco::VertexCollection>(regPSet.getParameter<edm::InputTag>("vertexCollection"));
77  m_maxNVertices = regPSet.getParameter<int>("maxNVertices");
78  }
79 
80  // RectangularEtaPhiTrackingRegion parameters:
81  m_ptMin = regPSet.getParameter<double>("ptMin");
82  m_originRadius = regPSet.getParameter<double>("originRadius");
83  m_zErrorBeamSpot = regPSet.getParameter<double>("zErrorBeamSpot");
84  m_deltaEta = regPSet.getParameter<double>("deltaEta");
85  m_deltaPhi = regPSet.getParameter<double>("deltaPhi");
86  m_precise = regPSet.getParameter<bool>("precise");
88  if (regPSet.exists("measurementTrackerName"))
89  {
90  // FIXME: when next time altering the configuration of this
91  // class, please change the types of the following parameters:
92  // - whereToUseMeasurementTracker to at least int32 or to a string
93  // corresponding to the UseMeasurementTracker enumeration
94  // - measurementTrackerName to InputTag
95  if (regPSet.exists("whereToUseMeasurementTracker"))
96  m_whereToUseMeasurementTracker = RectangularEtaPhiTrackingRegion::doubleToUseMeasurementTracker(regPSet.getParameter<double>("whereToUseMeasurementTracker"));
97  if(m_whereToUseMeasurementTracker != RectangularEtaPhiTrackingRegion::UseMeasurementTracker::kNever)
98  token_measurementTracker = iC.consumes<MeasurementTrackerEvent>(regPSet.getParameter<std::string>("measurementTrackerName"));
99  }
100  m_searchOpt = false;
101  if (regPSet.exists("searchOpt")) m_searchOpt = regPSet.getParameter<bool>("searchOpt");
102 
103  // mode-dependent z-halflength of tracking regions
104  if (m_mode == VERTICES_SIGMA) m_nSigmaZVertex = regPSet.getParameter<double>("nSigmaZVertex");
105  if (m_mode == VERTICES_FIXED) m_zErrorVetex = regPSet.getParameter<double>("zErrorVetex");
106  m_nSigmaZBeamSpot = -1.;
107  if (m_mode == BEAM_SPOT_SIGMA)
108  {
109  m_nSigmaZBeamSpot = regPSet.getParameter<double>("nSigmaZBeamSpot");
110  if (m_nSigmaZBeamSpot < 0.)
111  edm::LogError ("PointSeededTrackingRegionsProducer")<<"nSigmaZBeamSpot must be positive for BeamSpotSigma mode!";
112  }
113  }
114 
116 
117 
118  virtual std::vector<TrackingRegion* > regions(const edm::Event& e, const edm::EventSetup& es) const
119  {
120  std::vector<TrackingRegion* > 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  origins.push_back( std::make_pair( default_origin,
137  ));
138  } else if (m_mode == VERTICES_FIXED || m_mode == VERTICES_SIGMA) {
140  e.getByToken( token_vertex, vertices );
141  int n_vert = 0;
142  for (reco::VertexCollection::const_iterator iv = vertices->begin(), ev = vertices->end();
143  iv != ev && n_vert < m_maxNVertices; ++iv) {
144  if ( iv->isFake() || !iv->isValid() ) continue;
145 
146  origins.push_back( std::make_pair( GlobalPoint( iv->x(), iv->y(), iv->z() ),
147  (m_mode == VERTICES_FIXED) ? m_zErrorVetex : m_nSigmaZVertex*iv->zError()
148  ));
149  ++n_vert;
150  }
151  // no-vertex fall-back case:
152  if ( origins.empty() ) {
153  origins.push_back( std::make_pair( default_origin,
154  (m_nSigmaZBeamSpot > 0.) ? m_nSigmaZBeamSpot*bs->z0Error() : m_zErrorBeamSpot
155  ));
156  }
157  }
158 
163  measurementTracker = hmte.product();
164  }
165 
166  // create tracking regions (maximum MaxNRegions of them) in directions of the
167  // points of interest
168  size_t n_points = 1;
169  int n_regions = 0;
170  for(size_t i = 0; i < n_points && n_regions < m_maxNRegions; ++i ) {
171 
172  double x = TMath::Cos(phi_input);
173  double y = TMath::Sin(phi_input);
174  double theta = 2*TMath::ATan(TMath::Exp(-eta_input));
175  double z = (x*x+y*y)/TMath::Tan(theta);
176 
177  GlobalVector direction( x,y,z );
178 
179  for (size_t j=0; j<origins.size() && n_regions < m_maxNRegions; ++j) {
180 
181  result.push_back( new RectangularEtaPhiTrackingRegion(
182  direction, // GlobalVector
183  origins[j].first, // GlobalPoint
184  m_ptMin,
186  origins[j].second,
187  m_deltaEta,
188  m_deltaPhi,
190  m_precise,
191  measurementTracker,
193  ));
194  ++n_regions;
195  }
196  }
197  //std::cout<<"n_seeded_regions = "<<n_regions<<std::endl;
198  edm::LogInfo ("PointSeededTrackingRegionsProducer") << "produced "<<n_regions<<" regions";
199 
200  return result;
201  }
202 
203 private:
204 
206 
212 
213  float m_ptMin;
216  float m_deltaEta;
217  float m_deltaPhi;
218  bool m_precise;
222 
226 };
227 
228 #endif
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
PointSeededTrackingRegionsProducer(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
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
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
float float float z
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
U second(std::pair< T, U > const &p)
tuple result
Definition: query.py:137
int j
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< reco::VertexCollection > token_vertex
bool isValid() const
Definition: HandleBase.h:75
tuple conf
Definition: dbtoconf.py:185
static UseMeasurementTracker doubleToUseMeasurementTracker(double value)
T const * product() const
Definition: Handle.h:81
virtual std::vector< TrackingRegion * > regions(const edm::Event &e, const edm::EventSetup &es) const
bool isUninitialized() const
Definition: EDGetToken.h:71
Definition: DDAxes.h:10