test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 &, edm::ConsumesCollector &iC)
 Constructor. More...
 
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_
 
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,
edm::ConsumesCollector iC 
)

Constructor.

Definition at line 75 of file DTTimingExtractor.cc.

References edm::ParameterSet::getParameter(), MuonServiceProxy_cff::MuonServiceProxy, 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 = new MuonServiceProxy(serviceParameters);
89 
90  edm::ParameterSet matchParameters = iConfig.getParameter<edm::ParameterSet>("MatchParameters");
91 
92  theMatcher = new MuonSegmentMatcher(matchParameters, iC);
93 }
T getParameter(std::string const &) const
unsigned int theHitsMin_
MuonSegmentMatcher * theMatcher
MuonServiceProxy * theService
DTTimingExtractor::~DTTimingExtractor ( )

Destructor.

Definition at line 96 of file DTTimingExtractor.cc.

References theMatcher, and theService.

97 {
98  if (theService) delete theService;
99  if (theMatcher) delete theMatcher;
100 }
MuonSegmentMatcher * theMatcher
MuonServiceProxy * theService

Member Function Documentation

void DTTimingExtractor::fillTiming ( TimeMeasurementSequence tmSequence,
reco::TrackRef  muonTrack,
const edm::Event iEvent,
const edm::EventSetup iSetup 
)

Definition at line 109 of file DTTimingExtractor.cc.

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

Referenced by MuonTimingFiller::fillTiming().

