CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CosmicMuonSeedGenerator.cc
Go to the documentation of this file.
13 
15 
19 
25 
28 
32 
34 
35 #include <vector>
36 
39 
40 using namespace edm;
41 using namespace std;
42 
43 // Constructor
45  produces<TrajectorySeedCollection>();
46 
47  // enable the DT chamber
48  theEnableDTFlag = pset.getParameter<bool>("EnableDTMeasurement");
49  // enable the CSC chamber
50  theEnableCSCFlag = pset.getParameter<bool>("EnableCSCMeasurement");
51 
52  theDTRecSegmentLabel = pset.getParameter<InputTag>("DTRecSegmentLabel");
53 
54  theCSCRecSegmentLabel = pset.getParameter<InputTag>("CSCRecSegmentLabel");
55 
56  // the maximum number of TrajectorySeed
57  theMaxSeeds = pset.getParameter<int>("MaxSeeds");
58 
59  theMaxDTChi2 = pset.getParameter<double>("MaxDTChi2");
60  theMaxCSCChi2 = pset.getParameter<double>("MaxCSCChi2");
61 
62  theForcePointDownFlag = pset.existsAs<bool>("ForcePointDown") ? pset.getParameter<bool>("ForcePointDown") : true;
63 
64  // pre-determined parameters for seed pt calculation ( pt * dphi )
65  theParameters["topmb41"] = 0.87;
66  theParameters["bottommb41"] = 1.2;
67  theParameters["topmb42"] = 0.67;
68  theParameters["bottommb42"] = 0.98;
69  theParameters["topmb43"] = 0.34;
70  theParameters["bottommb43"] = 0.58;
71  theParameters["topmb31"] = 0.54;
72  theParameters["bottommb31"] = 0.77;
73  theParameters["topmb32"] = 0.35;
74  theParameters["bottommb32"] = 0.55;
75  theParameters["topmb21"] = 0.21;
76  theParameters["bottommb21"] = 0.31;
77 
78  edm::ConsumesCollector iC = consumesCollector();
79  muonMeasurements = new MuonDetLayerMeasurements(theDTRecSegmentLabel,
80  theCSCRecSegmentLabel,
81 
82  InputTag(),
83  InputTag(),
84  InputTag(),
85  iC,
86  theEnableDTFlag,
87  theEnableCSCFlag,
88  false,
89  false,
90  false);
91 }
92 
93 // Destructor
95  if (muonMeasurements)
96  delete muonMeasurements;
97 }
98 
99 // reconstruct muon's seeds
101  eSetup.get<IdealMagneticFieldRecord>().get(theField);
102 
103  auto output = std::make_unique<TrajectorySeedCollection>();
104 
106 
107  std::string category = "Muon|RecoMuon|CosmicMuonSeedGenerator";
108 
109  // Muon Geometry - DT, CSC and RPC
110  eSetup.get<MuonRecoGeometryRecord>().get(theMuonLayers);
111 
112  // get the DT layers
113  vector<const DetLayer*> dtLayers = theMuonLayers->allDTLayers();
114 
115  // get the CSC layers
116  vector<const DetLayer*> cscForwardLayers = theMuonLayers->forwardCSCLayers();
117  vector<const DetLayer*> cscBackwardLayers = theMuonLayers->backwardCSCLayers();
118 
119  muonMeasurements->setEvent(event);
120 
121  MuonRecHitContainer allHits;
122 
123  vector<MuonRecHitContainer> RHMBs;
124  vector<MuonRecHitContainer> RHMEFs;
125  vector<MuonRecHitContainer> RHMEBs;
126 
127  stable_sort(allHits.begin(), allHits.end(), DecreasingGlobalY());
128 
129  for (vector<const DetLayer*>::reverse_iterator icsclayer = cscForwardLayers.rbegin();
130  icsclayer != cscForwardLayers.rend() - 1;
131  ++icsclayer) {
132  MuonRecHitContainer RHMF = muonMeasurements->recHits(*icsclayer);
133  allHits.insert(allHits.end(), RHMF.begin(), RHMF.end());
134  }
135 
136  for (vector<const DetLayer*>::reverse_iterator icsclayer = cscBackwardLayers.rbegin();
137  icsclayer != cscBackwardLayers.rend() - 1;
138  ++icsclayer) {
139  MuonRecHitContainer RHMF = muonMeasurements->recHits(*icsclayer);
140  allHits.insert(allHits.end(), RHMF.begin(), RHMF.end());
141  }
142 
143  for (vector<const DetLayer*>::reverse_iterator idtlayer = dtLayers.rbegin(); idtlayer != dtLayers.rend();
144  ++idtlayer) {
145  MuonRecHitContainer RHMB = muonMeasurements->recHits(*idtlayer);
146  RHMBs.push_back(RHMB);
147 
148  if (idtlayer != dtLayers.rbegin())
149  allHits.insert(allHits.end(), RHMB.begin(), RHMB.end());
150  }
151 
152  // stable_sort(allHits.begin(),allHits.end(),DecreasingGlobalY());
153 
154  LogTrace(category) << "all RecHits: " << allHits.size();
155 
156  // CosmicMuonSeedGenerator::MuonRecHitPairVector mb41 = makeSegPairs(RHMBs[0], RHMBs[3], "mb41");
157  // createSeeds(seeds, mb41, eSetup);
158 
159  // CosmicMuonSeedGenerator::MuonRecHitPairVector mb43 = makeSegPairs(RHMBs[0],RHMBs[1], "mb43");
160  // createSeeds(seeds, mb43, eSetup);
161 
162  CosmicMuonSeedGenerator::MuonRecHitPairVector mb42 = makeSegPairs(RHMBs[0], RHMBs[2], "mb42");
163  createSeeds(seeds, mb42, eSetup);
164 
165  // CosmicMuonSeedGenerator::MuonRecHitPairVector mb32 = makeSegPairs(RHMBs[1], RHMBs[2], "mb32");
166  // createSeeds(seeds, mb32, eSetup);
167 
168  CosmicMuonSeedGenerator::MuonRecHitPairVector mb31 = makeSegPairs(RHMBs[1], RHMBs[3], "mb31");
169  createSeeds(seeds, mb31, eSetup);
170 
171  // CosmicMuonSeedGenerator::MuonRecHitPairVector mb21 = makeSegPairs(RHMBs[2], RHMBs[3], "mb21");
172  // createSeeds(seeds, mb21, eSetup);
173 
174  if (!allHits.empty()) {
175  MuonRecHitContainer goodhits = selectSegments(allHits);
176  LogTrace(category) << "good RecHits: " << goodhits.size();
177 
178  if (goodhits.empty()) {
179  LogTrace(category) << "No qualified Segments in Event! ";
180  LogTrace(category) << "Use 2D RecHit";
181 
182  createSeeds(seeds, allHits, eSetup);
183 
184  } else {
185  createSeeds(seeds, goodhits, eSetup);
186  }
187  }
188 
189  LogTrace(category) << "Seeds built: " << seeds.size();
190 
191  for (std::vector<TrajectorySeed>::iterator seed = seeds.begin(); seed != seeds.end(); ++seed) {
192  output->push_back(*seed);
193  }
194 
195  event.put(std::move(output));
196  seeds.clear();
197 }
198 
200  const std::string category = "Muon|RecoMuon|CosmicMuonSeedGenerator";
201 
202  // only use 4D segments
203  if (!hit->isValid())
204  return false;
205 
206  if (hit->dimension() < 4) {
207  LogTrace(category) << "dim < 4";
208  return false;
209  }
210 
211  if (hit->isDT() && (hit->chi2() > theMaxDTChi2)) {
212  LogTrace(category) << "DT chi2 too large";
213  return false;
214  } else if (hit->isCSC() && (hit->chi2() > theMaxCSCChi2)) {
215  LogTrace(category) << "CSC chi2 too large";
216  return false;
217  }
218  return true;
219 }
220 
223  const std::string category = "Muon|RecoMuon|CosmicMuonSeedGenerator";
224 
225  //Only select good quality Segments
226  for (MuonRecHitContainer::const_iterator hit = hits.begin(); hit != hits.end(); hit++) {
227  if (checkQuality(*hit))
228  result.push_back(*hit);
229  }
230 
231  if (result.size() < 2)
232  return result;
233 
234  MuonRecHitContainer result2;
235 
236  //avoid selecting Segments with similar direction
237  for (MuonRecHitContainer::iterator hit = result.begin(); hit != result.end(); hit++) {
238  if (*hit == nullptr)
239  continue;
240  if (!(*hit)->isValid())
241  continue;
242  bool good = true;
243  //UNUSED: GlobalVector dir1 = (*hit)->globalDirection();
244  //UNUSED: GlobalPoint pos1 = (*hit)->globalPosition();
245  for (MuonRecHitContainer::iterator hit2 = hit + 1; hit2 != result.end(); hit2++) {
246  if (*hit2 == nullptr)
247  continue;
248  if (!(*hit2)->isValid())
249  continue;
250 
251  //compare direction and position
252  //UNUSED: GlobalVector dir2 = (*hit2)->globalDirection();
253  //UNUSED: GlobalPoint pos2 = (*hit2)->globalPosition();
254  if (!areCorrelated((*hit), (*hit2)))
255  continue;
256 
257  if (!leftIsBetter((*hit), (*hit2))) {
258  good = false;
259  } else
260  (*hit2) = nullptr;
261  }
262 
263  if (good)
264  result2.push_back(*hit);
265  }
266 
267  result.clear();
268 
269  return result2;
270 }
271 
273  const MuonRecHitContainer& hits,
274  const edm::EventSetup& eSetup) const {
275  const std::string category = "Muon|RecoMuon|CosmicMuonSeedGenerator";
276 
277  if (hits.empty() || results.size() >= theMaxSeeds)
278  return;
279  for (MuonRecHitContainer::const_iterator ihit = hits.begin(); ihit != hits.end(); ihit++) {
280  const std::vector<TrajectorySeed>& sds = createSeed((*ihit), eSetup);
281  LogTrace(category) << "created seeds from rechit " << sds.size();
282  results.insert(results.end(), sds.begin(), sds.end());
283  if (results.size() >= theMaxSeeds)
284  break;
285  }
286  return;
287 }
288 
291  const edm::EventSetup& eSetup) const {
292  const std::string category = "Muon|RecoMuon|CosmicMuonSeedGenerator";
293 
294  if (hitpairs.empty() || results.size() >= theMaxSeeds)
295  return;
296  for (CosmicMuonSeedGenerator::MuonRecHitPairVector::const_iterator ihitpair = hitpairs.begin();
297  ihitpair != hitpairs.end();
298  ihitpair++) {
299  const std::vector<TrajectorySeed>& sds = createSeed((*ihitpair), eSetup);
300  LogTrace(category) << "created seeds from rechit " << sds.size();
301  results.insert(results.end(), sds.begin(), sds.end());
302  if (results.size() >= theMaxSeeds)
303  break;
304  }
305  return;
306 }
307 
308 std::vector<TrajectorySeed> CosmicMuonSeedGenerator::createSeed(const MuonRecHitPointer& hit,
309  const edm::EventSetup& eSetup) const {
310  std::vector<TrajectorySeed> result;
311 
312  const std::string category = "Muon|RecoMuon|CosmicMuonSeedGenerator";
313 
314  MuonPatternRecoDumper dumper;
315 
316  // set the pt by hand
317  double pt = 10.0;
318 
319  // AlgebraicVector4 t;
320  AlgebraicSymMatrix mat(5, 0);
321 
322  // Fill the LocalTrajectoryParameters
323  LocalPoint segPos = hit->localPosition();
324 
325  GlobalVector polar(GlobalVector::Spherical(hit->globalDirection().theta(), hit->globalDirection().phi(), 1.));
326  // Force all track downward for cosmic, not beam-halo
327  if (theForcePointDownFlag) {
328  if (hit->geographicalId().subdetId() == MuonSubdetId::DT && fabs(hit->globalDirection().eta()) < 4.0 &&
329  hit->globalDirection().phi() > 0)
330  polar = -polar;
331 
332  if (hit->geographicalId().subdetId() == MuonSubdetId::CSC && fabs(hit->globalDirection().eta()) > 2.3)
333  polar = -polar;
334  }
335 
336  polar *= fabs(pt) / polar.perp();
337 
338  LocalVector segDir = hit->det()->toLocal(polar);
339 
340  int charge = 1;
341  LocalTrajectoryParameters param(segPos, segDir, charge);
342 
343  charge = -1;
344  LocalTrajectoryParameters param2(segPos, segDir, charge);
345 
346  mat = hit->parametersError().similarityT(hit->projectionMatrix());
347 
348  float p_err = 0.2;
349  mat[0][0] = p_err;
350 
351  LocalTrajectoryError error(asSMatrix<5>(mat));
352 
353  // Create the TrajectoryStateOnSurface
354  TrajectoryStateOnSurface tsos(param, error, hit->det()->surface(), &*theField);
355  TrajectoryStateOnSurface tsos2(param2, error, hit->det()->surface(), &*theField);
356 
357  LogTrace(category) << "Trajectory State on Surface of Seed";
358  LogTrace(category) << "mom: " << tsos.globalMomentum() << " phi: " << tsos.globalMomentum().phi();
359  LogTrace(category) << "pos: " << tsos.globalPosition();
360  LogTrace(category) << "The RecSegment relies on: ";
361  LogTrace(category) << dumper.dumpMuonId(hit->geographicalId());
362 
364  container.push_back(hit->hit()->clone());
365 
366  result.push_back(tsosToSeed(tsos, hit->geographicalId().rawId(), container));
367  result.push_back(tsosToSeed(tsos2, hit->geographicalId().rawId(), container));
368 
369  return result;
370 }
371 
373  bool result = false;
374 
375  GlobalVector dir1 = lhs->globalDirection();
376  GlobalPoint pos1 = lhs->globalPosition();
377  GlobalVector dir2 = rhs->globalDirection();
378  GlobalPoint pos2 = rhs->globalPosition();
379 
380  GlobalVector dis = pos2 - pos1;
381 
382  if ((deltaR<double>(dir1.eta(), dir1.phi(), dir2.eta(), dir2.phi()) < 0.1 ||
383  deltaR<double>(dir1.eta(), dir1.phi(), -dir2.eta(), -dir2.phi()) < 0.1) &&
384  dis.mag() < 5.0)
385  result = true;
386 
387  if ((deltaR<double>(dir1.eta(), dir1.phi(), dir2.eta(), dir2.phi()) < 0.1 ||
388  deltaR<double>(dir1.eta(), dir1.phi(), -dir2.eta(), -dir2.phi()) < 0.1) &&
389  (deltaR<double>(dir1.eta(), dir1.phi(), dis.eta(), dis.phi()) < 0.1 ||
390  deltaR<double>(dir2.eta(), dir2.phi(), dis.eta(), dis.phi()) < 0.1))
391  result = true;
392 
393  if (fabs(dir1.eta()) > 4.0 || fabs(dir2.eta()) > 4.0) {
394  if ((fabs(dir1.theta() - dir2.theta()) < 0.07 || fabs(dir1.theta() + dir2.theta()) > 3.07) &&
395  (fabs(dir1.theta() - dis.theta()) < 0.07 || fabs(dir1.theta() - dis.theta()) < 0.07 ||
396  fabs(dir1.theta() + dis.theta()) > 3.07 || fabs(dir1.theta() + dis.theta()) > 3.07))
397 
398  result = true;
399  }
400 
401  return result;
402 }
403 
406  if ((lhs->degreesOfFreedom() > rhs->degreesOfFreedom()) ||
407  ((lhs->degreesOfFreedom() == rhs->degreesOfFreedom()) && (lhs)->chi2() < (rhs)->chi2()))
408  return true;
409  else
410  return false;
411 }
412 
416  std::string tag) const {
418  const std::string category = "Muon|RecoMuon|CosmicMuonSeedGenerator";
419 
420  if (hits1.empty() || hits2.empty())
421  return result;
422 
423  for (MuonRecHitContainer::const_iterator ihit1 = hits1.begin(); ihit1 != hits1.end(); ihit1++) {
424  if (!checkQuality(*ihit1))
425  continue;
426 
427  for (MuonRecHitContainer::const_iterator ihit2 = hits2.begin(); ihit2 != hits2.end(); ihit2++) {
428  if (!checkQuality(*ihit2))
429  continue;
430 
431  float dphi = deltaPhi((*ihit1)->globalPosition().barePhi(), (*ihit2)->globalPosition().barePhi());
432  if (dphi < 0.5) {
433  if ((*ihit1)->globalPosition().y() > 0.0 && ((*ihit1)->globalPosition().y() > (*ihit2)->globalPosition().y())) {
434  std::string tag2 = "top" + tag;
435 
436  result.push_back(MuonRecHitPair(*ihit1, *ihit2, tag2));
437  } else if ((*ihit1)->globalPosition().y() < 0.0 &&
438  ((*ihit1)->globalPosition().y() < (*ihit2)->globalPosition().y())) {
439  std::string tag2 = "bottom" + tag;
440  result.push_back(MuonRecHitPair(*ihit2, *ihit1, tag2));
441  }
442  }
443  }
444  }
445 
446  return result;
447 }
448 
450  const edm::EventSetup& eSetup) const {
451  std::vector<TrajectorySeed> result;
452 
453  const std::string category = "Muon|RecoMuon|CosmicMuonSeedGenerator";
454 
455  MuonPatternRecoDumper dumper;
456 
457  float dphi = deltaPhi((hitpair.first)->globalDirection().barePhi(), (hitpair.second)->globalDirection().barePhi());
458 
459  LogTrace(category) << "hitpair.type " << hitpair.type;
460 
461  map<string, float>::const_iterator iterPar = theParameters.find(hitpair.type);
462  if (iterPar == theParameters.end()) {
463  return result;
464  }
465 
466  // set the pt and charge by dphi
467  int charge = (dphi > 0) ? -1 : 1;
468 
469  double pt = 999.0;
470  float paraC = (iterPar->second);
471 
472  if (fabs(dphi) > 1e-5) {
473  pt = paraC / fabs(dphi);
474  }
475 
476  if (pt < 10.0) {
477  return result;
478  } //still use the old strategy for low pt
479 
480  AlgebraicVector t(4);
481  AlgebraicSymMatrix mat(5, 0);
482 
484  if (hit->dimension() < (hitpair.second)->dimension())
485  hit = hitpair.second;
486 
487  // Fill the LocalTrajectoryParameters
488  LocalPoint segPos = hit->localPosition();
489 
490  GlobalVector polar(GlobalVector::Spherical(hit->globalDirection().theta(), hit->globalDirection().phi(), 1.));
491  // Force all track downward for cosmic, not beam-halo
492  if (theForcePointDownFlag) {
493  if (hit->geographicalId().subdetId() == MuonSubdetId::DT && fabs(hit->globalDirection().eta()) < 4.0 &&
494  hit->globalDirection().phi() > 0)
495  polar = -polar;
496 
497  if (hit->geographicalId().subdetId() == MuonSubdetId::CSC && fabs(hit->globalDirection().eta()) > 2.3)
498  polar = -polar;
499  }
500 
501  polar *= fabs(pt) / polar.perp();
502 
503  LocalVector segDir = hit->det()->toLocal(polar);
504 
505  LocalTrajectoryParameters param(segPos, segDir, charge);
506 
507  mat = hit->parametersError().similarityT(hit->projectionMatrix());
508 
509  float p_err = 0.004 / paraC;
510  if (pt < 10.01)
511  p_err = 0.1;
512  mat[0][0] = p_err;
513 
514  LocalTrajectoryError error(asSMatrix<5>(mat));
515 
516  // Create the TrajectoryStateOnSurface
517  TrajectoryStateOnSurface tsos(param, error, hit->det()->surface(), &*theField);
518 
519  LogTrace(category) << "Trajectory State on Surface of Seed";
520  LogTrace(category) << "mom: " << tsos.globalMomentum() << " phi: " << tsos.globalMomentum().phi();
521  LogTrace(category) << "pos: " << tsos.globalPosition();
522  LogTrace(category) << "The RecSegment relies on: ";
523  LogTrace(category) << dumper.dumpMuonId(hit->geographicalId());
524 
526  container.push_back(hitpair.first->hit()->clone());
527  container.push_back(hitpair.second->hit()->clone());
528 
529  result.push_back(tsosToSeed(tsos, hit->geographicalId().rawId(), container));
530 
531  return result;
532 }
533 
536  return tsosToSeed(tsos, id, container);
537 }
538 
540  uint32_t id,
541  edm::OwnVector<TrackingRecHit>& container) const {
543  TrajectorySeed seed(seedTSOS, container, alongMomentum);
544  return seed;
545 }
T getParameter(std::string 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
T perp() const
Definition: PV3DBase.h:69
MuonTransientTrackingRecHit::MuonRecHitContainer MuonRecHitContainer
MuonTransientTrackingRecHit::MuonRecHitContainer selectSegments(const MuonTransientTrackingRecHit::MuonRecHitContainer &) const
select seed candidates from Segments in Event
MuonTransientTrackingRecHit::MuonRecHitPointer MuonRecHitPointer
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
MuonTransientTrackingRecHit::MuonRecHitPointer first
std::string dumpMuonId(const DetId &id) const
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
void push_back(D *&d)
Definition: OwnVector.h:326
bool areCorrelated(const MuonTransientTrackingRecHit::MuonRecHitPointer &, const MuonTransientTrackingRecHit::MuonRecHitPointer &) const
check if two rechits are correlated
T mag() const
Definition: PV3DBase.h:64
std::vector< TrajectorySeed > TrajectorySeedCollection
std::shared_ptr< MuonTransientTrackingRecHit > MuonRecHitPointer
bool leftIsBetter(const MuonTransientTrackingRecHit::MuonRecHitPointer &, const MuonTransientTrackingRecHit::MuonRecHitPointer &) const
compare quality of two rechits
std::vector< TrajectorySeed > createSeed(const MuonTransientTrackingRecHit::MuonRecHitPointer &, const edm::EventSetup &) const
create TrajectorySeed from MuonTransientTrackingRecHit
#define LogTrace(id)
std::vector< MuonRecHitPair > MuonRecHitPairVector
CLHEP::HepVector AlgebraicVector
bool checkQuality(const MuonTransientTrackingRecHit::MuonRecHitPointer &) const
determine if a MuonTransientTrackingRecHit is qualified to build seed
MuonTransientTrackingRecHit::MuonRecHitPointer MuonRecHitPointer
MuonTransientTrackingRecHit::MuonRecHitPointer second
T eta() const
Definition: PV3DBase.h:73
~CosmicMuonSeedGenerator() override
Destructor.
HLT enums.
T get() const
Definition: EventSetup.h:73
CLHEP::HepSymMatrix AlgebraicSymMatrix
MuonTransientTrackingRecHit::MuonRecHitContainer MuonRecHitContainer
void produce(edm::Event &, const edm::EventSetup &) override
reconstruct muon&#39;s seeds
std::vector< MuonRecHitPair > makeSegPairs(const MuonTransientTrackingRecHit::MuonRecHitContainer &, const MuonTransientTrackingRecHit::MuonRecHitContainer &, std::string) const
static constexpr int DT
Definition: MuonSubdetId.h:11
TrajectorySeed tsosToSeed(const TrajectoryStateOnSurface &, uint32_t) const
uint32_t dimension(pat::CandKinResolution::Parametrization parametrization)
Returns the number of free parameters in a parametrization (3 or 4)
static constexpr int CSC
Definition: MuonSubdetId.h:12
void createSeeds(TrajectorySeedCollection &results, const MuonTransientTrackingRecHit::MuonRecHitContainer &hits, const edm::EventSetup &eSetup) const
generate TrajectorySeeds and put them into results
std::vector< MuonRecHitPointer > MuonRecHitContainer
def move(src, dest)
Definition: eostools.py:511
CosmicMuonSeedGenerator(const edm::ParameterSet &)
Constructor.
Definition: event.py:1