CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TSGForOIFromL2.cc
Go to the documentation of this file.
1 
11 
12 #include <memory>
13 
15  : src_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("src"))),
16  maxSeeds_(iConfig.getParameter<uint32_t>("maxSeeds")),
17  maxHitlessSeeds_(iConfig.getParameter<uint32_t>("maxHitlessSeeds")),
18  maxHitSeeds_(iConfig.getParameter<uint32_t>("maxHitSeeds")),
19  numOfLayersToTry_(iConfig.getParameter<int32_t>("layersToTry")),
20  numOfHitsToTry_(iConfig.getParameter<int32_t>("hitsToTry")),
21  numL2ValidHitsCutAllEta_(iConfig.getParameter<uint32_t>("numL2ValidHitsCutAllEta")),
22  numL2ValidHitsCutAllEndcap_(iConfig.getParameter<uint32_t>("numL2ValidHitsCutAllEndcap")),
23  fixedErrorRescalingForHits_(iConfig.getParameter<double>("fixedErrorRescaleFactorForHits")),
24  fixedErrorRescalingForHitless_(iConfig.getParameter<double>("fixedErrorRescaleFactorForHitless")),
25  adjustErrorsDynamicallyForHits_(iConfig.getParameter<bool>("adjustErrorsDynamicallyForHits")),
26  adjustErrorsDynamicallyForHitless_(iConfig.getParameter<bool>("adjustErrorsDynamicallyForHitless")),
27  estimatorName_(iConfig.getParameter<std::string>("estimator")),
28  minEtaForTEC_(iConfig.getParameter<double>("minEtaForTEC")),
29  maxEtaForTOB_(iConfig.getParameter<double>("maxEtaForTOB")),
30  useHitLessSeeds_(iConfig.getParameter<bool>("UseHitLessSeeds")),
31  updator_(new KFUpdator()),
32  measurementTrackerTag_(
33  consumes<MeasurementTrackerEvent>(iConfig.getParameter<edm::InputTag>("MeasurementTrackerEvent"))),
34  pT1_(iConfig.getParameter<double>("pT1")),
35  pT2_(iConfig.getParameter<double>("pT2")),
36  pT3_(iConfig.getParameter<double>("pT3")),
37  eta1_(iConfig.getParameter<double>("eta1")),
38  eta2_(iConfig.getParameter<double>("eta2")),
39  eta3_(iConfig.getParameter<double>("eta3")),
40  eta4_(iConfig.getParameter<double>("eta4")),
41  eta5_(iConfig.getParameter<double>("eta5")),
42  eta6_(iConfig.getParameter<double>("eta6")),
43  eta7_(iConfig.getParameter<double>("eta7")),
44  SF1_(iConfig.getParameter<double>("SF1")),
45  SF2_(iConfig.getParameter<double>("SF2")),
46  SF3_(iConfig.getParameter<double>("SF3")),
47  SF4_(iConfig.getParameter<double>("SF4")),
48  SF5_(iConfig.getParameter<double>("SF5")),
49  SF6_(iConfig.getParameter<double>("SF6")),
50  tsosDiff1_(iConfig.getParameter<double>("tsosDiff1")),
51  tsosDiff2_(iConfig.getParameter<double>("tsosDiff2")),
52  propagatorName_(iConfig.getParameter<std::string>("propagatorName")),
53  theCategory_(std::string("Muon|RecoMuon|TSGForOIFromL2")),
54  estimatorToken_(esConsumes(edm::ESInputTag("", estimatorName_))),
55  magfieldToken_(esConsumes()),
56  propagatorToken_(esConsumes(edm::ESInputTag("", propagatorName_))),
57  tmpTkGeometryToken_(esConsumes()),
58  geometryToken_(esConsumes()),
59  sHPOppositeToken_(esConsumes(edm::ESInputTag("", "hltESPSteppingHelixPropagatorOpposite"))) {
60  produces<std::vector<TrajectorySeed> >();
61 }
62 
64 
65 //
66 // Produce seeds
67 //
69  // Initialize variables
70  unsigned int numSeedsMade = 0;
71  unsigned int layerCount = 0;
72  unsigned int hitlessSeedsMadeIP = 0;
73  unsigned int hitlessSeedsMadeMuS = 0;
74  unsigned int hitSeedsMade = 0;
75 
76  // Surface used to make a TSOS at the PCA to the beamline
78 
79  // Read ESHandles
80  edm::Handle<MeasurementTrackerEvent> measurementTrackerH;
82  const edm::ESHandle<MagneticField> magfieldH = iSetup.getHandle(magfieldToken_);
83  const edm::ESHandle<Propagator> propagatorAlongH = iSetup.getHandle(propagatorToken_);
84  const edm::ESHandle<Propagator>& propagatorOppositeH = propagatorAlongH;
85  const edm::ESHandle<TrackerGeometry> tmpTkGeometryH = iSetup.getHandle(tmpTkGeometryToken_);
87 
88  iEvent.getByToken(measurementTrackerTag_, measurementTrackerH);
89 
90  // Read L2 track collection
92  iEvent.getByToken(src_, l2TrackCol);
93 
94  // The product
95  std::unique_ptr<std::vector<TrajectorySeed> > result(new std::vector<TrajectorySeed>());
96 
97  // Get vector of Detector layers
98  std::vector<BarrelDetLayer const*> const& tob = measurementTrackerH->geometricSearchTracker()->tobLayers();
99  std::vector<ForwardDetLayer const*> const& tecPositive =
100  tmpTkGeometryH->isThere(GeomDetEnumerators::P2OTEC)
101  ? measurementTrackerH->geometricSearchTracker()->posTidLayers()
102  : measurementTrackerH->geometricSearchTracker()->posTecLayers();
103  std::vector<ForwardDetLayer const*> const& tecNegative =
104  tmpTkGeometryH->isThere(GeomDetEnumerators::P2OTEC)
105  ? measurementTrackerH->geometricSearchTracker()->negTidLayers()
106  : measurementTrackerH->geometricSearchTracker()->negTecLayers();
107 
108  // Get suitable propagators
109  std::unique_ptr<Propagator> propagatorAlong = SetPropagationDirection(*propagatorAlongH, alongMomentum);
110  std::unique_ptr<Propagator> propagatorOpposite = SetPropagationDirection(*propagatorOppositeH, oppositeToMomentum);
111 
112  // Stepping Helix Propagator for propogation from muon system to tracker
113  const edm::ESHandle<Propagator> SHPOpposite = iSetup.getHandle(sHPOppositeToken_);
114 
115  // Loop over the L2's and make seeds for all of them
116  LogTrace(theCategory_) << "TSGForOIFromL2::produce: Number of L2's: " << l2TrackCol->size();
117  for (unsigned int l2TrackColIndex(0); l2TrackColIndex != l2TrackCol->size(); ++l2TrackColIndex) {
118  const reco::TrackRef l2(l2TrackCol, l2TrackColIndex);
119 
120  // Container of Seeds
121  std::vector<TrajectorySeed> out;
122  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::produce: L2 muon pT, eta, phi --> " << l2->pt() << " , " << l2->eta()
123  << " , " << l2->phi() << std::endl;
124 
126 
127  dummyPlane->move(fts.position() - dummyPlane->position());
128  TrajectoryStateOnSurface tsosAtIP = TrajectoryStateOnSurface(fts, *dummyPlane);
129  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::produce: Created TSOSatIP: " << tsosAtIP << std::endl;
130 
131  // Get the TSOS on the innermost layer of the L2
132  TrajectoryStateOnSurface tsosAtMuonSystem =
133  trajectoryStateTransform::innerStateOnSurface(*l2, *geometryH, magfieldH.product());
134  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::produce: Created TSOSatMuonSystem: " << tsosAtMuonSystem
135  << std::endl;
136 
137  LogTrace("TSGForOIFromL2")
138  << "TSGForOIFromL2::produce: Check the error of the L2 parameter and use hit seeds if big errors" << std::endl;
139 
140  StateOnTrackerBound fromInside(propagatorAlong.get());
141  TrajectoryStateOnSurface outerTkStateInside = fromInside(fts);
142 
143  StateOnTrackerBound fromOutside(&*SHPOpposite);
144  TrajectoryStateOnSurface outerTkStateOutside = fromOutside(tsosAtMuonSystem);
145 
146  // Check if the two positions (using updated and not-updated TSOS) agree withing certain extent.
147  // If both TSOSs agree, use only the one at vertex, as it uses more information. If they do not agree, search for seeds based on both.
148  double L2muonEta = l2->eta();
149  double absL2muonEta = std::abs(L2muonEta);
150  bool useBoth = false;
151  if (outerTkStateInside.isValid() && outerTkStateOutside.isValid()) {
152  //following commented out variables dist1 (5 par compatibility of tsos at outertracker surface)
153  //dist2 (angle between two tsos) could further be explored in combination of L2 valid hits for seeding. So kept for
154  //future developers
155  //auto dist1 = match_Chi2(outerTkStateInside,outerTkStateOutside);//for future developers
156  //auto dist2 = deltaR(outerTkStateInside.globalMomentum(),outerTkStateOutside.globalMomentum());//for future developers
157  //if ((dist1 > tsosDiff1_ || dist2 > tsosDiff2_) && l2->numberOfValidHits() < 20) useBoth = true;//for future developers
158  if (l2->numberOfValidHits() < numL2ValidHitsCutAllEta_)
159  useBoth = true;
160  if (l2->numberOfValidHits() < numL2ValidHitsCutAllEndcap_ && absL2muonEta > eta7_)
161  useBoth = true;
162  if (absL2muonEta > eta1_ && absL2muonEta < eta1_)
163  useBoth = true;
164  }
165 
166  numSeedsMade = 0;
167  hitlessSeedsMadeIP = 0;
168  hitlessSeedsMadeMuS = 0;
169  hitSeedsMade = 0;
170 
171  // calculate scale factors
173  double errorSFHitless =
175 
176  // BARREL
177  if (absL2muonEta < maxEtaForTOB_) {
178  layerCount = 0;
179  for (auto it = tob.rbegin(); it != tob.rend(); ++it) {
180  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::produce: looping in TOB layer " << layerCount << std::endl;
181  if (useHitLessSeeds_ && hitlessSeedsMadeIP < maxHitlessSeeds_ && numSeedsMade < maxSeeds_)
183  tsosAtIP,
184  *(propagatorAlong.get()),
185  estimatorH,
186  errorSFHitless,
187  hitlessSeedsMadeIP,
188  numSeedsMade,
189  out);
190 
191  // Do not create hitbased seeds in barrel region
192  if (absL2muonEta > 1.0 && hitSeedsMade < maxHitSeeds_ && numSeedsMade < maxSeeds_)
193  makeSeedsFromHits(**it,
194  tsosAtIP,
195  *(propagatorAlong.get()),
196  estimatorH,
197  measurementTrackerH,
198  errorSFHits,
199  hitSeedsMade,
200  numSeedsMade,
201  layerCount,
202  out);
203 
204  if (useBoth) {
205  if (useHitLessSeeds_ && hitlessSeedsMadeMuS < maxHitlessSeeds_ && numSeedsMade < maxSeeds_)
207  outerTkStateOutside,
208  *(propagatorOpposite.get()),
209  estimatorH,
210  errorSFHitless,
211  hitlessSeedsMadeMuS,
212  numSeedsMade,
213  out);
214  }
215  }
216  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2:::produce: NumSeedsMade = " << numSeedsMade
217  << " , layerCount = " << layerCount << std::endl;
218  }
219 
220  // Reset number of seeds if in overlap region
221  if (absL2muonEta > minEtaForTEC_ && absL2muonEta < maxEtaForTOB_) {
222  numSeedsMade = 0;
223  hitlessSeedsMadeIP = 0;
224  hitlessSeedsMadeMuS = 0;
225  hitSeedsMade = 0;
226  }
227 
228  // ENDCAP+
229  if (L2muonEta > minEtaForTEC_) {
230  layerCount = 0;
231  for (auto it = tecPositive.rbegin(); it != tecPositive.rend(); ++it) {
232  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::produce: looping in TEC+ layer " << layerCount << std::endl;
233  if (useHitLessSeeds_ && hitlessSeedsMadeIP < maxHitlessSeeds_ && numSeedsMade < maxSeeds_)
235  tsosAtIP,
236  *(propagatorAlong.get()),
237  estimatorH,
238  errorSFHitless,
239  hitlessSeedsMadeIP,
240  numSeedsMade,
241  out);
242 
243  if (absL2muonEta > 1.0 && hitSeedsMade < maxHitSeeds_ && numSeedsMade < maxSeeds_)
244  makeSeedsFromHits(**it,
245  tsosAtIP,
246  *(propagatorAlong.get()),
247  estimatorH,
248  measurementTrackerH,
249  errorSFHits,
250  hitSeedsMade,
251  numSeedsMade,
252  layerCount,
253  out);
254 
255  if (useBoth) {
256  if (useHitLessSeeds_ && hitlessSeedsMadeMuS < maxHitlessSeeds_ && numSeedsMade < maxSeeds_)
258  outerTkStateOutside,
259  *(propagatorOpposite.get()),
260  estimatorH,
261  errorSFHitless,
262  hitlessSeedsMadeMuS,
263  numSeedsMade,
264  out);
265  }
266  }
267  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2:::produce: NumSeedsMade = " << numSeedsMade
268  << " , layerCount = " << layerCount << std::endl;
269  }
270 
271  // ENDCAP-
272  if (L2muonEta < -minEtaForTEC_) {
273  layerCount = 0;
274  for (auto it = tecNegative.rbegin(); it != tecNegative.rend(); ++it) {
275  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::produce: looping in TEC- layer " << layerCount << std::endl;
276  if (useHitLessSeeds_ && hitlessSeedsMadeIP < maxHitlessSeeds_ && numSeedsMade < maxSeeds_)
278  tsosAtIP,
279  *(propagatorAlong.get()),
280  estimatorH,
281  errorSFHitless,
282  hitlessSeedsMadeIP,
283  numSeedsMade,
284  out);
285 
286  if (absL2muonEta > 1.0 && hitSeedsMade < maxHitSeeds_ && numSeedsMade < maxSeeds_)
287  makeSeedsFromHits(**it,
288  tsosAtIP,
289  *(propagatorAlong.get()),
290  estimatorH,
291  measurementTrackerH,
292  errorSFHits,
293  hitSeedsMade,
294  numSeedsMade,
295  layerCount,
296  out);
297 
298  if (useBoth) {
299  if (useHitLessSeeds_ && hitlessSeedsMadeMuS < maxHitlessSeeds_ && numSeedsMade < maxSeeds_)
301  outerTkStateOutside,
302  *(propagatorOpposite.get()),
303  estimatorH,
304  errorSFHitless,
305  hitlessSeedsMadeMuS,
306  numSeedsMade,
307  out);
308  }
309  }
310  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2:::produce: NumSeedsMade = " << numSeedsMade
311  << " , layerCount = " << layerCount << std::endl;
312  }
313 
314  for (std::vector<TrajectorySeed>::iterator it = out.begin(); it != out.end(); ++it) {
315  result->push_back(*it);
316  }
317 
318  } // L2Collection
319 
320  edm::LogInfo(theCategory_) << "TSGForOIFromL2::produce: number of seeds made: " << result->size();
321 
322  iEvent.put(std::move(result));
323 }
324 
325 //
326 // Create seeds without hits on a given layer (TOB or TEC)
327 //
329  const TrajectoryStateOnSurface& tsos,
332  double errorSF,
333  unsigned int& hitlessSeedsMade,
334  unsigned int& numSeedsMade,
335  std::vector<TrajectorySeed>& out) const {
336  // create hitless seeds
337  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::makeSeedsWithoutHits: Start hitless" << std::endl;
338  std::vector<GeometricSearchDet::DetWithState> dets;
339  layer.compatibleDetsV(tsos, propagatorAlong, *estimator, dets);
340  if (!dets.empty()) {
341  auto const& detOnLayer = dets.front().first;
342  auto const& tsosOnLayer = dets.front().second;
343  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::makeSeedsWithoutHits: tsosOnLayer " << tsosOnLayer << std::endl;
344  if (!tsosOnLayer.isValid()) {
345  edm::LogInfo(theCategory_) << "ERROR!: Hitless TSOS is not valid!";
346  } else {
347  dets.front().second.rescaleError(errorSF);
348  PTrajectoryStateOnDet const& ptsod =
349  trajectoryStateTransform::persistentState(tsosOnLayer, detOnLayer->geographicalId().rawId());
351  out.push_back(TrajectorySeed(ptsod, rHC, oppositeToMomentum));
352  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::makeSeedsWithoutHits: TSOS (Hitless) done " << std::endl;
353  hitlessSeedsMade++;
354  numSeedsMade++;
355  }
356  }
357 }
358 
359 //
360 // Find hits on a given layer (TOB or TEC) and create seeds from updated TSOS with hit
361 //
363  const TrajectoryStateOnSurface& tsos,
367  double errorSF,
368  unsigned int& hitSeedsMade,
369  unsigned int& numSeedsMade,
370  unsigned int& layerCount,
371  std::vector<TrajectorySeed>& out) const {
372  if (layerCount > numOfLayersToTry_)
373  return;
374 
375  // Error Rescaling
376  TrajectoryStateOnSurface onLayer(tsos);
377  onLayer.rescaleError(errorSF);
378 
379  std::vector<GeometricSearchDet::DetWithState> dets;
380  layer.compatibleDetsV(onLayer, propagatorAlong, *estimator, dets);
381 
382  // Find Measurements on each DetWithState
383  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::makeSeedsFromHits: Find measurements on each detWithState "
384  << dets.size() << std::endl;
385  std::vector<TrajectoryMeasurement> meas;
386  for (std::vector<GeometricSearchDet::DetWithState>::iterator it = dets.begin(); it != dets.end(); ++it) {
387  MeasurementDetWithData det = measurementTracker->idToDet(it->first->geographicalId());
388  if (det.isNull())
389  continue;
390  if (!it->second.isValid())
391  continue; // Skip if TSOS is not valid
392 
393  std::vector<TrajectoryMeasurement> mymeas =
394  det.fastMeasurements(it->second, onLayer, propagatorAlong, *estimator); // Second TSOS is not used
395  for (std::vector<TrajectoryMeasurement>::const_iterator it2 = mymeas.begin(), ed2 = mymeas.end(); it2 != ed2;
396  ++it2) {
397  if (it2->recHit()->isValid())
398  meas.push_back(*it2); // Only save those which are valid
399  }
400  }
401 
402  // Update TSOS using TMs after sorting, then create Trajectory Seed and put into vector
403  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::makeSeedsFromHits: Update TSOS using TMs after sorting, then create "
404  "Trajectory Seed, number of TM = "
405  << meas.size() << std::endl;
406  std::sort(meas.begin(), meas.end(), TrajMeasLessEstim());
407 
408  unsigned int found = 0;
409  for (std::vector<TrajectoryMeasurement>::const_iterator it = meas.begin(); it != meas.end(); ++it) {
410  TrajectoryStateOnSurface updatedTSOS = updator_->update(it->forwardPredictedState(), *it->recHit());
411  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::makeSeedsFromHits: TSOS for TM " << found << std::endl;
412  if (not updatedTSOS.isValid())
413  continue;
414 
416  seedHits.push_back(*it->recHit()->hit());
417  PTrajectoryStateOnDet const& pstate =
418  trajectoryStateTransform::persistentState(updatedTSOS, it->recHit()->geographicalId().rawId());
419  LogTrace("TSGForOIFromL2") << "TSGForOIFromL2::makeSeedsFromHits: Number of seedHits: " << seedHits.size()
420  << std::endl;
421  TrajectorySeed seed(pstate, std::move(seedHits), oppositeToMomentum);
422  out.push_back(seed);
423  found++;
424  numSeedsMade++;
425  hitSeedsMade++;
426  if (found == numOfHitsToTry_)
427  break;
428  if (hitSeedsMade > maxHitSeeds_)
429  return;
430  }
431 
432  if (found)
433  layerCount++;
434 }
435 
436 //
437 // Calculate the dynamic error SF by analysing the L2
438 //
440  double theSF = 1.0;
441  // L2 direction vs pT blowup - as was previously done:
442  // Split into 4 pT ranges: <pT1_, pT1_<pT2_, pT2_<pT3_, <pT4_: 13,30,70
443  // Split into different eta ranges depending in pT
444  double abseta = std::abs(track->eta());
445  if (track->pt() <= pT1_)
446  theSF = SF1_;
447  else if (track->pt() > pT1_ && track->pt() <= pT2_) {
448  if (abseta <= eta3_)
449  theSF = SF3_;
450  else if (abseta > eta3_ && abseta <= eta6_)
451  theSF = SF2_;
452  else if (abseta > eta6_)
453  theSF = SF3_;
454  } else if (track->pt() > pT2_ && track->pt() <= pT3_) {
455  if (abseta <= eta1_)
456  theSF = SF6_;
457  else if (abseta > eta1_ && abseta <= eta2_)
458  theSF = SF4_;
459  else if (abseta > eta2_ && abseta <= eta3_)
460  theSF = SF6_;
461  else if (abseta > eta3_ && abseta <= eta4_)
462  theSF = SF1_;
463  else if (abseta > eta4_ && abseta <= eta5_)
464  theSF = SF1_;
465  else if (abseta > eta5_)
466  theSF = SF5_;
467  } else if (track->pt() > pT3_) {
468  if (abseta <= eta3_)
469  theSF = SF5_;
470  else if (abseta > eta3_ && abseta <= eta4_)
471  theSF = SF4_;
472  else if (abseta > eta4_ && abseta <= eta5_)
473  theSF = SF4_;
474  else if (abseta > eta5_)
475  theSF = SF5_;
476  }
477 
478  LogTrace(theCategory_) << "TSGForOIFromL2::calculateSFFromL2: SF has been calculated as: " << theSF;
479 
480  return theSF;
481 }
482 
483 //
484 // calculate Chi^2 of two trajectory states
485 //
487  if (!tsos1.isValid() || !tsos2.isValid())
488  return -1.;
489 
491  AlgebraicSymMatrix55 m(tsos1.localError().matrix() + tsos2.localError().matrix());
492 
493  bool ierr = !m.Invert();
494 
495  if (ierr) {
496  edm::LogInfo("TSGForOIFromL2") << "Error inverting covariance matrix";
497  return -1;
498  }
499 
500  double est = ROOT::Math::Similarity(v, m);
501 
502  return est;
503 }
504 
505 //
506 //
507 //
510  desc.add<edm::InputTag>("src", edm::InputTag("hltL2Muons", "UpdatedAtVtx"));
511  desc.add<int>("layersToTry", 2);
512  desc.add<double>("fixedErrorRescaleFactorForHitless", 2.0);
513  desc.add<int>("hitsToTry", 1);
514  desc.add<bool>("adjustErrorsDynamicallyForHits", false);
515  desc.add<bool>("adjustErrorsDynamicallyForHitless", true);
516  desc.add<edm::InputTag>("MeasurementTrackerEvent", edm::InputTag("hltSiStripClusters"));
517  desc.add<bool>("UseHitLessSeeds", true);
518  desc.add<std::string>("estimator", "hltESPChi2MeasurementEstimator100");
519  desc.add<double>("maxEtaForTOB", 1.8);
520  desc.add<double>("minEtaForTEC", 0.7);
521  desc.addUntracked<bool>("debug", false);
522  desc.add<double>("fixedErrorRescaleFactorForHits", 1.0);
523  desc.add<unsigned int>("maxSeeds", 20);
524  desc.add<unsigned int>("maxHitlessSeeds", 5);
525  desc.add<unsigned int>("maxHitSeeds", 1);
526  desc.add<unsigned int>("numL2ValidHitsCutAllEta", 20);
527  desc.add<unsigned int>("numL2ValidHitsCutAllEndcap", 30);
528  desc.add<double>("pT1", 13.0);
529  desc.add<double>("pT2", 30.0);
530  desc.add<double>("pT3", 70.0);
531  desc.add<double>("eta1", 0.2);
532  desc.add<double>("eta2", 0.3);
533  desc.add<double>("eta3", 1.0);
534  desc.add<double>("eta4", 1.2);
535  desc.add<double>("eta5", 1.6);
536  desc.add<double>("eta6", 1.4);
537  desc.add<double>("eta7", 2.1);
538  desc.add<double>("SF1", 3.0);
539  desc.add<double>("SF2", 4.0);
540  desc.add<double>("SF3", 5.0);
541  desc.add<double>("SF4", 7.0);
542  desc.add<double>("SF5", 10.0);
543  desc.add<double>("SF6", 2.0);
544  desc.add<double>("tsosDiff1", 0.2);
545  desc.add<double>("tsosDiff2", 0.02);
546  desc.add<std::string>("propagatorName", "PropagatorWithMaterialParabolicMf");
547  descriptions.add("TSGForOIFromL2", desc);
548 }
549 
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const bool useHitLessSeeds_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
TSGForOIFromL2(const edm::ParameterSet &iConfig)
const double eta2_
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const LocalTrajectoryParameters & localParameters() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
const unsigned int numOfHitsToTry_
How many hits to try per layer.
tuple propagatorAlong
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
tuple measurementTracker
const double minEtaForTEC_
Minimum eta value to activate searching in the TEC.
size_type size() const
Definition: OwnVector.h:300
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
const unsigned int maxHitSeeds_
Maximum number of hitbased seeds for each L2.
const double SF2_
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magfieldToken_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tmpTkGeometryToken_
const edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > geometryToken_
#define LogTrace(id)
constexpr std::array< uint8_t, layerIndexSize > layer
tuple result
Definition: mps_fire.py:311
const edm::ESGetToken< Chi2MeasurementEstimatorBase, TrackingComponentsRecord > estimatorToken_
const unsigned int maxHitlessSeeds_
Maximum number of hitless seeds for each L2.
AlgebraicVector5 vector() const
std::unique_ptr< Propagator > SetPropagationDirection(Propagator const &iprop, PropagationDirection dir)
void push_back(D *&d)
Definition: OwnVector.h:326
const double pT2_
const double SF6_
int iEvent
Definition: GenABIO.cc:224
const double eta7_
std::vector< TrajectoryMeasurement > fastMeasurements(const TrajectoryStateOnSurface &stateOnThisDet, const TrajectoryStateOnSurface &tsos2, const Propagator &prop, const MeasurementEstimator &est) const
static PlanePointer build(Args &&...args)
Definition: Plane.h:33
const edm::ESGetToken< Propagator, TrackingComponentsRecord > propagatorToken_
const unsigned int numL2ValidHitsCutAllEndcap_
def move
Definition: eostools.py:511
void produce(edm::StreamID sid, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
const unsigned int maxSeeds_
Maximum number of seeds for each L2.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ROOT::Math::SVector< double, 5 > AlgebraicVector5
const double eta1_
const std::unique_ptr< TrajectoryStateUpdator > updator_
KFUpdator defined in constructor.
const double SF3_
const double fixedErrorRescalingForHits_
Rescale L2 parameter uncertainties (fixed error vs pT, eta)
const double maxEtaForTOB_
Maximum eta value to activate searching in the TOB.
const AlgebraicSymMatrix55 & matrix() const
const LocalTrajectoryError & localError() const
const bool adjustErrorsDynamicallyForHits_
Whether or not to use an automatically calculated scale-factor value.
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const std::string theCategory_
const edm::EDGetTokenT< MeasurementTrackerEvent > measurementTrackerTag_
const unsigned int numL2ValidHitsCutAllEta_
L2 valid hit cuts to decide seed creation by both states.
~TSGForOIFromL2() override
Log< level::Info, false > LogInfo
const double eta6_
GlobalPoint position() const
double match_Chi2(const TrajectoryStateOnSurface &tsos1, const TrajectoryStateOnSurface &tsos2) const
Find compatability between two TSOSs.
const double eta5_
T const * product() const
Definition: ESHandle.h:86
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const double SF4_
void makeSeedsWithoutHits(const GeometricSearchDet &layer, const TrajectoryStateOnSurface &tsos, const Propagator &propagatorAlong, const edm::ESHandle< Chi2MeasurementEstimatorBase > &estimator, double errorSF, unsigned int &hitlessSeedsMade, unsigned int &numSeedsMade, std::vector< TrajectorySeed > &out) const
Create seeds without hits on a given layer (TOB or TEC)
const edm::EDGetTokenT< reco::TrackCollection > src_
Labels for input collections.
tuple propagatorOpposite
const double eta3_
double calculateSFFromL2(const reco::TrackRef track) const
Calculate the dynamic error SF by analysing the L2.
const double eta4_
const bool adjustErrorsDynamicallyForHitless_
const edm::ESGetToken< Propagator, TrackingComponentsRecord > sHPOppositeToken_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
Create L3MuonTrajectorySeeds from L2 Muons updated at vertex in an outside-in manner.
const double SF5_
void makeSeedsFromHits(const GeometricSearchDet &layer, const TrajectoryStateOnSurface &tsos, const Propagator &propagatorAlong, const edm::ESHandle< Chi2MeasurementEstimatorBase > &estimator, const edm::Handle< MeasurementTrackerEvent > &measurementTracker, double errorSF, unsigned int &hitSeedsMade, unsigned int &numSeedsMade, unsigned int &layerCount, std::vector< TrajectorySeed > &out) const
Find hits on a given layer (TOB or TEC) and create seeds from updated TSOS with hit.
const unsigned int numOfLayersToTry_
How many layers to try.
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const double SF1_
virtual void compatibleDetsV(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est, std::vector< DetWithState > &result) const
const double pT3_
const double pT1_
pT, eta ranges and scale factor values
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
const double fixedErrorRescalingForHitless_