CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
L2MuonSeedGeneratorFromL1TkMu Class Reference
Inheritance diagram for L2MuonSeedGeneratorFromL1TkMu:
edm::stream::EDProducer<>

Public Member Functions

 L2MuonSeedGeneratorFromL1TkMu (const edm::ParameterSet &)
 Constructor. More...
 
void produce (edm::Event &, const edm::EventSetup &) override
 
 ~L2MuonSeedGeneratorFromL1TkMu () override=default
 Destructor. More...
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

const TrajectorySeedassociateOfflineSeedToL1 (edm::Handle< edm::View< TrajectorySeed >> &, std::vector< int > &, TrajectoryStateOnSurface &, double)
 

Private Attributes

std::unique_ptr< MeasurementEstimatorestimator_
 
std::vector< double > etaBins_
 
const double l1MaxEta_
 
const double l1MinPt_
 
std::vector< double > matchingDR_
 
const double minPL1Tk_
 
const double minPtBarrel_
 
const double minPtEndcap_
 
const double minPtL1TkBarrel_
 
edm::EDGetTokenT< l1t::TrackerMuonCollectionmuCollToken_
 
edm::InputTag offlineSeedLabel_
 
edm::EDGetTokenT< edm::View< TrajectorySeed > > offlineSeedToken_
 
std::string propagatorName_
 
std::unique_ptr< MuonServiceProxyservice_
 the event setup proxy, it takes care the services update More...
 
edm::InputTag source_
 
const bool useOfflineSeed_
 
const bool useUnassociatedL1_
 

Static Private Attributes

static constexpr float distMB2_ {550.}
 
static constexpr float distME2_ {850.}
 
static constexpr float etaBoundary_ {1.1}
 
static constexpr float phiCorr0_ {1.464}
 
static constexpr float phiCorr1_ {1.7}
 
static constexpr float phiCorr2_ {144.}
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 53 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Constructor & Destructor Documentation

◆ L2MuonSeedGeneratorFromL1TkMu()

L2MuonSeedGeneratorFromL1TkMu::L2MuonSeedGeneratorFromL1TkMu ( const edm::ParameterSet iConfig)
explicit

Constructor.

Definition at line 105 of file L2MuonSeedGeneratorFromL1TkMu.cc.

References estimator_, etaBins_, Exception, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), matchingDR_, offlineSeedLabel_, offlineSeedToken_, service_, and useOfflineSeed_.