110 {
111 
112 // using reco::TrackCollection;
113 
114  if (debug)
115  std::cout << " *** Muon Timimng Extractor ***" << std::endl;
116 
117  theService->update(iSetup);
118 
119  const GlobalTrackingGeometry *theTrackingGeometry = &*theService->trackingGeometry();
120 
121  // get the DT geometry
122  edm::ESHandle<DTGeometry> theDTGeom;
123  iSetup.get<MuonGeometryRecord>().get(theDTGeom);
124 
126  iSetup.get<TrackingComponentsRecord>().get("SteppingHelixPropagatorAny", propagator);
127  const Propagator *propag = propagator.product();
128 
129  double invbeta=0;
130  double invbetaerr=0;
131  double totalWeightInvbeta=0;
132  double totalWeightVertex=0;
133  std::vector<TimeMeasurement> tms;
134 
135  math::XYZPoint pos=muonTrack->innerPosition();
136  math::XYZVector mom=muonTrack->innerMomentum();
137  GlobalPoint posp(pos.x(), pos.y(), pos.z());
138  GlobalVector momv(mom.x(), mom.y(), mom.z());
139  FreeTrajectoryState muonFTS(posp, momv, (TrackCharge)muonTrack->charge(), theService->magneticField().product());
140 
141  // get the DT segments that were used to construct the muon
142  std::vector<const DTRecSegment4D*> range = theMatcher->matchDT(*muonTrack,iEvent);
143 
144  if (debug)
145  std::cout << " The muon track matches " << range.size() << " segments." << std::endl;
146 
147 
148  // create a collection on TimeMeasurements for the track
149  for (std::vector<const DTRecSegment4D*>::iterator rechit = range.begin(); rechit!=range.end();++rechit) {
150 
151  // Create the ChamberId
152  DetId id = (*rechit)->geographicalId();
153  DTChamberId chamberId(id.rawId());
154  int station = chamberId.station();
155  if (debug) std::cout << "Matched DT segment in station " << station << std::endl;
156 
157  // use only segments with both phi and theta projections present (optional)
158  bool bothProjections = ( ((*rechit)->hasPhi()) && ((*rechit)->hasZed()) );
159 
160  if (requireBothProjections_ && !bothProjections) continue;
161 
162  // loop over (theta, phi) segments
163  for (int phi=0; phi<2; phi++) {
164 
165  if (dropTheta_ && !phi) continue;
166 
167  const DTRecSegment2D* segm;
168  if (phi) segm = dynamic_cast<const DTRecSegment2D*>((*rechit)->phiSegment());
169  else segm = dynamic_cast<const DTRecSegment2D*>((*rechit)->zSegment());
170 
171  if(segm == 0) continue;
172  if (!segm->specificRecHits().size()) continue;
173 
174  const GeomDet* geomDet = theTrackingGeometry->idToDet(segm->geographicalId());
175  const std::vector<DTRecHit1D> hits1d = segm->specificRecHits();
176 
177  // store all the hits from the segment
178  for (std::vector<DTRecHit1D>::const_iterator hiti=hits1d.begin(); hiti!=hits1d.end(); hiti++) {
179 
180  const GeomDet* dtcell = theTrackingGeometry->idToDet(hiti->geographicalId());
181  TimeMeasurement thisHit;
182 
183  std::pair< TrajectoryStateOnSurface, double> tsos;
184  tsos=propag->propagateWithPath(muonFTS,dtcell->surface());
185 
186  double dist;
187  double dist_straight = dtcell->toGlobal(hiti->localPosition()).mag();
188  if (tsos.first.isValid()) {
189  dist = tsos.second+posp.mag();
190 // std::cout << "Propagate distance: " << dist << " ( innermost: " << posp.mag() << ")" << std::endl;
191  } else {
192  dist = dist_straight;
193 // std::cout << "Geom distance: " << dist << std::endl;
194  }
195 
196  thisHit.driftCell = hiti->geographicalId();
197  if (hiti->lrSide()==DTEnums::Left) thisHit.isLeft=true; else thisHit.isLeft=false;
198  thisHit.isPhi = phi;
199  thisHit.posInLayer = geomDet->toLocal(dtcell->toGlobal(hiti->localPosition())).x();
200  thisHit.distIP = dist;
201  thisHit.station = station;
202  if (useSegmentT0_ && segm->ist0Valid()) thisHit.timeCorr=segm->t0();
203  else thisHit.timeCorr=0.;
204  thisHit.timeCorr += theTimeOffset_;
205 
206  // signal propagation along the wire correction for unmached theta or phi segment hits
207  if (doWireCorr_ && !bothProjections && tsos.first.isValid()) {
208  const DTLayer* layer = theDTGeom->layer(hiti->wireId());
209  float propgL = layer->toLocal( tsos.first.globalPosition() ).y();
210  float wirePropCorr = propgL/24.4*0.00543; // signal propagation speed along the wire
211  if (thisHit.isLeft) wirePropCorr=-wirePropCorr;
212  thisHit.posInLayer += wirePropCorr;
213  const DTSuperLayer *sl = layer->superLayer();
214  float tofCorr = sl->surface().position().mag()-tsos.first.globalPosition().mag();
215  tofCorr = (tofCorr/29.979)*0.00543;
216  if (thisHit.isLeft) tofCorr=-tofCorr;
217  thisHit.posInLayer += tofCorr;
218  } else {
219  // non straight-line path correction
220  float slCorr = (dist_straight-dist)/29.979*0.00543;
221  if (thisHit.isLeft) slCorr=-slCorr;
222  thisHit.posInLayer += slCorr;
223  }
224 
225  if (debug) std::cout << " dist: " << dist << " t0: " << thisHit.posInLayer << std::endl;
226 
227  tms.push_back(thisHit);
228  }
229  } // phi = (0,1)
230  } // rechit
231 
232  bool modified = false;
233  std::vector <double> dstnc, dsegm, dtraj, hitWeightVertex, hitWeightInvbeta, left;
234 
235  // Now loop over the measurements, calculate 1/beta and cut away outliers
236  do {
237 
238  modified = false;
239  dstnc.clear();
240  dsegm.clear();
241  dtraj.clear();
242  hitWeightVertex.clear();
243  hitWeightInvbeta.clear();
244  left.clear();
245 
246  std::vector <int> hit_idx;
247  totalWeightInvbeta=0;
248  totalWeightVertex=0;
249 
250  // Rebuild segments
251  for (int sta=1;sta<5;sta++)
252  for (int phi=0;phi<2;phi++) {
253  std::vector <TimeMeasurement> seg;
254  std::vector <int> seg_idx;
255  int tmpos=0;
256  for (std::vector<TimeMeasurement>::iterator tm=tms.begin(); tm!=tms.end(); ++tm) {
257  if ((tm->station==sta) && (tm->isPhi==phi)) {
258  seg.push_back(*tm);
259  seg_idx.push_back(tmpos);
260  }
261  tmpos++;
262  }
263 
264  unsigned int segsize = seg.size();
265  if (segsize<theHitsMin_) continue;
266 
267  double a=0, b=0;
268  std::vector <double> hitxl,hitxr,hityl,hityr;
269 
270  for (std::vector<TimeMeasurement>::iterator tm=seg.begin(); tm!=seg.end(); ++tm) {
271 
272  DetId id = tm->driftCell;
273  const GeomDet* dtcell = theTrackingGeometry->idToDet(id);
274  DTChamberId chamberId(id.rawId());
275  const GeomDet* dtcham = theTrackingGeometry->idToDet(chamberId);
276 
277  double celly=dtcham->toLocal(dtcell->position()).z();
278 
279  if (tm->isLeft) {
280  hitxl.push_back(celly);
281  hityl.push_back(tm->posInLayer);
282  } else {
283  hitxr.push_back(celly);
284  hityr.push_back(tm->posInLayer);
285  }
286  }
287 
288  if (!fitT0(a,b,hitxl,hityl,hitxr,hityr)) {
289  if (debug)
290  std::cout << " t0 = zero, Left hits: " << hitxl.size() << " Right hits: " << hitxr.size() << std::endl;
291  continue;
292  }
293 
294  // a segment must have at least one left and one right hit
295  if ((!hitxl.size()) || (!hityl.size())) continue;
296 
297  int segidx=0;
298  for (std::vector<TimeMeasurement>::const_iterator tm=seg.begin(); tm!=seg.end(); ++tm) {
299 
300  DetId id = tm->driftCell;
301  const GeomDet* dtcell = theTrackingGeometry->idToDet(id);
302  DTChamberId chamberId(id.rawId());
303  const GeomDet* dtcham = theTrackingGeometry->idToDet(chamberId);
304 
305  double layerZ = dtcham->toLocal(dtcell->position()).z();
306  double segmLocalPos = b+layerZ*a;
307  double hitLocalPos = tm->posInLayer;
308  int hitSide = -tm->isLeft*2+1;
309  double t0_segm = (-(hitSide*segmLocalPos)+(hitSide*hitLocalPos))/0.00543+tm->timeCorr;
310 
311  if (debug) std::cout << " Segm hit. dstnc: " << tm->distIP << " t0: " << t0_segm << std::endl;
312 
313  dstnc.push_back(tm->distIP);
314  dsegm.push_back(t0_segm);
315  left.push_back(hitSide);
316  hitWeightInvbeta.push_back(((double)seg.size()-2.)*tm->distIP*tm->distIP/((double)seg.size()*30.*30.*theError_*theError_));
317  hitWeightVertex.push_back(((double)seg.size()-2.)/((double)seg.size()*theError_*theError_));
318  hit_idx.push_back(seg_idx.at(segidx));
319  segidx++;
320  totalWeightInvbeta+=((double)seg.size()-2.)*tm->distIP*tm->distIP/((double)seg.size()*30.*30.*theError_*theError_);
321  totalWeightVertex+=((double)seg.size()-2.)/((double)seg.size()*theError_*theError_);
322  }
323  }
324 
325  if (totalWeightInvbeta==0) break;
326 
327  // calculate the value and error of 1/beta from the complete set of 1D hits
328  if (debug)
329  std::cout << " Points for global fit: " << dstnc.size() << std::endl;
330 
331  // inverse beta - weighted average of the contributions from individual hits
332  invbeta=0;
333  for (unsigned int i=0;i<dstnc.size();i++)
334  invbeta+=(1.+dsegm.at(i)/dstnc.at(i)*30.)*hitWeightInvbeta.at(i)/totalWeightInvbeta;
335 
336  double chimax=0.;
337  std::vector<TimeMeasurement>::iterator tmmax;
338 
339  // the dispersion of inverse beta
340  double diff;
341  for (unsigned int i=0;i<dstnc.size();i++) {
342  diff=(1.+dsegm.at(i)/dstnc.at(i)*30.)-invbeta;
343  diff=diff*diff*hitWeightInvbeta.at(i);
344  invbetaerr+=diff;
345  if (diff>chimax) {
346  tmmax=tms.begin()+hit_idx.at(i);
347  chimax=diff;
348  }
349  }
350 
351  invbetaerr=sqrt(invbetaerr/totalWeightInvbeta);
352 
353  // cut away the outliers
354  if (chimax>thePruneCut_) {
355  tms.erase(tmmax);
356  modified=true;
357  }
358 
359  if (debug)
360  std::cout << " Measured 1/beta: " << invbeta << " +/- " << invbetaerr << std::endl;
361 
362  } while (modified);
363 
364  for (unsigned int i=0;i<dstnc.size();i++) {
365  tmSequence.dstnc.push_back(dstnc.at(i));
366  tmSequence.local_t0.push_back(dsegm.at(i));
367  tmSequence.weightInvbeta.push_back(hitWeightInvbeta.at(i));
368  tmSequence.weightVertex.push_back(hitWeightVertex.at(i));
369  }
370 
371  tmSequence.totalWeightInvbeta=totalWeightInvbeta;
372  tmSequence.totalWeightVertex=totalWeightVertex;
373 
374 }
void update(const edm::EventSetup &setup)
update the services each event
int i
Definition: DBlmapReader.cc:9
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:52
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:67
edm::ESHandle< MagneticField > magneticField() const
get the magnetic field
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
int TrackCharge
Definition: TrackCharge.h:4
virtual const GeomDet * idToDet(DetId) const
T mag() const
Definition: PV3DBase.h:67
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:46
T sqrt(T t)
Definition: SSEVec.h:48
if(c.getParameter< edm::InputTag >("puppiValueMap").label().size()!=0)
std::vector< double > weightInvbeta
virtual std::pair< TrajectoryStateOnSurface, double > propagateWithPath(const FreeTrajectoryState &, const Surface &) const final
Definition: Propagator.cc:15
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
bool ist0Valid() const
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
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
double b
Definition: hdecay.h:120
MuonSegmentMatcher * theMatcher
const DTSuperLayer * superLayer() const
Definition: DTLayer.cc:54
edm::ESHandle< GlobalTrackingGeometry > trackingGeometry() const
get the tracking geometry
double a
Definition: hdecay.h:121
std::vector< const DTRecSegment4D * > matchDT(const reco::Track &muon, const edm::Event &event)
perform the matching
tuple cout
Definition: gather_cfg.py:121
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
std::vector< double > weightVertex
MuonServiceProxy * theService
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)
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 377 of file DTTimingExtractor.cc.

