CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
DTTimingExtractor Class Reference

#include <RecoMuon/MuonIdentification/src/DTTimingExtractor.cc>

Classes

class  TimeMeasurement
 

Public Member Functions

 DTTimingExtractor (const edm::ParameterSet &, MuonSegmentMatcher *segMatcher)
 Constructor. More...
 
void fillTiming (TimeMeasurementSequence &tmSequence, const std::vector< const DTRecSegment4D * > &segments, reco::TrackRef muonTrack, const edm::Event &iEvent, const edm::EventSetup &iSetup)
 
void fillTiming (TimeMeasurementSequence &tmSequence, reco::TrackRef muonTrack, const edm::Event &iEvent, const edm::EventSetup &iSetup)
 
 ~DTTimingExtractor ()
 Destructor. More...
 

Private Member Functions

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)
 

Private Attributes

bool debug
 
bool doWireCorr_
 
bool dropTheta_
 
edm::InputTag DTSegmentTags_
 
bool requireBothProjections_
 
double theError_
 
unsigned int theHitsMin_
 
MuonSegmentMatchertheMatcher
 
double thePruneCut_
 
std::unique_ptr< MuonServiceProxytheService
 
double theTimeOffset_
 
bool useSegmentT0_
 

Detailed Description

Extracts timing information associated to a muon track

Description: Produce timing information for a muon track using 1D DT hits from segments used to build the track

Definition at line 54 of file DTTimingExtractor.h.

Constructor & Destructor Documentation

DTTimingExtractor::DTTimingExtractor ( const edm::ParameterSet iConfig,
MuonSegmentMatcher segMatcher 
)

Constructor.

Definition at line 75 of file DTTimingExtractor.cc.

References edm::ParameterSet::getParameter(), theMatcher, and theService.

76  :
77  theHitsMin_(iConfig.getParameter<int>("HitsMin")),
78  thePruneCut_(iConfig.getParameter<double>("PruneCut")),
79  theTimeOffset_(iConfig.getParameter<double>("DTTimeOffset")),
80  theError_(iConfig.getParameter<double>("HitError")),
81  useSegmentT0_(iConfig.getParameter<bool>("UseSegmentT0")),
82  doWireCorr_(iConfig.getParameter<bool>("DoWireCorr")),
83  dropTheta_(iConfig.getParameter<bool>("DropTheta")),
84  requireBothProjections_(iConfig.getParameter<bool>("RequireBothProjections")),
85  debug(iConfig.getParameter<bool>("debug"))
86 {
87  edm::ParameterSet serviceParameters = iConfig.getParameter<edm::ParameterSet>("ServiceParameters");
88  theService = std::make_unique<MuonServiceProxy>(serviceParameters);
89  theMatcher = segMatcher;
90 }
T getParameter(std::string const &) const
std::unique_ptr< MuonServiceProxy > theService
unsigned int theHitsMin_
MuonSegmentMatcher * theMatcher
DTTimingExtractor::~DTTimingExtractor ( )

Destructor.

Definition at line 93 of file DTTimingExtractor.cc.

94 {
95 }

Member Function Documentation

void DTTimingExtractor::fillTiming ( TimeMeasurementSequence tmSequence,
const std::vector< const DTRecSegment4D * > &  segments,
reco::TrackRef  muonTrack,
const edm::Event iEvent,
const edm::EventSetup iSetup 
)

Definition at line 101 of file DTTimingExtractor.cc.