106  : source_(iConfig.getParameter<InputTag>("InputObjects")),
107  muCollToken_(consumes(source_)),
108  propagatorName_(iConfig.getParameter<string>("Propagator")),
109  l1MinPt_(iConfig.getParameter<double>("L1MinPt")),
110  l1MaxEta_(iConfig.getParameter<double>("L1MaxEta")),
111  minPtBarrel_(iConfig.getParameter<double>("SetMinPtBarrelTo")),
112  minPtEndcap_(iConfig.getParameter<double>("SetMinPtEndcapTo")),
113  minPL1Tk_(iConfig.getParameter<double>("MinPL1Tk")),
114  minPtL1TkBarrel_(iConfig.getParameter<double>("MinPtL1TkBarrel")),
115  useOfflineSeed_(iConfig.getUntrackedParameter<bool>("UseOfflineSeed", false)),
116  useUnassociatedL1_(iConfig.getParameter<bool>("UseUnassociatedL1")),
117  matchingDR_(iConfig.getParameter<std::vector<double>>("MatchDR")),
118  etaBins_(iConfig.getParameter<std::vector<double>>("EtaMatchingBins")) {
119  if (useOfflineSeed_) {
120  offlineSeedLabel_ = iConfig.getUntrackedParameter<InputTag>("OfflineSeedLabel");
121  offlineSeedToken_ = consumes<edm::View<TrajectorySeed>>(offlineSeedLabel_);
122 
123  // check that number of eta bins -1 matches number of dR cones
124  if (matchingDR_.size() != etaBins_.size() - 1) {
125  throw cms::Exception("Configuration") << "Size of MatchDR "
126  << "does not match number of eta bins." << endl;
127  }
128  }
129 
130  // service parameters
131  ParameterSet serviceParameters = iConfig.getParameter<ParameterSet>("ServiceParameters");
132 
133  // the services
134  service_ = std::make_unique<MuonServiceProxy>(serviceParameters, consumesCollector());
135 
136  // the estimator
137  estimator_ = std::make_unique<Chi2MeasurementEstimator>(10000.);
138 
139  produces<L2MuonTrajectorySeedCollection>();
140 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< edm::View< TrajectorySeed > > offlineSeedToken_
T getUntrackedParameter(std::string const &, T const &) const
std::unique_ptr< MuonServiceProxy > service_
the event setup proxy, it takes care the services update
std::unique_ptr< MeasurementEstimator > estimator_
edm::EDGetTokenT< l1t::TrackerMuonCollection > muCollToken_

◆ ~L2MuonSeedGeneratorFromL1TkMu()

L2MuonSeedGeneratorFromL1TkMu::~L2MuonSeedGeneratorFromL1TkMu ( )
overridedefault

Destructor.

Member Function Documentation

◆ associateOfflineSeedToL1()

const TrajectorySeed * L2MuonSeedGeneratorFromL1TkMu::associateOfflineSeedToL1 ( edm::Handle< edm::View< TrajectorySeed >> &  offseeds,
std::vector< int > &  offseedMap,
TrajectoryStateOnSurface newTsos,
double  dRcone 
)
private

Definition at line 394 of file L2MuonSeedGeneratorFromL1TkMu.cc.

References edm::View< T >::begin(), HLTMuonOfflineAnalyzer_cfi::deltaR2, JetPlusTrackCorrections_cff::dRcone, MuonPatternRecoDumper::dumpMuonId(), PV3DBase< T, PVType, FrameType >::eta(), TrajectoryStateOnSurface::globalPosition(), LogDebug, PV3DBase< T, PVType, FrameType >::phi(), Geom::Phi< T1, Range >::phi(), propagatorName_, service_, AlCaHLTBitMon_QueryRunRegistry::string, and TrajectoryStateOnSurface::surface().

Referenced by produce().

398  {
399  if (dRcone < 0.)
400  return nullptr;
401 
402  const std::string metlabel = "Muon|RecoMuon|L2MuonSeedGeneratorFromL1TkMu";
403  MuonPatternRecoDumper debugtmp;
404 
405  edm::View<TrajectorySeed>::const_iterator offseed, endOffseed = offseeds->end();
406  const TrajectorySeed *selOffseed = nullptr;
407  double bestDr2 = 99999.;
408  unsigned int nOffseed(0);
409  int lastOffseed(-1);
410 
411  for (offseed = offseeds->begin(); offseed != endOffseed; ++offseed, ++nOffseed) {
412  if (offseedMap[nOffseed] != 0)
413  continue;
414  GlobalPoint glbPos = service_->trackingGeometry()
415  ->idToDet(offseed->startingState().detId())
416  ->surface()
417  .toGlobal(offseed->startingState().parameters().position());
418  GlobalVector glbMom = service_->trackingGeometry()
419  ->idToDet(offseed->startingState().detId())
420  ->surface()
421  .toGlobal(offseed->startingState().parameters().momentum());
422 
423  // Preliminary check
424  double preDr2 = deltaR2(newTsos.globalPosition().eta(), newTsos.globalPosition().phi(), glbPos.eta(), glbPos.phi());
425  if (preDr2 > 1.0)
426  continue;
427 
428  const FreeTrajectoryState offseedFTS(
429  glbPos, glbMom, offseed->startingState().parameters().charge(), &*service_->magneticField());
430  TrajectoryStateOnSurface offseedTsos =
431  service_->propagator(propagatorName_)->propagate(offseedFTS, newTsos.surface());
432  LogDebug(metlabel) << "Offline seed info: Det and State" << std::endl;
433  LogDebug(metlabel) << debugtmp.dumpMuonId(offseed->startingState().detId()) << std::endl;
434  LogDebug(metlabel) << "pos: (r=" << offseedFTS.position().mag() << ", phi=" << offseedFTS.position().phi()
435  << ", eta=" << offseedFTS.position().eta() << ")" << std::endl;
436  LogDebug(metlabel) << "mom: (q*pt=" << offseedFTS.charge() * offseedFTS.momentum().perp()
437  << ", phi=" << offseedFTS.momentum().phi() << ", eta=" << offseedFTS.momentum().eta() << ")"
438  << std::endl
439  << std::endl;
440 
441  if (offseedTsos.isValid()) {
442  LogDebug(metlabel) << "Offline seed info after propagation to L1 layer:" << std::endl;
443  LogDebug(metlabel) << "pos: (r=" << offseedTsos.globalPosition().mag()
444  << ", phi=" << offseedTsos.globalPosition().phi()
445  << ", eta=" << offseedTsos.globalPosition().eta() << ")" << std::endl;
446  LogDebug(metlabel) << "mom: (q*pt=" << offseedTsos.charge() * offseedTsos.globalMomentum().perp()
447  << ", phi=" << offseedTsos.globalMomentum().phi()
448  << ", eta=" << offseedTsos.globalMomentum().eta() << ")" << std::endl
449  << std::endl;
450  double newDr2 = deltaR2(newTsos.globalPosition().eta(),
451  newTsos.globalPosition().phi(),
452  offseedTsos.globalPosition().eta(),
453  offseedTsos.globalPosition().phi());
454  LogDebug(metlabel) << " -- DR = " << newDr2 << std::endl;
455  if (newDr2 < bestDr2 && newDr2 < dRcone * dRcone) {
456  LogDebug(metlabel) << " --> OK! " << newDr2 << std::endl << std::endl;
457  selOffseed = &*offseed;
458  bestDr2 = newDr2;
459  offseedMap[nOffseed] = 1;
460  if (lastOffseed > -1)
461  offseedMap[lastOffseed] = 0;
462  lastOffseed = nOffseed;
463  } else {
464  LogDebug(metlabel) << " --> Rejected. " << newDr2 << std::endl << std::endl;
465  }
466  } else {
467  LogDebug(metlabel) << "Invalid offline seed TSOS after propagation!" << std::endl << std::endl;
468  }
469  }
470 
471  return selOffseed;
472 }
std::string dumpMuonId(const DetId &id) const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
T eta() const
Definition: PV3DBase.h:73
const SurfaceType & surface() const
GlobalPoint globalPosition() const
std::unique_ptr< MuonServiceProxy > service_
the event setup proxy, it takes care the services update
T1 phi() const
Definition: Phi.h:78
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
#define LogDebug(id)

◆ fillDescriptions()

void L2MuonSeedGeneratorFromL1TkMu::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 142 of file L2MuonSeedGeneratorFromL1TkMu.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::addUntracked(), submitPVResolutionJobs::desc, and ProducerED_cfi::InputTag.

142  {
144  desc.add<edm::InputTag>("InputObjects", edm::InputTag("hltGmtStage2Digis"));
145  desc.add<string>("Propagator", "");
146  desc.add<double>("L1MinPt", -1.);
147  desc.add<double>("L1MaxEta", 5.0);
148  desc.add<double>("SetMinPtBarrelTo", 3.5);
149  desc.add<double>("SetMinPtEndcapTo", 1.0);
150  desc.add<double>("MinPL1Tk", 3.5);
151  desc.add<double>("MinPtL1TkBarrel", 3.5);
152  desc.addUntracked<bool>("UseOfflineSeed", false);
153  desc.add<bool>("UseUnassociatedL1", true);
154  desc.add<std::vector<double>>("MatchDR", {0.3});
155  desc.add<std::vector<double>>("EtaMatchingBins", {0., 2.5});
156  desc.addUntracked<edm::InputTag>("OfflineSeedLabel", edm::InputTag(""));
157 
159  psd0.addUntracked<std::vector<std::string>>("Propagators", {"SteppingHelixPropagatorAny"});
160  psd0.add<bool>("RPCLayers", true);
161  psd0.addUntracked<bool>("UseMuonNavigation", true);
162  desc.add<edm::ParameterSetDescription>("ServiceParameters", psd0);
163  descriptions.add("L2MuonSeedGeneratorFromL1TkMu", desc);
164 }
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ produce()

void L2MuonSeedGeneratorFromL1TkMu::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 166 of file L2MuonSeedGeneratorFromL1TkMu.cc.

References funct::abs(), alongMomentum, associateOfflineSeedToL1(), Reference_intrackfit_cff::barrel, ALCARECOTkAlJpsiMuMu_cff::charge, TrajectoryStateOnSurface::charge(), GeometricSearchDet::compatibleDets(), funct::cos(), debug, dumpMFGeometry_cfg::delta, distMB2_, distME2_, JetPlusTrackCorrections_cff::dRcone, relativeConstraints::error, estimator_, PVValHelper::eta, PV3DBase< T, PVType, FrameType >::eta(), etaBins_, etaBoundary_, JetChargeProducer_cfi::exp, GeomDet::geographicalId(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), iEvent, TrajectoryStateOnSurface::isValid(), l1MaxEta_, l1MinPt_, LogDebug, M_PI, PV3DBase< T, PVType, FrameType >::mag(), matchingDR_, metname, minPL1Tk_, minPtBarrel_, minPtEndcap_, minPtL1TkBarrel_, eostools::move(), muCollToken_, offlineSeedToken_, convertSQLitetoXML_cfg::output, chargedHadronTrackResolutionFilter_cfi::p3, PV3DBase< T, PVType, FrameType >::perp(), trajectoryStateTransform::persistentState(), phi, PV3DBase< T, PVType, FrameType >::phi(), phiCorr0_, phiCorr1_, phiCorr2_, Geom::pi(), GeometricSearchDet::position(), propagatorName_, DiDispStaMuonMonitor_cfi::pt, edm::OwnVector< T, P >::push_back(), CosmicsPD_Skims::radius, DetId::rawId(), rpcPointValidation_cfi::recHit, TrajectorySeed::recHits(), reco::reduceRange(), service_, funct::sin(), TrajectorySeed::startingState(), AlCaHLTBitMon_QueryRunRegistry::string, GeometricSearchDet::surface(), nnet::tanh(), theta(), pfDeepBoostedJetPreprocessParams_cfi::upper_bound, useOfflineSeed_, useUnassociatedL1_, and PV3DBase< T, PVType, FrameType >::z().

166  {
167  const std::string metname = "Muon|RecoMuon|L2MuonSeedGeneratorFromL1TkMu";
169 
170  auto output = std::make_unique<L2MuonTrajectorySeedCollection>();
171 
172  auto const muColl = iEvent.getHandle(muCollToken_);
173  LogDebug(metname) << "Number of muons " << muColl->size() << endl;
174 
175  edm::Handle<edm::View<TrajectorySeed>> offlineSeedHandle;
176  vector<int> offlineSeedMap;
177  if (useOfflineSeed_) {
178  iEvent.getByToken(offlineSeedToken_, offlineSeedHandle);
179  offlineSeedMap = vector<int>(offlineSeedHandle->size(), 0);
180  }
181 
182  for (auto const &tkmu : *muColl) {
183  // L1 tracker track
184  auto const &it = tkmu.trkPtr();
185 
186  // propagate the L1 tracker track to GMT
187  auto p3 = it->momentum();
188 
189  float tk_p = p3.mag();
190  if (tk_p < minPL1Tk_)
191  continue;
192 
193  float tk_pt = p3.perp();
194  float tk_eta = p3.eta();
195  float tk_aeta = std::abs(tk_eta);
196 
197  bool barrel = tk_aeta < etaBoundary_;
198  if (barrel && tk_pt < minPtL1TkBarrel_)
199  continue;
200 
201  float tk_phi = p3.phi();
202  float tk_q = it->rInv() > 0 ? 1. : -1.;
203  float tk_z = it->POCA().z();
204 
205  float dzCorrPhi = 1.;
206  float deta = 0;
207  float etaProp = tk_aeta;
208 
209  if (barrel) {
210  etaProp = etaBoundary_;
211  deta = tk_z / distMB2_ / cosh(tk_aeta);
212  } else {
213  float delta = tk_z / distME2_; //roughly scales as distance to 2nd station
214  if (tk_eta > 0)
215  delta *= -1;
216  dzCorrPhi = 1. + delta;
217 
218  float zOzs = tk_z / distME2_;
219  if (tk_eta > 0)
220  deta = zOzs / (1. - zOzs);
221  else
222  deta = zOzs / (1. + zOzs);
223  deta = deta * tanh(tk_eta);
224  }
225  float resPhi = tk_phi - phiCorr0_ * tk_q * cosh(phiCorr1_) / cosh(etaProp) / tk_pt * dzCorrPhi - M_PI / phiCorr2_;
226  resPhi = reco::reduceRange(resPhi);
227 
228  float pt = tk_pt; //not corrected for eloss
229  float eta = tk_eta + deta;
230  float theta = 2 * atan(exp(-eta));
231  float phi = resPhi;
232  int charge = it->rInv() > 0 ? 1 : -1;
233 
234  if (pt < l1MinPt_ || std::abs(eta) > l1MaxEta_)
235  continue;
236 
237  LogDebug(metname) << "New L2 Muon Seed";
238  LogDebug(metname) << "Pt = " << pt << " GeV/c";
239  LogDebug(metname) << "eta = " << eta;
240  LogDebug(metname) << "theta = " << theta << " rad";
241  LogDebug(metname) << "phi = " << phi << " rad";
242  LogDebug(metname) << "charge = " << charge;
243  LogDebug(metname) << "In Barrel? = " << barrel;
244 
245  // Update the services
246  service_->update(iSetup);
247 
248  const DetLayer *detLayer = nullptr;
249  float radius = 0.;
250 
251  CLHEP::Hep3Vector vec(0., 1., 0.);
252  vec.setTheta(theta);
253  vec.setPhi(phi);
254 
255  DetId theid;
256  // Get the det layer on which the state should be put
257  if (barrel) {
258  LogDebug(metname) << "The seed is in the barrel";
259 
260  // MB2
261  theid = DTChamberId(0, 2, 0);
262  detLayer = service_->detLayerGeometry()->idToLayer(theid);
263  LogDebug(metname) << "L2 Layer: " << debug.dumpLayer(detLayer);
264 
265  const BoundSurface *sur = &(detLayer->surface());
266  const BoundCylinder *bc = dynamic_cast<const BoundCylinder *>(sur);
267 
268  radius = std::abs(bc->radius() / sin(theta));
269 
270  LogDebug(metname) << "radius " << radius;
271 
272  if (pt < minPtBarrel_)
273  pt = minPtBarrel_;
274  } else {
275  LogDebug(metname) << "The seed is in the endcap";
276 
277  // ME2
278  theid = theta < Geom::pi() / 2. ? CSCDetId(1, 2, 0, 0, 0) : CSCDetId(2, 2, 0, 0, 0);
279 
280  detLayer = service_->detLayerGeometry()->idToLayer(theid);
281  LogDebug(metname) << "L2 Layer: " << debug.dumpLayer(detLayer);
282 
283  radius = std::abs(detLayer->position().z() / cos(theta));
284 
285  if (pt < minPtEndcap_)
286  pt = minPtEndcap_;
287  }
288 
289  vec.setMag(radius);
290 
291  GlobalPoint pos(vec.x(), vec.y(), vec.z());
292 
293  GlobalVector mom(pt * cos(phi), pt * sin(phi), pt * cos(theta) / sin(theta));
294 
295  GlobalTrajectoryParameters param(pos, mom, charge, &*service_->magneticField());
297 
298  mat[0][0] = (0.25 / pt) * (0.25 / pt); // sigma^2(charge/abs_momentum)
299  if (!barrel)
300  mat[0][0] = (0.4 / pt) * (0.4 / pt);
301 
302  mat[1][1] = 0.05 * 0.05; // sigma^2(lambda)
303  mat[2][2] = 0.2 * 0.2; // sigma^2(phi)
304  mat[3][3] = 20. * 20.; // sigma^2(x_transverse))
305  mat[4][4] = 20. * 20.; // sigma^2(y_transverse))
306 
308 
309  const FreeTrajectoryState state(param, error);
310 
311  LogDebug(metname) << "Free trajectory State from the parameters";
312  LogDebug(metname) << debug.dumpFTS(state);
313 
314  // Propagate the state on the MB2/ME2 surface
315  TrajectoryStateOnSurface tsos = service_->propagator(propagatorName_)->propagate(state, detLayer->surface());
316 
317  LogDebug(metname) << "State after the propagation on the layer";
318  LogDebug(metname) << debug.dumpLayer(detLayer);
319  LogDebug(metname) << debug.dumpTSOS(tsos);
320 
321  double dRcone = matchingDR_[0];
322  if (std::abs(eta) < etaBins_.back()) {
323  std::vector<double>::iterator lowEdge = std::upper_bound(etaBins_.begin(), etaBins_.end(), std::abs(eta));
324  dRcone = matchingDR_.at(lowEdge - etaBins_.begin() - 1);
325  }
326 
327  if (tsos.isValid()) {
329 
330  if (useOfflineSeed_) {
331  // Get the compatible dets on the layer
332  std::vector<pair<const GeomDet *, TrajectoryStateOnSurface>> detsWithStates =
333  detLayer->compatibleDets(tsos, *service_->propagator(propagatorName_), *estimator_);
334 
335  if (detsWithStates.empty() && barrel) {
336  // Fallback solution using ME2, try again to propagate but using ME2 as reference
337  DetId fallback_id;
338  theta < Geom::pi() / 2. ? fallback_id = CSCDetId(1, 2, 0, 0, 0) : fallback_id = CSCDetId(2, 2, 0, 0, 0);
339  const DetLayer *ME2DetLayer = service_->detLayerGeometry()->idToLayer(fallback_id);
340 
341  tsos = service_->propagator(propagatorName_)->propagate(state, ME2DetLayer->surface());
342  detsWithStates = ME2DetLayer->compatibleDets(tsos, *service_->propagator(propagatorName_), *estimator_);
343  }
344 
345  if (!detsWithStates.empty()) {
346  TrajectoryStateOnSurface newTSOS = detsWithStates.front().second;
347  const GeomDet *newTSOSDet = detsWithStates.front().first;
348 
349  LogDebug(metname) << "Most compatible det";
350  LogDebug(metname) << debug.dumpMuonId(newTSOSDet->geographicalId());
351 
352  if (newTSOS.isValid()) {
353  LogDebug(metname) << "pos: (r=" << newTSOS.globalPosition().mag()
354  << ", phi=" << newTSOS.globalPosition().phi()
355  << ", eta=" << newTSOS.globalPosition().eta() << ")";
356  LogDebug(metname) << "mom: (q*pt=" << newTSOS.charge() * newTSOS.globalMomentum().perp()
357  << ", phi=" << newTSOS.globalMomentum().phi()
358  << ", eta=" << newTSOS.globalMomentum().eta() << ")";
359 
360  const TrajectorySeed *assoOffseed =
361  associateOfflineSeedToL1(offlineSeedHandle, offlineSeedMap, newTSOS, dRcone);
362 
363  if (assoOffseed != nullptr) {
364  PTrajectoryStateOnDet const &seedTSOS = assoOffseed->startingState();
365  for (auto const &recHit : assoOffseed->recHits()) {
366  container.push_back(recHit);
367  }
368  auto dummyRef = edm::Ref<MuonBxCollection>();
369  output->emplace_back(L2MuonTrajectorySeed(seedTSOS, container, alongMomentum, dummyRef));
370  } else {
371  if (useUnassociatedL1_) {
372  // convert the TSOS into a PTSOD
373  PTrajectoryStateOnDet const &seedTSOS =
375  auto dummyRef = edm::Ref<MuonBxCollection>();
376  output->emplace_back(L2MuonTrajectorySeed(seedTSOS, container, alongMomentum, dummyRef));
377  }
378  }
379  }
380  }
381  } else {
382  // convert the TSOS into a PTSOD
384  auto dummyRef = edm::Ref<MuonBxCollection>();
385  output->emplace_back(L2MuonTrajectorySeed(seedTSOS, container, alongMomentum, dummyRef));
386  }
387  }
388  }
389 
390  iEvent.put(std::move(output));
391 }
virtual const BoundSurface & surface() const =0
The surface of the GeometricSearchDet.
virtual const Surface::PositionType & position() const
Returns position of the surface.
T perp() const
Definition: PV3DBase.h:69
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
void tanh(data_T data[CONFIG_T::n_in], res_T res[CONFIG_T::n_in])
constexpr T reduceRange(T x)
Definition: deltaPhi.h:18
T z() const
Definition: PV3DBase.h:61
const std::string metname
RecHitRange recHits() const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T eta() const
Definition: PV3DBase.h:73
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
edm::EDGetTokenT< edm::View< TrajectorySeed > > offlineSeedToken_
const TrajectorySeed * associateOfflineSeedToL1(edm::Handle< edm::View< TrajectorySeed >> &, std::vector< int > &, TrajectoryStateOnSurface &, double)
void push_back(D *&d)
Definition: OwnVector.h:326
int iEvent
Definition: GenABIO.cc:224
GlobalPoint globalPosition() const
std::unique_ptr< MuonServiceProxy > service_
the event setup proxy, it takes care the services update
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::unique_ptr< MeasurementEstimator > estimator_
T mag() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
PTrajectoryStateOnDet const & startingState() const
#define M_PI
Definition: DetId.h:17
#define debug
Definition: HDRShower.cc:19
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
GlobalVector globalMomentum() const
edm::EDGetTokenT< l1t::TrackerMuonCollection > muCollToken_
constexpr double pi()
Definition: Pi.h:31
Geom::Theta< T > theta() const
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)

