CMS 3D CMS Logo

GlobalTrackQualityProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: GlobalTrackingTools
4 // Class: GlobalTrackQualityProducer
5 //
6 //
7 // Original Author: Adam Everett
8 //
9 //
10 
11 // system include files
12 #include <memory>
13 
14 // user include files
18 
21 
26 
28  : inputCollection_(iConfig.getParameter<edm::InputTag>("InputCollection")),
29  inputLinksCollection_(iConfig.getParameter<edm::InputTag>("InputLinksCollection")),
30  tTopoToken_(esConsumes()),
31  theService(nullptr),
32  theGlbRefitter(nullptr),
33  theGlbMatcher(nullptr) {
34  // service parameters
35  edm::ParameterSet serviceParameters = iConfig.getParameter<edm::ParameterSet>("ServiceParameters");
36  theService = new MuonServiceProxy(serviceParameters, consumesCollector());
37 
38  // TrackRefitter parameters
39  edm::ConsumesCollector iC = consumesCollector();
40  edm::ParameterSet refitterParameters = iConfig.getParameter<edm::ParameterSet>("RefitterParameters");
41  theGlbRefitter = new GlobalMuonRefitter(refitterParameters, theService, iC);
42 
43  edm::ParameterSet trackMatcherPSet = iConfig.getParameter<edm::ParameterSet>("GlobalMuonTrackMatcher");
44  theGlbMatcher = new GlobalMuonTrackMatcher(trackMatcherPSet, theService);
45 
46  double maxChi2 = iConfig.getParameter<double>("MaxChi2");
47  double nSigma = iConfig.getParameter<double>("nSigma");
49 
50  glbMuonsToken = consumes<reco::TrackCollection>(inputCollection_);
51  linkCollectionToken = consumes<reco::MuonTrackLinksCollection>(inputLinksCollection_);
52 
53  produces<edm::ValueMap<reco::MuonQuality>>();
54 }
55 
57  if (theService)
58  delete theService;
59  if (theGlbRefitter)
60  delete theGlbRefitter;
61  if (theGlbMatcher)
62  delete theGlbMatcher;
63  if (theEstimator)
64  delete theEstimator;
65 }
66 
68  const std::string theCategory = "Muon|RecoMuon|GlobalTrackQualityProducer";
69 
70  theService->update(iSetup);
71 
73 
75 
76  // Take the GLB muon container(s)
78  iEvent.getByToken(glbMuonsToken, glbMuons);
79 
81  iEvent.getByToken(linkCollectionToken, linkCollectionHandle);
82 
83  //Retrieve tracker topology from geometry
84  const TrackerTopology* tTopo = &iSetup.getData(tTopoToken_);
85 
86  // reserve some space
87  std::vector<reco::MuonQuality> valuesQual;
88  valuesQual.reserve(glbMuons->size());
89 
90  int trackIndex = 0;
91  for (reco::TrackCollection::const_iterator track = glbMuons->begin(); track != glbMuons->end();
92  ++track, ++trackIndex) {
93  reco::TrackRef glbRef(glbMuons, trackIndex);
94  reco::TrackRef staTrack = reco::TrackRef();
95 
96  std::vector<Trajectory> refitted = theGlbRefitter->refit(*track, 1, tTopo);
97 
98  LogTrace(theCategory) << "GLBQual N refitted " << refitted.size();
99 
100  std::pair<double, double> thisKink;
101  double relative_muon_chi2 = 0.0;
102  double relative_tracker_chi2 = 0.0;
103  double glbTrackProbability = 0.0;
104  if (!refitted.empty()) {
105  thisKink = kink(refitted.front());
106  std::pair<double, double> chi = newChi2(refitted.front());
107  relative_muon_chi2 = chi.second; //normalized inside to /sum(muHits.dimension)
108  relative_tracker_chi2 = chi.first; //normalized inside to /sum(tkHits.dimension)
109  glbTrackProbability = trackProbability(refitted.front());
110  }
111 
112  LogTrace(theCategory) << "GLBQual: Kink " << thisKink.first << " " << thisKink.second;
113  LogTrace(theCategory) << "GLBQual: Rel Chi2 " << relative_tracker_chi2 << " " << relative_muon_chi2;
114  LogTrace(theCategory) << "GLBQual: trackProbability " << glbTrackProbability;
115 
116  // Fill the STA-TK match information
117  float chi2, d, dist, Rpos;
118  chi2 = d = dist = Rpos = -1.0;
119  bool passTight = false;
121  if (linkCollectionHandle.isValid()) {
122  for (reco::MuonTrackLinksCollection::const_iterator links = linkCollectionHandle->begin();
123  links != linkCollectionHandle->end();
124  ++links) {
125  if (links->trackerTrack().isNull() || links->standAloneTrack().isNull() || links->globalTrack().isNull()) {
126  edm::LogWarning(theCategory) << "Global muon links to constituent tracks are invalid. There should be no "
127  "such object. Muon is skipped.";
128  continue;
129  }
130  if (links->globalTrack() == glbRef) {
131  staTrack = !links->standAloneTrack().isNull() ? links->standAloneTrack() : reco::TrackRef();
132  TrackCand staCand = TrackCand((Trajectory*)nullptr, links->standAloneTrack());
133  TrackCand tkCand = TrackCand((Trajectory*)nullptr, links->trackerTrack());
134  chi2 = theGlbMatcher->match(staCand, tkCand, 0, 0);
135  d = theGlbMatcher->match(staCand, tkCand, 1, 0);
136  Rpos = theGlbMatcher->match(staCand, tkCand, 2, 0);
137  dist = theGlbMatcher->match(staCand, tkCand, 3, 0);
138  passTight = theGlbMatcher->matchTight(staCand, tkCand);
139  }
140  }
141  }
142 
143  if (!staTrack.isNull())
144  LogTrace(theCategory) << "GLBQual: Used UpdatedAtVtx : "
145  << (iEvent.getStableProvenance(staTrack.id()).productInstanceName() ==
146  std::string("UpdatedAtVtx"));
147 
148  float maxFloat01 =
149  std::numeric_limits<float>::max() * 0.1; // a better solution would be to use float above .. m/be not
150  reco::MuonQuality muQual;
151  if (!staTrack.isNull())
152  muQual.updatedSta =
153  iEvent.getStableProvenance(staTrack.id()).productInstanceName() == std::string("UpdatedAtVtx");
154  muQual.trkKink = thisKink.first > maxFloat01 ? maxFloat01 : thisKink.first;
155  muQual.glbKink = thisKink.second > maxFloat01 ? maxFloat01 : thisKink.second;
156  muQual.trkRelChi2 = relative_tracker_chi2 > maxFloat01 ? maxFloat01 : relative_tracker_chi2;
157  muQual.staRelChi2 = relative_muon_chi2 > maxFloat01 ? maxFloat01 : relative_muon_chi2;
158  muQual.tightMatch = passTight;
159  muQual.chi2LocalPosition = dist;
160  muQual.chi2LocalMomentum = chi2;
161  muQual.localDistance = d;
162  muQual.globalDeltaEtaPhi = Rpos;
163  muQual.glbTrackProbability = glbTrackProbability;
164  valuesQual.push_back(muQual);
165  }
166 
167  /*
168  for(int i = 0; i < valuesTkRelChi2.size(); i++) {
169  LogTrace(theCategory)<<"value " << valuesTkRelChi2[i] ;
170  }
171  */
172 
173  // create and fill value maps
174  auto outQual = std::make_unique<edm::ValueMap<reco::MuonQuality>>();
175  edm::ValueMap<reco::MuonQuality>::Filler fillerQual(*outQual);
176  fillerQual.insert(glbMuons, valuesQual.begin(), valuesQual.end());
177  fillerQual.fill();
178 
179  // put value map into event
180  iEvent.put(std::move(outQual));
181 }
182 
183 std::pair<double, double> GlobalTrackQualityProducer::kink(Trajectory& muon) const {
184  const std::string theCategory = "Muon|RecoMuon|GlobalTrackQualityProducer";
185 
186  using namespace std;
187  using namespace edm;
188  using namespace reco;
189 
190  double result = 0.0;
191  double resultGlb = 0.0;
192 
194  typedef ConstRecHitPointer RecHit;
195  typedef std::vector<TrajectoryMeasurement>::const_iterator TMI;
196 
197  vector<TrajectoryMeasurement> meas = muon.measurements();
198 
199  for (TMI m = meas.begin(); m != meas.end(); m++) {
201 
202  //not used double estimate = 0.0;
203 
204  RecHit rhit = (*m).recHit();
205  bool ok = false;
206  if (rhit->isValid()) {
207  if (DetId::Tracker == rhit->geographicalId().det())
208  ok = true;
209  }
210 
211  //if ( !ok ) continue;
212 
213  const TrajectoryStateOnSurface& tsos = (*m).predictedState();
214 
215  if (tsos.isValid() && rhit->isValid() && rhit->hit()->isValid() &&
216  !edm::isNotFinite(rhit->localPositionError().xx()) //this is paranoia induced by reported case
217  && !edm::isNotFinite(rhit->localPositionError().xy()) //it's better to track down the origin of bad numbers
218  && !edm::isNotFinite(rhit->localPositionError().yy())) {
219  double phi1 = tsos.globalPosition().phi();
220  if (phi1 < 0)
221  phi1 = 2 * M_PI + phi1;
222 
223  double phi2 = rhit->globalPosition().phi();
224  if (phi2 < 0)
225  phi2 = 2 * M_PI + phi2;
226 
227  double diff = fabs(phi1 - phi2);
228  if (diff > M_PI)
229  diff = 2 * M_PI - diff;
230 
231  GlobalPoint hitPos = rhit->globalPosition();
232 
233  GlobalError hitErr = rhit->globalPositionError();
234  //LogDebug(theCategory)<<"hitPos " << hitPos;
235  double error = hitErr.phierr(hitPos); // error squared
236 
237  double s = (error > 0.0) ? (diff * diff) / error : (diff * diff);
238 
239  if (ok)
240  result += s;
241  resultGlb += s;
242  }
243  }
244 
245  return std::pair<double, double>(result, resultGlb);
246 }
247 
248 std::pair<double, double> GlobalTrackQualityProducer::newChi2(Trajectory& muon) const {
249  const std::string theCategory = "Muon|RecoMuon|GlobalTrackQualityProducer";
250 
251  using namespace std;
252  using namespace edm;
253  using namespace reco;
254 
255  double muChi2 = 0.0;
256  double tkChi2 = 0.0;
257  unsigned int muNdof = 0;
258  unsigned int tkNdof = 0;
259 
261  typedef ConstRecHitPointer RecHit;
262  typedef vector<TrajectoryMeasurement>::const_iterator TMI;
263 
264  vector<TrajectoryMeasurement> meas = muon.measurements();
265 
266  for (TMI m = meas.begin(); m != meas.end(); m++) {
268  const TrajectoryStateOnSurface& uptsos = (*m).updatedState();
269  // FIXME FIXME CLONE!!!
270  // TrackingRecHit::RecHitPointer preciseHit = hit->clone(uptsos);
271  const auto& preciseHit = hit;
272  double estimate = 0.0;
273  if (preciseHit->isValid() && uptsos.isValid()) {
274  estimate = theEstimator->estimate(uptsos, *preciseHit).second;
275  }
276 
277  //LogTrace(theCategory) << "estimate " << estimate << " TM.est " << m->estimate();
278  //UNUSED: double tkDiff = 0.0;
279  //UNUSED: double staDiff = 0.0;
280  if (hit->isValid() && (hit->geographicalId().det()) == DetId::Tracker) {
281  tkChi2 += estimate;
282  //UNUSED: tkDiff = estimate - m->estimate();
283  tkNdof += hit->dimension();
284  }
285  if (hit->isValid() && (hit->geographicalId().det()) == DetId::Muon) {
286  muChi2 += estimate;
287  //UNUSED staDiff = estimate - m->estimate();
288  muNdof += hit->dimension();
289  }
290  }
291 
292  //For tkNdof < 6, should a large number or something else
293  // be used instead of just tkChi2 directly?
294  if (tkNdof > 5) {
295  tkChi2 /= (tkNdof - 5.);
296  }
297 
298  //For muNdof < 6, should a large number or something else
299  // be used instead of just muChi2 directly?
300  if (muNdof > 5) {
301  muChi2 /= (muNdof - 5.);
302  }
303 
304  return std::pair<double, double>(tkChi2, muChi2);
305 }
306 
307 //
308 // calculate the tail probability (-ln(P)) of a fit
309 //
311  if (track.ndof() > 0 && track.chiSquared() > 0) {
312  return -LnChiSquaredProbability(track.chiSquared(), track.ndof());
313  } else {
314  return 0.0;
315  }
316 }
317 
320  {
322  psd1.setAllowAnything();
323  desc.add<edm::ParameterSetDescription>("ServiceParameters", psd1);
324  }
325  {
327  psd1.setAllowAnything();
328  desc.add<edm::ParameterSetDescription>("GlobalMuonTrackMatcher", psd1);
329  }
330  desc.add<edm::InputTag>("InputCollection", edm::InputTag("globalMuons"));
331  desc.add<edm::InputTag>("InputLinksCollection", edm::InputTag("globalMuons"));
332  desc.add<std::string>("BaseLabel", "GLB");
333  {
334  edm::ParameterSetDescription descGlbMuonRefitter;
335  descGlbMuonRefitter.setAllowAnything();
336  descGlbMuonRefitter.add<edm::InputTag>("DTRecSegmentLabel", edm::InputTag("dt1DRecHits"));
337  descGlbMuonRefitter.add<edm::InputTag>("CSCRecSegmentLabel", edm::InputTag("csc2DRecHits"));
338  descGlbMuonRefitter.add<edm::InputTag>("GEMRecHitLabel", edm::InputTag("gemRecHits"));
339  descGlbMuonRefitter.add<edm::InputTag>("ME0RecHitLabel", edm::InputTag("me0Segments"));
340  descGlbMuonRefitter.add<edm::InputTag>("RPCRecSegmentLabel", edm::InputTag("rpcRecHits"));
341 
342  descGlbMuonRefitter.add<std::string>("Fitter", "KFFitterForRefitInsideOut");
343  descGlbMuonRefitter.add<std::string>("Smoother", "KFSmootherForRefitInsideOut");
344  descGlbMuonRefitter.add<std::string>("Propagator", "SmartPropagatorAnyRK");
345  descGlbMuonRefitter.add<std::string>("TrackerRecHitBuilder", "WithAngleAndTemplate");
346  descGlbMuonRefitter.add<std::string>("MuonRecHitBuilder", "MuonRecHitBuilder");
347  descGlbMuonRefitter.add<bool>("DoPredictionsOnly", false);
348  descGlbMuonRefitter.add<std::string>("RefitDirection", "insideOut");
349  descGlbMuonRefitter.add<bool>("PropDirForCosmics", false);
350  descGlbMuonRefitter.add<bool>("RefitRPCHits", true);
351 
352  descGlbMuonRefitter.add<std::vector<int>>("DYTthrs", {10, 10});
353  descGlbMuonRefitter.add<int>("DYTselector", 1);
354  descGlbMuonRefitter.add<bool>("DYTupdator", false);
355  descGlbMuonRefitter.add<bool>("DYTuseAPE", false);
356  descGlbMuonRefitter.add<bool>("DYTuseThrsParametrization", true);
357  {
358  edm::ParameterSetDescription descDYTthrs;
359  descDYTthrs.add<std::vector<double>>("eta0p8", {1, -0.919853, 0.990742});
360  descDYTthrs.add<std::vector<double>>("eta1p2", {1, -0.897354, 0.987738});
361  descDYTthrs.add<std::vector<double>>("eta2p0", {4, -0.986855, 0.998516});
362  descDYTthrs.add<std::vector<double>>("eta2p2", {1, -0.940342, 0.992955});
363  descDYTthrs.add<std::vector<double>>("eta2p4", {1, -0.947633, 0.993762});
364  descGlbMuonRefitter.add<edm::ParameterSetDescription>("DYTthrsParameters", descDYTthrs);
365  }
366 
367  descGlbMuonRefitter.add<int>("SkipStation", -1);
368  descGlbMuonRefitter.add<int>("TrackerSkipSystem", -1);
369  descGlbMuonRefitter.add<int>("TrackerSkipSection", -1);
370  descGlbMuonRefitter.add<bool>("RefitFlag", true);
371 
372  desc.add<edm::ParameterSetDescription>("RefitterParameters", descGlbMuonRefitter);
373  }
374  desc.add<double>("nSigma", 3.0);
375  desc.add<double>("MaxChi2", 100000.0);
376 
377  descriptions.add("globalTrackQualityProducer", desc);
378 }
379 //#include "FWCore/Framework/interface/MakerMacros.h"
380 //DEFINE_FWK_MODULE(GlobalTrackQualityProducer);
float chi2LocalPosition
chi2 value for the STA-TK matching of local position
Definition: MuonQuality.h:19
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
float chi2LocalMomentum
chi2 value for the STA-TK matching of local momentum
Definition: MuonQuality.h:21
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::pair< const Trajectory *, reco::TrackRef > TrackCand
void setAllowAnything()
allow any parameter label/value pairs
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
bool updatedSta
bool returns true if standAloneMuon_updatedAtVtx was used in the fit
Definition: MuonQuality.h:9
float glbKink
value of the kink algorithm applied to the global track
Definition: MuonQuality.h:13
float LnChiSquaredProbability(double chiSquared, double nrDOF)
edm::EDGetTokenT< reco::TrackCollection > glbMuonsToken
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
float glbTrackProbability
the tail probability (-ln(P)) of the global fit
Definition: MuonQuality.h:29
void produce(edm::Event &, const edm::EventSetup &) override
void setServices(const edm::EventSetup &)
set the services needed by the TrackTransformer
float trkKink
value of the kink algorithm applied to the inner track stub
Definition: MuonQuality.h:11
float trkRelChi2
chi2 value for the inner track stub with respect to the global track
Definition: MuonQuality.h:15
T phierr(const GlobalPoint &aPoint) const
#define LogTrace(id)
TransientTrackingRecHit::ConstRecHitPointer ConstRecHitPointer
int iEvent
Definition: GenABIO.cc:224
bool tightMatch
if the STA-TK matching passed the tighter matching criteria
Definition: MuonQuality.h:27
GlobalPoint globalPosition() const
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
bool matchTight(const TrackCand &sta, const TrackCand &track) const
check if two tracks are compatible (less than Chi2Cut, DeltaDCut, DeltaRCut)
float globalDeltaEtaPhi
global delta-Eta-Phi of STA-TK matching
Definition: MuonQuality.h:25
virtual HitReturnType estimate(const TrajectoryStateOnSurface &ts, const TrackingRecHit &hit) const =0
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isNull() const
Checks for null.
Definition: Ref.h:235
double match(const TrackCand &sta, const TrackCand &track, int matchOption=0, int surfaceOption=1) const
d
Definition: ztail.py:151
#define M_PI
virtual std::pair< double, double > newChi2(Trajectory &muon) const
virtual void setEvent(const edm::Event &)
pass the Event to the algo at each event
float staRelChi2
chi2 value for the outer track stub with respect to the global track
Definition: MuonQuality.h:17
virtual double trackProbability(Trajectory &track) const
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< reco::MuonTrackLinksCollection > linkCollectionToken
bool isValid() const
Definition: HandleBase.h:70
fixed size matrix
HLT enums.
float localDistance
local distance seperation for STA-TK TSOS matching on same surface
Definition: MuonQuality.h:23
GlobalMuonTrackMatcher * theGlbMatcher
std::vector< Trajectory > refit(const reco::Track &globalTrack, const int theMuonHitsOption, const TrackerTopology *tTopo) const
build combined trajectory from sta Track and tracker RecHits
void update(const edm::EventSetup &setup, bool duringEvent=true)
update the services each event
const edm::EventSetup & eventSetup() const
virtual std::pair< double, double > kink(Trajectory &muon) const
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511
GlobalTrackQualityProducer(const edm::ParameterSet &iConfig)