CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GlobalTrackingRegionWithVerticesProducer.h
Go to the documentation of this file.
1 #ifndef RecoTracker_TkTrackingRegions_GlobalTrackingRegionWithVerticesProducer_H
2 #define RecoTracker_TkTrackingRegions_GlobalTrackingRegionWithVerticesProducer_H
3 
8 
14 
16 {
17 public:
18 
21  {
22  edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
23 
24  thePtMin = regionPSet.getParameter<double>("ptMin");
25  theOriginRadius = regionPSet.getParameter<double>("originRadius");
26  theNSigmaZ = regionPSet.getParameter<double>("nSigmaZ");
27  token_beamSpot = iC.consumes<reco::BeamSpot>(regionPSet.getParameter<edm::InputTag>("beamSpot"));
28  thePrecise = regionPSet.getParameter<bool>("precise");
29 
30  theSigmaZVertex = regionPSet.getParameter<double>("sigmaZVertex");
31  theFixedError = regionPSet.getParameter<double>("fixedError");
32 
33  theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices");
34  theUseFakeVertices = regionPSet.existsAs<bool>("useFakeVertices") ? regionPSet.getParameter<bool>("useFakeVertices") : false;
35  theUseFixedError = regionPSet.getParameter<bool>("useFixedError");
36  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("VertexCollection"));
37  }
38 
40 
41  virtual std::vector<TrackingRegion* > regions
42  (const edm::Event& ev, const edm::EventSetup&) const
43  {
44  std::vector<TrackingRegion* > result;
45 
46  GlobalPoint theOrigin;
48  ev.getByToken( token_beamSpot, bsHandle);
49  double bsSigmaZ;
50  if(bsHandle.isValid()) {
51  const reco::BeamSpot & bs = *bsHandle;
52  bsSigmaZ = theNSigmaZ*bs.sigmaZ();
53  theOrigin = GlobalPoint(bs.x0(), bs.y0(), bs.z0());
54  }else{
55  throw cms::Exception("Seeding") << "ERROR: input beamSpot is not valid in GlobalTrackingRegionWithVertices";
56  }
57 
59  {
61  ev.getByToken(token_vertex,vertexCollection);
62 
63  for(reco::VertexCollection::const_iterator iV=vertexCollection->begin(); iV != vertexCollection->end() ; iV++) {
64  if (!iV->isValid()) continue;
65  if (iV->isFake() && !(theUseFakeVertices && theUseFixedError)) continue;
66  GlobalPoint theOrigin_ = GlobalPoint(iV->x(),iV->y(),iV->z());
67  double theOriginHalfLength_ = (theUseFixedError ? theFixedError : (iV->zError())*theSigmaZVertex);
68  result.push_back( new GlobalTrackingRegion(thePtMin, theOrigin_, theOriginRadius, theOriginHalfLength_, thePrecise) );
69  }
70 
71  if (result.empty()) {
72  result.push_back( new GlobalTrackingRegion(thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise) );
73  }
74  }
75  else
76  {
77  result.push_back(
78  new GlobalTrackingRegion(thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise) );
79  }
80 
81  return result;
82  }
83 
84 private:
85  double thePtMin;
86  double theOriginRadius;
87  double theNSigmaZ;
89 
91  double theFixedError;
92  bool thePrecise;
93 
99 
100 
101 };
102 
103 #endif
virtual std::vector< TrackingRegion * > regions(const edm::Event &ev, const edm::EventSetup &) const
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:68
tuple cfg
Definition: looper.py:259
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:185
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
edm::EDGetTokenT< reco::VertexCollection > token_vertex
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool ev
tuple vertexCollection
tuple result
Definition: query.py:137
bool isValid() const
Definition: HandleBase.h:75
double sigmaZ() const
sigma z
Definition: BeamSpot.h:80
GlobalTrackingRegionWithVerticesProducer(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
double y0() const
y coordinate
Definition: BeamSpot.h:66
double x0() const
x coordinate
Definition: BeamSpot.h:64