CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TauRegionalPixelSeedGenerator.h
Go to the documentation of this file.
1 #ifndef TauRegionalPixelSeedGenerator_h
2 #define TauRegionalPixelSeedGenerator_h
3 
4 //
5 // Class: TauRegionalPixelSeedGenerator
6 
28 // Math
29 #include "Math/GenVector/VectorUtil.h"
30 #include "Math/GenVector/PxPyPzE4D.h"
31 
33 
36 
38 public:
41  edm::LogInfo("TauRegionalPixelSeedGenerator") << "Enter the TauRegionalPixelSeedGenerator";
42 
43  edm::ParameterSet regionPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
44 
45  m_ptMin = regionPSet.getParameter<double>("ptMin");
46  m_originRadius = regionPSet.getParameter<double>("originRadius");
47  m_halfLength = regionPSet.getParameter<double>("originHalfLength");
48  m_deltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
49  m_deltaPhi = regionPSet.getParameter<double>("deltaPhiRegion");
50  token_jet = iC.consumes<reco::CandidateView>(regionPSet.getParameter<edm::InputTag>("JetSrc"));
51  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("vertexSrc"));
52  if (regionPSet.exists("searchOpt")) {
53  m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
54  } else {
55  m_searchOpt = false;
56  }
58  regionPSet.getParameter<std::string>("howToUseMeasurementTracker"));
60  // temporary until everything migrated to InputTag
62  if (regionPSet.existsAs<edm::InputTag>("measurementTrackerName")) {
63  tag = regionPSet.getParameter<edm::InputTag>("measurementTrackerName");
64  } else {
65  tag = edm::InputTag(regionPSet.getParameter<std::string>("measurementTrackerName"));
66  }
67 
69  }
70  }
71 
73 
74  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
76 
77  desc.add<double>("ptMin", 5.0);
78  desc.add<double>("originHalfLength", 0.2);
79  desc.add<double>("originRadius", 0.2);
80  desc.add<double>("deltaEtaRegion", 0.1);
81  desc.add<double>("deltaPhiRegion", 0.1);
82  desc.add<edm::InputTag>("JetSrc", edm::InputTag("icone5Tau1"));
83  desc.add<edm::InputTag>("vertexSrc", edm::InputTag("pixelVertices"));
84  desc.add<bool>("searchOpt", false);
85 
86  desc.add<std::string>("howToUseMeasurementTracker", "ForSiStrips");
87 
88  // allow both InputTag and string for the moment, use InputTag as the default
90  "measurementTrackerName", edm::InputTag("MeasurementTrackerEvent"), true) xor
91  edm::ParameterDescription<std::string>("measurementTrackerName", "MeasurementTrackerEvent", true));
92 
93  // Only for backwards-compatibility
95  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
96 
97  descriptions.add("tauRegionalPixelSeedTrackingRegions", descRegion);
98  }
99 
100  std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& e, const edm::EventSetup& es) const override {
101  std::vector<std::unique_ptr<TrackingRegion> > result;
102 
103  // double originZ;
104  double deltaZVertex, deltaRho;
105  GlobalPoint vertex;
106  // get the primary vertex
108  e.getByToken(token_vertex, h_vertices);
109  const reco::VertexCollection& vertices = *h_vertices;
110  if (not vertices.empty()) {
111  // originZ = vertices.front().z();
112  GlobalPoint myTmp(vertices.at(0).position().x(), vertices.at(0).position().y(), vertices.at(0).position().z());
113  vertex = myTmp;
114  deltaZVertex = m_halfLength;
115  deltaRho = m_originRadius;
116  } else {
117  // originZ = 0.;
118  GlobalPoint myTmp(0., 0., 0.);
119  vertex = myTmp;
120  deltaRho = 1.;
121  deltaZVertex = 15.;
122  }
123 
124  // get the jet direction
126  e.getByToken(token_jet, h_jets);
127 
132  measurementTracker = hmte.product();
133  }
134 
135  const auto& bfield = es.getData(token_bfield);
136  const auto& msmaker = es.getData(token_msmaker);
137 
138  for (const reco::Candidate& myJet : *h_jets) {
139  GlobalVector jetVector(myJet.momentum().x(), myJet.momentum().y(), myJet.momentum().z());
140  // GlobalPoint vertex(0, 0, originZ);
141  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(jetVector,
142  vertex,
143  m_ptMin,
144  deltaRho,
145  deltaZVertex,
146  m_deltaEta,
147  m_deltaPhi,
148  bfield,
149  &msmaker,
150  true,
152  measurementTracker,
153  m_searchOpt));
154  }
155 
156  return result;
157  }
158 
159 private:
160  float m_ptMin;
163  float m_deltaEta;
164  float m_deltaPhi;
172 };
173 
174 #endif
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
tuple measurementTracker
constexpr bool isUninitialized() const noexcept
Definition: EDGetToken.h:99
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
tuple result
Definition: mps_fire.py:311
bool getData(T &iHolder) const
Definition: EventSetup.h:122
TauRegionalPixelSeedGenerator(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > token_bfield
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< reco::VertexCollection > token_vertex
Log< level::Info, false > LogInfo
T const * product() const
Definition: Handle.h:70
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > token_msmaker
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::EDGetTokenT< reco::CandidateView > token_jet