CMS 3D CMS Logo

L2MuonSeedGeneratorFromL1TkMu.cc
Go to the documentation of this file.
1 /* \class L2MuonSeedGeneratorFromL1TkMu
2  *
3  * L2 muon seed generator:
4  * Transform the L1TkMuon informations in seeds
5  * for the L2 muon reconstruction
6  * (mimicking L2MuonSeedGeneratorFromL1T)
7  *
8  * Author: H. Kwon
9  * Modified by M. Oh
10  */
11 
12 // Framework
23 
24 // Data Formats
36 
37 #include "CLHEP/Vector/ThreeVector.h"
46 
49 
50 using namespace std;
51 using namespace edm;
52 using namespace l1t;
53 
55 public:
58 
60  ~L2MuonSeedGeneratorFromL1TkMu() override = default;
61 
62  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
63  void produce(edm::Event &, const edm::EventSetup &) override;
64 
65 private:
68 
71 
73 
74  const double l1MinPt_;
75  const double l1MaxEta_;
76  const double minPtBarrel_;
77  const double minPtEndcap_;
78  const double minPL1Tk_;
79  const double minPtL1TkBarrel_;
80  const bool useOfflineSeed_;
81  const bool useUnassociatedL1_;
82  std::vector<double> matchingDR_;
83  std::vector<double> etaBins_;
84 
85  // parameters used in propagating L1 tracker track
86  // to the second muon station, numbers are
87  // taken from L1TkMuonProducer::propagateToGMT
88  static constexpr float etaBoundary_{1.1};
89  static constexpr float distMB2_{550.};
90  static constexpr float distME2_{850.};
91  static constexpr float phiCorr0_{1.464};
92  static constexpr float phiCorr1_{1.7};
93  static constexpr float phiCorr2_{144.};
94 
96  std::unique_ptr<MuonServiceProxy> service_;
97  std::unique_ptr<MeasurementEstimator> estimator_;
98 
99  const TrajectorySeed *associateOfflineSeedToL1(edm::Handle<edm::View<TrajectorySeed>> &,
100  std::vector<int> &,
102  double);
103 };
104 
105 // constructors
107  : source_(iConfig.getParameter<InputTag>("InputObjects")),
108  muCollToken_(consumes<l1t::TkMuonCollection>(source_)),
109  propagatorName_(iConfig.getParameter<string>("Propagator")),
110  l1MinPt_(iConfig.getParameter<double>("L1MinPt")),
111  l1MaxEta_(iConfig.getParameter<double>("L1MaxEta")),
112  minPtBarrel_(iConfig.getParameter<double>("SetMinPtBarrelTo")),
113  minPtEndcap_(iConfig.getParameter<double>("SetMinPtEndcapTo")),
114  minPL1Tk_(iConfig.getParameter<double>("MinPL1Tk")),
115  minPtL1TkBarrel_(iConfig.getParameter<double>("MinPtL1TkBarrel")),
116  useOfflineSeed_(iConfig.getUntrackedParameter<bool>("UseOfflineSeed", false)),
117  useUnassociatedL1_(iConfig.getParameter<bool>("UseUnassociatedL1")),
118  matchingDR_(iConfig.getParameter<std::vector<double>>("MatchDR")),
119  etaBins_(iConfig.getParameter<std::vector<double>>("EtaMatchingBins")) {
120  if (useOfflineSeed_) {
121  offlineSeedLabel_ = iConfig.getUntrackedParameter<InputTag>("OfflineSeedLabel");
122  offlineSeedToken_ = consumes<edm::View<TrajectorySeed>>(offlineSeedLabel_);
123 
124  // check that number of eta bins -1 matches number of dR cones
125  if (matchingDR_.size() != etaBins_.size() - 1) {
126  throw cms::Exception("Configuration") << "Size of MatchDR "
127  << "does not match number of eta bins." << endl;
128  }
129  }
130 
131  // service parameters
132  ParameterSet serviceParameters = iConfig.getParameter<ParameterSet>("ServiceParameters");
133 
134  // the services
135  service_ = std::make_unique<MuonServiceProxy>(serviceParameters, consumesCollector());
136 
137  // the estimator
138  estimator_ = std::make_unique<Chi2MeasurementEstimator>(10000.);
139 
140  produces<L2MuonTrajectorySeedCollection>();
141 }
142 
145  desc.add<edm::InputTag>("InputObjects", edm::InputTag("hltGmtStage2Digis"));
146  desc.add<string>("Propagator", "");
147  desc.add<double>("L1MinPt", -1.);
148  desc.add<double>("L1MaxEta", 5.0);
149  desc.add<double>("SetMinPtBarrelTo", 3.5);
150  desc.add<double>("SetMinPtEndcapTo", 1.0);
151  desc.add<double>("MinPL1Tk", 3.5);
152  desc.add<double>("MinPtL1TkBarrel", 3.5);
153  desc.addUntracked<bool>("UseOfflineSeed", false);
154  desc.add<bool>("UseUnassociatedL1", true);
155  desc.add<std::vector<double>>("MatchDR", {0.3});
156  desc.add<std::vector<double>>("EtaMatchingBins", {0., 2.5});
157  desc.addUntracked<edm::InputTag>("OfflineSeedLabel", edm::InputTag(""));
158 
160  psd0.addUntracked<std::vector<std::string>>("Propagators", {"SteppingHelixPropagatorAny"});
161  psd0.add<bool>("RPCLayers", true);
162  psd0.addUntracked<bool>("UseMuonNavigation", true);
163  desc.add<edm::ParameterSetDescription>("ServiceParameters", psd0);
164  descriptions.add("L2MuonSeedGeneratorFromL1TkMu", desc);
165 }
166 
168  const std::string metname = "Muon|RecoMuon|L2MuonSeedGeneratorFromL1TkMu";
170 
171  auto output = std::make_unique<L2MuonTrajectorySeedCollection>();
172 
174  iEvent.getByToken(muCollToken_, muColl);
175  LogDebug(metname) << "Number of muons " << muColl->size() << endl;
176 
177  edm::Handle<edm::View<TrajectorySeed>> offlineSeedHandle;
178  vector<int> offlineSeedMap;
179  if (useOfflineSeed_) {
180  iEvent.getByToken(offlineSeedToken_, offlineSeedHandle);
181  offlineSeedMap = vector<int>(offlineSeedHandle->size(), 0);
182  }
183 
184  for (auto const &tkmu : *muColl) {
185  // L1 tracker track
186  auto const &it = tkmu.trkPtr();
187 
188  // propagate the L1 tracker track to GMT
189  auto p3 = it->momentum();
190 
191  float tk_p = p3.mag();
192  if (tk_p < minPL1Tk_)
193  continue;
194 
195  float tk_pt = p3.perp();
196  float tk_eta = p3.eta();
197  float tk_aeta = std::abs(tk_eta);
198 
199  bool barrel = tk_aeta < etaBoundary_;
200  if (barrel && tk_pt < minPtL1TkBarrel_)
201  continue;
202 
203  float tk_phi = p3.phi();
204  float tk_q = it->rInv() > 0 ? 1. : -1.;
205  float tk_z = it->POCA().z();
206 
207  float dzCorrPhi = 1.;
208  float deta = 0;
209  float etaProp = tk_aeta;
210 
211  if (barrel) {
212  etaProp = etaBoundary_;
213  deta = tk_z / distMB2_ / cosh(tk_aeta);
214  } else {
215  float delta = tk_z / distME2_; //roughly scales as distance to 2nd station
216  if (tk_eta > 0)
217  delta *= -1;
218  dzCorrPhi = 1. + delta;
219 
220  float zOzs = tk_z / distME2_;
221  if (tk_eta > 0)
222  deta = zOzs / (1. - zOzs);
223  else
224  deta = zOzs / (1. + zOzs);
225  deta = deta * tanh(tk_eta);
226  }
227  float resPhi = tk_phi - phiCorr0_ * tk_q * cosh(phiCorr1_) / cosh(etaProp) / tk_pt * dzCorrPhi - M_PI / phiCorr2_;
228  resPhi = reco::reduceRange(resPhi);
229 
230  float pt = tk_pt; //not corrected for eloss
231  float eta = tk_eta + deta;
232  float theta = 2 * atan(exp(-eta));
233  float phi = resPhi;
234  int charge = it->rInv() > 0 ? 1 : -1;
235 
236  if (pt < l1MinPt_ || std::abs(eta) > l1MaxEta_)
237  continue;
238 
239  LogDebug(metname) << "New L2 Muon Seed";
240  LogDebug(metname) << "Pt = " << pt << " GeV/c";
241  LogDebug(metname) << "eta = " << eta;
242  LogDebug(metname) << "theta = " << theta << " rad";
243  LogDebug(metname) << "phi = " << phi << " rad";
244  LogDebug(metname) << "charge = " << charge;
245  LogDebug(metname) << "In Barrel? = " << barrel;
246 
247  // Update the services
248  service_->update(iSetup);
249 
250  const DetLayer *detLayer = nullptr;
251  float radius = 0.;
252 
253  CLHEP::Hep3Vector vec(0., 1., 0.);
254  vec.setTheta(theta);
255  vec.setPhi(phi);
256 
257  DetId theid;
258  // Get the det layer on which the state should be put
259  if (barrel) {
260  LogDebug(metname) << "The seed is in the barrel";
261 
262  // MB2
263  theid = DTChamberId(0, 2, 0);
264  detLayer = service_->detLayerGeometry()->idToLayer(theid);
265  LogDebug(metname) << "L2 Layer: " << debug.dumpLayer(detLayer);
266 
267  const BoundSurface *sur = &(detLayer->surface());
268  const BoundCylinder *bc = dynamic_cast<const BoundCylinder *>(sur);
269 
270  radius = std::abs(bc->radius() / sin(theta));
271 
272  LogDebug(metname) << "radius " << radius;
273 
274  if (pt < minPtBarrel_)
275  pt = minPtBarrel_;
276  } else {
277  LogDebug(metname) << "The seed is in the endcap";
278 
279  // ME2
280  theid = theta < Geom::pi() / 2. ? CSCDetId(1, 2, 0, 0, 0) : CSCDetId(2, 2, 0, 0, 0);
281 
282  detLayer = service_->detLayerGeometry()->idToLayer(theid);
283  LogDebug(metname) << "L2 Layer: " << debug.dumpLayer(detLayer);
284 
285  radius = std::abs(detLayer->position().z() / cos(theta));
286 
287  if (pt < minPtEndcap_)
288  pt = minPtEndcap_;
289  }
290 
291  vec.setMag(radius);
292 
293  GlobalPoint pos(vec.x(), vec.y(), vec.z());
294 
295  GlobalVector mom(pt * cos(phi), pt * sin(phi), pt * cos(theta) / sin(theta));
296 
297  GlobalTrajectoryParameters param(pos, mom, charge, &*service_->magneticField());
299 
300  mat[0][0] = (0.25 / pt) * (0.25 / pt); // sigma^2(charge/abs_momentum)
301  if (!barrel)
302  mat[0][0] = (0.4 / pt) * (0.4 / pt);
303 
304  mat[1][1] = 0.05 * 0.05; // sigma^2(lambda)
305  mat[2][2] = 0.2 * 0.2; // sigma^2(phi)
306  mat[3][3] = 20. * 20.; // sigma^2(x_transverse))
307  mat[4][4] = 20. * 20.; // sigma^2(y_transverse))
308 
310 
311  const FreeTrajectoryState state(param, error);
312 
313  LogDebug(metname) << "Free trajectory State from the parameters";
314  LogDebug(metname) << debug.dumpFTS(state);
315 
316  // Propagate the state on the MB2/ME2 surface
317  TrajectoryStateOnSurface tsos = service_->propagator(propagatorName_)->propagate(state, detLayer->surface());
318 
319  LogDebug(metname) << "State after the propagation on the layer";
320  LogDebug(metname) << debug.dumpLayer(detLayer);
321  LogDebug(metname) << debug.dumpTSOS(tsos);
322 
323  double dRcone = matchingDR_[0];
324  if (std::abs(eta) < etaBins_.back()) {
325  std::vector<double>::iterator lowEdge = std::upper_bound(etaBins_.begin(), etaBins_.end(), std::abs(eta));
326  dRcone = matchingDR_.at(lowEdge - etaBins_.begin() - 1);
327  }
328 
329  if (tsos.isValid()) {
331 
332  if (useOfflineSeed_) {
333  // Get the compatible dets on the layer
334  std::vector<pair<const GeomDet *, TrajectoryStateOnSurface>> detsWithStates =
335  detLayer->compatibleDets(tsos, *service_->propagator(propagatorName_), *estimator_);
336 
337  if (detsWithStates.empty() && barrel) {
338  // Fallback solution using ME2, try again to propagate but using ME2 as reference
339  DetId fallback_id;
340  theta < Geom::pi() / 2. ? fallback_id = CSCDetId(1, 2, 0, 0, 0) : fallback_id = CSCDetId(2, 2, 0, 0, 0);
341  const DetLayer *ME2DetLayer = service_->detLayerGeometry()->idToLayer(fallback_id);
342 
343  tsos = service_->propagator(propagatorName_)->propagate(state, ME2DetLayer->surface());
344  detsWithStates = ME2DetLayer->compatibleDets(tsos, *service_->propagator(propagatorName_), *estimator_);
345  }
346 
347  if (!detsWithStates.empty()) {
348  TrajectoryStateOnSurface newTSOS = detsWithStates.front().second;
349  const GeomDet *newTSOSDet = detsWithStates.front().first;
350 
351  LogDebug(metname) << "Most compatible det";
352  LogDebug(metname) << debug.dumpMuonId(newTSOSDet->geographicalId());
353 
354  if (newTSOS.isValid()) {
355  LogDebug(metname) << "pos: (r=" << newTSOS.globalPosition().mag()
356  << ", phi=" << newTSOS.globalPosition().phi()
357  << ", eta=" << newTSOS.globalPosition().eta() << ")";
358  LogDebug(metname) << "mom: (q*pt=" << newTSOS.charge() * newTSOS.globalMomentum().perp()
359  << ", phi=" << newTSOS.globalMomentum().phi()
360  << ", eta=" << newTSOS.globalMomentum().eta() << ")";
361 
362  const TrajectorySeed *assoOffseed =
363  associateOfflineSeedToL1(offlineSeedHandle, offlineSeedMap, newTSOS, dRcone);
364 
365  if (assoOffseed != nullptr) {
366  PTrajectoryStateOnDet const &seedTSOS = assoOffseed->startingState();
367  for (auto const &recHit : assoOffseed->recHits()) {
368  container.push_back(recHit);
369  }
370  auto dummyRef = edm::Ref<MuonBxCollection>();
371  output->emplace_back(L2MuonTrajectorySeed(seedTSOS, container, alongMomentum, dummyRef));
372  } else {
373  if (useUnassociatedL1_) {
374  // convert the TSOS into a PTSOD
375  PTrajectoryStateOnDet const &seedTSOS =
377  auto dummyRef = edm::Ref<MuonBxCollection>();
378  output->emplace_back(L2MuonTrajectorySeed(seedTSOS, container, alongMomentum, dummyRef));
379  }
380  }
381  }
382  }
383  } else {
384  // convert the TSOS into a PTSOD
386  auto dummyRef = edm::Ref<MuonBxCollection>();
387  output->emplace_back(L2MuonTrajectorySeed(seedTSOS, container, alongMomentum, dummyRef));
388  }
389  }
390  }
391 
392  iEvent.put(std::move(output));
393 }
394 
395 // FIXME: does not resolve ambiguities yet!
398  std::vector<int> &offseedMap,
399  TrajectoryStateOnSurface &newTsos,
400  double dRcone) {
401  if (dRcone < 0.)
402  return nullptr;
403 
404  const std::string metlabel = "Muon|RecoMuon|L2MuonSeedGeneratorFromL1TkMu";
405  MuonPatternRecoDumper debugtmp;
406 
407  edm::View<TrajectorySeed>::const_iterator offseed, endOffseed = offseeds->end();
408  const TrajectorySeed *selOffseed = nullptr;
409  double bestDr2 = 99999.;
410  unsigned int nOffseed(0);
411  int lastOffseed(-1);
412 
413  for (offseed = offseeds->begin(); offseed != endOffseed; ++offseed, ++nOffseed) {
414  if (offseedMap[nOffseed] != 0)
415  continue;
416  GlobalPoint glbPos = service_->trackingGeometry()
417  ->idToDet(offseed->startingState().detId())
418  ->surface()
419  .toGlobal(offseed->startingState().parameters().position());
420  GlobalVector glbMom = service_->trackingGeometry()
421  ->idToDet(offseed->startingState().detId())
422  ->surface()
423  .toGlobal(offseed->startingState().parameters().momentum());
424 
425  // Preliminary check
426  double preDr2 = deltaR2(newTsos.globalPosition().eta(), newTsos.globalPosition().phi(), glbPos.eta(), glbPos.phi());
427  if (preDr2 > 1.0)
428  continue;
429 
430  const FreeTrajectoryState offseedFTS(
431  glbPos, glbMom, offseed->startingState().parameters().charge(), &*service_->magneticField());
432  TrajectoryStateOnSurface offseedTsos =
433  service_->propagator(propagatorName_)->propagate(offseedFTS, newTsos.surface());
434  LogDebug(metlabel) << "Offline seed info: Det and State" << std::endl;
435  LogDebug(metlabel) << debugtmp.dumpMuonId(offseed->startingState().detId()) << std::endl;
436  LogDebug(metlabel) << "pos: (r=" << offseedFTS.position().mag() << ", phi=" << offseedFTS.position().phi()
437  << ", eta=" << offseedFTS.position().eta() << ")" << std::endl;
438  LogDebug(metlabel) << "mom: (q*pt=" << offseedFTS.charge() * offseedFTS.momentum().perp()
439  << ", phi=" << offseedFTS.momentum().phi() << ", eta=" << offseedFTS.momentum().eta() << ")"
440  << std::endl
441  << std::endl;
442 
443  if (offseedTsos.isValid()) {
444  LogDebug(metlabel) << "Offline seed info after propagation to L1 layer:" << std::endl;
445  LogDebug(metlabel) << "pos: (r=" << offseedTsos.globalPosition().mag()
446  << ", phi=" << offseedTsos.globalPosition().phi()
447  << ", eta=" << offseedTsos.globalPosition().eta() << ")" << std::endl;
448  LogDebug(metlabel) << "mom: (q*pt=" << offseedTsos.charge() * offseedTsos.globalMomentum().perp()
449  << ", phi=" << offseedTsos.globalMomentum().phi()
450  << ", eta=" << offseedTsos.globalMomentum().eta() << ")" << std::endl
451  << std::endl;
452  double newDr2 = deltaR2(newTsos.globalPosition().eta(),
453  newTsos.globalPosition().phi(),
454  offseedTsos.globalPosition().eta(),
455  offseedTsos.globalPosition().phi());
456  LogDebug(metlabel) << " -- DR = " << newDr2 << std::endl;
457  if (newDr2 < bestDr2 && newDr2 < dRcone * dRcone) {
458  LogDebug(metlabel) << " --> OK! " << newDr2 << std::endl << std::endl;
459  selOffseed = &*offseed;
460  bestDr2 = newDr2;
461  offseedMap[nOffseed] = 1;
462  if (lastOffseed > -1)
463  offseedMap[lastOffseed] = 0;
464  lastOffseed = nOffseed;
465  } else {
466  LogDebug(metlabel) << " --> Rejected. " << newDr2 << std::endl << std::endl;
467  }
468  } else {
469  LogDebug(metlabel) << "Invalid offline seed TSOS after propagation!" << std::endl << std::endl;
470  }
471  }
472 
473  return selOffseed;
474 }
475 
476 //define this as a plug-in
Vector3DBase
Definition: Vector3DBase.h:8
ConfigurationDescriptions.h
Chi2MeasurementEstimator.h
pfDeepBoostedJetPreprocessParams_cfi.upper_bound
upper_bound
Definition: pfDeepBoostedJetPreprocessParams_cfi.py:16
edm::View::begin
const_iterator begin() const
TrajectoryStateOnSurface.h
Handle.h
L2MuonSeedGeneratorFromL1TkMu::distMB2_
static constexpr float distMB2_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:89
electrons_cff.bool
bool
Definition: electrons_cff.py:366
FreeTrajectoryState.h
Reference_intrackfit_cff.barrel
list barrel
Definition: Reference_intrackfit_cff.py:37
GlobalTrajectoryParameters.h
MessageLogger.h
MuonPatternRecoDumper.h
funct::false
false
Definition: Factorize.h:29
GeomDet
Definition: GeomDet.h:27
DetLayer
Definition: DetLayer.h:21
MuonPatternRecoDumper::dumpMuonId
std::string dumpMuonId(const DetId &id) const
Definition: MuonPatternRecoDumper.cc:66
L2MuonSeedGeneratorFromL1TkMu::source_
edm::InputTag source_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:66
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:72
edm::EDGetTokenT< l1t::TkMuonCollection >
edm
HLT enums.
Definition: AlignableModifier.h:19
trajectoryStateTransform::persistentState
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
Definition: TrajectoryStateTransform.cc:14
TrajectoryStateOnSurface::globalPosition
GlobalPoint globalPosition() const
Definition: TrajectoryStateOnSurface.h:65
Muon.h
pos
Definition: PixelAliasList.h:18
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89285
L2MuonSeedGeneratorFromL1TkMu::minPL1Tk_
const double minPL1Tk_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:78
L2MuonSeedGeneratorFromL1TkMu::l1MaxEta_
const double l1MaxEta_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:75
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
TrajectoryStateOnSurface::charge
TrackCharge charge() const
Definition: TrajectoryStateOnSurface.h:68
L2MuonSeedGeneratorFromL1TkMu::minPtL1TkBarrel_
const double minPtL1TkBarrel_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:79
EDProducer.h
Surface
Definition: Surface.h:36
L2MuonSeedGeneratorFromL1TkMu::service_
std::unique_ptr< MuonServiceProxy > service_
the event setup proxy, it takes care the services update
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:96
GeometricSearchDet::compatibleDets
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
Definition: GeometricSearchDet.cc:35
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
BoundCylinder.h
L2MuonSeedGeneratorFromL1TkMu::useUnassociatedL1_
const bool useUnassociatedL1_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:81
L2MuonSeedGeneratorFromL1TkMu::muCollToken_
edm::EDGetTokenT< l1t::TkMuonCollection > muCollToken_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:67
GeometricSearchDet::position
virtual const Surface::PositionType & position() const
Returns position of the surface.
Definition: GeometricSearchDet.h:31
edm::Handle
Definition: AssociativeIterator.h:50
relativeConstraints.error
error
Definition: relativeConstraints.py:53
rpcPointValidation_cfi.recHit
recHit
Definition: rpcPointValidation_cfi.py:7
CSCDetId.h
JetPlusTrackCorrections_cff.dRcone
dRcone
Definition: JetPlusTrackCorrections_cff.py:38
edm::Ref< MuonBxCollection >
TkMuon.h
L2MuonSeedGeneratorFromL1TkMu::etaBoundary_
static constexpr float etaBoundary_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:88
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
L2MuonSeedGeneratorFromL1TkMu::L2MuonSeedGeneratorFromL1TkMu
L2MuonSeedGeneratorFromL1TkMu(const edm::ParameterSet &)
Constructor.
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:106
DetId
Definition: DetId.h:17
L2MuonSeedGeneratorFromL1TkMu::matchingDR_
std::vector< double > matchingDR_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:82
L2MuonSeedGeneratorFromL1TkMu::associateOfflineSeedToL1
const TrajectorySeed * associateOfflineSeedToL1(edm::Handle< edm::View< TrajectorySeed >> &, std::vector< int > &, TrajectoryStateOnSurface &, double)
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:396
TrajectoryStateOnSurface
Definition: TrajectoryStateOnSurface.h:16
L2MuonTrajectorySeedCollection.h
MakerMacros.h
debug
#define debug
Definition: HDRShower.cc:19
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
CurvilinearTrajectoryError
Definition: CurvilinearTrajectoryError.h:27
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
PTrajectoryStateOnDet.h
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
MuonPatternRecoDumper
Definition: MuonPatternRecoDumper.h:18
CurvilinearTrajectoryError.h
PVValHelper::eta
Definition: PVValidationHelpers.h:70
fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Geom::pi
constexpr double pi()
Definition: Pi.h:31
L2MuonSeedGeneratorFromL1TkMu::minPtEndcap_
const double minPtEndcap_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:77
theta
Geom::Theta< T > theta() const
Definition: Basic3DVectorLD.h:150
GlobalTrajectoryParameters
Definition: GlobalTrajectoryParameters.h:15
Point3DBase< float, GlobalTag >
L2MuonSeedGeneratorFromL1TkMu
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:54
DTChamberId.h
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
L2MuonSeedGeneratorFromL1TkMu::propagatorName_
std::string propagatorName_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:72
L2MuonSeedGeneratorFromL1TkMu::distME2_
static constexpr float distME2_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:90
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
edm::View
Definition: CaloClusterFwd.h:14
L2MuonSeedGeneratorFromL1TkMu::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:167
GeomDet::geographicalId
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
edm::ParameterSetDescription::addUntracked
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:100
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
edm::ParameterSet
Definition: ParameterSet.h:47
Geom::Phi::phi
T1 phi() const
Definition: Phi.h:78
Event.h
L2MuonSeedGeneratorFromL1TkMu::estimator_
std::unique_ptr< MeasurementEstimator > estimator_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:97
deltaR.h
l1t
delete x;
Definition: CaloConfig.h:22
dumpMFGeometry_cfg.delta
delta
Definition: dumpMFGeometry_cfg.py:25
CSCDetId
Definition: CSCDetId.h:26
L2MuonSeedGeneratorFromL1TkMu::phiCorr0_
static constexpr float phiCorr0_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:91
PV3DBase::eta
T eta() const
Definition: PV3DBase.h:73
iEvent
int iEvent
Definition: GenABIO.cc:224
L2MuonSeedGeneratorFromL1TkMu::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:143
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
edm::stream::EDProducer
Definition: EDProducer.h:38
edm::EventSetup
Definition: EventSetup.h:58
L2MuonSeedGeneratorFromL1TkMu::etaBins_
std::vector< double > etaBins_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:83
TrajectorySeed::recHits
RecHitRange recHits() const
Definition: TrajectorySeed.h:52
InputTag.h
DDAxes::phi
PV3DBase::mag
T mag() const
Definition: PV3DBase.h:64
L2MuonTrajectorySeed
Definition: L2MuonTrajectorySeed.h:17
L2MuonSeedGeneratorFromL1TkMu::l1MinPt_
const double l1MinPt_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:74
L2MuonSeedGeneratorFromL1TkMu::minPtBarrel_
const double minPtBarrel_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:76
TrajectoryStateOnSurface::globalMomentum
GlobalVector globalMomentum() const
Definition: TrajectoryStateOnSurface.h:66
L2MuonSeedGeneratorFromL1TkMu::useOfflineSeed_
const bool useOfflineSeed_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:80
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
DetId::rawId
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
RunInfoPI::state
state
Definition: RunInfoPayloadInspectoHelper.h:16
FreeTrajectoryState
Definition: FreeTrajectoryState.h:27
L2MuonSeedGeneratorFromL1TkMu::phiCorr1_
static constexpr float phiCorr1_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:92
HLTMuonOfflineAnalyzer_cfi.deltaR2
deltaR2
Definition: HLTMuonOfflineAnalyzer_cfi.py:105
Frameworkfwd.h
L2MuonSeedGeneratorFromL1TkMu::offlineSeedLabel_
edm::InputTag offlineSeedLabel_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:69
Exception
Definition: hltDiff.cc:245
MuonServiceProxy.h
CosmicsPD_Skims.radius
radius
Definition: CosmicsPD_Skims.py:135
GeometricSearchDet::surface
virtual const BoundSurface & surface() const =0
The surface of the GeometricSearchDet.
DetLayer.h
TrajectorySeed
Definition: TrajectorySeed.h:18
EventSetup.h
BoundCylinder
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
p3
double p3[4]
Definition: TauolaWrapper.h:91
TrajectoryStateTransform.h
GeomDetEnumerators.h
edm::OwnVector::push_back
void push_back(D *&d)
Definition: OwnVector.h:326
TrajectoryStateOnSurface::surface
const SurfaceType & surface() const
Definition: TrajectoryStateOnSurface.h:78
edm::View::const_iterator
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
TrajectorySeed::startingState
PTrajectoryStateOnDet const & startingState() const
Definition: TrajectorySeed.h:55
L2MuonSeedGeneratorFromL1TkMu::offlineSeedToken_
edm::EDGetTokenT< edm::View< TrajectorySeed > > offlineSeedToken_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:70
l1t::TkMuonCollection
std::vector< TkMuon > TkMuonCollection
Definition: TkMuonFwd.h:16
ConsumesCollector.h
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
DTChamberId
Definition: DTChamberId.h:14
PTrajectoryStateOnDet
Definition: PTrajectoryStateOnDet.h:10
ParameterSet.h
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
L2MuonTrajectorySeed.h
edm::Event
Definition: Event.h:73
AlgebraicSymMatrix55
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
Definition: AlgebraicROOTObjects.h:23
TkMuonFwd.h
PV3DBase::perp
T perp() const
Definition: PV3DBase.h:69
edm::InputTag
Definition: InputTag.h:15
alongMomentum
Definition: PropagationDirection.h:4
PV3DBase::phi
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
TrajectoryStateOnSurface::isValid
bool isValid() const
Definition: TrajectoryStateOnSurface.h:54
edm::OwnVector< TrackingRecHit >
L2MuonSeedGeneratorFromL1TkMu::phiCorr2_
static constexpr float phiCorr2_
Definition: L2MuonSeedGeneratorFromL1TkMu.cc:93
metname
const std::string metname
Definition: MuonSeedOrcaPatternRecognition.cc:40
reco::reduceRange
constexpr T reduceRange(T x)
Definition: deltaPhi.h:18