CMS 3D CMS Logo

DTTimingExtractor.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: MuonIdentification
4 // Class: DTTimingExtractor
5 //
11 //
12 // Original Author: Traczyk Piotr
13 // Created: Thu Oct 11 15:01:28 CEST 2007
14 //
15 //
16 
18 
19 // user include files
22 
24 
26 
28 
37 
41 
49 
54 
55 // system include files
56 #include <memory>
57 #include <vector>
58 #include <string>
59 #include <iostream>
60 
61 namespace edm {
62  class ParameterSet;
63  class EventSetup;
64  class InputTag;
65 } // namespace edm
66 
67 class MuonServiceProxy;
68 
69 //
70 // constructors and destructor
71 //
73  MuonSegmentMatcher* segMatcher,
75  : theHitsMin_(iConfig.getParameter<int>("HitsMin")),
76  thePruneCut_(iConfig.getParameter<double>("PruneCut")),
77  theTimeOffset_(iConfig.getParameter<double>("DTTimeOffset")),
78  theError_(iConfig.getParameter<double>("HitError")),
79  useSegmentT0_(iConfig.getParameter<bool>("UseSegmentT0")),
80  doWireCorr_(iConfig.getParameter<bool>("DoWireCorr")),
81  dropTheta_(iConfig.getParameter<bool>("DropTheta")),
82  requireBothProjections_(iConfig.getParameter<bool>("RequireBothProjections")),
83  debug(iConfig.getParameter<bool>("debug")),
84  theDTGeomToken(iC.esConsumes()),
85  thePropagatorToken(iC.esConsumes(edm::ESInputTag("", "SteppingHelixPropagatorAny"))) {
86  edm::ParameterSet serviceParameters = iConfig.getParameter<edm::ParameterSet>("ServiceParameters");
87  theService = std::make_unique<MuonServiceProxy>(serviceParameters, edm::ConsumesCollector(iC));
88  theMatcher = segMatcher;
89 }
90 
92 
93 //
94 // member functions
95 //
97  const std::vector<const DTRecSegment4D*>& segments,
98  reco::TrackRef muonTrack,
99  const edm::Event& iEvent,
100  const edm::EventSetup& iSetup) {
101  if (debug)
102  std::cout << " *** DT Timimng Extractor ***" << std::endl;
103 
104  theService->update(iSetup);
105 
106  const GlobalTrackingGeometry* theTrackingGeometry = &*theService->trackingGeometry();
107 
108  // get the DT geometry
109  DTGeometry const& theDTGeom = iSetup.getData(theDTGeomToken);
110 
111  // get the propagator
112  const Propagator* propag = &iSetup.getData(thePropagatorToken);
113 
114  math::XYZPoint pos = muonTrack->innerPosition();
115  math::XYZVector mom = muonTrack->innerMomentum();
116  GlobalPoint posp(pos.x(), pos.y(), pos.z());
117  GlobalVector momv(mom.x(), mom.y(), mom.z());
118  FreeTrajectoryState muonFTS(posp, momv, (TrackCharge)muonTrack->charge(), theService->magneticField().product());
119 
120  // create a collection on TimeMeasurements for the track
121  std::vector<TimeMeasurement> tms;
122  for (const auto& rechit : segments) {
123  // Create the ChamberId
124  DetId id = rechit->geographicalId();
125  DTChamberId chamberId(id.rawId());
126  int station = chamberId.station();
127  if (debug)
128  std::cout << "Matched DT segment in station " << station << std::endl;
129 
130  // use only segments with both phi and theta projections present (optional)
131  bool bothProjections = ((rechit->hasPhi()) && (rechit->hasZed()));
132  if (requireBothProjections_ && !bothProjections)
133  continue;
134 
135  // loop over (theta, phi) segments
136  for (int phi = 0; phi < 2; phi++) {
137  if (dropTheta_ && !phi)
138  continue;
139 
140  const DTRecSegment2D* segm;
141  if (phi) {
142  segm = dynamic_cast<const DTRecSegment2D*>(rechit->phiSegment());
143  if (debug)
144  std::cout << " *** Segment t0: " << segm->t0() << std::endl;
145  } else
146  segm = dynamic_cast<const DTRecSegment2D*>(rechit->zSegment());
147 
148  if (segm == nullptr)
149  continue;
150  if (segm->specificRecHits().empty())
151  continue;
152 
153  const GeomDet* geomDet = theTrackingGeometry->idToDet(segm->geographicalId());
154  const std::vector<DTRecHit1D>& hits1d(segm->specificRecHits());
155 
156  // store all the hits from the segment
157  for (const auto& hiti : hits1d) {
158  const GeomDet* dtcell = theTrackingGeometry->idToDet(hiti.geographicalId());
159  TimeMeasurement thisHit;
160 
161  std::pair<TrajectoryStateOnSurface, double> tsos;
162  tsos = propag->propagateWithPath(muonFTS, dtcell->surface());
163 
164  double dist;
165  double dist_straight = dtcell->toGlobal(hiti.localPosition()).mag();
166  if (tsos.first.isValid()) {
167  dist = tsos.second + posp.mag();
168  // std::cout << "Propagate distance: " << dist << " ( innermost: " << posp.mag() << ")" << std::endl;
169  } else {
170  dist = dist_straight;
171  // std::cout << "Geom distance: " << dist << std::endl;
172  }
173 
174  thisHit.driftCell = hiti.geographicalId();
175  if (hiti.lrSide() == DTEnums::Left)
176  thisHit.isLeft = true;
177  else
178  thisHit.isLeft = false;
179  thisHit.isPhi = phi;
180  thisHit.posInLayer = geomDet->toLocal(dtcell->toGlobal(hiti.localPosition())).x();
181  thisHit.distIP = dist;
182  thisHit.station = station;
183  if (useSegmentT0_ && segm->ist0Valid())
184  thisHit.timeCorr = segm->t0();
185  else
186  thisHit.timeCorr = 0.;
187  thisHit.timeCorr += theTimeOffset_;
188 
189  // signal propagation along the wire correction for unmached theta or phi segment hits
190  if (doWireCorr_ && !bothProjections && tsos.first.isValid()) {
191  const DTLayer* layer = theDTGeom.layer(hiti.wireId());
192  float propgL = layer->toLocal(tsos.first.globalPosition()).y();
193  float wirePropCorr = propgL / 24.4 * 0.00543; // signal propagation speed along the wire
194  if (thisHit.isLeft)
195  wirePropCorr = -wirePropCorr;
196  thisHit.posInLayer += wirePropCorr;
197  const DTSuperLayer* sl = layer->superLayer();
198  float tofCorr = sl->surface().position().mag() - tsos.first.globalPosition().mag();
199  tofCorr = (tofCorr / 29.979) * 0.00543;
200  if (thisHit.isLeft)
201  tofCorr = -tofCorr;
202  thisHit.posInLayer += tofCorr;
203  } else {
204  // non straight-line path correction
205  float slCorr = (dist_straight - dist) / 29.979 * 0.00543;
206  if (thisHit.isLeft)
207  slCorr = -slCorr;
208  thisHit.posInLayer += slCorr;
209  }
210 
211  if (debug)
212  std::cout << " dist: " << dist << " t0: " << thisHit.posInLayer << std::endl;
213 
214  tms.push_back(thisHit);
215  }
216  } // phi = (0,1)
217  } // rechit
218 
219  bool modified = false;
220  std::vector<double> dstnc, local_t0, hitWeightTimeVtx, hitWeightInvbeta, left;
221  double totalWeightInvbeta = 0;
222  double totalWeightTimeVtx = 0;
223 
224  // Now loop over the measurements, calculate 1/beta and time at vertex and cut away outliers
225  do {
226  modified = false;
227  dstnc.clear();
228  local_t0.clear();
229  hitWeightTimeVtx.clear();
230  hitWeightInvbeta.clear();
231  left.clear();
232 
233  std::vector<int> hit_idx;
234  totalWeightInvbeta = 0;
235  totalWeightTimeVtx = 0;
236 
237  // Rebuild segments
238  for (int sta = 1; sta < 5; sta++)
239  for (int phi = 0; phi < 2; phi++) {
240  std::vector<TimeMeasurement> seg;
241  std::vector<int> seg_idx;
242  int tmpos = 0;
243  for (auto& tm : tms) {
244  if ((tm.station == sta) && (tm.isPhi == phi)) {
245  seg.push_back(tm);
246  seg_idx.push_back(tmpos);
247  }
248  tmpos++;
249  }
250 
251  unsigned int segsize = seg.size();
252  if (segsize < theHitsMin_)
253  continue;
254 
255  double a = 0, b = 0;
256  std::vector<double> hitxl, hitxr, hityl, hityr;
257 
258  for (auto& tm : seg) {
259  DetId id = tm.driftCell;
260  const GeomDet* dtcell = theTrackingGeometry->idToDet(id);
261  DTChamberId chamberId(id.rawId());
262  const GeomDet* dtcham = theTrackingGeometry->idToDet(chamberId);
263 
264  double celly = dtcham->toLocal(dtcell->position()).z();
265 
266  if (tm.isLeft) {
267  hitxl.push_back(celly);
268  hityl.push_back(tm.posInLayer);
269  } else {
270  hitxr.push_back(celly);
271  hityr.push_back(tm.posInLayer);
272  }
273  }
274 
275  if (!fitT0(a, b, hitxl, hityl, hitxr, hityr)) {
276  if (debug)
277  std::cout << " t0 = zero, Left hits: " << hitxl.size() << " Right hits: " << hitxr.size() << std::endl;
278  continue;
279  }
280 
281  // a segment must have at least one left and one right hit
282  if ((hitxl.empty()) || (hityl.empty()))
283  continue;
284 
285  int segidx = 0;
286  for (const auto& tm : seg) {
287  DetId id = tm.driftCell;
288  const GeomDet* dtcell = theTrackingGeometry->idToDet(id);
289  DTChamberId chamberId(id.rawId());
290  const GeomDet* dtcham = theTrackingGeometry->idToDet(chamberId);
291 
292  double layerZ = dtcham->toLocal(dtcell->position()).z();
293  double segmLocalPos = b + layerZ * a;
294  double hitLocalPos = tm.posInLayer;
295  int hitSide = -tm.isLeft * 2 + 1;
296  double t0_segm = (-(hitSide * segmLocalPos) + (hitSide * hitLocalPos)) / 0.00543 + tm.timeCorr;
297 
298  if (debug)
299  std::cout << " Segm hit. dstnc: " << tm.distIP << " t0: " << t0_segm << std::endl;
300 
301  dstnc.push_back(tm.distIP);
302  local_t0.push_back(t0_segm);
303  left.push_back(hitSide);
304  hitWeightInvbeta.push_back(((double)seg.size() - 2.) * tm.distIP * tm.distIP /
305  ((double)seg.size() * 30. * 30. * theError_ * theError_));
306  hitWeightTimeVtx.push_back(((double)seg.size() - 2.) / ((double)seg.size() * theError_ * theError_));
307  hit_idx.push_back(seg_idx.at(segidx));
308  segidx++;
309  totalWeightInvbeta += ((double)seg.size() - 2.) * tm.distIP * tm.distIP /
310  ((double)seg.size() * 30. * 30. * theError_ * theError_);
311  totalWeightTimeVtx += ((double)seg.size() - 2.) / ((double)seg.size() * theError_ * theError_);
312  }
313  }
314 
315  if (totalWeightInvbeta == 0)
316  break;
317 
318  // calculate the value and error of 1/beta and timeVtx from the complete set of 1D hits
319  if (debug)
320  std::cout << " Points for global fit: " << dstnc.size() << std::endl;
321 
322  double invbeta = 0, invbetaErr = 0;
323  double timeVtx = 0, timeVtxErr = 0;
324 
325  for (unsigned int i = 0; i < dstnc.size(); i++) {
326  invbeta += (1. + local_t0.at(i) / dstnc.at(i) * 30.) * hitWeightInvbeta.at(i) / totalWeightInvbeta;
327  timeVtx += local_t0.at(i) * hitWeightTimeVtx.at(i) / totalWeightTimeVtx;
328  }
329 
330  double chimax = 0.;
331  std::vector<TimeMeasurement>::iterator tmmax;
332 
333  // Calculate the inv beta and time at vertex dispersion
334  double diff_ibeta, diff_tvtx;
335  for (unsigned int i = 0; i < dstnc.size(); i++) {
336  diff_ibeta = (1. + local_t0.at(i) / dstnc.at(i) * 30.) - invbeta;
337  diff_ibeta = diff_ibeta * diff_ibeta * hitWeightInvbeta.at(i);
338  diff_tvtx = local_t0.at(i) - timeVtx;
339  diff_tvtx = diff_tvtx * diff_tvtx * hitWeightTimeVtx.at(i);
340  invbetaErr += diff_ibeta;
341  timeVtxErr += diff_tvtx;
342 
343  // decide if we cut away time at vertex outliers or inverse beta outliers
344  // currently not configurable.
345  if (diff_tvtx > chimax) {
346  tmmax = tms.begin() + hit_idx.at(i);
347  chimax = diff_tvtx;
348  }
349  }
350 
351  // cut away the outliers
352  if (chimax > thePruneCut_) {
353  tms.erase(tmmax);
354  modified = true;
355  }
356 
357  if (debug) {
358  double cf = 1. / (dstnc.size() - 1);
359  invbetaErr = sqrt(invbetaErr / totalWeightInvbeta * cf);
360  timeVtxErr = sqrt(timeVtxErr / totalWeightTimeVtx * cf);
361  std::cout << " Measured 1/beta: " << invbeta << " +/- " << invbetaErr << std::endl;
362  std::cout << " Measured time: " << timeVtx << " +/- " << timeVtxErr << std::endl;
363  }
364 
365  } while (modified);
366 
367  for (unsigned int i = 0; i < dstnc.size(); i++) {
368  tmSequence.dstnc.push_back(dstnc.at(i));
369  tmSequence.local_t0.push_back(local_t0.at(i));
370  tmSequence.weightInvbeta.push_back(hitWeightInvbeta.at(i));
371  tmSequence.weightTimeVtx.push_back(hitWeightTimeVtx.at(i));
372  }
373 
374  tmSequence.totalWeightInvbeta = totalWeightInvbeta;
375  tmSequence.totalWeightTimeVtx = totalWeightTimeVtx;
376 }
377 
378 // ------------ method called to produce the data ------------
380  reco::TrackRef muonTrack,
381  const edm::Event& iEvent,
382  const edm::EventSetup& iSetup) {
383  // get the DT segments that were used to construct the muon
384  std::vector<const DTRecSegment4D*> range = theMatcher->matchDT(*muonTrack, iEvent);
385 
386  if (debug)
387  std::cout << " The muon track matches " << range.size() << " segments." << std::endl;
388  fillTiming(tmSequence, range, muonTrack, iEvent, iSetup);
389 }
390 
392  double& b,
393  const std::vector<double>& xl,
394  const std::vector<double>& yl,
395  const std::vector<double>& xr,
396  const std::vector<double>& yr) {
397  double sx = 0, sy = 0, sxy = 0, sxx = 0, ssx = 0, ssy = 0, s = 0, ss = 0;
398 
399  for (unsigned int i = 0; i < xl.size(); i++) {
400  sx += xl[i];
401  sy += yl[i];
402  sxy += xl[i] * yl[i];
403  sxx += xl[i] * xl[i];
404  s++;
405  ssx += xl[i];
406  ssy += yl[i];
407  ss++;
408  }
409 
410  for (unsigned int i = 0; i < xr.size(); i++) {
411  sx += xr[i];
412  sy += yr[i];
413  sxy += xr[i] * yr[i];
414  sxx += xr[i] * xr[i];
415  s++;
416  ssx -= xr[i];
417  ssy -= yr[i];
418  ss--;
419  }
420 
421  double delta = ss * ss * sxx + s * sx * sx + s * ssx * ssx - s * s * sxx - 2 * ss * sx * ssx;
422 
423  double t0_corr = 0.;
424 
425  if (delta) {
426  a = (ssy * s * ssx + sxy * ss * ss + sy * sx * s - sy * ss * ssx - ssy * sx * ss - sxy * s * s) / delta;
427  b = (ssx * sy * ssx + sxx * ssy * ss + sx * sxy * s - sxx * sy * s - ssx * sxy * ss - sx * ssy * ssx) / delta;
428  t0_corr = (ssx * s * sxy + sxx * ss * sy + sx * sx * ssy - sxx * s * ssy - sx * ss * sxy - ssx * sx * sy) / delta;
429  }
430 
431  // convert drift distance to time
432  t0_corr /= -0.00543;
433 
434  return t0_corr;
435 }
436 
437 //define this as a plug-in
438 //DEFINE_FWK_MODULE(DTTimingExtractor);
int station() const
Return the station number.
Definition: DTChamberId.h:42
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
bool ist0Valid() const
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::vector< double > local_t0
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
double t0() const
Get the segment t0 (if recomputed, 0 is returned otherwise)
DTTimingExtractor(const edm::ParameterSet &, MuonSegmentMatcher *segMatcher, edm::ConsumesCollector &)
Constructor.
const GeomDet * idToDet(DetId) const override
int TrackCharge
Definition: TrackCharge.h:4
T x() const
Definition: PV3DBase.h:59
void fillTiming(TimeMeasurementSequence &tmSequence, const std::vector< const DTRecSegment4D *> &segments, reco::TrackRef muonTrack, const edm::Event &iEvent, const edm::EventSetup &iSetup)
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
T mag() const
Definition: PV3DBase.h:64
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
std::vector< double > weightInvbeta
~DTTimingExtractor()
Destructor.
virtual std::pair< TrajectoryStateOnSurface, double > propagateWithPath(const FreeTrajectoryState &, const Surface &) const final
Definition: Propagator.cc:10
std::unique_ptr< MuonServiceProxy > theService
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
unsigned int theHitsMin_
Definition: DetId.h:17
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
#define debug
Definition: HDRShower.cc:19
const PositionType & position() const
DetId geographicalId() const
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< double > weightTimeVtx
double b
Definition: hdecay.h:120
MuonSegmentMatcher * theMatcher
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:43
HLT enums.
double a
Definition: hdecay.h:121
std::vector< const DTRecSegment4D * > matchDT(const reco::Track &muon, const edm::Event &event)
perform the matching
constexpr float xr
std::vector< double > dstnc
edm::ESGetToken< DTGeometry, MuonGeometryRecord > theDTGeomToken
double fitT0(double &a, double &b, const std::vector< double > &xl, const std::vector< double > &yl, const std::vector< double > &xr, const std::vector< double > &yr)
const DTLayer * layer(const DTLayerId &id) const
Return a layer given its id.
Definition: DTGeometry.cc:96
edm::ESGetToken< Propagator, TrackingComponentsRecord > thePropagatorToken