Member Data Documentation

◆ distMB2_

constexpr float L2MuonSeedGeneratorFromL1TkMu::distMB2_ {550.}
staticprivate

Definition at line 88 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ distME2_

constexpr float L2MuonSeedGeneratorFromL1TkMu::distME2_ {850.}
staticprivate

Definition at line 89 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ estimator_

std::unique_ptr<MeasurementEstimator> L2MuonSeedGeneratorFromL1TkMu::estimator_
private

Definition at line 96 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by L2MuonSeedGeneratorFromL1TkMu(), and produce().

◆ etaBins_

std::vector<double> L2MuonSeedGeneratorFromL1TkMu::etaBins_
private

Definition at line 82 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by L2MuonSeedGeneratorFromL1TkMu(), and produce().

◆ etaBoundary_

constexpr float L2MuonSeedGeneratorFromL1TkMu::etaBoundary_ {1.1}
staticprivate

Definition at line 87 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ l1MaxEta_

const double L2MuonSeedGeneratorFromL1TkMu::l1MaxEta_
private

Definition at line 74 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ l1MinPt_

const double L2MuonSeedGeneratorFromL1TkMu::l1MinPt_
private

Definition at line 73 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ matchingDR_

std::vector<double> L2MuonSeedGeneratorFromL1TkMu::matchingDR_
private

