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  }
56  if (regionPSet.exists("measurementTrackerName")){
57  // FIXME: when next time altering the configuration of this
58  // class, please change the types of the following parameters:
59  // - howToUseMeasurementTracker to at least int32 or to a string
60  // corresponding to the UseMeasurementTracker enumeration
61  // - measurementTrackerName to InputTag
62  if (regionPSet.exists("howToUseMeasurementTracker")){
64  }
66  token_measurementTracker = iC.consumes<MeasurementTrackerEvent>(regionPSet.getParameter<std::string>("measurementTrackerName"));
67  }
68  }
69  }
70 
72 
73 
74  virtual std::vector<TrackingRegion* > regions(const edm::Event& e, const edm::EventSetup& es) const {
75  std::vector<TrackingRegion* > result;
76 
77  // double originZ;
78  double deltaZVertex, deltaRho;
79  GlobalPoint vertex;
80  // get the primary vertex
82  e.getByToken(token_vertex, h_vertices);
83  const reco::VertexCollection & vertices = * h_vertices;
84  if (not vertices.empty()) {
85 // originZ = vertices.front().z();
86  GlobalPoint myTmp(vertices.at(0).position().x(),vertices.at(0).position().y(), vertices.at(0).position().z());
87  vertex = myTmp;
88  deltaZVertex = m_halfLength;
89  deltaRho = m_originRadius;
90  } else {
91  // originZ = 0.;
92  GlobalPoint myTmp(0.,0.,0.);
93  vertex = myTmp;
94  deltaRho = 1.;
95  deltaZVertex = 15.;
96  }
97 
98  // get the jet direction
100  e.getByToken(token_jet, h_jets);
101 
102  const MeasurementTrackerEvent *measurementTracker = nullptr;
106  measurementTracker = hmte.product();
107  }
108 
109  for(const reco::Candidate& myJet: *h_jets)
110  {
111  GlobalVector jetVector(myJet.momentum().x(),myJet.momentum().y(),myJet.momentum().z());
112 // GlobalPoint vertex(0, 0, originZ);
114  vertex,
115  m_ptMin,
116  deltaRho,
117  deltaZVertex,
118  m_deltaEta,
119  m_deltaPhi,
121  true,
122  measurementTracker,
123  m_searchOpt);
124  result.push_back(etaphiRegion);
125  }
126 
127  return result;
128  }
129 
130  private:
132 
133  float m_ptMin;
136  float m_deltaEta;
137  float m_deltaPhi;
143 };
144 
145 #endif
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
TauRegionalPixelSeedGenerator(const edm::ParameterSet &conf_, edm::ConsumesCollector &&iC)
virtual std::vector< TrackingRegion * > regions(const edm::Event &e, const edm::EventSetup &es) const
tuple result
Definition: query.py:137
edm::EDGetTokenT< reco::VertexCollection > token_vertex
static UseMeasurementTracker doubleToUseMeasurementTracker(double value)
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
T const * product() const
Definition: Handle.h:81
bool isUninitialized() const
Definition: EDGetToken.h:71
edm::EDGetTokenT< reco::CandidateView > token_jet