References a, b, gather_cfg::cout, debug, DTTimingExtractor::TimeMeasurement::distIP, doWireCorr_, DTTimingExtractor::TimeMeasurement::driftCell, dropTheta_, TimeMeasurementSequence::dstnc, fitT0(), TrackingRecHit::geographicalId(), edm::EventSetup::get(), mps_fire::i, GlobalTrackingGeometry::idToDet(), DTTimingExtractor::TimeMeasurement::isLeft, DTTimingExtractor::TimeMeasurement::isPhi, DTRecSegment2D::ist0Valid(), DTGeometry::layer(), DTEnums::Left, TimeMeasurementSequence::local_t0, PV3DBase< T, PVType, FrameType >::mag(), mag(), phi, DTTimingExtractor::TimeMeasurement::posInLayer, GloballyPositioned< T >::position(), GeomDet::position(), edm::ESHandle< T >::product(), Propagator::propagateWithPath(), PhotonConversionTrajectorySeedProducerFromQuadruplets_cfi::propagator, requireBothProjections_, DTRecSegment2D::specificRecHits(), mathSSE::sqrt(), DTChamberId::station(), relativeConstraints::station, DTTimingExtractor::TimeMeasurement::station, DTLayer::superLayer(), GeomDet::surface(), DTRecSegment2D::t0(), theError_, theHitsMin_, thePruneCut_, theService, theTimeOffset_, MuonServiceProxy::theTrackingGeometry, DTTimingExtractor::TimeMeasurement::timeCorr, GeomDet::toGlobal(), GeomDet::toLocal(), TimeMeasurementSequence::totalWeightInvbeta, TimeMeasurementSequence::totalWeightTimeVtx, useSegmentT0_, TimeMeasurementSequence::weightInvbeta, TimeMeasurementSequence::weightTimeVtx, PV3DBase< T, PVType, FrameType >::x(), y, and z.

Referenced by fillTiming().

