CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
L2MuonSeedGenerator Class Reference

#include <L2MuonSeedGenerator.h>

Inheritance diagram for L2MuonSeedGenerator:
edm::stream::EDProducer<>

Public Member Functions

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

Private Member Functions

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

Private Attributes

edm::EDGetTokenT< L1MuGMTReadoutCollectiongmtToken_
 
edm::EDGetTokenT< l1extra::L1MuonParticleCollectionmuCollToken_
 
edm::EDGetTokenT< edm::View< TrajectorySeed > > offlineSeedToken_
 
MeasurementEstimatortheEstimator
 
edm::InputTag theL1GMTReadoutCollection
 
const double theL1MaxEta
 
const double theL1MinPt
 
const unsigned theL1MinQuality
 
edm::InputTag theOfflineSeedLabel
 
std::string thePropagatorName
 
MuonServiceProxytheService
 the event setup proxy, it takes care the services update More...
 
edm::InputTag theSource
 
const bool useOfflineSeed
 
const bool useUnassociatedL1
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

L2 muon seed generator: Transform the L1 informations in seeds for the L2 muon reconstruction

Author
A.Everett, R.Bellan, J. Alcaraz

ORCA's author: N. Neumeister

L2 muon seed generator: Transform the L1 informations in seeds for the L2 muon reconstruction

Author
A.Everett, R.Bellan

ORCA's author: N. Neumeister

Definition at line 54 of file L2MuonSeedGenerator.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 43 of file L2MuonSeedGenerator.cc.

References Chi2MeasurementEstimator_cfi::Chi2MeasurementEstimator, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), gmtToken_, muCollToken_, MuonServiceProxy_cff::MuonServiceProxy, offlineSeedToken_, theEstimator, theL1GMTReadoutCollection, theOfflineSeedLabel, theService, theSource, and useOfflineSeed.

