CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // $Id: GlobalTrackQualityProducer.cc,v 1.7 2011/10/28 22:10:33 slava77 Exp $
9 //
10 //
11 
12 // system include files
13 #include <memory>
14 
15 // user include files
18 
21 
27 
29  inputCollection_(iConfig.getParameter<edm::InputTag>("InputCollection")),inputLinksCollection_(iConfig.getParameter<edm::InputTag>("InputLinksCollection")),theService(0),theGlbRefitter(0),theGlbMatcher(0)
30 {
31  // service parameters
32  edm::ParameterSet serviceParameters = iConfig.getParameter<edm::ParameterSet>("ServiceParameters");
33  theService = new MuonServiceProxy(serviceParameters);
34 
35  // TrackRefitter parameters
36  edm::ParameterSet refitterParameters = iConfig.getParameter<edm::ParameterSet>("RefitterParameters");
37  theGlbRefitter = new GlobalMuonRefitter(refitterParameters, theService);
38 
39  edm::ParameterSet trackMatcherPSet = iConfig.getParameter<edm::ParameterSet>("GlobalMuonTrackMatcher");
40  theGlbMatcher = new GlobalMuonTrackMatcher(trackMatcherPSet,theService);
41 
42  double maxChi2 = iConfig.getParameter<double>("MaxChi2");
43  double nSigma = iConfig.getParameter<double>("nSigma");
44  theEstimator = new Chi2MeasurementEstimator(maxChi2,nSigma);
45 
46  produces<edm::ValueMap<reco::MuonQuality> >();
47 }
48 
50  if (theService) delete theService;
51  if (theGlbRefitter) delete theGlbRefitter;
52  if (theGlbMatcher) delete theGlbMatcher;
53  if (theEstimator) delete theEstimator;
54 }
55 
56 void
58 {
59  const std::string theCategory = "Muon|RecoMuon|GlobalTrackQualityProducer";
60 
61  theService->update(iSetup);
62 
63  theGlbRefitter->setEvent(iEvent);
64 
66 
67  // Take the GLB muon container(s)
69  iEvent.getByLabel(inputCollection_,glbMuons);
70 
72  iEvent.getByLabel(inputLinksCollection_, linkCollectionHandle);
73 
74  // reserve some space
75  std::vector<reco::MuonQuality> valuesQual;
76  valuesQual.reserve(glbMuons->size());
77 
78  int trackIndex = 0;
79  for (reco::TrackCollection::const_iterator track = glbMuons->begin(); track!=glbMuons->end(); ++track , ++trackIndex) {
80  reco::TrackRef glbRef(glbMuons,trackIndex);
81  reco::TrackRef staTrack = reco::TrackRef();
82 
83  std::vector<Trajectory> refitted=theGlbRefitter->refit(*track,1);
84 
85  LogTrace(theCategory)<<"GLBQual N refitted " << refitted.size();
86 
87  std::pair<double,double> thisKink;
88  double relative_muon_chi2 = 0.0;
89  double relative_tracker_chi2 = 0.0;
90  double glbTrackProbability = 0.0;
91  if(refitted.size()>0) {
92  thisKink = kink(refitted.front()) ;
93  std::pair<double,double> chi = newChi2(refitted.front());
94  relative_muon_chi2 = chi.second; //normalized inside to /sum(muHits.dimension)
95  relative_tracker_chi2 = chi.first; //normalized inside to /sum(tkHits.dimension)
96  glbTrackProbability = trackProbability(refitted.front());
97  }
98 
99  LogTrace(theCategory)<<"GLBQual: Kink " << thisKink.first << " " << thisKink.second;
100  LogTrace(theCategory)<<"GLBQual: Rel Chi2 " << relative_tracker_chi2 << " " << relative_muon_chi2;
101  LogTrace(theCategory)<<"GLBQual: trackProbability " << glbTrackProbability;
102 
103  // Fill the STA-TK match information
104  float chi2, d, dist, Rpos;
105  chi2 = d = dist = Rpos = -1.0;
106  bool passTight = false;
107  typedef MuonTrajectoryBuilder::TrackCand TrackCand;
108  if ( linkCollectionHandle.isValid() ) {
109  for ( reco::MuonTrackLinksCollection::const_iterator links = linkCollectionHandle->begin();
110  links != linkCollectionHandle->end(); ++links )
111  {
112  if ( links->trackerTrack().isNull() ||
113  links->standAloneTrack().isNull() ||
114  links->globalTrack().isNull() )
115  {
116  edm::LogWarning(theCategory) << "Global muon links to constituent tracks are invalid. There should be no such object. Muon is skipped.";
117  continue;
118  }
119  if (links->globalTrack() == glbRef) {
120  staTrack = !links->standAloneTrack().isNull() ? links->standAloneTrack() : reco::TrackRef();
121  TrackCand staCand = TrackCand((Trajectory*)(0),links->standAloneTrack());
122  TrackCand tkCand = TrackCand((Trajectory*)(0),links->trackerTrack());
123  chi2 = theGlbMatcher->match(staCand,tkCand,0,0);
124  d = theGlbMatcher->match(staCand,tkCand,1,0);
125  Rpos = theGlbMatcher->match(staCand,tkCand,2,0);
126  dist = theGlbMatcher->match(staCand,tkCand,3,0);
127  passTight = theGlbMatcher->matchTight(staCand,tkCand);
128  }
129  }
130  }
131 
132  if(!staTrack.isNull()) LogTrace(theCategory)<<"GLBQual: Used UpdatedAtVtx : " << (iEvent.getProvenance(staTrack.id()).productInstanceName() == std::string("UpdatedAtVtx"));
133 
134  float maxFloat01 = std::numeric_limits<float>::max()*0.1; // a better solution would be to use float above .. m/be not
135  reco::MuonQuality muQual;
136  if(!staTrack.isNull()) muQual.updatedSta = iEvent.getProvenance(staTrack.id()).productInstanceName() == std::string("UpdatedAtVtx");
137  muQual.trkKink = thisKink.first > maxFloat01 ? maxFloat01 : thisKink.first;
138  muQual.glbKink = thisKink.second > maxFloat01 ? maxFloat01 : thisKink.second;
139  muQual.trkRelChi2 = relative_tracker_chi2 > maxFloat01 ? maxFloat01 : relative_tracker_chi2;
140  muQual.staRelChi2 = relative_muon_chi2 > maxFloat01 ? maxFloat01 : relative_muon_chi2;
141  muQual.tightMatch = passTight;
142  muQual.chi2LocalPosition = dist;
143  muQual.chi2LocalMomentum = chi2;
144  muQual.localDistance = d;
145  muQual.globalDeltaEtaPhi = Rpos;
146  muQual.glbTrackProbability = glbTrackProbability;
147  valuesQual.push_back(muQual);
148  }
149 
150  /*
151  for(int i = 0; i < valuesTkRelChi2.size(); i++) {
152  LogTrace(theCategory)<<"value " << valuesTkRelChi2[i] ;
153  }
154  */
155 
156  // create and fill value maps
157  std::auto_ptr<edm::ValueMap<reco::MuonQuality> > outQual(new edm::ValueMap<reco::MuonQuality>());
158  edm::ValueMap<reco::MuonQuality>::Filler fillerQual(*outQual);
159  fillerQual.insert(glbMuons, valuesQual.begin(), valuesQual.end());
160  fillerQual.fill();
161 
162  // put value map into event
163  iEvent.put(outQual);
164 }
165 
166 std::pair<double,double> GlobalTrackQualityProducer::kink(Trajectory& muon) const {
167 
168  const std::string theCategory = "Muon|RecoMuon|GlobalTrackQualityProducer";
169 
170  using namespace std;
171  using namespace edm;
172  using namespace reco;
173 
174  double result = 0.0;
175  double resultGlb = 0.0;
176 
177 
179  typedef ConstRecHitPointer RecHit;
180  typedef std::vector<TrajectoryMeasurement>::const_iterator TMI;
181 
182 
183  vector<TrajectoryMeasurement> meas = muon.measurements();
184 
185  for ( TMI m = meas.begin(); m != meas.end(); m++ ) {
187 
188  //not used double estimate = 0.0;
189 
190  RecHit rhit = (*m).recHit();
191  bool ok = false;
192  if ( rhit->isValid() ) {
193  if(DetId::Tracker == rhit->geographicalId().det()) ok = true;
194  }
195 
196  //if ( !ok ) continue;
197 
198  const TrajectoryStateOnSurface& tsos = (*m).predictedState();
199 
200 
201  if ( tsos.isValid() && rhit->isValid() && rhit->hit()->isValid()
202  && !std::isinf(rhit->localPositionError().xx()) //this is paranoia induced by reported case
203  && !std::isinf(rhit->localPositionError().xy()) //it's better to track down the origin of bad numbers
204  && !std::isinf(rhit->localPositionError().yy())
205  && !std::isnan(rhit->localPositionError().xx()) //this is paranoia induced by reported case
206  && !std::isnan(rhit->localPositionError().xy()) //it's better to track down the origin of bad numbers
207  && !std::isnan(rhit->localPositionError().yy())
208  ) {
209 
210  double phi1 = tsos.globalPosition().phi();
211  if ( phi1 < 0 ) phi1 = 2*M_PI + phi1;
212 
213  double phi2 = rhit->globalPosition().phi();
214  if ( phi2 < 0 ) phi2 = 2*M_PI + phi2;
215 
216  double diff = fabs(phi1 - phi2);
217  if ( diff > M_PI ) diff = 2*M_PI - diff;
218 
219  GlobalPoint hitPos = rhit->globalPosition();
220 
221  GlobalError hitErr = rhit->globalPositionError();
222  //LogDebug(theCategory)<<"hitPos " << hitPos;
223  double error = hitErr.phierr(hitPos); // error squared
224 
225  double s = ( error > 0.0 ) ? (diff*diff)/error : (diff*diff);
226 
227  if(ok) result += s;
228  resultGlb += s;
229  }
230 
231  }
232 
233 
234  return std::pair<double,double>(result,resultGlb);
235 
236 }
237 
238 std::pair<double,double> GlobalTrackQualityProducer::newChi2(Trajectory& muon) const {
239  const std::string theCategory = "Muon|RecoMuon|GlobalTrackQualityProducer";
240 
241  using namespace std;
242  using namespace edm;
243  using namespace reco;
244 
245  double muChi2 = 0.0;
246  double tkChi2 = 0.0;
247  unsigned int muNdof = 0;
248  unsigned int tkNdof = 0;
249 
251  typedef ConstRecHitPointer RecHit;
252  typedef vector<TrajectoryMeasurement>::const_iterator TMI;
253 
254  vector<TrajectoryMeasurement> meas = muon.measurements();
255 
256  for ( TMI m = meas.begin(); m != meas.end(); m++ ) {
258  const TrajectoryStateOnSurface& uptsos = (*m).updatedState();
259  TransientTrackingRecHit::RecHitPointer preciseHit = hit->clone(uptsos);
260  double estimate = 0.0;
261  if (preciseHit->isValid() && uptsos.isValid()) {
262  estimate = theEstimator->estimate(uptsos, *preciseHit ).second;
263  }
264 
265  //LogTrace(theCategory) << "estimate " << estimate << " TM.est " << m->estimate();
266  //UNUSED: double tkDiff = 0.0;
267  //UNUSED: double staDiff = 0.0;
268  if ( hit->isValid() && (hit->geographicalId().det()) == DetId::Tracker ) {
269  tkChi2 += estimate;
270  //UNUSED: tkDiff = estimate - m->estimate();
271  tkNdof += hit->dimension();
272  }
273  if ( hit->isValid() && (hit->geographicalId().det()) == DetId::Muon ) {
274  muChi2 += estimate;
275  //UNUSED staDiff = estimate - m->estimate();
276  muNdof += hit->dimension();
277  }
278  }
279 
280  if (tkNdof < 6 ) tkChi2 = tkChi2; // or should I set it to a large number ?
281  else tkChi2 /= (tkNdof-5.);
282 
283  if (muNdof < 6 ) muChi2 = muChi2; // or should I set it to a large number ?
284  else muChi2 /= (muNdof-5.);
285 
286  return std::pair<double,double>(tkChi2,muChi2);
287 
288 }
289 
290 //
291 // calculate the tail probability (-ln(P)) of a fit
292 //
293 double
295 
296  if ( track.ndof() > 0 && track.chiSquared() > 0 ) {
297  return -LnChiSquaredProbability(track.chiSquared(), track.ndof());
298  } else {
299  return 0.0;
300  }
301 
302 }
303 
304 //#include "FWCore/Framework/interface/MakerMacros.h"
305 //DEFINE_FWK_MODULE(GlobalTrackQualityProducer);
void update(const edm::EventSetup &setup)
update the services each event
float chi2LocalPosition
chi2 value for the STA-TK matching of local position
Definition: MuonQuality.h:19
T getParameter(std::string const &) const
const edm::EventSetup & eventSetup() const
get the whole EventSetup
float chi2LocalMomentum
chi2 value for the STA-TK matching of local momentum
Definition: MuonQuality.h:21
std::pair< const Trajectory *, reco::TrackRef > TrackCand
virtual HitReturnType estimate(const TrajectoryStateOnSurface &ts, const TransientTrackingRecHit &hit) const =0
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)
bool matchTight(const TrackCand &sta, const TrackCand &track) const
check if two tracks are compatible (less than Chi2Cut, DeltaDCut, DeltaRCut)
Geom::Phi< T > phi() const
Definition: PV3DBase.h:68
float glbTrackProbability
the tail probability (-ln(P)) of the global fit
Definition: MuonQuality.h:29
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
GlobalPoint globalPosition() const
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
virtual std::pair< double, double > newChi2(Trajectory &muon) const
virtual double trackProbability(Trajectory &track) const
T phierr(const GlobalPoint &aPoint) const
TransientTrackingRecHit::ConstRecHitPointer ConstRecHitPointer
virtual std::pair< double, double > kink(Trajectory &muon) const
DataContainer const & measurements() const
Definition: Trajectory.h:203
int iEvent
Definition: GenABIO.cc:243
bool isNull() const
Checks for null.
Definition: Ref.h:247
bool tightMatch
if the STA-TK matching passed the tighter matching criteria
Definition: MuonQuality.h:27
const T & max(const T &a, const T &b)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
bool isnan(float x)
Definition: math.h:13
tuple result
Definition: query.py:137
float globalDeltaEtaPhi
global delta-Eta-Phi of STA-TK matching
Definition: MuonQuality.h:25
double match(const TrackCand &sta, const TrackCand &track, int matchOption=0, int surfaceOption=1) const
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
#define LogTrace(id)
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
int ndof(bool bon=true) const
Definition: Trajectory.cc:74
#define M_PI
Definition: BFit3D.cc:3
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:14
std::vector< Trajectory > refit(const reco::Track &globalTrack, const int theMuonHitsOption) const
build combined trajectory from sta Track and tracker RecHits
float localDistance
local distance seperation for STA-TK TSOS matching on same surface
Definition: MuonQuality.h:23
GlobalMuonTrackMatcher * theGlbMatcher
virtual void produce(edm::Event &, const edm::EventSetup &)
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
Provenance getProvenance(BranchID const &theID) const
Definition: Event.cc:60
double chiSquared() const
Definition: Trajectory.h:242
GlobalTrackQualityProducer(const edm::ParameterSet &iConfig)