105 {
106  if (debug)
107  std::cout << " *** DT Timimng Extractor ***" << std::endl;
108 
109  theService->update(iSetup);
110 
111  const GlobalTrackingGeometry *theTrackingGeometry = &*theService->trackingGeometry();
112 
113  // get the DT geometry
114  edm::ESHandle<DTGeometry> theDTGeom;
115  iSetup.get<MuonGeometryRecord>().get(theDTGeom);
116 
117  // get the propagator
119  iSetup.get<TrackingComponentsRecord>().get("SteppingHelixPropagatorAny", propagator);
120  const Propagator *propag = propagator.product();
121 
122  math::XYZPoint pos=muonTrack->innerPosition();
123  math::XYZVector mom=muonTrack->innerMomentum();
124  GlobalPoint posp(pos.x(), pos.y(), pos.z());
125  GlobalVector momv(mom.x(), mom.y(), mom.z());
126  FreeTrajectoryState muonFTS(posp, momv, (TrackCharge)muonTrack->charge(), theService->magneticField().product());
127 
128  // create a collection on TimeMeasurements for the track
129  std::vector<TimeMeasurement> tms;
130  for (const auto& rechit : segments ) {
131 
132  // Create the ChamberId
133  DetId id = rechit->geographicalId();
134  DTChamberId chamberId(id.rawId());
135  int station = chamberId.station();
136  if (debug) std::cout << "Matched DT segment in station " << station << std::endl;
137 
138  // use only segments with both phi and theta projections present (optional)
139  bool bothProjections = ( (rechit->hasPhi()) && (rechit->hasZed()) );
140  if (requireBothProjections_ && !bothProjections) continue;
141 
142  // loop over (theta, phi) segments
143  for (int phi=0; phi<2; phi++) {
144 
145  if (dropTheta_ && !phi) continue;
146 
147  const DTRecSegment2D* segm;
148  if (phi) {
149  segm = dynamic_cast<const DTRecSegment2D*>(rechit->phiSegment());
150  if (debug) std::cout << " *** Segment t0: " << segm->t0() << std::endl;
151  }
152  else segm = dynamic_cast<const DTRecSegment2D*>(rechit->zSegment());
153 
154  if(segm == nullptr) continue;
155  if (segm->specificRecHits().empty()) continue;
156 
157  const GeomDet* geomDet = theTrackingGeometry->idToDet(segm->geographicalId());
158  const std::vector<DTRecHit1D>& hits1d(segm->specificRecHits());
159 
160  // store all the hits from the segment
161  for (const auto& hiti : hits1d) {
162 
163  const GeomDet* dtcell = theTrackingGeometry->idToDet(hiti.geographicalId());
164  TimeMeasurement thisHit;
165 
166  std::pair< TrajectoryStateOnSurface, double> tsos;
167  tsos=propag->propagateWithPath(muonFTS,dtcell->surface());
168 
169  double dist;
170  double dist_straight = dtcell->toGlobal(hiti.localPosition()).mag();
171  if (tsos.first.isValid()) {
172  dist = tsos.second+posp.mag();
173 // std::cout << "Propagate distance: " << dist << " ( innermost: " << posp.mag() << ")" << std::endl;
174  } else {
175  dist = dist_straight;
176 // std::cout << "Geom distance: " << dist << std::endl;
177  }
178 
179  thisHit.driftCell = hiti.geographicalId();
180  if (hiti.lrSide()==DTEnums::Left) thisHit.isLeft=true; else thisHit.isLeft=false;
181  thisHit.isPhi = phi;
182  thisHit.posInLayer = geomDet->toLocal(dtcell->toGlobal(hiti.localPosition())).x();
183  thisHit.distIP = dist;
184  thisHit.station = station;
185  if (useSegmentT0_ && segm->ist0Valid()) thisHit.timeCorr=segm->t0();
186  else 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) wirePropCorr=-wirePropCorr;
195  thisHit.posInLayer += wirePropCorr;
196  const DTSuperLayer *sl = layer->superLayer();
197  float tofCorr = sl->surface().position().mag()-tsos.first.globalPosition().mag();
198  tofCorr = (tofCorr/29.979)*0.00543;
199  if (thisHit.isLeft) tofCorr=-tofCorr;
200  thisHit.posInLayer += tofCorr;
201  } else {
202  // non straight-line path correction
203  float slCorr = (dist_straight-dist)/29.979*0.00543;
204  if (thisHit.isLeft) slCorr=-slCorr;
205  thisHit.posInLayer += slCorr;
206  }
207 
208  if (debug) std::cout << " dist: " << dist << " t0: " << thisHit.posInLayer << std::endl;
209 
210  tms.push_back(thisHit);
211  }
212  } // phi = (0,1)
213  } // rechit
214 
215  bool modified = false;
216  std::vector <double> dstnc, local_t0, hitWeightTimeVtx, hitWeightInvbeta, left;
217  double totalWeightInvbeta=0;
218  double totalWeightTimeVtx=0;
219 
220  // Now loop over the measurements, calculate 1/beta and time at vertex and cut away outliers
221  do {
222 
223  modified = false;
224  dstnc.clear();
225  local_t0.clear();
226  hitWeightTimeVtx.clear();
227  hitWeightInvbeta.clear();
228  left.clear();
229 
230  std::vector <int> hit_idx;
231  totalWeightInvbeta=0;
232  totalWeightTimeVtx=0;
233 
234  // Rebuild segments
235  for (int sta=1;sta<5;sta++)
236  for (int phi=0;phi<2;phi++) {
237  std::vector <TimeMeasurement> seg;
238  std::vector <int> seg_idx;
239  int tmpos=0;
240  for (auto& tm : tms) {
241  if ((tm.station==sta) && (tm.isPhi==phi)) {
242  seg.push_back(tm);
243  seg_idx.push_back(tmpos);
244  }
245  tmpos++;
246  }
247 
248  unsigned int segsize = seg.size();
249  if (segsize<theHitsMin_) continue;
250 
251  double a=0, b=0;
252  std::vector <double> hitxl,hitxr,hityl,hityr;
253 
254  for (auto& tm : seg) {
255 
256  DetId id = tm.driftCell;
257  const GeomDet* dtcell = theTrackingGeometry->idToDet(id);
258  DTChamberId chamberId(id.rawId());
259  const GeomDet* dtcham = theTrackingGeometry->idToDet(chamberId);
260 
261  double celly=dtcham->toLocal(dtcell->position()).z();
262 
263  if (tm.isLeft) {
264  hitxl.push_back(celly);
265  hityl.push_back(tm.posInLayer);
266  } else {
267  hitxr.push_back(celly);
268  hityr.push_back(tm.posInLayer);
269  }
270  }
271 
272  if (!fitT0(a,b,hitxl,hityl,hitxr,hityr)) {
273  if (debug)
274  std::cout << " t0 = zero, Left hits: " << hitxl.size() << " Right hits: " << hitxr.size() << std::endl;
275  continue;
276  }
277 
278  // a segment must have at least one left and one right hit
279  if ((hitxl.empty()) || (hityl.empty())) continue;
280 
281  int segidx=0;
282  for (const auto& tm : seg) {
283 
284  DetId id = tm.driftCell;
285  const GeomDet* dtcell = theTrackingGeometry->idToDet(id);
286  DTChamberId chamberId(id.rawId());
287  const GeomDet* dtcham = theTrackingGeometry->idToDet(chamberId);
288 
289  double layerZ = dtcham->toLocal(dtcell->position()).z();
290  double segmLocalPos = b+layerZ*a;
291  double hitLocalPos = tm.posInLayer;
292  int hitSide = -tm.isLeft*2+1;
293  double t0_segm = (-(hitSide*segmLocalPos)+(hitSide*hitLocalPos))/0.00543+tm.timeCorr;
294 
295  if (debug) std::cout << " Segm hit. dstnc: " << tm.distIP << " t0: " << t0_segm << std::endl;
296 
297  dstnc.push_back(tm.distIP);
298  local_t0.push_back(t0_segm);
299  left.push_back(hitSide);
300  hitWeightInvbeta.push_back(((double)seg.size()-2.)*tm.distIP*tm.distIP/((double)seg.size()*30.*30.*theError_*theError_));
301  hitWeightTimeVtx.push_back(((double)seg.size()-2.)/((double)seg.size()*theError_*theError_));
302  hit_idx.push_back(seg_idx.at(segidx));
303  segidx++;
304  totalWeightInvbeta+=((double)seg.size()-2.)*tm.distIP*tm.distIP/((double)seg.size()*30.*30.*theError_*theError_);
305  totalWeightTimeVtx+=((double)seg.size()-2.)/((double)seg.size()*theError_*theError_);
306  }
307  }
308 
309  if (totalWeightInvbeta==0) break;
310 
311  // calculate the value and error of 1/beta and timeVtx from the complete set of 1D hits
312  if (debug)
313  std::cout << " Points for global fit: " << dstnc.size() << std::endl;
314 
315  double invbeta=0,invbetaErr=0;
316  double timeVtx=0,timeVtxErr=0;
317 
318  for (unsigned int i=0;i<dstnc.size();i++) {
319  invbeta+=(1.+local_t0.at(i)/dstnc.at(i)*30.)*hitWeightInvbeta.at(i)/totalWeightInvbeta;
320  timeVtx+=local_t0.at(i)*hitWeightTimeVtx.at(i)/totalWeightTimeVtx;
321  }
322 
323  double chimax=0.;
324  std::vector<TimeMeasurement>::iterator tmmax;
325 
326  // Calculate the inv beta and time at vertex dispersion
327  double diff_ibeta,diff_tvtx;
328  for (unsigned int i=0;i<dstnc.size();i++) {
329  diff_ibeta=(1.+local_t0.at(i)/dstnc.at(i)*30.)-invbeta;
330  diff_ibeta=diff_ibeta*diff_ibeta*hitWeightInvbeta.at(i);
331  diff_tvtx=local_t0.at(i)-timeVtx;
332  diff_tvtx=diff_tvtx*diff_tvtx*hitWeightTimeVtx.at(i);
333  invbetaErr+=diff_ibeta;
334  timeVtxErr+=diff_tvtx;
335 
336  // decide if we cut away time at vertex outliers or inverse beta outliers
337  // currently not configurable.
338  if (diff_tvtx>chimax) {
339  tmmax=tms.begin()+hit_idx.at(i);
340  chimax=diff_tvtx;
341  }
342  }
343 
344  // cut away the outliers
345  if (chimax>thePruneCut_) {
346  tms.erase(tmmax);
347  modified=true;
348  }
349 
350  if (debug) {
351  double cf = 1./(dstnc.size()-1);
352  invbetaErr=sqrt(invbetaErr/totalWeightInvbeta*cf);
353  timeVtxErr=sqrt(timeVtxErr/totalWeightTimeVtx*cf);
354  std::cout << " Measured 1/beta: " << invbeta << " +/- " << invbetaErr << std::endl;
355  std::cout << " Measured time: " << timeVtx << " +/- " << timeVtxErr << std::endl;
356  }
357 
358  } while (modified);
359 
360  for (unsigned int i=0;i<dstnc.size();i++) {
361  tmSequence.dstnc.push_back(dstnc.at(i));
362  tmSequence.local_t0.push_back(local_t0.at(i));
363  tmSequence.weightInvbeta.push_back(hitWeightInvbeta.at(i));
364  tmSequence.weightTimeVtx.push_back(hitWeightTimeVtx.at(i));
365  }
366 
367  tmSequence.totalWeightInvbeta=totalWeightInvbeta;
368  tmSequence.totalWeightTimeVtx=totalWeightTimeVtx;
369 }
std::vector< double > local_t0
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:69
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
int TrackCharge
Definition: TrackCharge.h:4
T mag() const
Definition: PV3DBase.h:67
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:48
T sqrt(T t)
Definition: SSEVec.h:18
virtual std::pair< TrajectoryStateOnSurface, double > propagateWithPath(const FreeTrajectoryState &, const Surface &) const final
Definition: Propagator.cc:15
std::vector< double > weightInvbeta
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
bool ist0Valid() const
std::unique_ptr< MuonServiceProxy > theService
unsigned int theHitsMin_
Definition: DetId.h:18
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< double > weightTimeVtx
double b
Definition: hdecay.h:120
const DTSuperLayer * superLayer() const
Definition: DTLayer.cc:54
double a
Definition: hdecay.h:121
T get() const
Definition: EventSetup.h:63
const GeomDet * idToDet(DetId) const override
const DTLayer * layer(const DTLayerId &id) const
Return a layer given its id.
Definition: DTGeometry.cc:109
DetId geographicalId() const
double t0() const
Get the segment t0 (if recomputed, 0 is returned otherwise)
std::vector< double > dstnc
T x() const
Definition: PV3DBase.h:62
const PositionType & position() const
T const * product() const
Definition: ESHandle.h:86
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)
void DTTimingExtractor::fillTiming ( TimeMeasurementSequence tmSequence,
reco::TrackRef  muonTrack,
const edm::Event iEvent,
const edm::EventSetup iSetup 
)