44  : theSource(iConfig.getParameter<InputTag>("InputObjects")),
45  theL1GMTReadoutCollection(iConfig.getParameter<InputTag>("GMTReadoutCollection")),
46  thePropagatorName(iConfig.getParameter<string>("Propagator")),
47  theL1MinPt(iConfig.getParameter<double>("L1MinPt")),
48  theL1MaxEta(iConfig.getParameter<double>("L1MaxEta")),
49  theL1MinQuality(iConfig.getParameter<unsigned int>("L1MinQuality")),
50  useOfflineSeed(iConfig.getUntrackedParameter<bool>("UseOfflineSeed", false)),
51  useUnassociatedL1(iConfig.existsAs<bool>("UseUnassociatedL1") ? iConfig.getParameter<bool>("UseUnassociatedL1")
52  : true) {
53  gmtToken_ = consumes<L1MuGMTReadoutCollection>(theL1GMTReadoutCollection);
54  muCollToken_ = consumes<L1MuonParticleCollection>(theSource);
55 
56  if (useOfflineSeed) {
57  theOfflineSeedLabel = iConfig.getUntrackedParameter<InputTag>("OfflineSeedLabel");
58  offlineSeedToken_ = consumes<edm::View<TrajectorySeed> >(theOfflineSeedLabel);
59  }
60 
61  // service parameters
62  ParameterSet serviceParameters = iConfig.getParameter<ParameterSet>("ServiceParameters");
63 
64  // the services
65  theService = new MuonServiceProxy(serviceParameters);
66 
67  // the estimator
69 
70  produces<L2MuonTrajectorySeedCollection>();
71 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:160
edm::EDGetTokenT< l1extra::L1MuonParticleCollection > muCollToken_
edm::InputTag theL1GMTReadoutCollection
edm::EDGetTokenT< edm::View< TrajectorySeed > > offlineSeedToken_
const unsigned theL1MinQuality
MeasurementEstimator * theEstimator
MuonServiceProxy * theService
the event setup proxy, it takes care the services update
edm::EDGetTokenT< L1MuGMTReadoutCollection > gmtToken_
edm::InputTag theOfflineSeedLabel
L2MuonSeedGenerator::~L2MuonSeedGenerator ( )
override

Destructor.

Definition at line 74 of file L2MuonSeedGenerator.cc.

References theEstimator, and theService.

74  {
75  if (theService)
76  delete theService;
77  if (theEstimator)
78  delete theEstimator;
79 }
MeasurementEstimator * theEstimator
MuonServiceProxy * theService
the event setup proxy, it takes care the services update

Member Function Documentation

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

Definition at line 325 of file L2MuonSeedGenerator.cc.

References edm::View< T >::begin(), PbPb_ZMuSkimMuonDPG_cff::deltaR, MuonPatternRecoDumper::dumpMuonId(), PV3DBase< T, PVType, FrameType >::eta(), TrajectoryStateOnSurface::globalPosition(), LogDebug, PV3DBase< T, PVType, FrameType >::phi(), Geom::Phi< T1, Range >::phi(), AlCaHLTBitMon_QueryRunRegistry::string, TrajectoryStateOnSurface::surface(), thePropagatorName, and theService.

Referenced by produce().

327  {
328  const std::string metlabel = "Muon|RecoMuon|L2MuonSeedGenerator";
329  MuonPatternRecoDumper debugtmp;
330 
331  edm::View<TrajectorySeed>::const_iterator offseed, endOffseed = offseeds->end();
332  const TrajectorySeed* selOffseed = nullptr;
333  double bestDr = 99999.;
334  unsigned int nOffseed(0);
335  int lastOffseed(-1);
336 
337  for (offseed = offseeds->begin(); offseed != endOffseed; ++offseed, ++nOffseed) {
338  if (offseedMap[nOffseed] != 0)
339  continue;
340  GlobalPoint glbPos = theService->trackingGeometry()
341  ->idToDet(offseed->startingState().detId())
342  ->surface()
343  .toGlobal(offseed->startingState().parameters().position());
344  GlobalVector glbMom = theService->trackingGeometry()
345  ->idToDet(offseed->startingState().detId())
346  ->surface()
347  .toGlobal(offseed->startingState().parameters().momentum());
348 
349  // Preliminary check
350  double preDr = deltaR(newTsos.globalPosition().eta(), newTsos.globalPosition().phi(), glbPos.eta(), glbPos.phi());
351  if (preDr > 1.0)
352  continue;
353 
354  const FreeTrajectoryState offseedFTS(
355  glbPos, glbMom, offseed->startingState().parameters().charge(), &*theService->magneticField());
356  TrajectoryStateOnSurface offseedTsos =
357  theService->propagator(thePropagatorName)->propagate(offseedFTS, newTsos.surface());
358  LogDebug(metlabel) << "Offline seed info: Det and State" << std::endl;
359  LogDebug(metlabel) << debugtmp.dumpMuonId(offseed->startingState().detId()) << std::endl;
360  //LogDebug(metlabel) << "(x, y, z) = (" << newTSOS.globalPosition().x() << ", "
361  // << newTSOS.globalPosition().y() << ", " << newTSOS.globalPosition().z() << ")" << std::endl;
362  LogDebug(metlabel) << "pos: (r=" << offseedFTS.position().mag() << ", phi=" << offseedFTS.position().phi()
363  << ", eta=" << offseedFTS.position().eta() << ")" << std::endl;
364  LogDebug(metlabel) << "mom: (q*pt=" << offseedFTS.charge() * offseedFTS.momentum().perp()
365  << ", phi=" << offseedFTS.momentum().phi() << ", eta=" << offseedFTS.momentum().eta() << ")"
366  << std::endl
367  << std::endl;
368  //LogDebug(metlabel) << debugtmp.dumpFTS(offseedFTS) << std::endl;
369 
370  if (offseedTsos.isValid()) {
371  LogDebug(metlabel) << "Offline seed info after propagation to L1 layer:" << std::endl;
372  //LogDebug(metlabel) << "(x, y, z) = (" << offseedTsos.globalPosition().x() << ", "
373  // << offseedTsos.globalPosition().y() << ", " << offseedTsos.globalPosition().z() << ")" << std::endl;
374  LogDebug(metlabel) << "pos: (r=" << offseedTsos.globalPosition().mag()
375  << ", phi=" << offseedTsos.globalPosition().phi()
376  << ", eta=" << offseedTsos.globalPosition().eta() << ")" << std::endl;
377  LogDebug(metlabel) << "mom: (q*pt=" << offseedTsos.charge() * offseedTsos.globalMomentum().perp()
378  << ", phi=" << offseedTsos.globalMomentum().phi()
379  << ", eta=" << offseedTsos.globalMomentum().eta() << ")" << std::endl
380  << std::endl;
381  //LogDebug(metlabel) << debugtmp.dumpTSOS(offseedTsos) << std::endl;
382  double newDr = deltaR(newTsos.globalPosition().eta(),
383  newTsos.globalPosition().phi(),
384  offseedTsos.globalPosition().eta(),
385  offseedTsos.globalPosition().phi());
386  LogDebug(metlabel) << " -- DR = " << newDr << std::endl;
387  if (newDr < 0.3 && newDr < bestDr) {
388  LogDebug(metlabel) << " --> OK! " << newDr << std::endl << std::endl;
389  selOffseed = &*offseed;
390  bestDr = newDr;
391  offseedMap[nOffseed] = 1;
392  if (lastOffseed > -1)
393  offseedMap[lastOffseed] = 0;
394  lastOffseed = nOffseed;
395  } else {
396  LogDebug(metlabel) << " --> Rejected. " << newDr << std::endl << std::endl;
397  }
398  } else {
399  LogDebug(metlabel) << "Invalid offline seed TSOS after propagation!" << std::endl << std::endl;
400  }
401  }
402 
403  return selOffseed;
404 }
#define LogDebug(id)
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
GlobalPoint globalPosition() const
T1 phi() const
Definition: Phi.h:78
std::string dumpMuonId(const DetId &id) const
const SurfaceType & surface() const
T eta() const
Definition: PV3DBase.h:73
MuonServiceProxy * theService
the event setup proxy, it takes care the services update
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
void L2MuonSeedGenerator::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 81 of file L2MuonSeedGenerator.cc.

References alongMomentum, associateOfflineSeedToL1(), Reference_intrackfit_cff::barrel, ALCARECOTkAlJpsiMuMu_cff::charge, TrajectoryStateOnSurface::charge(), L1MuGMTCand::charge_valid(), GeometricSearchDet::compatibleDets(), funct::cos(), debug, MuonPatternRecoDumper::dumpFTS(), MuonPatternRecoDumper::dumpLayer(), MuonPatternRecoDumper::dumpMuonId(), L1MuGMTCand::empty(), relativeConstraints::error, PVValHelper::eta, PV3DBase< T, PVType, FrameType >::eta(), JetChargeProducer_cfi::exp, GeomDet::geographicalId(), L1MuGMTReadoutRecord::getBrlRPCCands(), edm::Event::getByToken(), L1MuGMTReadoutRecord::getCSCCands(), L1MuGMTReadoutRecord::getDTBXCands(), L1MuGMTExtendedCand::getDTCSCIndex(), L1MuGMTReadoutRecord::getFwdRPCCands(), L1MuGMTExtendedCand::getRPCIndex(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), gmtToken_, triggerObjects_cff::id, training_settings::idx, L1MuGMTExtendedCand::isFwd(), L1MuGMTExtendedCand::isRPC(), TrajectoryStateOnSurface::isValid(), LogDebug, LogTrace, PV3DBase< T, PVType, FrameType >::mag(), metname, eostools::move(), muCollToken_, offlineSeedToken_, convertSQLitetoXML_cfg::output, PV3DBase< T, PVType, FrameType >::perp(), trajectoryStateTransform::persistentState(), phi, PV3DBase< T, PVType, FrameType >::phi(), Geom::pi(), GeometricSearchDet::position(), DiDispStaMuonMonitor_cfi::pt, edm::OwnVector< T, P >::push_back(), edm::Event::put(), qcdUeDQM_cfi::quality, L1MuGMTCand::quality(), CosmicsPD_Skims::radius, DetId::rawId(), TrajectorySeed::recHits(), funct::sin(), TrajectorySeed::startingState(), AlCaHLTBitMon_QueryRunRegistry::string, GeometricSearchDet::surface(), theEstimator, theL1MaxEta, theL1MinPt, theL1MinQuality, thePropagatorName, theService, theta(), useOfflineSeed, useUnassociatedL1, and PV3DBase< T, PVType, FrameType >::z().

81  {
82  const std::string metname = "Muon|RecoMuon|L2MuonSeedGenerator";
84 
85  auto output = std::make_unique<L2MuonTrajectorySeedCollection>();
86 
87  // Muon particles and GMT readout collection
89  iEvent.getByToken(gmtToken_, gmtrc_handle);
90  L1MuGMTReadoutRecord const& gmtrr = gmtrc_handle.product()->getRecord(0);
91 
93  iEvent.getByToken(muCollToken_, muColl);
94  LogTrace(metname) << "Number of muons " << muColl->size() << endl;
95 
96  edm::Handle<edm::View<TrajectorySeed> > offlineSeedHandle;
97  vector<int> offlineSeedMap;
98  if (useOfflineSeed) {
99  iEvent.getByToken(offlineSeedToken_, offlineSeedHandle);
100  LogTrace(metname) << "Number of offline seeds " << offlineSeedHandle->size() << endl;
101  offlineSeedMap = vector<int>(offlineSeedHandle->size(), 0);
102  }
103 
104  L1MuonParticleCollection::const_iterator it;
105  L1MuonParticleRef::key_type l1ParticleIndex = 0;
106 
107  for (it = muColl->begin(); it != muColl->end(); ++it, ++l1ParticleIndex) {
108  const L1MuGMTExtendedCand muonCand = (*it).gmtMuonCand();
109  unsigned int quality = 0;
110  bool valid_charge = false;
111  ;
112 
113  if (muonCand.empty()) {
114  LogWarning(metname) << "L2MuonSeedGenerator: WARNING, no L1MuGMTCand! " << endl;
115  LogWarning(metname) << "L2MuonSeedGenerator: this should make sense only within MC tests" << endl;
116  // FIXME! Temporary to handle the MC input
117  quality = 7;
118  valid_charge = true;
119  } else {
120  quality = muonCand.quality();
121  valid_charge = muonCand.charge_valid();
122  }
123 
124  float pt = (*it).pt();
125  float eta = (*it).eta();
126  float theta = 2 * atan(exp(-eta));
127  float phi = (*it).phi();
128  int charge = (*it).charge();
129  // Set charge=0 for the time being if the valid charge bit is zero
130  if (!valid_charge)
131  charge = 0;
132  bool barrel = !(*it).isForward();
133 
134  // Get a better eta and charge from regional information
135  // Phi has the same resolution in GMT than regionally, is not it?
136  if (!(muonCand.empty())) {
137  int idx = -1;
138  vector<L1MuRegionalCand> rmc;
139  if (!muonCand.isRPC()) {
140  idx = muonCand.getDTCSCIndex();
141  if (muonCand.isFwd())
142  rmc = gmtrr.getCSCCands();
143  else
144  rmc = gmtrr.getDTBXCands();
145  } else {
146  idx = muonCand.getRPCIndex();
147  if (muonCand.isFwd())
148  rmc = gmtrr.getFwdRPCCands();
149  else
150  rmc = gmtrr.getBrlRPCCands();
151  }
152  if (idx >= 0) {
153  eta = rmc[idx].etaValue();
154  //phi = rmc[idx].phiValue();
155  // Use this charge if the valid charge bit is zero
156  if (!valid_charge)
157  charge = rmc[idx].chargeValue();
158  }
159  }
160 
161  if (pt < theL1MinPt || fabs(eta) > theL1MaxEta)
162  continue;
163 
164  LogTrace(metname) << "New L2 Muon Seed";
165  LogTrace(metname) << "Pt = " << pt << " GeV/c";
166  LogTrace(metname) << "eta = " << eta;
167  LogTrace(metname) << "theta = " << theta << " rad";
168  LogTrace(metname) << "phi = " << phi << " rad";
169  LogTrace(metname) << "charge = " << charge;
170  LogTrace(metname) << "In Barrel? = " << barrel;
171 
172  if (quality <= theL1MinQuality)
173  continue;
174  LogTrace(metname) << "quality = " << quality;
175 
176  // Update the services
177  theService->update(iSetup);
178 
179  const DetLayer* detLayer = nullptr;
180  float radius = 0.;
181 
182  CLHEP::Hep3Vector vec(0., 1., 0.);
183  vec.setTheta(theta);
184  vec.setPhi(phi);
185 
186  // Get the det layer on which the state should be put
187  if (barrel) {
188  LogTrace(metname) << "The seed is in the barrel";
189 
190  // MB2
191  DetId id = DTChamberId(0, 2, 0);
192  detLayer = theService->detLayerGeometry()->idToLayer(id);
193  LogTrace(metname) << "L2 Layer: " << debug.dumpLayer(detLayer);
194 
195  const BoundSurface* sur = &(detLayer->surface());
196  const BoundCylinder* bc = dynamic_cast<const BoundCylinder*>(sur);
197 
198  radius = fabs(bc->radius() / sin(theta));
199 
200  LogTrace(metname) << "radius " << radius;
201 
202  if (pt < 3.5)
203  pt = 3.5;
204  } else {
205  LogTrace(metname) << "The seed is in the endcap";
206 
207  DetId id;
208  // ME2
209  if (theta < Geom::pi() / 2.)
210  id = CSCDetId(1, 2, 0, 0, 0);
211  else
212  id = CSCDetId(2, 2, 0, 0, 0);
213 
214  detLayer = theService->detLayerGeometry()->idToLayer(id);
215  LogTrace(metname) << "L2 Layer: " << debug.dumpLayer(detLayer);
216 
217  radius = fabs(detLayer->position().z() / cos(theta));
218 
219  if (pt < 1.0)
220  pt = 1.0;
221  }
222 
223  vec.setMag(radius);
224 
225  GlobalPoint pos(vec.x(), vec.y(), vec.z());
226 
227  GlobalVector mom(pt * cos(phi), pt * sin(phi), pt * cos(theta) / sin(theta));
228 
229  GlobalTrajectoryParameters param(pos, mom, charge, &*theService->magneticField());
231 
232  mat[0][0] = (0.25 / pt) * (0.25 / pt); // sigma^2(charge/abs_momentum)
233  if (!barrel)
234  mat[0][0] = (0.4 / pt) * (0.4 / pt);
235 
236  //Assign q/pt = 0 +- 1/pt if charge has been declared invalid
237  if (!valid_charge)
238  mat[0][0] = (1. / pt) * (1. / pt);
239 
240  mat[1][1] = 0.05 * 0.05; // sigma^2(lambda)
241  mat[2][2] = 0.2 * 0.2; // sigma^2(phi)
242  mat[3][3] = 20. * 20.; // sigma^2(x_transverse))
243  mat[4][4] = 20. * 20.; // sigma^2(y_transverse))
244 
246 
247  const FreeTrajectoryState state(param, error);
248 
249  LogTrace(metname) << "Free trajectory State from the parameters";
250  LogTrace(metname) << debug.dumpFTS(state);
251 
252  // Propagate the state on the MB2/ME2 surface
253  TrajectoryStateOnSurface tsos = theService->propagator(thePropagatorName)->propagate(state, detLayer->surface());
254 
255  LogTrace(metname) << "State after the propagation on the layer";
256  LogTrace(metname) << debug.dumpLayer(detLayer);
257  LogTrace(metname) << debug.dumpFTS(state);
258 
259  if (tsos.isValid()) {
260  // Get the compatible dets on the layer
261  std::vector<pair<const GeomDet*, TrajectoryStateOnSurface> > detsWithStates =
262  detLayer->compatibleDets(tsos, *theService->propagator(thePropagatorName), *theEstimator);
263  if (!detsWithStates.empty()) {
264  TrajectoryStateOnSurface newTSOS = detsWithStates.front().second;
265  const GeomDet* newTSOSDet = detsWithStates.front().first;
266 
267  LogTrace(metname) << "Most compatible det";
268  LogTrace(metname) << debug.dumpMuonId(newTSOSDet->geographicalId());
269 
270  LogDebug(metname) << "L1 info: Det and State:";
271  LogDebug(metname) << debug.dumpMuonId(newTSOSDet->geographicalId());
272 
273  if (newTSOS.isValid()) {
274  //LogDebug(metname) << "(x, y, z) = (" << newTSOS.globalPosition().x() << ", "
275  // << newTSOS.globalPosition().y() << ", " << newTSOS.globalPosition().z() << ")";
276  LogDebug(metname) << "pos: (r=" << newTSOS.globalPosition().mag()
277  << ", phi=" << newTSOS.globalPosition().phi() << ", eta=" << newTSOS.globalPosition().eta()
278  << ")";
279  LogDebug(metname) << "mom: (q*pt=" << newTSOS.charge() * newTSOS.globalMomentum().perp()
280  << ", phi=" << newTSOS.globalMomentum().phi() << ", eta=" << newTSOS.globalMomentum().eta()
281  << ")";
282 
283  //LogDebug(metname) << "State on it";
284  //LogDebug(metname) << debug.dumpTSOS(newTSOS);
285 
286  //PTrajectoryStateOnDet seedTSOS;
288 
289  if (useOfflineSeed) {
290  const TrajectorySeed* assoOffseed = associateOfflineSeedToL1(offlineSeedHandle, offlineSeedMap, newTSOS);
291 
292  if (assoOffseed != nullptr) {
293  PTrajectoryStateOnDet const& seedTSOS = assoOffseed->startingState();
294  TrajectorySeed::const_iterator tsci = assoOffseed->recHits().first, tscie = assoOffseed->recHits().second;
295  for (; tsci != tscie; ++tsci) {
296  container.push_back(*tsci);
297  }
298  output->push_back(
299  L2MuonTrajectorySeed(seedTSOS, container, alongMomentum, L1MuonParticleRef(muColl, l1ParticleIndex)));
300  } else {
301  if (useUnassociatedL1) {
302  // convert the TSOS into a PTSOD
303  PTrajectoryStateOnDet const& seedTSOS =
305  output->push_back(L2MuonTrajectorySeed(
306  seedTSOS, container, alongMomentum, L1MuonParticleRef(muColl, l1ParticleIndex)));
307  }
308  }
309  } else {
310  // convert the TSOS into a PTSOD
311  PTrajectoryStateOnDet const& seedTSOS =
313  output->push_back(
314  L2MuonTrajectorySeed(seedTSOS, container, alongMomentum, L1MuonParticleRef(muColl, l1ParticleIndex)));
315  }
316  }
317  }
318  }
319  }
320 
321  iEvent.put(std::move(output));
322 }
#define LogDebug(id)
std::remove_cv< typename std::remove_reference< argument_type >::type >::type key_type
Definition: Ref.h:164
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
T perp() const
Definition: PV3DBase.h:69
edm::EDGetTokenT< l1extra::L1MuonParticleCollection > muCollToken_
edm::EDGetTokenT< edm::View< TrajectorySeed > > offlineSeedToken_
std::string dumpLayer(const DetLayer *layer) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
const std::string metname
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
Geom::Theta< T > theta() const
GlobalPoint globalPosition() const
std::vector< L1MuRegionalCand > getBrlRPCCands() const
get barrel RPC candidates vector
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
std::string dumpMuonId(const DetId &id) const
std::string dumpFTS(const FreeTrajectoryState &fts) const
std::vector< L1MuRegionalCand > getFwdRPCCands() const
get forward RPC candidates vector
const TrajectorySeed * associateOfflineSeedToL1(edm::Handle< edm::View< TrajectorySeed > > &, std::vector< int > &, TrajectoryStateOnSurface &)
const unsigned theL1MinQuality
std::vector< L1MuRegionalCand > getCSCCands() const
get CSC candidates vector
void push_back(D *&d)
Definition: OwnVector.h:326
bool empty() const
is it an empty muon candidate?
Definition: L1MuGMTCand.h:61
T mag() const
Definition: PV3DBase.h:64
bool isRPC() const
get RPC bit (true=RPC, false = DT/CSC or matched)
recHitContainer::const_iterator const_iterator
bool charge_valid() const
is the charge valid ?
Definition: L1MuGMTCand.h:135
T z() const
Definition: PV3DBase.h:61
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
unsigned getRPCIndex() const
get index of contributing RPC muon
virtual const BoundSurface & surface() const =0
The surface of the GeometricSearchDet.
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
#define LogTrace(id)
unsigned int quality() const
get quality
Definition: L1MuGMTCand.h:90
Definition: DetId.h:17
PTrajectoryStateOnDet const & startingState() const
#define debug
Definition: HDRShower.cc:19
unsigned getDTCSCIndex() const
get index of contributing DT/CSC muon
MeasurementEstimator * theEstimator
virtual const Surface::PositionType & position() const
Returns position of the surface.
std::vector< L1MuRegionalCand > getDTBXCands() const
get DT candidates vector
range recHits() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
T eta() const
Definition: PV3DBase.h:73
MuonServiceProxy * theService
the event setup proxy, it takes care the services update
GlobalVector globalMomentum() const
edm::EDGetTokenT< L1MuGMTReadoutCollection > gmtToken_
constexpr double pi()
Definition: Pi.h:31
bool isFwd() const
get forward bit (true=forward, false=barrel)
edm::Ref< L1MuonParticleCollection > L1MuonParticleRef
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

edm::EDGetTokenT<L1MuGMTReadoutCollection> L2MuonSeedGenerator::gmtToken_
private

Definition at line 70 of file L2MuonSeedGenerator.h.

Referenced by L2MuonSeedGenerator(), and produce().

edm::EDGetTokenT<l1extra::L1MuonParticleCollection> L2MuonSeedGenerator::muCollToken_
private

Definition at line 71 of file L2MuonSeedGenerator.h.

Referenced by L2MuonSeedGenerator(), and produce().

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

Definition at line 72 of file L2MuonSeedGenerator.h.

Referenced by L2MuonSeedGenerator(), and produce().

MeasurementEstimator* L2MuonSeedGenerator::theEstimator
private

Definition at line 83 of file L2MuonSeedGenerator.h.

Referenced by L2MuonSeedGenerator(), produce(), and ~L2MuonSeedGenerator().

edm::InputTag L2MuonSeedGenerator::theL1GMTReadoutCollection
private

Definition at line 66 of file L2MuonSeedGenerator.h.

Referenced by L2MuonSeedGenerator().

const double L2MuonSeedGenerator::theL1MaxEta
private

Definition at line 75 of file L2MuonSeedGenerator.h.

Referenced by produce().

const double L2MuonSeedGenerator::theL1MinPt
private

Definition at line 74 of file L2MuonSeedGenerator.h.

Referenced by produce().

const unsigned L2MuonSeedGenerator::theL1MinQuality
private

Definition at line 76 of file L2MuonSeedGenerator.h.

Referenced by produce().

edm::InputTag L2MuonSeedGenerator::theOfflineSeedLabel
private

Definition at line 67 of file L2MuonSeedGenerator.h.

Referenced by L2MuonSeedGenerator().

std::string L2MuonSeedGenerator::thePropagatorName
private

Definition at line 68 of file L2MuonSeedGenerator.h.

Referenced by associateOfflineSeedToL1(), and produce().

MuonServiceProxy* L2MuonSeedGenerator::theService
private

the event setup proxy, it takes care the services update

Definition at line 81 of file L2MuonSeedGenerator.h.

Referenced by associateOfflineSeedToL1(), L2MuonSeedGenerator(), produce(), and ~L2MuonSeedGenerator().

edm::InputTag L2MuonSeedGenerator::theSource
private

Definition at line 65 of file L2MuonSeedGenerator.h.

Referenced by L2MuonSeedGenerator().

const bool L2MuonSeedGenerator::useOfflineSeed
private

Definition at line 77 of file L2MuonSeedGenerator.h.

Referenced by L2MuonSeedGenerator(), and produce().

const bool L2MuonSeedGenerator::useUnassociatedL1
private

Definition at line 78 of file L2MuonSeedGenerator.h.

Referenced by produce().