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  theUseFixedError = regionPSet.getParameter<bool>("useFixedError");
35  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("VertexCollection"));
36  }
37 
39 
40  virtual std::vector<TrackingRegion* > regions
41  (const edm::Event& ev, const edm::EventSetup&) const
42  {
43  std::vector<TrackingRegion* > result;
44 
45  GlobalPoint theOrigin;
47  ev.getByToken( token_beamSpot, bsHandle);
48  double bsSigmaZ;
49  if(bsHandle.isValid()) {
50  const reco::BeamSpot & bs = *bsHandle;
51  bsSigmaZ = theNSigmaZ*bs.sigmaZ();
52  theOrigin = GlobalPoint(bs.x0(), bs.y0(), bs.z0());
53  }else{
54  throw cms::Exception("Seeding") << "ERROR: input beamSpot is not valid in GlobalTrackingRegionWithVertices";
55  }
56 
58  {
60  ev.getByToken(token_vertex,vertexCollection);
61 
62  for(reco::VertexCollection::const_iterator iV=vertexCollection->begin(); iV != vertexCollection->end() ; iV++) {
63  if (iV->isFake() || !iV->isValid()) continue;
64  GlobalPoint theOrigin_ = GlobalPoint(iV->x(),iV->y(),iV->z());
65  double theOriginHalfLength_ = (theUseFixedError ? theFixedError : (iV->zError())*theSigmaZVertex);
66  result.push_back( new GlobalTrackingRegion(thePtMin, theOrigin_, theOriginRadius, theOriginHalfLength_, thePrecise) );
67  }
68 
69  if (result.empty()) {
70  result.push_back( new GlobalTrackingRegion(thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise) );
71  }
72  }
73  else
74  {
75  result.push_back(
76  new GlobalTrackingRegion(thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise) );
77  }
78 
79  return result;
80  }
81 
82 private:
83  double thePtMin;
84  double theOriginRadius;
85  double theNSigmaZ;
87 
89  double theFixedError;
90  bool thePrecise;
91 
96 
97 
98 };
99 
100 #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
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
edm::EDGetTokenT< reco::VertexCollection > token_vertex
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
tuple vertexCollection
tuple result
Definition: query.py:137
bool isValid() const
Definition: HandleBase.h:76
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