Definition at line 373 of file DTTimingExtractor.cc.

References gather_cfg::cout, debug, fillTiming(), MuonSegmentMatcher::matchDT(), and theMatcher.

375 {
376  // get the DT segments that were used to construct the muon
377  std::vector<const DTRecSegment4D*> range = theMatcher->matchDT(*muonTrack,iEvent);
378 
379  if (debug)
380  std::cout << " The muon track matches " << range.size() << " segments." << std::endl;
381  fillTiming(tmSequence,range,muonTrack,iEvent,iSetup);
382 }
void fillTiming(TimeMeasurementSequence &tmSequence, const std::vector< const DTRecSegment4D * > &segments, reco::TrackRef muonTrack, const edm::Event &iEvent, const edm::EventSetup &iSetup)
MuonSegmentMatcher * theMatcher
std::vector< const DTRecSegment4D * > matchDT(const reco::Track &muon, const edm::Event &event)
perform the matching
double DTTimingExtractor::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 
)
private

Definition at line 385 of file DTTimingExtractor.cc.

References delta, mps_fire::i, alignCSCRings::s, fftjetcommon_cfi::sx, and fftjetcommon_cfi::sy.

Referenced by fillTiming().

385  {
386 
387  double sx=0,sy=0,sxy=0,sxx=0,ssx=0,ssy=0,s=0,ss=0;
388 
389  for (unsigned int i=0; i<xl.size(); i++) {
390  sx+=xl[i];
391  sy+=yl[i];
392  sxy+=xl[i]*yl[i];
393  sxx+=xl[i]*xl[i];
394  s++;
395  ssx+=xl[i];
396  ssy+=yl[i];
397  ss++;
398  }
399 
400  for (unsigned int i=0; i<xr.size(); i++) {
401  sx+=xr[i];
402  sy+=yr[i];
403  sxy+=xr[i]*yr[i];
404  sxx+=xr[i]*xr[i];
405  s++;
406  ssx-=xr[i];
407  ssy-=yr[i];
408  ss--;
409  }
410 
411  double delta = ss*ss*sxx+s*sx*sx+s*ssx*ssx-s*s*sxx-2*ss*sx*ssx;
412 
413  double t0_corr=0.;
414 
415  if (delta) {
416  a=(ssy*s*ssx+sxy*ss*ss+sy*sx*s-sy*ss*ssx-ssy*sx*ss-sxy*s*s)/delta;
417  b=(ssx*sy*ssx+sxx*ssy*ss+sx*sxy*s-sxx*sy*s-ssx*sxy*ss-sx*ssy*ssx)/delta;
418  t0_corr=(ssx*s*sxy+sxx*ss*sy+sx*sx*ssy-sxx*s*ssy-sx*ss*sxy-ssx*sx*sy)/delta;
419  }
420 
421  // convert drift distance to time
422  t0_corr/=-0.00543;
423 
424  return t0_corr;
425 }
dbl * delta
Definition: mlp_gen.cc:36
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121

