CMS 3D CMS Logo

TauRegionalPixelSeedGenerator.h
Go to the documentation of this file.
1 #ifndef TauRegionalPixelSeedGenerator_h
2 #define TauRegionalPixelSeedGenerator_h
3 
4 //
5 // Class: TauRegionalPixelSeedGenerator
6 
24 // Math
25 #include "Math/GenVector/VectorUtil.h"
26 #include "Math/GenVector/PxPyPzE4D.h"
27 
29 
32 
34 public:
36  edm::LogInfo("TauRegionalPixelSeedGenerator") << "Enter the TauRegionalPixelSeedGenerator";
37 
38  edm::ParameterSet regionPSet = conf_.getParameter<edm::ParameterSet>("RegionPSet");
39 
40  m_ptMin = regionPSet.getParameter<double>("ptMin");
41  m_originRadius = regionPSet.getParameter<double>("originRadius");
42  m_halfLength = regionPSet.getParameter<double>("originHalfLength");
43  m_deltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
44  m_deltaPhi = regionPSet.getParameter<double>("deltaPhiRegion");
45  token_jet = iC.consumes<reco::CandidateView>(regionPSet.getParameter<edm::InputTag>("JetSrc"));
46  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("vertexSrc"));
47  if (regionPSet.exists("searchOpt")) {
48  m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
49  } else {
50  m_searchOpt = false;
51  }
53  regionPSet.getParameter<std::string>("howToUseMeasurementTracker"));
55  // temporary until everything migrated to InputTag
57  if (regionPSet.existsAs<edm::InputTag>("measurementTrackerName")) {
58  tag = regionPSet.getParameter<edm::InputTag>("measurementTrackerName");
59  } else {
60  tag = edm::InputTag(regionPSet.getParameter<std::string>("measurementTrackerName"));
61  }
62 
64  }
65  }
66 
68 
69  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
71 
72  desc.add<double>("ptMin", 5.0);
73  desc.add<double>("originHalfLength", 0.2);
74  desc.add<double>("originRadius", 0.2);
75  desc.add<double>("deltaEtaRegion", 0.1);
76  desc.add<double>("deltaPhiRegion", 0.1);
77  desc.add<edm::InputTag>("JetSrc", edm::InputTag("icone5Tau1"));
78  desc.add<edm::InputTag>("vertexSrc", edm::InputTag("pixelVertices"));
79  desc.add<bool>("searchOpt", false);
80 
81  desc.add<std::string>("howToUseMeasurementTracker", "ForSiStrips");
82 
83  // allow both InputTag and string for the moment, use InputTag as the default
85  "measurementTrackerName", edm::InputTag("MeasurementTrackerEvent"), true) xor
86  edm::ParameterDescription<std::string>("measurementTrackerName", "MeasurementTrackerEvent", true));
87 
88  // Only for backwards-compatibility
90  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
91 
92  descriptions.add("tauRegionalPixelSeedTrackingRegions", descRegion);
93  }
94 
95  std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& e, const edm::EventSetup& es) const override {
96  std::vector<std::unique_ptr<TrackingRegion> > result;
97 
98  // double originZ;
99  double deltaZVertex, deltaRho;
101  // get the primary vertex
103  e.getByToken(token_vertex, h_vertices);
104  const reco::VertexCollection& vertices = *h_vertices;
105  if (not vertices.empty()) {
106  // originZ = vertices.front().z();
107  GlobalPoint myTmp(vertices.at(0).position().x(), vertices.at(0).position().y(), vertices.at(0).position().z());
108  vertex = myTmp;
109  deltaZVertex = m_halfLength;
110  deltaRho = m_originRadius;
111  } else {
112  // originZ = 0.;
113  GlobalPoint myTmp(0., 0., 0.);
114  vertex = myTmp;
115  deltaRho = 1.;
116  deltaZVertex = 15.;
117  }
118 
119  // get the jet direction
121  e.getByToken(token_jet, h_jets);
122 
127  measurementTracker = hmte.product();
128  }
129 
130  for (const reco::Candidate& myJet : *h_jets) {
131  GlobalVector jetVector(myJet.momentum().x(), myJet.momentum().y(), myJet.momentum().z());
132  // GlobalPoint vertex(0, 0, originZ);
133  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(jetVector,
134  vertex,
135  m_ptMin,
136  deltaRho,
137  deltaZVertex,
138  m_deltaEta,
139  m_deltaPhi,
141  true,
142  measurementTracker,
143  m_searchOpt));
144  }
145 
146  return result;
147  }
148 
149 private:
151 
152  float m_ptMin;
155  float m_deltaEta;
156  float m_deltaPhi;
162 };
163 
164 #endif
T getParameter(std::string const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:160
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
bool exists(std::string const &parameterName) const
checks if a parameter exists
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override
TauRegionalPixelSeedGenerator(const edm::ParameterSet &conf_, edm::ConsumesCollector &&iC)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< reco::VertexCollection > token_vertex
T const * product() const
Definition: Handle.h:69
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
bool isUninitialized() const
Definition: EDGetToken.h:70
edm::EDGetTokenT< reco::CandidateView > token_jet