CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Private Attributes
CandidateSeededTrackingRegionsProducer Class Reference

#include <CandidateSeededTrackingRegionsProducer.h>

Inheritance diagram for CandidateSeededTrackingRegionsProducer:
TrackingRegionProducer

Public Types

enum  Mode { BEAM_SPOT_FIXED, BEAM_SPOT_SIGMA, VERTICES_FIXED, VERTICES_SIGMA }
 

Public Member Functions

 CandidateSeededTrackingRegionsProducer (const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
 
virtual std::vector
< std::unique_ptr
< TrackingRegion > > 
regions (const edm::Event &e, const edm::EventSetup &es) const override
 
virtual ~CandidateSeededTrackingRegionsProducer ()
 
- Public Member Functions inherited from TrackingRegionProducer
virtual ~TrackingRegionProducer ()
 

Private Attributes

float m_deltaEta
 
float m_deltaPhi
 
int m_maxNRegions
 
int m_maxNVertices
 
Mode m_mode
 
float m_nSigmaZBeamSpot
 
float m_nSigmaZVertex
 
float m_originRadius
 
bool m_precise
 
float m_ptMin
 
bool m_searchOpt
 
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
 
float m_zErrorBeamSpot
 
float m_zErrorVetex
 
edm::EDGetTokenT< reco::BeamSpottoken_beamSpot
 
edm::EDGetTokenT
< reco::CandidateView
token_input
 
edm::EDGetTokenT
< MeasurementTrackerEvent
token_measurementTracker
 
edm::EDGetTokenT
< reco::VertexCollection
token_vertex
 

Detailed Description

class CandidateSeededTrackingRegionsProducer

eta-phi TrackingRegions producer in directions defined by Candidate-based objects of interest from a collection defined by the "input" parameter.

Four operational modes are supported ("mode" parameter):

BeamSpotFixed: origin is defined by the beam spot z-half-length is defined by a fixed zErrorBeamSpot parameter BeamSpotSigma: origin is defined by the beam spot z-half-length is defined by nSigmaZBeamSpot * beamSpot.sigmaZ VerticesFixed: origins are defined by vertices from VertexCollection (use maximum MaxNVertices of them) z-half-length is defined by a fixed zErrorVetex parameter VerticesSigma: origins are defined by vertices from VertexCollection (use maximum MaxNVertices of them) z-half-length is defined by nSigmaZVertex * vetex.zError

If, while using one of the "Vertices" modes, there's no vertices in an event, we fall back into either BeamSpotSigma or BeamSpotFixed mode, depending on the positiveness of nSigmaZBeamSpot.

Author
Vadim Khotilovich

Definition at line 44 of file CandidateSeededTrackingRegionsProducer.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

CandidateSeededTrackingRegionsProducer::CandidateSeededTrackingRegionsProducer ( const edm::ParameterSet conf,
edm::ConsumesCollector &&  iC 
)
inlineexplicit

Definition at line 50 of file CandidateSeededTrackingRegionsProducer.h.

References BEAM_SPOT_FIXED, BEAM_SPOT_SIGMA, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), RectangularEtaPhiTrackingRegion::kNever, m_deltaEta, m_deltaPhi, m_maxNRegions, m_maxNVertices, m_mode, m_nSigmaZBeamSpot, m_nSigmaZVertex, m_originRadius, m_precise, m_ptMin, m_searchOpt, m_whereToUseMeasurementTracker, m_zErrorBeamSpot, m_zErrorVetex, AlCaHLTBitMon_QueryRunRegistry::string, RectangularEtaPhiTrackingRegion::stringToUseMeasurementTracker(), token_beamSpot, token_input, token_measurementTracker, token_vertex, VERTICES_FIXED, and VERTICES_SIGMA.

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
65  m_maxNRegions = regPSet.getParameter<int>("maxNRegions");
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) {
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  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
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
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
virtual CandidateSeededTrackingRegionsProducer::~CandidateSeededTrackingRegionsProducer ( )
inlinevirtual

Definition at line 100 of file CandidateSeededTrackingRegionsProducer.h.

100 {}

Member Function Documentation

virtual std::vector<std::unique_ptr<TrackingRegion> > CandidateSeededTrackingRegionsProducer::regions ( const edm::Event e,
const edm::EventSetup es 
) const
inlineoverridevirtual

Implements TrackingRegionProducer.

Definition at line 103 of file CandidateSeededTrackingRegionsProducer.h.

References BEAM_SPOT_FIXED, BEAM_SPOT_SIGMA, plotBeamSpotDB::first, edm::Event::getByToken(), i, edm::EDGetTokenT< T >::isUninitialized(), edm::HandleBase::isValid(), j, m_deltaEta, m_deltaPhi, m_maxNRegions, m_maxNVertices, m_mode, m_nSigmaZBeamSpot, m_nSigmaZVertex, m_originRadius, m_precise, m_ptMin, m_searchOpt, m_whereToUseMeasurementTracker, m_zErrorBeamSpot, m_zErrorVetex, HLT_FULL_cff::measurementTracker, edm::Handle< T >::product(), query::result, edm::second(), token_beamSpot, token_input, token_measurementTracker, token_vertex, findQualityFiles::v, HLT_FULL_cff::vertices, VERTICES_FIXED, VERTICES_SIGMA, x, y, and z.

104  {
105  std::vector<std::unique_ptr<TrackingRegion> > result;
106 
107  // pick up the candidate objects of interest
109  e.getByToken( token_input, objects );
110  size_t n_objects = objects->size();
111  if (n_objects == 0) return result;
112 
113  // always need the beam spot (as a fall back strategy for vertex modes)
115  e.getByToken( token_beamSpot, bs );
116  if( !bs.isValid() ) return result;
117 
118  // this is a default origin for all modes
119  GlobalPoint default_origin( bs->x0(), bs->y0(), bs->z0() );
120 
121  // vector of origin & halfLength pairs:
122  std::vector< std::pair< GlobalPoint, float > > origins;
123 
124  // fill the origins and halfLengths depending on the mode
126  {
127  origins.push_back( std::make_pair(
128  default_origin,
130  ));
131  }
132  else if (m_mode == VERTICES_FIXED || m_mode == VERTICES_SIGMA)
133  {
135  e.getByToken( token_vertex, vertices );
136  int n_vert = 0;
137  for (reco::VertexCollection::const_iterator v = vertices->begin(); v != vertices->end() && n_vert < m_maxNVertices; ++v)
138  {
139  if ( v->isFake() || !v->isValid() ) continue;
140 
141  origins.push_back( std::make_pair(
142  GlobalPoint( v->x(), v->y(), v->z() ),
144  ));
145  ++n_vert;
146  }
147  // no-vertex fall-back case:
148  if (origins.empty())
149  {
150  origins.push_back( std::make_pair(
151  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  // objects of interest (we expect that the collection was sorted in decreasing pt order)
166  int n_regions = 0;
167  for(size_t i = 0; i < n_objects && n_regions < m_maxNRegions; ++i )
168  {
169  const reco::Candidate & object = (*objects)[i];
170  GlobalVector direction( object.momentum().x(), object.momentum().y(), object.momentum().z() );
171 
172  for (size_t j=0; j<origins.size() && n_regions < m_maxNRegions; ++j)
173  {
174  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(
175  direction,
176  origins[j].first,
177  m_ptMin,
179  origins[j].second,
180  m_deltaEta,
181  m_deltaPhi,
183  m_precise,
184  measurementTracker,
186  ));
187  ++n_regions;
188  }
189  }
190  //std::cout<<"n_seeded_regions = "<<n_regions<<std::endl;
191  edm::LogInfo ("CandidateSeededTrackingRegionsProducer") << "produced "<<n_regions<<" regions";
192 
193  return result;
194  }
int i
Definition: DBlmapReader.cc:9
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
tuple measurementTracker
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
edm::EDGetTokenT< reco::VertexCollection > token_vertex
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 isValid() const
Definition: HandleBase.h:75
T const * product() const
Definition: Handle.h:81
bool isUninitialized() const
Definition: EDGetToken.h:73

Member Data Documentation

float CandidateSeededTrackingRegionsProducer::m_deltaEta
private
float CandidateSeededTrackingRegionsProducer::m_deltaPhi
private
int CandidateSeededTrackingRegionsProducer::m_maxNRegions
private
int CandidateSeededTrackingRegionsProducer::m_maxNVertices
private
Mode CandidateSeededTrackingRegionsProducer::m_mode
private
float CandidateSeededTrackingRegionsProducer::m_nSigmaZBeamSpot
private
float CandidateSeededTrackingRegionsProducer::m_nSigmaZVertex
private
float CandidateSeededTrackingRegionsProducer::m_originRadius
private
bool CandidateSeededTrackingRegionsProducer::m_precise
private
float CandidateSeededTrackingRegionsProducer::m_ptMin
private
bool CandidateSeededTrackingRegionsProducer::m_searchOpt
private
RectangularEtaPhiTrackingRegion::UseMeasurementTracker CandidateSeededTrackingRegionsProducer::m_whereToUseMeasurementTracker
private
float CandidateSeededTrackingRegionsProducer::m_zErrorBeamSpot
private
float CandidateSeededTrackingRegionsProducer::m_zErrorVetex
private
edm::EDGetTokenT<reco::BeamSpot> CandidateSeededTrackingRegionsProducer::token_beamSpot
private
edm::EDGetTokenT<reco::CandidateView> CandidateSeededTrackingRegionsProducer::token_input
private
edm::EDGetTokenT<MeasurementTrackerEvent> CandidateSeededTrackingRegionsProducer::token_measurementTracker
private
edm::EDGetTokenT<reco::VertexCollection> CandidateSeededTrackingRegionsProducer::token_vertex
private