References delta, i, alignCSCRings::s, and contentValuesCheck::ss.

Referenced by fillTiming().

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

Member Data Documentation

bool DTTimingExtractor::debug
private

Definition at line 90 of file DTTimingExtractor.h.

Referenced by fillTiming().

bool DTTimingExtractor::doWireCorr_
private

Definition at line 87 of file DTTimingExtractor.h.

Referenced by fillTiming().

bool DTTimingExtractor::dropTheta_
private

Definition at line 88 of file DTTimingExtractor.h.

Referenced by fillTiming().

edm::InputTag DTTimingExtractor::DTSegmentTags_
private

Definition at line 81 of file DTTimingExtractor.h.

bool DTTimingExtractor::requireBothProjections_
private

Definition at line 89 of file DTTimingExtractor.h.

Referenced by fillTiming().

double DTTimingExtractor::theError_
private

Definition at line 85 of file DTTimingExtractor.h.

Referenced by fillTiming().

unsigned int DTTimingExtractor::theHitsMin_
private

Definition at line 82 of file DTTimingExtractor.h.

Referenced by fillTiming().

MuonSegmentMatcher* DTTimingExtractor::theMatcher
private

Definition at line 94 of file DTTimingExtractor.h.

Referenced by DTTimingExtractor(), fillTiming(), and ~DTTimingExtractor().

double DTTimingExtractor::thePruneCut_
private

Definition at line 83 of file DTTimingExtractor.h.

Referenced by fillTiming().

MuonServiceProxy* DTTimingExtractor::theService
private

Definition at line 92 of file DTTimingExtractor.h.

Referenced by DTTimingExtractor(), fillTiming(), and ~DTTimingExtractor().

double DTTimingExtractor::theTimeOffset_
private

Definition at line 84 of file DTTimingExtractor.h.

Referenced by fillTiming().

bool DTTimingExtractor::useSegmentT0_
private

Definition at line 86 of file DTTimingExtractor.h.

Referenced by fillTiming().