CMS 3D CMS Logo

Functions
CandidateSimMuonMatcher.cc File Reference
#include "L1Trigger/L1TMuonOverlapPhase1/interface/Tools/CandidateSimMuonMatcher.h"
#include "L1Trigger/L1TMuon/interface/MicroGMTConfiguration.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "DataFormats/GeometrySurface/interface/BoundCylinder.h"
#include "DataFormats/GeometrySurface/interface/SimpleCylinderBounds.h"
#include "DataFormats/GeometrySurface/interface/BoundDisk.h"
#include "DataFormats/GeometrySurface/interface/SimpleDiskBounds.h"
#include "DataFormats/MuonDetId/interface/RPCDetId.h"
#include "TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h"
#include "boost/dynamic_bitset.hpp"
#include "TFile.h"
#include "TH1D.h"

Go to the source code of this file.

Functions

double foldPhi (double phi)
 
double hwGmtPhiToGlobalPhi (int phi)
 
float normal_pdf (float x, float m, float s)
 
bool simTrackIsMuonInBx0 (const SimTrack &simTrack)
 
bool simTrackIsMuonInOmtf (const SimTrack &simTrack)
 
bool simTrackIsMuonInOmtfBx0 (const SimTrack &simTrack)
 
bool trackingParticleIsMuonInBx0 (const TrackingParticle &trackingParticle)
 
bool trackingParticleIsMuonInOmtfBx0 (const TrackingParticle &trackingParticle)
 
bool trackingParticleIsMuonInOmtfEvent0 (const TrackingParticle &trackingParticle)
 

Function Documentation

◆ foldPhi()

double foldPhi ( double  phi)

Definition at line 33 of file CandidateSimMuonMatcher.cc.

References M_PI, and PVValHelper::phi.

Referenced by CandidateSimMuonMatcher::match(), and CandidateSimMuonMatcher::matchSimple().

33  {
34  if (phi > M_PI)
35  return (phi - 2 * M_PI);
36  else if (phi < -M_PI)
37  return (phi + 2 * M_PI);
38 
39  return phi;
40 }
#define M_PI

◆ hwGmtPhiToGlobalPhi()

double hwGmtPhiToGlobalPhi ( int  phi)

Definition at line 28 of file CandidateSimMuonMatcher.cc.

References M_PI, and PVValHelper::phi.

Referenced by CandidateSimMuonMatcher::ghostBust(), CandidateSimMuonMatcher::match(), and CandidateSimMuonMatcher::matchSimple().

28  {
29  double phiGmtUnit = 2. * M_PI / 576.;
30  return phi * phiGmtUnit;
31 }
#define M_PI

◆ normal_pdf()

float normal_pdf ( float  x,
float  m,
float  s 
)

Definition at line 403 of file CandidateSimMuonMatcher.cc.

References a, JetChargeProducer_cfi::exp, visualization-live-secondInstance_cfg::m, alignCSCRings::s, and x.

Referenced by CandidateSimMuonMatcher::match().

403  {
404  static const float inv_sqrt_2pi = 0.3989422804014327;
405  float a = (x - m) / s;
406 
407  return inv_sqrt_2pi / s * std::exp(-0.5 * a * a);
408 }
double a
Definition: hdecay.h:121
float x

◆ simTrackIsMuonInBx0()

bool simTrackIsMuonInBx0 ( const SimTrack simTrack)

Definition at line 136 of file CandidateSimMuonMatcher.cc.

References funct::abs(), and cscDigiValidation_cfi::simTrack.

Referenced by CandidateSimMuonMatcher::observeEventEnd().

