CMS 3D CMS Logo

L3MumuTrackingRegion.h
Go to the documentation of this file.
1 #ifndef HLTrigger_btau_L3MumuTrackingRegion_H
2 #define HLTrigger_btau_L3MumuTrackingRegion_H
3 
6 
20 
22 public:
25  edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
26 
27  theVertexTag = regionPSet.getParameter<edm::InputTag>("vertexSrc");
28  theVertex = (theVertexTag.label().length() > 1);
29  theInputTrkTag = regionPSet.getParameter<edm::InputTag>("TrkSrc");
30  useVtxTks = regionPSet.getParameter<bool>("UseVtxTks");
31 
32  if (theVertex)
34  if (!(theVertex && useVtxTks))
36 
37  thePtMin = regionPSet.getParameter<double>("ptMin");
38  theOriginRadius = regionPSet.getParameter<double>("originRadius");
39  theOriginHalfLength = regionPSet.getParameter<double>("originHalfLength");
40  theOriginZPos = regionPSet.getParameter<double>("vertexZDefault");
41 
42  theDeltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
43  theDeltaPhi = regionPSet.getParameter<double>("deltaPhiRegion");
44  if (regionPSet.exists("searchOpt")) {
45  m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
46  } else {
47  m_searchOpt = false;
48  }
50  regionPSet.getParameter<std::string>("howToUseMeasurementTracker"));
53  iC.consumes<MeasurementTrackerEvent>(regionPSet.getParameter<edm::InputTag>("measurementTrackerName"));
54  }
55  }
56 
57  ~L3MumuTrackingRegion() override = default;
58 
59  std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& ev,
60  const edm::EventSetup& es) const override {
61  std::vector<std::unique_ptr<TrackingRegion> > result;
62 
66  ev.getByToken(theMeasurementTrackerToken, hmte);
67  measurementTracker = hmte.product();
68  }
69  const auto& field = es.getData(theFieldToken);
70  const auto& msmaker = es.getData(theMSMakerToken);
71 
72  // optional constraint for vertex
73  // get highest Pt pixel vertex (if existing)
74  double deltaZVertex = theOriginHalfLength;
75  double originz = theOriginZPos;
76  if (theVertex) {
78  ev.getByToken(theVertexToken, vertices);
79  const reco::VertexCollection vertCollection = *(vertices.product());
80  reco::VertexCollection::const_iterator ci = vertCollection.begin();
81  if (!vertCollection.empty()) {
82  originz = ci->z();
83  } else {
84  originz = theOriginZPos;
85  deltaZVertex = 15.;
86  }
87  if (useVtxTks) {
88  for (ci = vertCollection.begin(); ci != vertCollection.end(); ci++)
89  for (reco::Vertex::trackRef_iterator trackIt = ci->tracks_begin(); trackIt != ci->tracks_end(); trackIt++) {
90  reco::TrackRef iTrk = (*trackIt).castTo<reco::TrackRef>();
91  GlobalVector dirVector((iTrk)->px(), (iTrk)->py(), (iTrk)->pz());
92  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(dirVector,
93  GlobalPoint(0, 0, float(ci->z())),
94  thePtMin,
96  deltaZVertex,
99  field,
100  &msmaker,
101  true,
104  m_searchOpt));
105  }
106  return result;
107  }
108  }
109 
112  ev.getByToken(theInputTrkToken, trks);
113  for (reco::TrackCollection::const_iterator iTrk = trks->begin(); iTrk != trks->end(); iTrk++) {
114  GlobalVector dirVector((iTrk)->px(), (iTrk)->py(), (iTrk)->pz());
115  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(dirVector,
116  GlobalPoint(0, 0, float(originz)),
117  thePtMin,
119  deltaZVertex,
120  theDeltaEta,
121  theDeltaPhi,
122  field,
123  &msmaker,
124  true,
127  m_searchOpt));
128  }
129  return result;
130  }
131 
132 private:
134  bool theVertex;
138 
139  bool useVtxTks;
140 
141  double thePtMin;
145 
146  double theDeltaEta;
147  double theDeltaPhi;
153 };
154 
155 #endif
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &ev, const edm::EventSetup &es) const override
edm::EDGetTokenT< reco::TrackCollection > theInputTrkToken
bool exists(std::string const &parameterName) const
checks if a parameter exists
T const * product() const
Definition: Handle.h:70
constexpr bool isUninitialized() const noexcept
Definition: EDGetToken.h:104
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
std::string const & label() const
Definition: InputTag.h:36
edm::EDGetTokenT< reco::VertexCollection > theVertexToken
bool getData(T &iHolder) const
Definition: EventSetup.h:122
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > theMSMakerToken
edm::EDGetTokenT< MeasurementTrackerEvent > theMeasurementTrackerToken
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > theFieldToken
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
L3MumuTrackingRegion(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
~L3MumuTrackingRegion() override=default
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:38