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 
7 
25 // Math
26 #include "Math/GenVector/VectorUtil.h"
27 #include "Math/GenVector/PxPyPzE4D.h"
28 
30 
33 
34 
36  public:
37 
40  edm::LogInfo ("TauRegionalPixelSeedGenerator") << "Enter the TauRegionalPixelSeedGenerator";
41 
42  edm::ParameterSet regionPSet = conf_.getParameter<edm::ParameterSet>("RegionPSet");
43 
44  m_ptMin = regionPSet.getParameter<double>("ptMin");
45  m_originRadius = regionPSet.getParameter<double>("originRadius");
46  m_halfLength = regionPSet.getParameter<double>("originHalfLength");
47  m_deltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
48  m_deltaPhi = regionPSet.getParameter<double>("deltaPhiRegion");
49  token_jet = iC.consumes<reco::CandidateView>(regionPSet.getParameter<edm::InputTag>("JetSrc"));
50  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("vertexSrc"));
51  if (regionPSet.exists("searchOpt")){
52  m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
53  }
54  else{
55  m_searchOpt = false;
56  }
59  // temporary until everything migrated to InputTag
61  if(regionPSet.existsAs<edm::InputTag>("measurementTrackerName")) {
62  tag = regionPSet.getParameter<edm::InputTag>("measurementTrackerName");
63  }
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
89  desc.addNode(edm::ParameterDescription<edm::InputTag>("measurementTrackerName", edm::InputTag("MeasurementTrackerEvent"), true) xor
90  edm::ParameterDescription<std::string>("measurementTrackerName", "MeasurementTrackerEvent", true));
91 
92  // Only for backwards-compatibility
94  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
95 
96  descriptions.add("tauRegionalPixelSeedTrackingRegions", descRegion);
97  }
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 
128  const MeasurementTrackerEvent *measurementTracker = nullptr;
132  measurementTracker = hmte.product();
133  }
134 
135  for(const reco::Candidate& myJet: *h_jets)
136  {
137  GlobalVector jetVector(myJet.momentum().x(),myJet.momentum().y(),myJet.momentum().z());
138 // GlobalPoint vertex(0, 0, originZ);
139  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>( jetVector,
140  vertex,
141  m_ptMin,
142  deltaRho,
143  deltaZVertex,
144  m_deltaEta,
145  m_deltaPhi,
147  true,
148  measurementTracker,
149  m_searchOpt));
150  }
151 
152  return result;
153  }
154 
155  private:
157 
158  float m_ptMin;
161  float m_deltaEta;
162  float m_deltaPhi;
168 };
169 
170 #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:186
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
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:81
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:73
edm::EDGetTokenT< reco::CandidateView > token_jet