136  {
137  if (std::abs(simTrack.type()) == 13 ||
138  std::abs(simTrack.type()) == 1000015) { // 1000015 is stau, todo use other selection (e.g. pt>20) if needed
139  //only muons
140  if (simTrack.eventId().bunchCrossing() == 0)
141  return true;
142  }
143  return false;
144 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ simTrackIsMuonInOmtf()

bool simTrackIsMuonInOmtf ( const SimTrack simTrack)

Definition at line 100 of file CandidateSimMuonMatcher.cc.

References funct::abs(), LogTrace, and cscDigiValidation_cfi::simTrack.

Referenced by simTrackIsMuonInOmtfBx0().

100  {
101  if (std::abs(simTrack.type()) == 13 ||
102  std::abs(simTrack.type()) == 1000015) { // 1000015 is stau, todo use other selection (e.g. pt>20) if needed
103  //only muons
104  } else
105  return false;
106 
107  //in the overlap, the propagation of muons with pt less then ~3.2 fails - the actual threshold depends slightly on eta,
108  if (simTrack.momentum().pt() < 2.5)
109  return false;
110 
111  LogTrace("l1tOmtfEventPrint") << "simTrackIsMuonInOmtf: simTrack type " << std::setw(3) << simTrack.type() << " pt "
112  << std::setw(9) << simTrack.momentum().pt() << " eta " << std::setw(9)
113  << simTrack.momentum().eta() << " phi " << std::setw(9) << simTrack.momentum().phi()
114  << std::endl;
115 
116  //some margin for matching must be used on top of actual OMTF region,
117  //i.e. (0.82-1.24)=>(0.72-1.3),
118  //otherwise many candidates are marked as ghosts
119  if ((std::abs(simTrack.momentum().eta()) >= 0.72) && (std::abs(simTrack.momentum().eta()) <= 1.3)) {
120  LogTrace("l1tOmtfEventPrint") << "simTrackIsMuonInOmtf: is in OMTF";
121  } else {
122  LogTrace("l1tOmtfEventPrint") << "simTrackIsMuonInOmtf: not in OMTF";
123  return false;
124  }
125 
126  return true;
127 }
#define LogTrace(id)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ simTrackIsMuonInOmtfBx0()

bool simTrackIsMuonInOmtfBx0 ( const SimTrack simTrack)

Definition at line 129 of file CandidateSimMuonMatcher.cc.

References cscDigiValidation_cfi::simTrack, and simTrackIsMuonInOmtf().

Referenced by CandidateSimMuonMatcher::observeEventEnd().

129  {
130  if (simTrack.eventId().bunchCrossing() != 0)
131  return false;
132 
134 }
bool simTrackIsMuonInOmtf(const SimTrack &simTrack)

◆ trackingParticleIsMuonInBx0()

bool trackingParticleIsMuonInBx0 ( const TrackingParticle trackingParticle)

Definition at line 146 of file CandidateSimMuonMatcher.cc.

References funct::abs(), EncodedEventId::bunchCrossing(), TrackingParticle::eventId(), and TrackingParticle::pdgId().

Referenced by CandidateSimMuonMatcher::observeEventEnd().

146  {
147  if (std::abs(trackingParticle.pdgId()) == 13 ||
148  std::abs(trackingParticle.pdgId()) ==
149  1000015) { // 1000015 is stau, todo use other selection (e.g. pt>20) if needed
150  //only muons
151  if (trackingParticle.eventId().bunchCrossing() == 0)
152  return true;
153  }
154  return false;
155 }
int pdgId() const
PDG ID.
EncodedEventId eventId() const
Signal source, crossing number.
int bunchCrossing() const
get the detector field from this detid
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ trackingParticleIsMuonInOmtfBx0()

bool trackingParticleIsMuonInOmtfBx0 ( const TrackingParticle trackingParticle)

Definition at line 157 of file CandidateSimMuonMatcher.cc.

References funct::abs(), EncodedEventId::bunchCrossing(), EncodedEventId::event(), TrackingParticle::eventId(), TrackingParticle::g4Tracks(), edm::Ref< C, T, F >::isNonnull(), LogTrace, TrackingParticle::momentum(), TrackingParticle::parentVertex(), TrackingParticle::pdgId(), and TrackingParticle::pt().

Referenced by trackingParticleIsMuonInOmtfEvent0().

157  {
158  if (trackingParticle.eventId().bunchCrossing() != 0)
159  return false;
160 
161  if (std::abs(trackingParticle.pdgId()) == 13 || std::abs(trackingParticle.pdgId()) == 1000015) {
162  // 1000015 is stau, todo use other selection (e.g. pt>20) if needed
163  //only muons
164  } else
165  return false;
166 
167  //in the overlap, the propagation of muons with pt less then ~3.2 fails, the actual threshold depends slightly on eta,
168  if (trackingParticle.pt() < 2.5)
169  return false;
170 
171  if (trackingParticle.parentVertex().isNonnull())
172  LogTrace("l1tOmtfEventPrint") << "trackingParticleIsMuonInOmtfBx0, pdgId " << std::setw(3)
173  << trackingParticle.pdgId() << " pt " << std::setw(9) << trackingParticle.pt()
174  << " eta " << std::setw(9) << trackingParticle.momentum().eta() << " phi "
175  << std::setw(9) << trackingParticle.momentum().phi() << " event "
176  << trackingParticle.eventId().event() << " trackId "
177  << trackingParticle.g4Tracks().at(0).trackId() << " parentVertex Rho "
178  << trackingParticle.parentVertex()->position().Rho() << " eta "
179  << trackingParticle.parentVertex()->position().eta() << " phi "
180  << trackingParticle.parentVertex()->position().phi() << std::endl;
181  else
182  LogTrace("l1tOmtfEventPrint") << "trackingParticleIsMuonInOmtfBx0, pdgId " << std::setw(3)
183  << trackingParticle.pdgId() << " pt " << std::setw(9) << trackingParticle.pt()
184  << " eta " << std::setw(9) << trackingParticle.momentum().eta() << " phi "
185  << std::setw(9) << trackingParticle.momentum().phi() << " trackId "
186  << trackingParticle.g4Tracks().at(0).trackId();
187 
188  //some margin for matching must be used on top of actual OMTF region,
189  //i.e. (0.82-1.24)=>(0.72-1.3),
190  //otherwise many candidates are marked as ghosts
191  if ((std::abs(trackingParticle.momentum().eta()) >= 0.72) && (std::abs(trackingParticle.momentum().eta()) <= 1.3)) {
192  } else
193  return false;
194 
195  return true;
196 }
int pdgId() const
PDG ID.
EncodedEventId eventId() const
Signal source, crossing number.
Vector momentum() const
spatial momentum vector
int event() const
get the contents of the subdetector field (should be protected?)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
const std::vector< SimTrack > & g4Tracks() const
#define LogTrace(id)
int bunchCrossing() const
get the detector field from this detid
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const TrackingVertexRef & parentVertex() const
double pt() const
Transverse momentum. Note this is taken from the first SimTrack only.

◆ trackingParticleIsMuonInOmtfEvent0()

bool trackingParticleIsMuonInOmtfEvent0 ( const TrackingParticle trackingParticle)

Definition at line 198 of file CandidateSimMuonMatcher.cc.

References EncodedEventId::event(), TrackingParticle::eventId(), and trackingParticleIsMuonInOmtfBx0().

198  {
199  if (trackingParticle.eventId().event() != 0)
200  return false;
201 
202  return trackingParticleIsMuonInOmtfBx0(trackingParticle);
203 }
EncodedEventId eventId() const
Signal source, crossing number.
int event() const
get the contents of the subdetector field (should be protected?)
bool trackingParticleIsMuonInOmtfBx0(const TrackingParticle &trackingParticle)