Definition at line 81 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by L2MuonSeedGeneratorFromL1TkMu(), and produce().

◆ minPL1Tk_

const double L2MuonSeedGeneratorFromL1TkMu::minPL1Tk_
private

Definition at line 77 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ minPtBarrel_

const double L2MuonSeedGeneratorFromL1TkMu::minPtBarrel_
private

Definition at line 75 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ minPtEndcap_

const double L2MuonSeedGeneratorFromL1TkMu::minPtEndcap_
private

Definition at line 76 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ minPtL1TkBarrel_

const double L2MuonSeedGeneratorFromL1TkMu::minPtL1TkBarrel_
private

Definition at line 78 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ muCollToken_

edm::EDGetTokenT<l1t::TrackerMuonCollection> L2MuonSeedGeneratorFromL1TkMu::muCollToken_
private

Definition at line 66 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ offlineSeedLabel_

edm::InputTag L2MuonSeedGeneratorFromL1TkMu::offlineSeedLabel_
private

Definition at line 68 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by L2MuonSeedGeneratorFromL1TkMu().

◆ offlineSeedToken_

edm::EDGetTokenT<edm::View<TrajectorySeed> > L2MuonSeedGeneratorFromL1TkMu::offlineSeedToken_
private

Definition at line 69 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by L2MuonSeedGeneratorFromL1TkMu(), and produce().

