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 {
48 public:
49 
51 
54  {
55  edm::ParameterSet regPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
56 
57  // operation mode
58  std::string modeString = regPSet.getParameter<std::string>("mode");
59  if (modeString == "BeamSpotFixed") m_mode = BEAM_SPOT_FIXED;
60  else if (modeString == "BeamSpotSigma") m_mode = BEAM_SPOT_SIGMA;
61  else if (modeString == "VerticesFixed") m_mode = VERTICES_FIXED;
62  else if (modeString == "VerticesSigma") m_mode = VERTICES_SIGMA;
63  else edm::LogError ("PointSeededTrackingRegionsProducer")<<"Unknown mode string: "<<modeString;
64 
65  // basic inputs
66  edm::ParameterSet point_input = regPSet.getParameter<edm::ParameterSet>("point_input");
67  eta_input = point_input.getParameter<double>("eta");
68  phi_input = point_input.getParameter<double>("phi");
69  m_maxNRegions = regPSet.getParameter<int>("maxNRegions");
70  token_beamSpot = iC.consumes<reco::BeamSpot>(regPSet.getParameter<edm::InputTag>("beamSpot"));
71  m_maxNVertices = 1;
73  {
74  token_vertex = iC.consumes<reco::VertexCollection>(regPSet.getParameter<edm::InputTag>("vertexCollection"));
75  m_maxNVertices = regPSet.getParameter<int>("maxNVertices");
76  }
77 
78  // RectangularEtaPhiTrackingRegion parameters:
79  m_ptMin = regPSet.getParameter<double>("ptMin");
80  m_originRadius = regPSet.getParameter<double>("originRadius");
81  m_zErrorBeamSpot = regPSet.getParameter<double>("zErrorBeamSpot");
82  m_deltaEta = regPSet.getParameter<double>("deltaEta");
83  m_deltaPhi = regPSet.getParameter<double>("deltaPhi");
84  m_precise = regPSet.getParameter<bool>("precise");
86  if (regPSet.exists("measurementTrackerName"))
87  {
88  // FIXME: when next time altering the configuration of this
89  // class, please change the types of the following parameters:
90  // - whereToUseMeasurementTracker to at least int32 or to a string
91  // corresponding to the UseMeasurementTracker enumeration
92  // - measurementTrackerName to InputTag
93  if (regPSet.exists("whereToUseMeasurementTracker"))
94  m_whereToUseMeasurementTracker = RectangularEtaPhiTrackingRegion::doubleToUseMeasurementTracker(regPSet.getParameter<double>("whereToUseMeasurementTracker"));
95  if(m_whereToUseMeasurementTracker != RectangularEtaPhiTrackingRegion::UseMeasurementTracker::kNever)
96  token_measurementTracker = iC.consumes<MeasurementTrackerEvent>(regPSet.getParameter<std::string>("measurementTrackerName"));
97  }
98  m_searchOpt = false;
99  if (regPSet.exists("searchOpt")) m_searchOpt = regPSet.getParameter<bool>("searchOpt");
100 
101  // mode-dependent z-halflength of tracking regions
102  if (m_mode == VERTICES_SIGMA) m_nSigmaZVertex = regPSet.getParameter<double>("nSigmaZVertex");
103  if (m_mode == VERTICES_FIXED) m_zErrorVetex = regPSet.getParameter<double>("zErrorVetex");
104  m_nSigmaZBeamSpot = -1.;
105  if (m_mode == BEAM_SPOT_SIGMA)
106  {
107  m_nSigmaZBeamSpot = regPSet.getParameter<double>("nSigmaZBeamSpot");
108  if (m_nSigmaZBeamSpot < 0.)
109  edm::LogError ("PointSeededTrackingRegionsProducer")<<"nSigmaZBeamSpot must be positive for BeamSpotSigma mode!";
110  }
111  }
112 
114 
115 
116  virtual std::vector<TrackingRegion* > regions(const edm::Event& e, const edm::EventSetup& es) const
117  {
118  std::vector<TrackingRegion* > result;
119 
120  // always need the beam spot (as a fall back strategy for vertex modes)
122  e.getByToken( token_beamSpot, bs );
123  if( !bs.isValid() ) return result;
124 
125  // this is a default origin for all modes
126  GlobalPoint default_origin( bs->x0(), bs->y0(), bs->z0() );
127 
128  // vector of origin & halfLength pairs:
129  std::vector< std::pair< GlobalPoint, float > > origins;
130 
131  // fill the origins and halfLengths depending on the mode
133  origins.push_back( std::make_pair( default_origin,
135  ));
136  } else if (m_mode == VERTICES_FIXED || m_mode == VERTICES_SIGMA) {
138  e.getByToken( token_vertex, vertices );
139  int n_vert = 0;
140  for (reco::VertexCollection::const_iterator iv = vertices->begin(), ev = vertices->end();
141  iv != ev && n_vert < m_maxNVertices; ++iv) {
142  if ( iv->isFake() || !iv->isValid() ) continue;
143 
144  origins.push_back( std::make_pair( GlobalPoint( iv->x(), iv->y(), iv->z() ),
145  (m_mode == VERTICES_FIXED) ? m_zErrorVetex : m_nSigmaZVertex*iv->zError()
146  ));
147  ++n_vert;
148  }
149  // no-vertex fall-back case:
150  if ( origins.empty() ) {
151  origins.push_back( std::make_pair( default_origin,
152  (m_nSigmaZBeamSpot > 0.) ? m_nSigmaZBeamSpot*bs->z0Error() : m_zErrorBeamSpot
153  ));
154  }
155  }
156 
161  measurementTracker = hmte.product();
162  }
163 
164  // create tracking regions (maximum MaxNRegions of them) in directions of the
165  // points of interest
166  size_t n_points = 1;
167  int n_regions = 0;
168  for(size_t i = 0; i < n_points && n_regions < m_maxNRegions; ++i ) {
169 
170  double x = TMath::Cos(phi_input);
171  double y = TMath::Sin(phi_input);
172  double theta = 2*TMath::ATan(TMath::Exp(-eta_input));
173  double z = (x*x+y*y)/TMath::Tan(theta);
174 
175  GlobalVector direction( x,y,z );
176 
177  for (size_t j=0; j<origins.size() && n_regions < m_maxNRegions; ++j) {
178 
179  result.push_back( new RectangularEtaPhiTrackingRegion(
180  direction, // GlobalVector
181  origins[j].first, // GlobalPoint
182  m_ptMin,
184  origins[j].second,
185  m_deltaEta,
186  m_deltaPhi,
188  m_precise,
189  measurementTracker,
191  ));
192  ++n_regions;
193  }
194  }
195  //std::cout<<"n_seeded_regions = "<<n_regions<<std::endl;
196  edm::LogInfo ("PointSeededTrackingRegionsProducer") << "produced "<<n_regions<<" regions";
197 
198  return result;
199  }
200 
201 private:
202 
204 
210 
211  float m_ptMin;
214  float m_deltaEta;
215  float m_deltaPhi;
216  bool m_precise;
220 
224 };
225 
226 #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:457
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
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