Member Data Documentation

bool DTTimingExtractor::debug
private

Definition at line 96 of file DTTimingExtractor.h.

Referenced by fillTiming().

bool DTTimingExtractor::doWireCorr_
private

Definition at line 93 of file DTTimingExtractor.h.

Referenced by fillTiming().

bool DTTimingExtractor::dropTheta_
private

Definition at line 94 of file DTTimingExtractor.h.

Referenced by fillTiming().

edm::InputTag DTTimingExtractor::DTSegmentTags_
private

Definition at line 87 of file DTTimingExtractor.h.

bool DTTimingExtractor::requireBothProjections_
private

Definition at line 95 of file DTTimingExtractor.h.

Referenced by fillTiming().

double DTTimingExtractor::theError_
private

Definition at line 91 of file DTTimingExtractor.h.

Referenced by fillTiming().

unsigned int DTTimingExtractor::theHitsMin_
private

Definition at line 88 of file DTTimingExtractor.h.

Referenced by fillTiming().

MuonSegmentMatcher* DTTimingExtractor::theMatcher
private

Definition at line 99 of file DTTimingExtractor.h.

Referenced by DTTimingExtractor(), and fillTiming().

double DTTimingExtractor::thePruneCut_
private

Definition at line 89 of file DTTimingExtractor.h.

Referenced by fillTiming().

std::unique_ptr<MuonServiceProxy> DTTimingExtractor::theService
private

Definition at line 98 of file DTTimingExtractor.h.

Referenced by DTTimingExtractor(), and fillTiming().

double DTTimingExtractor::theTimeOffset_
private

Definition at line 90 of file DTTimingExtractor.h.

Referenced by fillTiming().

bool DTTimingExtractor::useSegmentT0_
private

Definition at line 92 of file DTTimingExtractor.h.

Referenced by fillTiming().