◆ phiCorr0_

constexpr float L2MuonSeedGeneratorFromL1TkMu::phiCorr0_ {1.464}
staticprivate

Definition at line 90 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ phiCorr1_

constexpr float L2MuonSeedGeneratorFromL1TkMu::phiCorr1_ {1.7}
staticprivate

Definition at line 91 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ phiCorr2_

constexpr float L2MuonSeedGeneratorFromL1TkMu::phiCorr2_ {144.}
staticprivate

Definition at line 92 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().

◆ propagatorName_

std::string L2MuonSeedGeneratorFromL1TkMu::propagatorName_
private

Definition at line 71 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by associateOfflineSeedToL1(), and produce().

◆ service_

std::unique_ptr<MuonServiceProxy> L2MuonSeedGeneratorFromL1TkMu::service_
private

the event setup proxy, it takes care the services update

Definition at line 95 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by associateOfflineSeedToL1(), L2MuonSeedGeneratorFromL1TkMu(), and produce().

◆ source_

edm::InputTag L2MuonSeedGeneratorFromL1TkMu::source_
private

◆ useOfflineSeed_

const bool L2MuonSeedGeneratorFromL1TkMu::useOfflineSeed_
private

Definition at line 79 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by L2MuonSeedGeneratorFromL1TkMu(), and produce().

◆ useUnassociatedL1_

const bool L2MuonSeedGeneratorFromL1TkMu::useUnassociatedL1_
private

Definition at line 80 of file L2MuonSeedGeneratorFromL1TkMu.cc.

Referenced by produce().