CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TauRegionalPixelSeedGenerator.h
Go to the documentation of this file.
1 #ifndef TauRegionalPixelSeedGenerator_h
2 #define TauRegionalPixelSeedGenerator_h
3 
4 //
5 // Class: TauRegionalPixelSeedGenerator
6 
7 
23 // Math
24 #include "Math/GenVector/VectorUtil.h"
25 #include "Math/GenVector/PxPyPzE4D.h"
26 
28 
31 
32 
34  public:
35 
38  edm::LogInfo ("TauRegionalPixelSeedGenerator") << "Enter the TauRegionalPixelSeedGenerator";
39 
40  edm::ParameterSet regionPSet = conf_.getParameter<edm::ParameterSet>("RegionPSet");
41 
42  m_ptMin = regionPSet.getParameter<double>("ptMin");
43  m_originRadius = regionPSet.getParameter<double>("originRadius");
44  m_halfLength = regionPSet.getParameter<double>("originHalfLength");
45  m_deltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
46  m_deltaPhi = regionPSet.getParameter<double>("deltaPhiRegion");
47  token_jet = iC.consumes<reco::CandidateView>(regionPSet.getParameter<edm::InputTag>("JetSrc"));
48  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("vertexSrc"));
49  if (regionPSet.exists("searchOpt")){
50  m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
51  }
52  else{
53  m_searchOpt = false;
54  }
57  token_measurementTracker = iC.consumes<MeasurementTrackerEvent>(regionPSet.getParameter<std::string>("measurementTrackerName"));
58  }
59  }
60 
62 
63 
64  virtual std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& e, const edm::EventSetup& es) const override {
65  std::vector<std::unique_ptr<TrackingRegion> > result;
66 
67  // double originZ;
68  double deltaZVertex, deltaRho;
69  GlobalPoint vertex;
70  // get the primary vertex
72  e.getByToken(token_vertex, h_vertices);
73  const reco::VertexCollection & vertices = * h_vertices;
74  if (not vertices.empty()) {
75 // originZ = vertices.front().z();
76  GlobalPoint myTmp(vertices.at(0).position().x(),vertices.at(0).position().y(), vertices.at(0).position().z());
77  vertex = myTmp;
78  deltaZVertex = m_halfLength;
79  deltaRho = m_originRadius;
80  } else {
81  // originZ = 0.;
82  GlobalPoint myTmp(0.,0.,0.);
83  vertex = myTmp;
84  deltaRho = 1.;
85  deltaZVertex = 15.;
86  }
87 
88  // get the jet direction
90  e.getByToken(token_jet, h_jets);
91 
96  measurementTracker = hmte.product();
97  }
98 
99  for(const reco::Candidate& myJet: *h_jets)
100  {
101  GlobalVector jetVector(myJet.momentum().x(),myJet.momentum().y(),myJet.momentum().z());
102 // GlobalPoint vertex(0, 0, originZ);
103  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>( jetVector,
104  vertex,
105  m_ptMin,
106  deltaRho,
107  deltaZVertex,
108  m_deltaEta,
109  m_deltaPhi,
111  true,
112  measurementTracker,
113  m_searchOpt));
114  }
115 
116  return result;
117  }
118 
119  private:
121 
122  float m_ptMin;
125  float m_deltaEta;
126  float m_deltaPhi;
132 };
133 
134 #endif
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
tuple result
Definition: mps_fire.py:83
TauRegionalPixelSeedGenerator(const edm::ParameterSet &conf_, edm::ConsumesCollector &&iC)
edm::EDGetTokenT< reco::VertexCollection > token_vertex
T const * product() const
Definition: Handle.h:81
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
virtual std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override
bool isUninitialized() const
Definition: EDGetToken.h:73
edm::EDGetTokenT< reco::CandidateView > token_jet